2015-02-05 08:05:38 -06:00
|
|
|
/*
|
2016-09-12 05:31:59 -05:00
|
|
|
Copyright 2015-2016 Dr. Blatt - HPC-Simulation-Software & Services.
|
|
|
|
Coypright 2015 NTNU
|
|
|
|
Copyright 2015-2016 Statoil AS
|
|
|
|
Copyright 2015 IRIS AS
|
2015-02-05 08:05:38 -06:00
|
|
|
|
|
|
|
This file is part of the Open Porous Media project (OPM).
|
|
|
|
|
|
|
|
OPM is free software: you can redistribute it and/or modify
|
|
|
|
it under the terms of the GNU General Public License as published by
|
|
|
|
the Free Software Foundation, either version 3 of the License, or
|
|
|
|
(at your option) any later version.
|
|
|
|
|
|
|
|
OPM is distributed in the hope that it will be useful,
|
|
|
|
but WITHOUT ANY WARRANTY; without even the implied warranty of
|
|
|
|
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
|
|
|
GNU General Public License for more details.
|
|
|
|
|
|
|
|
You should have received a copy of the GNU General Public License
|
|
|
|
along with OPM. If not, see <http://www.gnu.org/licenses/>.
|
|
|
|
*/
|
2015-06-16 09:17:56 -05:00
|
|
|
#ifndef OPM_REDISTRIBUTEDATAHANDLES_HEADER
|
|
|
|
#define OPM_REDISTRIBUTEDATAHANDLES_HEADER
|
2015-02-05 08:05:38 -06:00
|
|
|
|
2016-09-12 05:31:59 -05:00
|
|
|
#include <unordered_set>
|
|
|
|
#include <string>
|
2017-01-25 13:36:37 -06:00
|
|
|
#include <type_traits>
|
|
|
|
#include <iterator>
|
2016-09-12 05:31:59 -05:00
|
|
|
|
2015-06-16 09:17:56 -05:00
|
|
|
#include <opm/core/simulator/BlackoilState.hpp>
|
2015-02-05 08:05:38 -06:00
|
|
|
|
|
|
|
#include <opm/autodiff/BlackoilPropsAdFromDeck.hpp>
|
2015-06-16 09:17:56 -05:00
|
|
|
#include <opm/autodiff/ExtractParallelGridInformationToISTL.hpp>
|
2017-01-26 02:25:58 -06:00
|
|
|
#include <opm/autodiff/createGlobalCellArray.hpp>
|
2015-06-16 09:17:56 -05:00
|
|
|
|
|
|
|
#include<boost/any.hpp>
|
2015-02-05 08:05:38 -06:00
|
|
|
|
|
|
|
namespace Opm
|
|
|
|
{
|
2015-05-07 04:21:09 -05:00
|
|
|
|
2015-06-16 09:17:56 -05:00
|
|
|
template <class Grid>
|
2016-09-12 05:31:59 -05:00
|
|
|
inline std::unordered_set<std::string>
|
|
|
|
distributeGridAndData( Grid& ,
|
2016-10-14 02:23:26 -05:00
|
|
|
const Opm::Deck& ,
|
|
|
|
const EclipseState& ,
|
2017-10-29 15:06:19 -05:00
|
|
|
const Schedule&,
|
|
|
|
BlackoilState&,
|
2016-09-12 05:31:59 -05:00
|
|
|
BlackoilPropsAdFromDeck& ,
|
|
|
|
DerivedGeology&,
|
|
|
|
std::shared_ptr<BlackoilPropsAdFromDeck::MaterialLawManager>&,
|
|
|
|
std::vector<double>&,
|
|
|
|
boost::any& ,
|
|
|
|
const bool )
|
2015-06-16 09:17:56 -05:00
|
|
|
{
|
2016-09-12 05:31:59 -05:00
|
|
|
return std::unordered_set<std::string>();
|
2015-06-16 09:17:56 -05:00
|
|
|
}
|
|
|
|
|
2017-01-25 13:36:37 -06:00
|
|
|
/// \brief A handle that copies a fixed number data per index.
|
|
|
|
///
|
|
|
|
/// It works on Iterators to allow for communicating C arrays.
|
|
|
|
/// \tparam Iter1 Constant random access iterator type.
|
|
|
|
/// \tparam Iter1 Mutable random access iterator type.
|
|
|
|
template<class Iter1, class Iter2=Iter1>
|
|
|
|
class FixedSizeIterCopyHandle
|
|
|
|
{
|
|
|
|
typedef typename std::iterator_traits<Iter1>::value_type DataType2;
|
|
|
|
|
|
|
|
public:
|
|
|
|
typedef typename std::iterator_traits<Iter1>::value_type DataType;
|
|
|
|
|
|
|
|
/// \brief Constructor.
|
|
|
|
/// \param send_begin The begin iterator for sending.
|
|
|
|
/// \param receive_begin The begin iterator for receiving.
|
|
|
|
FixedSizeIterCopyHandle(const Iter1& send_begin,
|
|
|
|
const Iter2& receive_begin,
|
|
|
|
std::size_t size = 1)
|
|
|
|
: send_(send_begin), receive_(receive_begin), size_(size)
|
|
|
|
{
|
|
|
|
static_assert(std::is_same<DataType,DataType2>::value,
|
|
|
|
"Iter1 and Iter2 have to have the same value_type!");
|
|
|
|
}
|
|
|
|
|
|
|
|
template<class Buffer>
|
|
|
|
void gather(Buffer& buffer, std::size_t i)
|
|
|
|
{
|
|
|
|
for(auto index = i*size(i), end = (i+1)*size(i);
|
|
|
|
index < end; ++index)
|
|
|
|
{
|
|
|
|
buffer.write(send_[index]);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
template<class Buffer>
|
2017-06-16 06:52:51 -05:00
|
|
|
void scatter(Buffer& buffer, std::size_t i, std::size_t s OPM_OPTIM_UNUSED)
|
2017-01-25 13:36:37 -06:00
|
|
|
{
|
|
|
|
assert(s==size(i));
|
2017-04-12 04:13:10 -05:00
|
|
|
static_cast<void>(s);
|
2017-01-25 13:36:37 -06:00
|
|
|
|
|
|
|
for(auto index = i*size(i), end = (i+1)*size(i);
|
|
|
|
index < end; ++index)
|
|
|
|
{
|
|
|
|
buffer.read(receive_[index]);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
bool fixedsize()
|
|
|
|
{
|
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
|
|
|
std::size_t size(std::size_t)
|
|
|
|
{
|
|
|
|
return size_;
|
|
|
|
}
|
|
|
|
private:
|
|
|
|
Iter1 send_;
|
|
|
|
Iter2 receive_;
|
|
|
|
std::size_t size_;
|
|
|
|
};
|
|
|
|
|
|
|
|
|
2017-11-20 08:15:33 -06:00
|
|
|
/// \brief Data handle for gathering the rank that owns a cell
|
|
|
|
template<class Mapper>
|
|
|
|
class CellOwnerDataHandle
|
|
|
|
{
|
|
|
|
public:
|
|
|
|
using DataType = int;
|
|
|
|
|
2017-11-22 07:18:32 -06:00
|
|
|
CellOwnerDataHandle(const Mapper& globalMapper, std::vector<int>& globalData)
|
|
|
|
: globalMapper_(globalMapper), globalData_(globalData)
|
2017-11-20 08:15:33 -06:00
|
|
|
{
|
|
|
|
int argc = 0;
|
|
|
|
char** argv = nullptr;
|
|
|
|
my_rank_ = Dune::MPIHelper::instance(argc,argv).rank();
|
|
|
|
}
|
|
|
|
bool fixedsize(int /*dim*/, int /*codim*/)
|
|
|
|
{
|
|
|
|
return true;
|
|
|
|
}
|
|
|
|
template<class T>
|
|
|
|
std::size_t size(const T& e)
|
|
|
|
{
|
|
|
|
if ( T::codimension == 0)
|
|
|
|
{
|
|
|
|
return 1;
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
OPM_THROW(std::logic_error, "Data handle can only be used for elements");
|
|
|
|
}
|
|
|
|
}
|
|
|
|
template<class B, class T>
|
|
|
|
void gather(B& buffer, const T& e)
|
|
|
|
{
|
|
|
|
buffer.write(my_rank_);
|
|
|
|
}
|
|
|
|
template<class B, class T>
|
|
|
|
void scatter(B& buffer, const T& e, std::size_t /* size */)
|
|
|
|
{
|
2017-11-22 07:18:32 -06:00
|
|
|
const auto& index = globalMapper_.index(e);
|
2017-11-20 08:15:33 -06:00
|
|
|
buffer.read(globalData_[index]);
|
|
|
|
}
|
|
|
|
bool contains(int dim, int codim)
|
|
|
|
{
|
|
|
|
return codim==0;
|
|
|
|
}
|
|
|
|
|
|
|
|
private:
|
|
|
|
int my_rank_;
|
|
|
|
const Mapper& globalMapper_;
|
|
|
|
std::vector<int>& globalData_;
|
|
|
|
};
|
|
|
|
|
2017-11-21 08:06:36 -06:00
|
|
|
#if HAVE_OPM_GRID && HAVE_MPI
|
2016-02-18 09:52:54 -06:00
|
|
|
/// \brief a data handle to distribute the threshold pressures
|
|
|
|
class ThresholdPressureDataHandle
|
|
|
|
{
|
|
|
|
public:
|
|
|
|
/// \brief type of the data we send
|
|
|
|
typedef double DataType;
|
|
|
|
/// \brief Constructor
|
|
|
|
/// \param sendGrid The grid that the data is attached to when sending.
|
|
|
|
/// \param recvGrid The grid that the data is attached to when receiving.
|
|
|
|
/// \param sendPressures The container where we will retrieve the values to be sent.
|
|
|
|
/// \param numFaces Number of faces of the distributed grid.
|
|
|
|
ThresholdPressureDataHandle(const Dune::CpGrid& sendGrid,
|
|
|
|
const Dune::CpGrid& recvGrid,
|
|
|
|
const std::vector<double>& sendPressures,
|
|
|
|
std::vector<double>& recvPressures)
|
|
|
|
: sendGrid_(sendGrid), recvGrid_(recvGrid), sendPressures_(sendPressures),
|
|
|
|
recvPressures_(recvPressures)
|
|
|
|
{}
|
|
|
|
|
|
|
|
bool fixedsize(int /*dim*/, int /*codim*/)
|
|
|
|
{
|
|
|
|
return false;
|
|
|
|
}
|
|
|
|
template<class T>
|
|
|
|
std::size_t size(const T& e)
|
|
|
|
{
|
|
|
|
if ( T::codimension == 0)
|
|
|
|
{
|
|
|
|
return sendGrid_.numCellFaces(e.index());
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
OPM_THROW(std::logic_error, "Data handle can only be used for elements");
|
|
|
|
}
|
|
|
|
}
|
|
|
|
template<class B, class T>
|
|
|
|
void gather(B& buffer, const T& e)
|
|
|
|
{
|
|
|
|
assert( T::codimension == 0);
|
|
|
|
for ( int i=0; i< sendGrid_.numCellFaces(e.index()); ++i )
|
|
|
|
{
|
|
|
|
buffer.write(sendPressures_[sendGrid_.cellFace(e.index(), i)]);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
template<class B, class T>
|
|
|
|
void scatter(B& buffer, const T& e, std::size_t /* size */)
|
|
|
|
{
|
|
|
|
assert( T::codimension == 0);
|
|
|
|
for ( int i=0; i< recvGrid_.numCellFaces(e.index()); ++i )
|
|
|
|
{
|
|
|
|
double val;
|
|
|
|
buffer.read(val);
|
|
|
|
recvPressures_[recvGrid_.cellFace(e.index(), i)]=val;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
bool contains(int dim, int codim)
|
|
|
|
{
|
|
|
|
return dim==3 && codim==0;
|
|
|
|
}
|
|
|
|
|
|
|
|
private:
|
|
|
|
/// \brief The grid that the data we send is associated with.
|
|
|
|
const Dune::CpGrid& sendGrid_;
|
|
|
|
/// \brief The grid that the data we receive is associated with.
|
|
|
|
const Dune::CpGrid& recvGrid_;
|
|
|
|
/// \brief The data to send.
|
|
|
|
const std::vector<double>& sendPressures_;
|
|
|
|
/// \brief The data to receive.
|
|
|
|
std::vector<double>& recvPressures_;
|
|
|
|
};
|
|
|
|
|
2015-05-07 04:21:09 -05:00
|
|
|
/// \brief a data handle to distribute Derived Geology
|
|
|
|
class GeologyDataHandle
|
|
|
|
{
|
|
|
|
public:
|
|
|
|
/// \brief type of the data we send
|
|
|
|
typedef double DataType;
|
|
|
|
/// \brief Constructor
|
|
|
|
/// \param sendGrid The grid that the data is attached to when sending.
|
|
|
|
/// \param recvGrid The grid that the data is attached to when receiving.
|
|
|
|
/// \param sendGeology The state where we will retieve the values to be sent.
|
|
|
|
/// \param recvGeology The state where we will store the received values.
|
|
|
|
GeologyDataHandle(const Dune::CpGrid& sendGrid,
|
|
|
|
const Dune::CpGrid& recvGrid,
|
|
|
|
const DerivedGeology& sendGeology,
|
|
|
|
DerivedGeology& recvGeology)
|
|
|
|
: sendGrid_(sendGrid), recvGrid_(recvGrid), sendGeology_(sendGeology),
|
|
|
|
recvGeology_(recvGeology)
|
|
|
|
{}
|
|
|
|
|
|
|
|
bool fixedsize(int /*dim*/, int /*codim*/)
|
|
|
|
{
|
|
|
|
return false;
|
|
|
|
}
|
|
|
|
template<class T>
|
|
|
|
std::size_t size(const T& e)
|
|
|
|
{
|
|
|
|
if ( T::codimension == 0)
|
|
|
|
{
|
|
|
|
return 1 + sendGrid_.numCellFaces(e.index());
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
OPM_THROW(std::logic_error, "Data handle can only be used for elements");
|
|
|
|
}
|
|
|
|
}
|
|
|
|
template<class B, class T>
|
|
|
|
void gather(B& buffer, const T& e)
|
|
|
|
{
|
|
|
|
assert( T::codimension == 0);
|
|
|
|
buffer.write(sendGeology_.poreVolume()[e.index()]);
|
|
|
|
for ( int i=0; i< sendGrid_.numCellFaces(e.index()); ++i )
|
|
|
|
{
|
|
|
|
buffer.write(sendGeology_.transmissibility()[sendGrid_.cellFace(e.index(), i)]);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
template<class B, class T>
|
2015-05-19 12:40:59 -05:00
|
|
|
void scatter(B& buffer, const T& e, std::size_t /* size */)
|
2015-05-07 04:21:09 -05:00
|
|
|
{
|
|
|
|
assert( T::codimension == 0);
|
|
|
|
double val;
|
|
|
|
buffer.read(val);
|
|
|
|
recvGeology_.poreVolume()[e.index()]=val;
|
|
|
|
for ( int i=0; i< recvGrid_.numCellFaces(e.index()); ++i )
|
|
|
|
{
|
|
|
|
buffer.read(val);
|
|
|
|
recvGeology_.transmissibility()[recvGrid_.cellFace(e.index(), i)]=val;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
bool contains(int dim, int codim)
|
|
|
|
{
|
|
|
|
return dim==3 && codim==0;
|
|
|
|
}
|
|
|
|
private:
|
|
|
|
/// \brief The grid that the data we send is associated with.
|
|
|
|
const Dune::CpGrid& sendGrid_;
|
|
|
|
/// \brief The grid that the data we receive is associated with.
|
|
|
|
const Dune::CpGrid& recvGrid_;
|
|
|
|
/// \brief The data to send.
|
|
|
|
const DerivedGeology& sendGeology_;
|
|
|
|
/// \brief The data to receive.
|
|
|
|
DerivedGeology& recvGeology_;
|
|
|
|
};
|
|
|
|
|
2015-02-05 08:05:38 -06:00
|
|
|
/// \brief a data handle to distribute the BlackoilState
|
|
|
|
class BlackoilStateDataHandle
|
|
|
|
{
|
|
|
|
public:
|
|
|
|
/// \brief The data that we send.
|
|
|
|
typedef double DataType;
|
|
|
|
/// \brief Constructor.
|
|
|
|
/// \param sendGrid The grid that the data is attached to when sending.
|
|
|
|
/// \param recvGrid The grid that the data is attached to when receiving.
|
|
|
|
/// \param sendState The state where we will retieve the values to be sent.
|
2015-05-07 04:21:09 -05:00
|
|
|
/// \param recvState The state where we will store the received values.
|
2015-02-05 08:05:38 -06:00
|
|
|
BlackoilStateDataHandle(const Dune::CpGrid& sendGrid,
|
|
|
|
const Dune::CpGrid& recvGrid,
|
|
|
|
const BlackoilState& sendState,
|
|
|
|
BlackoilState& recvState)
|
|
|
|
: sendGrid_(sendGrid), recvGrid_(recvGrid), sendState_(sendState), recvState_(recvState)
|
2016-05-25 07:47:24 -05:00
|
|
|
{
|
|
|
|
// construction does not resize surfacevol and hydroCarbonState. Do it manually.
|
2016-05-25 08:07:29 -05:00
|
|
|
recvState.surfacevol().resize(recvGrid.numCells()*sendState.numPhases(),
|
2016-05-25 07:47:24 -05:00
|
|
|
std::numeric_limits<double>::max());
|
2016-05-25 08:07:29 -05:00
|
|
|
recvState.hydroCarbonState().resize(recvGrid.numCells());
|
2016-05-25 07:47:24 -05:00
|
|
|
}
|
2015-02-05 08:05:38 -06:00
|
|
|
|
|
|
|
bool fixedsize(int /*dim*/, int /*codim*/)
|
|
|
|
{
|
|
|
|
return false;
|
|
|
|
}
|
|
|
|
|
|
|
|
template<class T>
|
|
|
|
std::size_t size(const T& e)
|
|
|
|
{
|
|
|
|
if ( T::codimension == 0)
|
|
|
|
{
|
2016-05-24 14:52:09 -05:00
|
|
|
return 2 * sendState_.numPhases() + 5 + 2*sendGrid_.numCellFaces(e.index());
|
2015-02-05 08:05:38 -06:00
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
OPM_THROW(std::logic_error, "Data handle can only be used for elements");
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
template<class B, class T>
|
|
|
|
void gather(B& buffer, const T& e)
|
|
|
|
{
|
|
|
|
assert( T::codimension == 0);
|
|
|
|
|
2016-03-15 04:40:48 -05:00
|
|
|
for ( size_t i=0; i<sendState_.numPhases(); ++i )
|
2015-02-05 08:05:38 -06:00
|
|
|
{
|
|
|
|
buffer.write(sendState_.surfacevol()[e.index()*sendState_.numPhases()+i]);
|
|
|
|
}
|
|
|
|
buffer.write(sendState_.gasoilratio()[e.index()]);
|
|
|
|
buffer.write(sendState_.rv()[e.index()]);
|
|
|
|
buffer.write(sendState_.pressure()[e.index()]);
|
|
|
|
buffer.write(sendState_.temperature()[e.index()]);
|
2016-05-24 14:52:09 -05:00
|
|
|
//We can only send one type with this buffer. Ergo we convert the enum to a double.
|
|
|
|
double hydroCarbonState_ = sendState_.hydroCarbonState()[e.index()];
|
|
|
|
buffer.write(hydroCarbonState_);
|
|
|
|
|
2016-03-15 04:40:48 -05:00
|
|
|
for ( size_t i=0; i<sendState_.numPhases(); ++i )
|
2015-05-05 12:14:27 -05:00
|
|
|
{
|
|
|
|
buffer.write(sendState_.saturation()[e.index()*sendState_.numPhases()+i]);
|
|
|
|
}
|
2015-02-05 08:05:38 -06:00
|
|
|
for ( int i=0; i<sendGrid_.numCellFaces(e.index()); ++i )
|
|
|
|
{
|
|
|
|
buffer.write(sendState_.facepressure()[sendGrid_.cellFace(e.index(), i)]);
|
|
|
|
}
|
|
|
|
for ( int i=0; i<sendGrid_.numCellFaces(e.index()); ++i )
|
|
|
|
{
|
2015-03-02 04:04:54 -06:00
|
|
|
buffer.write(sendState_.faceflux()[sendGrid_.cellFace(e.index(), i)]);
|
2015-02-05 08:05:38 -06:00
|
|
|
}
|
|
|
|
}
|
|
|
|
template<class B, class T>
|
2015-09-02 06:02:27 -05:00
|
|
|
void scatter(B& buffer, const T& e, std::size_t size_arg)
|
2015-02-05 08:05:38 -06:00
|
|
|
{
|
|
|
|
assert( T::codimension == 0);
|
2016-05-25 09:36:33 -05:00
|
|
|
assert( size_arg == 2 * recvState_.numPhases() + 5 +2*recvGrid_.numCellFaces(e.index()));
|
2015-09-02 06:02:27 -05:00
|
|
|
static_cast<void>(size_arg);
|
2015-02-05 08:05:38 -06:00
|
|
|
|
2015-02-17 06:44:52 -06:00
|
|
|
double val;
|
2016-03-15 04:40:48 -05:00
|
|
|
for ( size_t i=0; i<recvState_.numPhases(); ++i )
|
2015-02-05 08:05:38 -06:00
|
|
|
{
|
|
|
|
buffer.read(val);
|
2015-05-05 12:14:27 -05:00
|
|
|
recvState_.surfacevol()[e.index()*sendState_.numPhases()+i]=val;
|
2015-02-05 08:05:38 -06:00
|
|
|
}
|
|
|
|
buffer.read(val);
|
|
|
|
recvState_.gasoilratio()[e.index()]=val;
|
|
|
|
buffer.read(val);
|
|
|
|
recvState_.rv()[e.index()]=val;
|
|
|
|
buffer.read(val);
|
|
|
|
recvState_.pressure()[e.index()]=val;
|
|
|
|
buffer.read(val);
|
|
|
|
recvState_.temperature()[e.index()]=val;
|
2016-05-24 14:52:09 -05:00
|
|
|
//We can only send one type with this buffer. Ergo we convert the enum to a double.
|
|
|
|
buffer.read(val);
|
|
|
|
recvState_.hydroCarbonState()[e.index()]=static_cast<HydroCarbonState>(val);
|
|
|
|
|
2016-03-15 04:40:48 -05:00
|
|
|
for ( size_t i=0; i<recvState_.numPhases(); ++i )
|
2015-05-05 12:14:27 -05:00
|
|
|
{
|
|
|
|
buffer.read(val);
|
|
|
|
recvState_.saturation()[e.index()*sendState_.numPhases()+i]=val;
|
|
|
|
}
|
2015-02-05 08:05:38 -06:00
|
|
|
for ( int i=0; i<recvGrid_.numCellFaces(e.index()); ++i )
|
|
|
|
{
|
|
|
|
buffer.read(val);
|
|
|
|
recvState_.facepressure()[recvGrid_.cellFace(e.index(), i)]=val;
|
|
|
|
}
|
|
|
|
for ( int i=0; i<recvGrid_.numCellFaces(e.index()); ++i )
|
|
|
|
{
|
|
|
|
buffer.read(val);
|
|
|
|
recvState_.faceflux()[recvGrid_.cellFace(e.index(), i)]=val;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
bool contains(int dim, int codim)
|
|
|
|
{
|
|
|
|
return dim==3 && codim==0;
|
|
|
|
}
|
|
|
|
private:
|
|
|
|
/// \brief The grid that the data is attached to when sending
|
|
|
|
const Dune::CpGrid& sendGrid_;
|
|
|
|
/// \brief The grid that the data is attached to when receiving
|
|
|
|
const Dune::CpGrid& recvGrid_;
|
|
|
|
/// \brief The state where we will retieve the values to be sent.
|
|
|
|
const BlackoilState& sendState_;
|
|
|
|
// \brief The state where we will store the received values.
|
|
|
|
BlackoilState& recvState_;
|
|
|
|
};
|
|
|
|
|
2015-10-28 10:23:07 -05:00
|
|
|
/// \brief A DUNE data handle for sending the blackoil properties
|
2015-02-05 08:05:38 -06:00
|
|
|
class BlackoilPropsDataHandle
|
|
|
|
{
|
|
|
|
public:
|
|
|
|
/// \brief The data that we send.
|
|
|
|
typedef double DataType;
|
|
|
|
/// \brief Constructor.
|
|
|
|
/// \param sendProps The properties where we will retieve the values to be sent.
|
|
|
|
/// \parame recvProps The properties where we will store the received values.
|
2015-02-19 02:35:18 -06:00
|
|
|
BlackoilPropsDataHandle(const BlackoilPropsAdFromDeck& sendProps,
|
2015-02-05 08:05:38 -06:00
|
|
|
BlackoilPropsAdFromDeck& recvProps)
|
2015-02-19 02:35:18 -06:00
|
|
|
: sendProps_(sendProps), recvProps_(recvProps),
|
2015-10-28 10:23:07 -05:00
|
|
|
size_(11) // full permeability tensor 9 + porosity 1 + pvt region index
|
2015-02-05 08:05:38 -06:00
|
|
|
{
|
|
|
|
// satOilMax might be non empty. In this case we will need to send it, too.
|
|
|
|
if ( sendProps.satOilMax_.size()>0 )
|
|
|
|
{
|
2016-02-18 09:52:54 -06:00
|
|
|
|
2015-10-28 10:23:07 -05:00
|
|
|
// satOilMax has to have the same size as the cellPvtRegionIdx_
|
2015-09-01 13:20:16 -05:00
|
|
|
recvProps_.satOilMax_.resize(recvProps_.cellPvtRegionIdx_.size(),
|
2015-02-05 08:05:38 -06:00
|
|
|
-std::numeric_limits<double>::max());
|
|
|
|
++size_;
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
bool fixedsize(int /*dim*/, int /*codim*/)
|
|
|
|
{
|
|
|
|
return true;
|
|
|
|
}
|
|
|
|
|
|
|
|
template<class T>
|
|
|
|
std::size_t size(const T&)
|
|
|
|
{
|
|
|
|
if ( T::codimension == 0)
|
|
|
|
{
|
|
|
|
return size_;
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
OPM_THROW(std::logic_error, "Data handle can only be used for elements");
|
|
|
|
}
|
|
|
|
}
|
|
|
|
template<class B, class T>
|
|
|
|
void gather(B& buffer, const T& e)
|
|
|
|
{
|
|
|
|
assert( T::codimension == 0);
|
|
|
|
|
|
|
|
buffer.write(sendProps_.cellPvtRegionIndex()[e.index()]);
|
2015-10-26 04:39:59 -05:00
|
|
|
for( std::size_t i = 0; i < 9; ++i )
|
|
|
|
{
|
|
|
|
buffer.write(sendProps_.rock_.permeability_[e.index()*9+i]);
|
|
|
|
}
|
|
|
|
buffer.write(sendProps_.rock_.porosity_[e.index()]);
|
|
|
|
if ( size_ > 11 ) {
|
2015-02-26 06:25:14 -06:00
|
|
|
buffer.write(sendProps_.satOilMax_[e.index()]);
|
2015-02-16 04:42:42 -06:00
|
|
|
}
|
2015-02-05 08:05:38 -06:00
|
|
|
}
|
|
|
|
template<class B, class T>
|
2015-09-02 06:02:27 -05:00
|
|
|
void scatter(B& buffer, const T& e, std::size_t size_arg)
|
2015-02-05 08:05:38 -06:00
|
|
|
{
|
|
|
|
assert( T::codimension == 0);
|
2015-09-02 06:02:27 -05:00
|
|
|
assert( size_arg==size_ ); (void) size_arg;
|
2015-02-05 08:05:38 -06:00
|
|
|
double val;
|
|
|
|
buffer.read(val);
|
|
|
|
recvProps_.cellPvtRegionIdx_[e.index()]=val;
|
2015-10-26 04:39:59 -05:00
|
|
|
for( std::size_t i = 0; i < 9; ++i )
|
|
|
|
{
|
|
|
|
buffer.read(val);
|
|
|
|
recvProps_.rock_.permeability_[e.index()*9+i]
|
|
|
|
= val;
|
|
|
|
}
|
|
|
|
buffer.read(val);
|
|
|
|
recvProps_.rock_.porosity_[e.index()]=val;
|
|
|
|
if ( size_ > 11 ) {
|
2015-02-26 06:25:14 -06:00
|
|
|
buffer.read(val);
|
|
|
|
recvProps_.satOilMax_[e.index()]=val;
|
|
|
|
}
|
2015-02-05 08:05:38 -06:00
|
|
|
}
|
|
|
|
bool contains(int dim, int codim)
|
|
|
|
{
|
|
|
|
return dim==3 && codim==0;
|
|
|
|
}
|
|
|
|
private:
|
|
|
|
/// \brief The properties where we will retieve the values to be sent.
|
|
|
|
const BlackoilPropsAdFromDeck& sendProps_;
|
2015-10-28 10:23:07 -05:00
|
|
|
/// \brief The properties where we will store the received values.
|
2015-02-05 08:05:38 -06:00
|
|
|
BlackoilPropsAdFromDeck& recvProps_;
|
2015-10-28 10:23:07 -05:00
|
|
|
/// \brief The number of entries to send.
|
|
|
|
///
|
|
|
|
/// full permeability tensor 9 + porosity 1 + pvt region index and
|
|
|
|
/// in some case satOilMax
|
2015-02-05 08:05:38 -06:00
|
|
|
std::size_t size_;
|
|
|
|
};
|
|
|
|
|
2015-06-16 09:17:56 -05:00
|
|
|
inline
|
2016-09-12 05:31:59 -05:00
|
|
|
std::unordered_set<std::string>
|
|
|
|
distributeGridAndData( Dune::CpGrid& grid,
|
2016-10-14 02:23:26 -05:00
|
|
|
const Opm::Deck& deck,
|
|
|
|
const EclipseState& eclipseState,
|
2017-10-29 15:06:19 -05:00
|
|
|
const Schedule& schedule,
|
2016-09-12 05:31:59 -05:00
|
|
|
BlackoilState& state,
|
|
|
|
BlackoilPropsAdFromDeck& properties,
|
|
|
|
DerivedGeology& geology,
|
|
|
|
std::shared_ptr<BlackoilPropsAdFromDeck::MaterialLawManager>& material_law_manager,
|
|
|
|
std::vector<double>& threshold_pressures,
|
|
|
|
boost::any& parallelInformation,
|
|
|
|
const bool useLocalPerm)
|
2015-06-16 09:17:56 -05:00
|
|
|
{
|
|
|
|
Dune::CpGrid global_grid ( grid );
|
|
|
|
global_grid.switchToGlobalView();
|
|
|
|
|
|
|
|
// distribute the grid and switch to the distributed view
|
2016-09-12 05:31:59 -05:00
|
|
|
using std::get;
|
2017-10-29 15:06:19 -05:00
|
|
|
auto wells = schedule.getWells();
|
|
|
|
auto my_defunct_wells = get<1>(grid.loadBalance(&wells, geology.transmissibility().data()));
|
2015-06-16 09:17:56 -05:00
|
|
|
grid.switchToDistributedView();
|
2015-10-12 08:13:42 -05:00
|
|
|
std::vector<int> compressedToCartesianIdx;
|
|
|
|
Opm::createGlobalCellArray(grid, compressedToCartesianIdx);
|
|
|
|
typedef BlackoilPropsAdFromDeck::MaterialLawManager MaterialLawManager;
|
|
|
|
auto distributed_material_law_manager = std::make_shared<MaterialLawManager>();
|
|
|
|
distributed_material_law_manager->initFromDeck(deck, eclipseState, compressedToCartesianIdx);
|
|
|
|
// copy the values from the global to the local MaterialLawManager
|
|
|
|
// We should actually communicate these to be future proof. But that is
|
|
|
|
// really, really cumbersome for the underlying vector<shared_ptr>
|
|
|
|
// where the classes pointed to even have more shared_ptr stored in them.
|
|
|
|
typedef Dune::CpGrid::ParallelIndexSet IndexSet;
|
|
|
|
const IndexSet& local_indices = grid.getCellIndexSet();
|
2015-10-28 10:33:58 -05:00
|
|
|
for ( auto index : local_indices )
|
2015-10-12 08:13:42 -05:00
|
|
|
{
|
2016-07-06 06:00:35 -05:00
|
|
|
distributed_material_law_manager->materialLawParamsPointerReferenceHack(index.local()) =
|
|
|
|
material_law_manager->materialLawParamsPointerReferenceHack(index.global());
|
2015-06-16 09:17:56 -05:00
|
|
|
|
2016-07-06 06:00:35 -05:00
|
|
|
distributed_material_law_manager->oilWaterScaledEpsInfoDrainagePointerReferenceHack(index.local()) =
|
|
|
|
material_law_manager->oilWaterScaledEpsInfoDrainagePointerReferenceHack(index.global());
|
2015-10-12 08:13:42 -05:00
|
|
|
}
|
|
|
|
BlackoilPropsAdFromDeck distributed_props(properties,
|
|
|
|
distributed_material_law_manager,
|
|
|
|
grid.numCells());
|
2016-03-15 04:40:48 -05:00
|
|
|
BlackoilState distributed_state(grid.numCells(), grid.numFaces(), state.numPhases());
|
2015-07-09 06:45:15 -05:00
|
|
|
BlackoilStateDataHandle state_handle(global_grid, grid,
|
|
|
|
state, distributed_state);
|
|
|
|
BlackoilPropsDataHandle props_handle(properties,
|
2015-10-12 04:55:54 -05:00
|
|
|
distributed_props);
|
2015-06-16 09:17:56 -05:00
|
|
|
grid.scatterData(state_handle);
|
|
|
|
grid.scatterData(props_handle);
|
|
|
|
// Create a distributed Geology. Some values will be updated using communication
|
|
|
|
// below
|
2015-10-12 04:55:54 -05:00
|
|
|
DerivedGeology distributed_geology(grid,
|
|
|
|
distributed_props, eclipseState,
|
|
|
|
useLocalPerm, geology.gravity());
|
2015-07-09 06:45:15 -05:00
|
|
|
GeologyDataHandle geo_handle(global_grid, grid,
|
2015-10-12 04:55:54 -05:00
|
|
|
geology, distributed_geology);
|
2015-06-16 09:17:56 -05:00
|
|
|
grid.scatterData(geo_handle);
|
|
|
|
|
2016-02-18 09:52:54 -06:00
|
|
|
std::vector<double> distributed_pressures;
|
|
|
|
|
|
|
|
if( !threshold_pressures.empty() ) // Might be empty if not specified
|
|
|
|
{
|
|
|
|
if( threshold_pressures.size() !=
|
|
|
|
static_cast<std::size_t>(UgGridHelpers::numFaces(global_grid)) )
|
|
|
|
{
|
|
|
|
OPM_THROW(std::runtime_error, "NNCs not yet supported for parallel runs. "
|
|
|
|
<< UgGridHelpers::numFaces(grid) << " faces but " <<
|
2016-02-18 15:05:01 -06:00
|
|
|
threshold_pressures.size()<<" threshold pressure values");
|
2016-02-18 09:52:54 -06:00
|
|
|
}
|
|
|
|
distributed_pressures.resize(UgGridHelpers::numFaces(grid));
|
|
|
|
ThresholdPressureDataHandle press_handle(global_grid, grid,
|
|
|
|
threshold_pressures,
|
|
|
|
distributed_pressures);
|
|
|
|
grid.scatterData(press_handle);
|
|
|
|
}
|
|
|
|
|
2015-06-16 09:17:56 -05:00
|
|
|
// copy states
|
2015-10-12 04:55:54 -05:00
|
|
|
properties = distributed_props;
|
|
|
|
geology = distributed_geology;
|
|
|
|
state = distributed_state;
|
2015-10-12 08:13:42 -05:00
|
|
|
material_law_manager = distributed_material_law_manager;
|
2016-02-18 09:52:54 -06:00
|
|
|
threshold_pressures = distributed_pressures;
|
2015-07-09 06:45:15 -05:00
|
|
|
extractParallelGridInformationToISTL(grid, parallelInformation);
|
2016-09-12 05:31:59 -05:00
|
|
|
|
|
|
|
return my_defunct_wells;
|
2015-06-16 09:17:56 -05:00
|
|
|
}
|
|
|
|
#endif
|
|
|
|
|
2015-02-05 08:05:38 -06:00
|
|
|
} // end namespace Opm
|
2015-06-16 09:17:56 -05:00
|
|
|
#endif
|