mirror of
https://github.com/OPM/opm-simulators.git
synced 2024-11-25 10:40:21 -06:00
changed: avoid usage of equilGrid on non-root processes
This commit is contained in:
parent
40f91b0e32
commit
7efbee9fb0
@ -189,33 +189,34 @@ public:
|
||||
{
|
||||
std::set<int> send, recv;
|
||||
typedef typename Vanguard::EquilGrid::LeafGridView EquilGridView;
|
||||
const EquilGridView equilGridView = vanguard.equilGrid().leafGridView();
|
||||
|
||||
#if DUNE_VERSION_NEWER(DUNE_GRID, 2,6)
|
||||
typedef Dune::MultipleCodimMultipleGeomTypeMapper<EquilGridView> EquilElementMapper;
|
||||
EquilElementMapper equilElemMapper(equilGridView, Dune::mcmgElementLayout());
|
||||
#else
|
||||
typedef Dune::MultipleCodimMultipleGeomTypeMapper<EquilGridView, Dune::MCMGElementLayout> EquilElementMapper;
|
||||
EquilElementMapper equilElemMapper(equilGridView);
|
||||
#endif
|
||||
|
||||
// We need a mapping from local to global grid, here we
|
||||
// use equilGrid which represents a view on the global grid
|
||||
const size_t globalSize = vanguard.equilGrid().leafGridView().size(0);
|
||||
// reserve memory
|
||||
auto logSize = vanguard.grid().logicalCartesianSize();
|
||||
const size_t globalSize = logSize[0]*logSize[1]*logSize[2];
|
||||
globalCartesianIndex_.resize(globalSize, -1);
|
||||
|
||||
// loop over all elements (global grid) and store Cartesian index
|
||||
auto elemIt = vanguard.equilGrid().leafGridView().template begin<0>();
|
||||
const auto& elemEndIt = vanguard.equilGrid().leafGridView().template end<0>();
|
||||
for (; elemIt != elemEndIt; ++elemIt) {
|
||||
int elemIdx = equilElemMapper.index(*elemIt);
|
||||
int cartElemIdx = vanguard.equilCartesianIndexMapper().cartesianIndex(elemIdx);
|
||||
globalCartesianIndex_[elemIdx] = cartElemIdx;
|
||||
}
|
||||
|
||||
// the I/O rank receives from all other ranks
|
||||
if (isIORank()) {
|
||||
const EquilGridView equilGridView = vanguard.equilGrid().leafGridView();
|
||||
|
||||
#if DUNE_VERSION_NEWER(DUNE_GRID, 2,6)
|
||||
typedef Dune::MultipleCodimMultipleGeomTypeMapper<EquilGridView> EquilElementMapper;
|
||||
EquilElementMapper equilElemMapper(equilGridView, Dune::mcmgElementLayout());
|
||||
#else
|
||||
typedef Dune::MultipleCodimMultipleGeomTypeMapper<EquilGridView, Dune::MCMGElementLayout> EquilElementMapper;
|
||||
EquilElementMapper equilElemMapper(equilGridView);
|
||||
#endif
|
||||
|
||||
// loop over all elements (global grid) and store Cartesian index
|
||||
auto elemIt = vanguard.equilGrid().leafGridView().template begin<0>();
|
||||
const auto& elemEndIt = vanguard.equilGrid().leafGridView().template end<0>();
|
||||
for (; elemIt != elemEndIt; ++elemIt) {
|
||||
int elemIdx = equilElemMapper.index(*elemIt);
|
||||
int cartElemIdx = vanguard.equilCartesianIndexMapper().cartesianIndex(elemIdx);
|
||||
globalCartesianIndex_[elemIdx] = cartElemIdx;
|
||||
}
|
||||
|
||||
for (int i = 0; i < comm.size(); ++i) {
|
||||
if (i != ioRank)
|
||||
recv.insert(i);
|
||||
|
@ -119,7 +119,12 @@ public:
|
||||
* same as the grid which is used for the actual simulation.
|
||||
*/
|
||||
const EquilGrid& equilGrid() const
|
||||
{ return *equilGrid_; }
|
||||
{
|
||||
int mpiRank;
|
||||
MPI_Comm_rank(MPI_COMM_WORLD, &mpiRank);
|
||||
assert(mpiRank == 0);
|
||||
return *equilGrid_;
|
||||
}
|
||||
|
||||
/*!
|
||||
* \brief Indicates that the initial condition has been computed and the memory used
|
||||
@ -232,7 +237,12 @@ public:
|
||||
* \brief Returns mapper from compressed to cartesian indices for the EQUIL grid
|
||||
*/
|
||||
const CartesianIndexMapper& equilCartesianIndexMapper() const
|
||||
{ return *equilCartesianIndexMapper_; }
|
||||
{
|
||||
int mpiRank;
|
||||
MPI_Comm_rank(MPI_COMM_WORLD, &mpiRank);
|
||||
assert(mpiRank == 0);
|
||||
return *equilCartesianIndexMapper_;
|
||||
}
|
||||
|
||||
std::unordered_set<std::string> defunctWellNames() const
|
||||
{ return defunctWellNames_; }
|
||||
|
Loading…
Reference in New Issue
Block a user