48#ifndef OPM_CPGRIDDATA_HEADER
49#define OPM_CPGRIDDATA_HEADER
52#include <dune/common/parallel/mpihelper.hh>
54#include <dune/istl/owneroverlapcopy.hh>
57#include <dune/common/parallel/communication.hh>
58#include <dune/common/parallel/variablesizecommunicator.hh>
59#include <dune/grid/common/gridenums.hh>
62#include <opm/input/eclipse/EclipseState/Grid/NNC.hpp>
67#include "Entity2IndexDataHandle.hpp"
68#include "CpGridDataTraits.hpp"
71#include "Geometry.hpp"
74#include <initializer_list>
91class LevelGlobalIdSet;
92class PartitionTypeIndicator;
93template<
int,
int>
class Geometry;
94template<
int>
class Entity;
95template<
int>
class EntityRep;
100 const std::array<int, 3>&,
103void refinePatch_and_check(
const std::array<int,3>&,
104 const std::array<int,3>&,
105 const std::array<int,3>&);
108 const std::vector<std::array<int,3>>&,
109 const std::vector<std::array<int,3>>&,
110 const std::vector<std::array<int,3>>&,
111 const std::vector<std::string>&);
116void fieldProp_check(
const Dune::CpGrid& grid, Opm::EclipseGrid eclGrid,
const std::string& deck_string);
118void testInactiveCellsLgrs(
const std::string&,
119 const std::vector<std::array<int,3>>&,
120 const std::vector<std::array<int,3>>&,
121 const std::vector<std::array<int,3>>&,
122 const std::vector<std::string>&);
130template<
class T,
int i>
struct Mover;
148 const std::array<int, 3>&,
151 void ::refinePatch_and_check(
const std::array<int,3>&,
152 const std::array<int,3>&,
153 const std::array<int,3>&);
157 const std::vector<std::array<int,3>>&,
158 const std::vector<std::array<int,3>>&,
159 const std::vector<std::array<int,3>>&,
160 const std::vector<std::string>&);
166 void ::fieldProp_check(
const Dune::CpGrid& grid, Opm::EclipseGrid eclGrid,
const std::string& deck_string);
168 void ::testInactiveCellsLgrs(
const std::string&,
169 const std::vector<std::array<int,3>>&,
170 const std::vector<std::array<int,3>>&,
171 const std::vector<std::array<int,3>>&,
172 const std::vector<std::string>&);
179#ifndef MAX_DATA_COMMUNICATED_PER_ENTITY
201 explicit CpGridData(MPIHelper::MPICommunicator comm, std::vector<std::shared_ptr<CpGridData>>& data);
206 explicit CpGridData(std::vector<std::shared_ptr<CpGridData>>& data);
214 int size(
int codim)
const;
217 int size (GeometryType type)
const
220 return size(3 - type.dim());
231 void readEclipseFormat(
const std::string& filename,
bool periodic_extension,
bool turn_normals =
false);
242 void processEclipseFormat(
const Opm::Deck& deck,
bool periodic_extension,
bool turn_normals =
false,
bool clip_z =
false,
243 const std::vector<double>& poreVolume = std::vector<double>());
258 std::vector<std::size_t>
260 bool periodic_extension,
bool turn_normals =
false,
bool clip_z =
false,
261 bool pinchActive =
true);
277 Opm::EclipseState* ecl_state,
279 std::array<std::set<std::pair<int, int>>, 2>& nnc,
280 bool remove_ij_boundary,
bool turn_normals,
bool pinchActive,
281 double tolerance_unique_points);
290 void getIJK(
int c, std::array<int,3>& ijk)
const
292 ijk =
getIJK(global_cell_[c], logical_cartesian_size_);
295 int cellFace(
int cell,
int local_index)
const
300 int faceToCellSize(
int face)
const {
302 return face_to_cell_[faceEntity].
size();
322 std::array<int,3>
getIJK(
int idx_in_parent_cell,
const std::array<int,3>& cells_per_dim)
const
326 assert(cells_per_dim[0]);
327 assert(cells_per_dim[1]);
328 assert(cells_per_dim[2]);
330 std::array<int,3> ijk = {0,0,0};
331 ijk[0] = idx_in_parent_cell % cells_per_dim[0]; idx_in_parent_cell /= cells_per_dim[0];
332 ijk[1] = idx_in_parent_cell % cells_per_dim[1];
333 ijk[2] = idx_in_parent_cell /cells_per_dim[1];
342 bool disjointPatches(
const std::vector<std::array<int,3>>& startIJK_vec,
const std::vector<std::array<int,3>>& endIJK_vec)
const;
352 getPatchesCells(
const std::vector<std::array<int,3>>& startIJK_vec,
const std::vector<std::array<int,3>>& endIJK_vec)
const;
356 bool hasNNCs(
const std::vector<int>& cellIndices)
const;
364 void validStartEndIJKs(
const std::vector<std::array<int,3>>& startIJK_vec,
const std::vector<std::array<int,3>>& endIJK_vec)
const;
374 bool patchesShareFace(
const std::vector<std::array<int,3>>& startIJK_vec,
const std::vector<std::array<int,3>>& endIJK_vec)
const;
376 int sharedFaceTag(
const std::vector<std::array<int,3>>& startIJK_2Patches,
const std::vector<std::array<int,3>>& endIJK_2Patches)
const;
424 bool compatibleSubdivisions(
const std::vector<std::array<int,3>>& cells_per_dim_vec,
425 const std::vector<std::array<int,3>>& startIJK_vec,
426 const std::vector<std::array<int,3>>& endIJK_vec)
const;
428 std::array<Dune::FieldVector<double,3>,8> getReferenceRefinedCorners(
int idx_in_parent_cell,
const std::array<int,3>& cells_per_dim)
const;
437 std::array<int,3> getPatchDim(
const std::array<int,3>& startIJK,
const std::array<int,3>& endIJK)
const;
446 std::vector<int> getPatchCorners(
const std::array<int,3>& startIJK,
const std::array<int,3>& endIJK)
const;
455 std::vector<int> getPatchFaces(
const std::array<int,3>& startIJK,
const std::array<int,3>& endIJK)
const;
464 std::vector<int> getPatchCells(
const std::array<int,3>& startIJK,
const std::array<int,3>& endIJK)
const;
473 std::vector<int> getPatchBoundaryCorners(
const std::array<int,3>& startIJK,
const std::array<int,3>& endIJK)
const;
482 std::array<std::vector<int>,6> getBoundaryPatchFaces(
const std::array<int,3>& startIJK,
const std::array<int,3>& endIJK)
const;
485 std::array<std::vector<double>,3> getWidthsLengthsHeights(
const std::array<int,3>& startIJK,
const std::array<int,3>& endIJK)
const;
501 Geometry<3,3> cellifyPatch(
const std::array<int,3>& startIJK,
const std::array<int,3>& endIJK,
503 std::array<int,8>& cellifiedPatch_to_point,
504 std::array<int,8>& allcorners_cellifiedPatch)
const;
510 std::array<double,3> getAverageArr(
const std::vector<std::array<double,3>>& vec)
const;
522 if ((level_data_ptr_ ->
size() >1) && (level_ == 0) && (!child_to_parent_cells_.empty())) {
523 return level_data_ptr_->size() -1;
528 const std::vector<std::shared_ptr<Dune::cpgrid::CpGridData>>&
levelData()
const
530 if (level_data_ptr_->empty()) {
531 OPM_THROW(std::logic_error,
"Level data has not been initialized\n");
533 return *level_data_ptr_;
544 return parent_to_children_cells_[elemIdx];
547 const std::vector<std::tuple<int,std::vector<int>>>& getParentToChildren()
const {
548 return parent_to_children_cells_;
571 std::tuple< const std::shared_ptr<CpGridData>,
572 const std::vector<std::array<int,2>>,
573 const std::vector<std::tuple<int,std::vector<int>>>,
574 const std::tuple<int, std::vector<int>>,
575 const std::vector<std::array<int,2>>,
576 const std::vector<std::array<int,2>>>
577 refineSingleCell(
const std::array<int,3>& cells_per_dim,
const int& parent_idx)
const;
591 std::tuple< std::shared_ptr<CpGridData>,
592 const std::vector<std::array<int,2>>,
593 const std::vector<std::tuple<int,std::vector<int>>>,
594 const std::vector<std::tuple<int,std::vector<int>>>,
595 const std::vector<std::tuple<int,std::vector<int>>>,
596 const std::vector<std::array<int,2>>,
597 const std::vector<std::array<int,2>>>
598 refinePatch(
const std::array<int,3>& cells_per_dim,
const std::array<int,3>& startIJK,
const std::array<int,3>& endIJK)
const;
606 std::array<double,3> computeEclCentroid(
const int idx)
const;
614 std::array<double,3> computeEclCentroid(
const Entity<0>& elem)
const;
617 void computeUniqueBoundaryIds();
624 return use_unique_boundary_ids_;
631 use_unique_boundary_ids_ = uids;
632 if (use_unique_boundary_ids_ && unique_boundary_ids_.empty()) {
633 computeUniqueBoundaryIds();
655 return *local_id_set_;
661 return *global_id_set_;
669 return logical_cartesian_size_;
677 const std::vector<int>& cell_part);
684 template<
class DataHandle>
685 void communicate(DataHandle& data, InterfaceType iftype, CommunicationDirection dir);
687 void computeCellPartitionType();
689 void computePointPartitionType();
691 void computeCommunicationInterfaces(
int noexistingPoints);
697 using CollectiveCommunication = CpGridDataTraits::CollectiveCommunication;
703 using Communicator = CpGridDataTraits::Communicator;
706 using InterfaceMap = CpGridDataTraits::InterfaceMap;
709 using CommunicationType = CpGridDataTraits::CommunicationType;
712 using ParallelIndexSet = CpGridDataTraits::ParallelIndexSet;
715 using RemoteIndices = CpGridDataTraits::RemoteIndices;
720 CommunicationType& cellCommunication()
728 const CommunicationType& cellCommunication()
const
733 ParallelIndexSet& cellIndexSet()
735 return cellCommunication().indexSet();
738 const ParallelIndexSet& cellIndexSet()
const
740 return cellCommunication().indexSet();
743 RemoteIndices& cellRemoteIndices()
745 return cellCommunication().remoteIndices();
748 const RemoteIndices& cellRemoteIndices()
const
750 return cellCommunication().remoteIndices();
757 return aquifer_cells_;
763 void populateGlobalCellIndexSet();
772 template<
class DataHandle>
773 void gatherData(DataHandle& data,
CpGridData* global_view,
783 template<
int codim,
class DataHandle>
784 void gatherCodimData(DataHandle& data,
CpGridData* global_data,
793 template<
class DataHandle>
794 void scatterData(DataHandle& data,
const CpGridData* global_data,
795 const CpGridData* distributed_data,
const InterfaceMap& cell_inf,
796 const InterfaceMap& point_inf);
805 template<
int codim,
class DataHandle>
806 void scatterCodimData(DataHandle& data,
CpGridData* global_data,
817 template<
int codim,
class DataHandle>
819 const Interface& interface);
829 template<
int codim,
class DataHandle>
831 const InterfaceMap& interface);
835 void computeGeometry(
const CpGrid& grid,
837 const std::vector<int>& globalAquiferCells,
840 std::vector<int>& aquiferCells,
842 const std::vector< std::array<int,8> >& cell2Points);
858 std::vector< std::array<int,8> > cell_to_point_;
865 std::array<int, 3> logical_cartesian_size_{};
872 std::vector<int> global_cell_;
878 typedef FieldVector<double, 3> PointType;
884 std::unique_ptr<cpgrid::IndexSet> index_set_;
886 std::shared_ptr<const cpgrid::IdSet> local_id_set_;
888 std::shared_ptr<LevelGlobalIdSet> global_id_set_;
890 std::shared_ptr<PartitionTypeIndicator> partition_type_indicator_;
892 std::vector<int> mark_;
896 std::vector<std::shared_ptr<CpGridData>>* level_data_ptr_;
899 std::vector<int> level_to_leaf_cells_;
// {level LGR, {child0, child1, ...}}
901 std::vector<std::tuple<int,std::vector<int>>> parent_to_children_cells_;
903 std::array<int,3> cells_per_dim_;
// {level, cell index in that level}
906 std::vector<std::array<int,2>> leaf_to_level_cells_;
908 std::vector<std::array<int,2>> corner_history_;
// {level parent cell, parent cell index}
911 std::vector<std::array<int,2>> child_to_parent_cells_;
914 std::vector<int> cell_to_idxInParentCell_;
922 bool use_unique_boundary_ids_;
929 std::vector<double> zcorn;
932 std::vector<int> aquifer_cells_;
937 CommunicationType cell_comm_;
940 std::tuple<Interface,Interface,Interface,Interface,Interface> cell_interfaces_;
945 std::tuple<InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap>
949 std::tuple<InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap>
956 const EntityVariable<Geometry<3 - codim, 3>, codim>& geomVector()
const
962 template<
int>
friend class Entity;
963 template<
int>
friend class EntityRep;
964 friend class Intersection;
965 friend class PartitionTypeIndicator;
979T& getInterface(InterfaceType iftype,
980 std::tuple<T,T,T,T,T>& interfaces)
985 return std::get<0>(interfaces);
987 return std::get<1>(interfaces);
989 return std::get<2>(interfaces);
991 return std::get<3>(interfaces);
993 return std::get<4>(interfaces);
995 OPM_THROW(std::runtime_error,
"Invalid Interface type was used during communication");
1000template<
int codim,
class DataHandle>
1001void CpGridData::communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data, CommunicationDirection dir,
1002 const Interface& interface)
1004 this->
template communicateCodim<codim>(data, dir, interface.interfaces());
1007template<
int codim,
class DataHandle>
1008void CpGridData::communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data_wrapper, CommunicationDirection dir,
1009 const InterfaceMap& interface)
1011 Communicator comm(ccobj_, interface);
1013 if(dir==ForwardCommunication)
1014 comm.forward(data_wrapper);
1016 comm.backward(data_wrapper);
1020template<
class DataHandle>
1022 CommunicationDirection dir)
1025 if(data.contains(3,0))
1028 communicateCodim<0>(data_wrapper, dir, getInterface(iftype, cell_interfaces_));
1030 if(data.contains(3,3))
1033 communicateCodim<3>(data_wrapper, dir, getInterface(iftype, point_interfaces_));
1045#include <opm/grid/cpgrid/Entity.hpp>
1046#include <opm/grid/cpgrid/Indexsets.hpp>
1060 data=buffer_[index_++];
1062 void write(
const T& data)
1064 buffer_[index_++]=data;
1070 void resize(std::size_t size)
1072 buffer_.resize(size);
1076 std::vector<T> buffer_;
1077 typename std::vector<T>::size_type index_;
1079template<
class DataHandle,
int codim>
1084template<
class DataHandle>
1087 explicit BaseMover(DataHandle& data)
1091 void moveData(
const E& from,
const E& to)
1093 std::size_t size=data_.size(from);
1094 buffer.resize(size);
1095 data_.gather(buffer, from);
1097 data_.scatter(buffer, to, size);
1100 MoveBuffer<typename DataHandle::DataType> buffer;
1104template<
class DataHandle>
1105struct Mover<DataHandle,0> :
public BaseMover<DataHandle>
1107 Mover(DataHandle& data, CpGridData* gatherView,
1108 CpGridData* scatterView)
1109 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
1112 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
1114 Entity<0> from_entity=Entity<0>(*gatherView_, from_cell_index,
true);
1115 Entity<0> to_entity=Entity<0>(*scatterView_, to_cell_index,
true);
1116 this->moveData(from_entity, to_entity);
1118 CpGridData* gatherView_;
1119 CpGridData* scatterView_;
1122template<
class DataHandle>
1123struct Mover<DataHandle,1> :
public BaseMover<DataHandle>
1125 Mover(DataHandle& data, CpGridData* gatherView,
1126 CpGridData* scatterView)
1127 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
1130 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
1132 typedef typename OrientedEntityTable<0,1>::row_type row_type;
1133 EntityRep<0> from_cell=EntityRep<0>(from_cell_index,
true);
1134 EntityRep<0> to_cell=EntityRep<0>(to_cell_index,
true);
1135 const OrientedEntityTable<0,1>& table = gatherView_->cell_to_face_;
1136 row_type from_faces=table.operator[](from_cell);
1137 row_type to_faces=scatterView_->cell_to_face_[to_cell];
1139 for(
int i=0; i<from_faces.size(); ++i)
1140 this->moveData(from_faces[i], to_faces[i]);
1142 CpGridData *gatherView_;
1143 CpGridData *scatterView_;
1146template<
class DataHandle>
1147struct Mover<DataHandle,3> :
public BaseMover<DataHandle>
1149 Mover(DataHandle& data, CpGridData* gatherView,
1150 CpGridData* scatterView)
1151 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
1153 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
1155 const std::array<int,8>& from_cell_points=
1156 gatherView_->cell_to_point_[from_cell_index];
1157 const std::array<int,8>& to_cell_points=
1158 scatterView_->cell_to_point_[to_cell_index];
1159 for(std::size_t i=0; i<8; ++i)
1161 this->moveData(Entity<3>(*gatherView_, from_cell_points[i],
true),
1162 Entity<3>(*scatterView_, to_cell_points[i],
true));
1165 CpGridData* gatherView_;
1166 CpGridData* scatterView_;
1171template<
class DataHandle>
1172void CpGridData::scatterData(DataHandle& data,
const CpGridData* global_data,
1173 const CpGridData* distributed_data,
const InterfaceMap& cell_inf,
1174 const InterfaceMap& point_inf)
1177 if(data.contains(3,0))
1179 Entity2IndexDataHandle<DataHandle, 0> data_wrapper(*global_data, *distributed_data, data);
1180 communicateCodim<0>(data_wrapper, ForwardCommunication, cell_inf);
1182 if(data.contains(3,3))
1184 Entity2IndexDataHandle<DataHandle, 3> data_wrapper(*global_data, *distributed_data, data);
1185 communicateCodim<3>(data_wrapper, ForwardCommunication, point_inf);
1190template<
int codim,
class DataHandle>
1191void CpGridData::scatterCodimData(DataHandle& data, CpGridData* global_data,
1192 CpGridData* distributed_data)
1194 CpGridData *gather_view, *scatter_view;
1195 gather_view=global_data;
1196 scatter_view=distributed_data;
1198 mover::Mover<DataHandle,codim> mover(data, gather_view, scatter_view);
1201 for(
auto index=distributed_data->cellIndexSet().begin(),
1202 end = distributed_data->cellIndexSet().end();
1203 index!=end; ++index)
1205 std::size_t from=index->global();
1206 std::size_t to=index->local();
1214template<
int codim,
class T,
class F>
1215void visitInterior(CpGridData& distributed_data, T begin, T endit, F& func)
1217 for(T it=begin; it!=endit; ++it)
1219 Entity<codim> entity(distributed_data, it-begin,
true);
1220 PartitionType pt = entity.partitionType();
1221 if(pt==Dune::InteriorEntity)
1228template<
class DataHandle>
1229struct GlobalIndexSizeGatherer
1231 GlobalIndexSizeGatherer(DataHandle& data_,
1232 std::vector<int>& ownedGlobalIndices_,
1233 std::vector<int>& ownedSizes_)
1234 : data(data_), ownedGlobalIndices(ownedGlobalIndices_), ownedSizes(ownedSizes_)
1237 template<
class T,
class E>
1238 void operator()(T& i, E& entity)
1240 ownedGlobalIndices.push_back(i);
1241 ownedSizes.push_back(data.size(entity));
1244 std::vector<int>& ownedGlobalIndices;
1245 std::vector<int>& ownedSizes;
1248template<
class DataHandle>
1251 DataGatherer(mover::MoveBuffer<typename DataHandle::DataType>& buffer_,
1253 : buffer(buffer_), data(data_)
1256 template<
class T,
class E>
1257 void operator()(T& , E& entity)
1259 data.gather(buffer, entity);
1261 mover::MoveBuffer<typename DataHandle::DataType>& buffer;
1267template<
class DataHandle>
1268void CpGridData::gatherData(DataHandle& data, CpGridData* global_data,
1269 CpGridData* distributed_data)
1272 if(data.contains(3,0))
1273 gatherCodimData<0>(data, global_data, distributed_data);
1274 if(data.contains(3,3))
1275 gatherCodimData<3>(data, global_data, distributed_data);
1279template<
int codim,
class DataHandle>
1280void CpGridData::gatherCodimData(DataHandle& data, CpGridData* global_data,
1281 CpGridData* distributed_data)
1285 const std::vector<int>& mapping =
1286 distributed_data->global_id_set_->getMapping<codim>();
1290 std::vector<int> owned_global_indices;
1291 std::vector<int> owned_sizes;
1292 owned_global_indices.reserve(mapping.size());
1293 owned_sizes.reserve(mapping.size());
1295 GlobalIndexSizeGatherer<DataHandle> gisg(data, owned_global_indices, owned_sizes);
1296 visitInterior<codim>(*distributed_data, mapping.begin(), mapping.end(), gisg);
1299 int no_indices=owned_sizes.size();
1302 if ( owned_global_indices.empty() )
1303 owned_global_indices.resize(1);
1304 if ( owned_sizes.empty() )
1305 owned_sizes.resize(1);
1306 std::vector<int> no_indices_to_recv(distributed_data->ccobj_.size());
1307 distributed_data->ccobj_.allgather(&no_indices, 1, &(no_indices_to_recv[0]));
1311 std::vector<int> displ(distributed_data->ccobj_.size()+1, 0);
1312 std::transform(displ.begin(), displ.end()-1, no_indices_to_recv.begin(), displ.begin()+1,
1314 int global_size=displ[displ.size()-1];
1315 std::vector<int> global_indices(global_size);
1316 std::vector<int> global_sizes(global_size);
1317 MPI_Allgatherv(&(owned_global_indices[0]), no_indices, MPITraits<int>::getType(),
1318 &(global_indices[0]), &(no_indices_to_recv[0]), &(displ[0]),
1319 MPITraits<int>::getType(),
1320 distributed_data->ccobj_);
1321 MPI_Allgatherv(&(owned_sizes[0]), no_indices, MPITraits<int>::getType(),
1322 &(global_sizes[0]), &(no_indices_to_recv[0]), &(displ[0]),
1323 MPITraits<int>::getType(),
1324 distributed_data->ccobj_);
1325 std::vector<int>().swap(owned_global_indices);
1327 std::vector<int> no_data_send(distributed_data->ccobj_.size());
1328 for(
typename std::vector<int>::iterator begin=no_data_send.begin(),
1329 i=begin, end=no_data_send.end(); i!=end; ++i)
1330 *i = std::accumulate(global_sizes.begin()+displ[i-begin],
1331 global_sizes.begin()+displ[i-begin+1], std::size_t());
1333 std::vector<int>().swap(owned_sizes);
1336 std::transform(displ.begin(), displ.end()-1, no_data_send.begin(), displ.begin()+1,
1337 std::plus<std::size_t>());
1339 int no_data_recv = displ[displ.size()-1];
1342 mover::MoveBuffer<typename DataHandle::DataType> local_data_buffer, global_data_buffer;
1343 if ( no_data_send[distributed_data->ccobj_.rank()] )
1345 local_data_buffer.resize(no_data_send[distributed_data->ccobj_.rank()]);
1349 local_data_buffer.resize(1);
1351 global_data_buffer.resize(no_data_recv);
1353 DataGatherer<DataHandle> gatherer(local_data_buffer, data);
1354 visitInterior<codim>(*distributed_data, mapping.begin(), mapping.end(), gatherer);
1355 MPI_Allgatherv(&(local_data_buffer.buffer_[0]), no_data_send[distributed_data->ccobj_.rank()],
1356 MPITraits<typename DataHandle::DataType>::getType(),
1357 &(global_data_buffer.buffer_[0]), &(no_data_send[0]), &(displ[0]),
1358 MPITraits<typename DataHandle::DataType>::getType(),
1359 distributed_data->ccobj_);
1360 Entity2IndexDataHandle<DataHandle, codim> edata(*global_data, data);
1362 for(
int i=0; i< codim; ++i)
1363 offset+=global_data->size(i);
1365 typename std::vector<int>::const_iterator s=global_sizes.begin();
1366 for(
typename std::vector<int>::const_iterator i=global_indices.begin(),
1367 end=global_indices.end();
1370 edata.scatter(global_data_buffer, *i-offset, *s);
[ provides Dune::Grid ]
Definition CpGrid.hpp:198
Struct that hods all the data needed to represent a Cpgrid.
Definition CpGridData.hpp:138
void postAdapt()
Clean up refinement/coarsening markers - set every element to the mark 0 which represents 'doing noth...
Definition CpGridData.cpp:2735
const cpgrid::LevelGlobalIdSet & globalIdSet() const
Get the global index set.
Definition CpGridData.hpp:659
@ MAX_DATA_PER_CELL
The maximum data items allowed per cell (DUNE < 2.5.2)
Definition CpGridData.hpp:186
std::vector< int > getPatchesCells(const std::vector< std::array< int, 3 > > &startIJK_vec, const std::vector< std::array< int, 3 > > &endIJK_vec) const
Compute cell indices of selected patches of cells (Cartesian grid required).
Definition CpGridData.cpp:2183
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition CpGridData.hpp:217
const std::array< int, 3 > & logicalCartesianSize() const
The logical cartesian size of the grid.
Definition CpGridData.hpp:667
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir)
communicate objects for all codims on a given level
Definition CpGridData.hpp:1021
bool uniqueBoundaryIds() const
Is the grid currently using unique boundary ids?
Definition CpGridData.hpp:622
const std::tuple< int, std::vector< int > > & getChildrenLevelAndIndexList(int elemIdx) const
Retrieves the level and child indices of a given parent cell.
Definition CpGridData.hpp:543
std::array< int, 3 > getIJK(int idx_in_parent_cell, const std::array< int, 3 > &cells_per_dim) const
Extract Cartesian index triplet (i,j,k) given an index between 0 and NXxNYxNZ -1 where NX,...
Definition CpGridData.hpp:322
bool patchesShareFace(const std::vector< std::array< int, 3 > > &startIJK_vec, const std::vector< std::array< int, 3 > > &endIJK_vec) const
Determine if a finite amount of patches (of cells) share a face.
Definition CpGridData.cpp:2020
int size(int codim) const
number of leaf entities per codim in this process
Definition CpGridData.cpp:147
const std::vector< int > & globalCell() const
Return global_cell_ of any level grid, or the leaf grid view (in presence of refinement).
Definition CpGridData.hpp:311
void readEclipseFormat(const std::string &filename, bool periodic_extension, bool turn_normals=false)
Read the Eclipse grid format ('grdecl').
CpGridDataTraits::Communication Communication
The type of the collective communication.
Definition CpGridData.hpp:696
~CpGridData()
Destructor.
Definition CpGridData.cpp:100
void getIJK(int c, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell.
Definition CpGridData.hpp:290
const IndexSet & indexSet() const
Get the index set.
Definition CpGridData.hpp:647
std::tuple< std::shared_ptr< CpGridData >, const std::vector< std::array< int, 2 > >, const std::vector< std::tuple< int, std::vector< int > > >, const std::vector< std::tuple< int, std::vector< int > > >, const std::vector< std::tuple< int, std::vector< int > > >, const std::vector< std::array< int, 2 > >, const std::vector< std::array< int, 2 > > > refinePatch(const std::array< int, 3 > &cells_per_dim, const std::array< int, 3 > &startIJK, const std::array< int, 3 > &endIJK) const
Refine a (connected block-shaped) patch of cells.
Definition CpGridData.cpp:2475
int getMark(const cpgrid::Entity< 0 > &element) const
Return refinement mark for entity.
Definition CpGridData.cpp:2709
CpGridDataTraits::MPICommunicator MPICommunicator
The type of the mpi communicator.
Definition CpGridData.hpp:694
void processEclipseFormat(const grdecl &input_data, std::array< std::set< std::pair< int, int > >, 2 > &nnc, bool remove_ij_boundary, bool turn_normals, bool pinchActive, double tolerance_unique_points)
Read the Eclipse grid format ('grdecl').
Definition processEclipseFormat.cpp:409
std::tuple< const std::shared_ptr< CpGridData >, const std::vector< std::array< int, 2 > >, const std::vector< std::tuple< int, std::vector< int > > >, const std::tuple< int, std::vector< int > >, const std::vector< std::array< int, 2 > >, const std::vector< std::array< int, 2 > > > refineSingleCell(const std::array< int, 3 > &cells_per_dim, const int &parent_idx) const
Refine a single cell and return a shared pointer of CpGridData type.
Definition CpGridData.cpp:2344
void distributeGlobalGrid(CpGrid &grid, const CpGridData &view_data, const std::vector< int > &cell_part)
Redistribute a global grid.
Definition CpGridData.cpp:1512
const std::vector< std::shared_ptr< Dune::cpgrid::CpGridData > > & levelData() const
Add doc/or remove method and replace it with better approach.
Definition CpGridData.hpp:528
const std::vector< double > & zcornData() const
Return the internalized zcorn copy from the grid processing, if no cells were adjusted during the min...
Definition CpGridData.hpp:640
void setUniqueBoundaryIds(bool uids)
Set whether we want to have unique boundary ids.
Definition CpGridData.hpp:629
bool preAdapt()
Set mightVanish flags for elements that will be refined in the next adapt() call Need to be called af...
Definition CpGridData.cpp:2714
int getGridIdx() const
Add doc/or remove method and replace it with better approach.
Definition CpGridData.hpp:514
bool hasNNCs(const std::vector< int > &cellIndices) const
Check all cells selected for refinement have no NNCs (no neighbor connections).
Definition CpGridData.cpp:2194
bool mark(int refCount, const cpgrid::Entity< 0 > &element)
Mark entity for refinement or coarsening.
Definition CpGridData.cpp:2692
const cpgrid::IdSet & localIdSet() const
Get the local index set.
Definition CpGridData.hpp:653
bool adapt()
TO DO: Documentation. Triggers the grid refinement process - Currently, returns preAdapt()
Definition CpGridData.cpp:2730
void checkCuboidShape(const std::vector< int > &cellIdx_vec) const
Check that every cell to be refined has cuboid shape.
Definition CpGridData.cpp:1824
bool disjointPatches(const std::vector< std::array< int, 3 > > &startIJK_vec, const std::vector< std::array< int, 3 > > &endIJK_vec) const
Determine if a finite amount of patches (of cells) are disjoint, namely, they do not share any corner...
Definition CpGridData.cpp:1971
const std::vector< int > & sortedNumAquiferCells() const
Get sorted active cell indices of numerical aquifer.
Definition CpGridData.hpp:755
void validStartEndIJKs(const std::vector< std::array< int, 3 > > &startIJK_vec, const std::vector< std::array< int, 3 > > &endIJK_vec) const
Check startIJK and endIJK of each patch of cells to be refined are valid, i.e.
Definition CpGridData.cpp:2218
Definition DefaultGeometryPolicy.hpp:53
const EntityVariable< cpgrid::Geometry< 3 - codim, 3 >, codim > & geomVector() const
Definition DefaultGeometryPolicy.hpp:86
Wrapper that turns a data handle suitable for dune-grid into one based on integers instead of entitie...
Definition Entity2IndexDataHandle.hpp:56
Represents an entity of a given codim, with positive or negative orientation.
Definition EntityRep.hpp:99
A class design to hold a variable with a value for each entity of the given codimension,...
Definition EntityRep.hpp:272
This class encapsulates geometry for vertices, intersections, and cells.
Definition Geometry.hpp:75
The global id set for Dune.
Definition Indexsets.hpp:450
Only needs to provide interface for doing nothing.
Definition Iterators.hpp:108
Definition Indexsets.hpp:199
Definition Indexsets.hpp:56
Definition Indexsets.hpp:350
Represents the topological relationships between sets of entities, for example cells and faces.
Definition OrientedEntityTable.hpp:139
int size() const
Returns the number of rows in the table.
Definition SparseTable.hpp:121
A class design to hold a variable with a value for each entity of the given codimension,...
Definition EntityRep.hpp:307
A SparseTable stores a table with rows of varying size as efficiently as possible.
Definition SparseTable.hpp:55
The namespace Dune is the main namespace for all Dune code.
Definition CartesianIndexMapper.hpp:10
Holds the implementation of the CpGrid as a pimple.
Definition CellQuadrature.cpp:68
Low-level corner-point processing routines and supporting data structures.
MPIHelper::MPICommunicator MPICommunicator
The type of the collective communication.
Definition CpGridDataTraits.hpp:56
AttributeSet
The type of the set of the attributes.
Definition CpGridDataTraits.hpp:66
Definition CpGridData.hpp:130
Raw corner-point specification of a particular geological model.
Definition preprocess.h:56