47 #ifndef OPM_CPGRIDDATA_HEADER
48 #define OPM_CPGRIDDATA_HEADER
51 #include <opm/grid/utility/platform_dependent/disable_warnings.h>
53 #include <dune/common/version.hh>
54 #include <dune/common/parallel/mpihelper.hh>
56 #include <dune/istl/owneroverlapcopy.hh>
59 #if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 7)
60 #include <dune/common/parallel/communication.hh>
62 #include <dune/common/parallel/collectivecommunication.hh>
64 #include <dune/common/parallel/indexset.hh>
65 #include <dune/common/parallel/interface.hh>
66 #include <dune/common/parallel/plocalindex.hh>
67 #if DUNE_VERSION_NEWER(DUNE_GRID, 2, 7)
68 #include <dune/common/parallel/variablesizecommunicator.hh>
70 #include <opm/grid/utility/VariableSizeCommunicator.hpp>
72 #include <dune/grid/common/gridenums.hh>
74 #include <opm/grid/utility/platform_dependent/reenable_warnings.h>
77 #include <opm/parser/eclipse/EclipseState/Grid/NNC.hpp>
85 #include "OrientedEntityTable.hpp"
86 #include "DefaultGeometryPolicy.hpp"
89 #include <opm/grid/utility/OpmParserIncludes.hpp>
91 #include "Entity2IndexDataHandle.hpp"
92 #include "DataHandleWrappers.hpp"
93 #include "GlobalIdMapping.hpp"
107 class LevelGlobalIdSet;
108 class PartitionTypeIndicator;
109 template<
int,
int>
class Geometry;
110 template<
int>
class Entity;
111 template<
int>
class EntityRep;
115 template<
class T,
int i>
struct Mover;
133 #ifndef MAX_DATA_COMMUNICATED_PER_ENTITY
156 explicit CpGridData(MPIHelper::MPICommunicator comm);
163 int size(
int codim)
const;
166 int size (GeometryType type)
const
169 return size(3 - type.dim());
189 void readEclipseFormat(
const std::string& filename,
bool periodic_extension,
bool turn_normals =
false);
200 void processEclipseFormat(
const Opm::Deck& deck,
bool periodic_extension,
bool turn_normals =
false,
bool clip_z =
false,
201 const std::vector<double>& poreVolume = std::vector<double>());
214 std::vector<std::size_t>
processEclipseFormat(
const Opm::EclipseGrid* ecl_grid, Opm::EclipseState* ecl_state,
215 bool periodic_extension,
bool turn_normals =
false,
bool clip_z =
false,
bool pinchActive =
true);
228 std::array<std::set<std::pair<int, int>>, 2>& nnc,
229 bool remove_ij_boundary,
bool turn_normals,
bool pinchActive);
238 void getIJK(
int c, std::array<int,3>& ijk)
const
240 int gc = global_cell_[c];
241 ijk[0] = gc % logical_cartesian_size_[0]; gc /= logical_cartesian_size_[0];
242 ijk[1] = gc % logical_cartesian_size_[1];
243 ijk[2] = gc / logical_cartesian_size_[1];
247 void computeUniqueBoundaryIds();
254 return use_unique_boundary_ids_;
261 use_unique_boundary_ids_ = uids;
262 if (use_unique_boundary_ids_ && unique_boundary_ids_.empty()) {
263 computeUniqueBoundaryIds();
286 return logical_cartesian_size_;
294 const std::vector<int>& cell_part);
301 template<
class DataHandle>
302 void communicate(DataHandle& data, InterfaceType iftype, CommunicationDirection dir);
305 #if DUNE_VERSION_NEWER(DUNE_GRID, 2, 7)
307 using Communicator = VariableSizeCommunicator<>;
310 using Communicator = Opm::VariableSizeCommunicator<>;
314 using InterfaceMap = Communicator::InterfaceMap;
317 using CommunicationType = Dune::OwnerOverlapCopyCommunication<int,int>;
320 using ParallelIndexSet = CommunicationType::ParallelIndexSet;
323 using RemoteIndices = Dune::RemoteIndices<ParallelIndexSet>;
328 CommunicationType& cellCommunication()
336 const CommunicationType& cellCommunication()
const
341 ParallelIndexSet& cellIndexSet()
343 return cellCommunication().indexSet();
346 const ParallelIndexSet& cellIndexSet()
const
348 return cellCommunication().indexSet();
351 RemoteIndices& cellRemoteIndices()
353 return cellCommunication().remoteIndices();
356 const RemoteIndices& cellRemoteIndices()
const
358 return cellCommunication().remoteIndices();
362 #ifdef HAVE_DUNE_ISTL
364 typedef Dune::OwnerOverlapCopyAttributeSet::AttributeSet
AttributeSet;
373 void populateGlobalCellIndexSet();
382 template<
class DataHandle>
383 void gatherData(DataHandle& data,
CpGridData* global_view,
393 template<
int codim,
class DataHandle>
394 void gatherCodimData(DataHandle& data,
CpGridData* global_data,
403 template<
class DataHandle>
404 void scatterData(DataHandle& data,
CpGridData* global_data,
405 CpGridData* distributed_data,
const InterfaceMap& cell_inf,
406 const InterfaceMap& point_inf);
415 template<
int codim,
class DataHandle>
416 void scatterCodimData(DataHandle& data,
CpGridData* global_data,
427 template<
int codim,
class DataHandle>
429 const Interface& interface);
439 template<
int codim,
class DataHandle>
441 const InterfaceMap& interface);
445 void computeGeometry(
CpGrid& grid,
450 const std::vector< std::array<int,8> >& cell2Points);
466 std::vector< std::array<int,8> > cell_to_point_;
473 std::array<int, 3> logical_cartesian_size_;
480 std::vector<int> global_cell_;
486 typedef FieldVector<double, 3> PointType;
501 typedef MPIHelper::MPICommunicator MPICommunicator;
502 typedef Dune::CollectiveCommunication<MPICommunicator> CollectiveCommunication;
504 CollectiveCommunication ccobj_;
507 bool use_unique_boundary_ids_;
514 std::vector<double> zcorn;
519 CommunicationType cell_comm_;
522 std::tuple<Interface,Interface,Interface,Interface,Interface> cell_interfaces_;
531 std::tuple<InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap>
544 template<
int>
friend class Entity;
545 template<
int>
friend class EntityRep;
546 template<
int>
friend class EntityPointer;
547 friend class Intersection;
548 friend class PartitionTypeIndicator;
560 T& getInterface(InterfaceType iftype,
561 std::tuple<T,T,T,T,T>& interfaces)
566 return std::get<0>(interfaces);
568 return std::get<1>(interfaces);
570 return std::get<2>(interfaces);
572 return std::get<3>(interfaces);
574 return std::get<4>(interfaces);
576 OPM_THROW(std::runtime_error,
"Invalid Interface type was used during communication");
581 template<
int codim,
class DataHandle>
582 void CpGridData::communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data, CommunicationDirection dir,
583 const Interface& interface)
585 this->
template communicateCodim<codim>(data, dir, interface.interfaces());
588 template<
int codim,
class DataHandle>
589 void CpGridData::communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data_wrapper, CommunicationDirection dir,
590 const InterfaceMap& interface)
592 Communicator comm(ccobj_, interface);
594 if(dir==ForwardCommunication)
595 comm.forward(data_wrapper);
597 comm.backward(data_wrapper);
601 template<
class DataHandle>
603 CommunicationDirection dir)
606 if(data.contains(3,0))
609 communicateCodim<0>(data_wrapper, dir, getInterface(iftype, cell_interfaces_));
611 if(data.contains(3,3))
614 communicateCodim<3>(data_wrapper, dir, getInterface(iftype, point_interfaces_));
626 #include <opm/grid/cpgrid/Entity.hpp>
627 #include <opm/grid/cpgrid/Indexsets.hpp>
641 data=buffer_[index_++];
643 void write(
const T& data)
645 buffer_[index_++]=data;
651 void resize(std::size_t size)
653 buffer_.resize(size);
657 std::vector<T> buffer_;
658 typename std::vector<T>::size_type index_;
660 template<
class DataHandle,
int codim>
665 template<
class DataHandle>
668 BaseMover(DataHandle& data)
672 void moveData(
const E& from,
const E& to)
674 std::size_t size=data_.size(from);
676 data_.gather(buffer, from);
678 data_.scatter(buffer, to, size);
681 MoveBuffer<typename DataHandle::DataType> buffer;
685 template<
class DataHandle>
686 struct Mover<DataHandle,0> :
public BaseMover<DataHandle>
688 Mover<DataHandle,0>(DataHandle& data, CpGridData* gatherView,
689 CpGridData* scatterView)
690 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
693 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
695 Entity<0> from_entity=Entity<0>(*gatherView_, from_cell_index,
true);
696 Entity<0> to_entity=Entity<0>(*scatterView_, to_cell_index,
true);
697 this->moveData(from_entity, to_entity);
699 CpGridData* gatherView_;
700 CpGridData* scatterView_;
703 template<
class DataHandle>
704 struct Mover<DataHandle,1> :
public BaseMover<DataHandle>
706 Mover<DataHandle,1>(DataHandle& data, CpGridData* gatherView,
707 CpGridData* scatterView)
708 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
711 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
713 typedef typename OrientedEntityTable<0,1>::row_type row_type;
714 EntityRep<0> from_cell=EntityRep<0>(from_cell_index,
true);
715 EntityRep<0> to_cell=EntityRep<0>(to_cell_index,
true);
716 OrientedEntityTable<0,1>& table = gatherView_->cell_to_face_;
717 row_type from_faces=table.operator[](from_cell);
718 row_type to_faces=scatterView_->cell_to_face_[to_cell];
720 for(
int i=0; i<from_faces.size(); ++i)
721 this->moveData(from_faces[i], to_faces[i]);
723 CpGridData *gatherView_;
724 CpGridData *scatterView_;
727 template<
class DataHandle>
728 struct Mover<DataHandle,3> :
public BaseMover<DataHandle>
730 Mover<DataHandle,3>(DataHandle& data, CpGridData* gatherView,
731 CpGridData* scatterView)
732 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
734 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
736 const std::array<int,8>& from_cell_points=
737 gatherView_->cell_to_point_[from_cell_index];
738 const std::array<int,8>& to_cell_points=
739 scatterView_->cell_to_point_[to_cell_index];
740 for(std::size_t i=0; i<8; ++i)
742 this->moveData(Entity<3>(*gatherView_, from_cell_points[i],
true),
743 Entity<3>(*scatterView_, to_cell_points[i],
true));
746 CpGridData* gatherView_;
747 CpGridData* scatterView_;
752 template<
class DataHandle>
753 void CpGridData::scatterData(DataHandle& data, CpGridData* global_data,
754 CpGridData* distributed_data,
const InterfaceMap& cell_inf,
755 const InterfaceMap& point_inf)
758 if(data.contains(3,0))
760 Entity2IndexDataHandle<DataHandle, 0> data_wrapper(*global_data, *distributed_data, data);
761 communicateCodim<0>(data_wrapper, ForwardCommunication, cell_inf);
763 if(data.contains(3,3))
765 Entity2IndexDataHandle<DataHandle, 3> data_wrapper(*global_data, *distributed_data, data);
766 communicateCodim<3>(data_wrapper, ForwardCommunication, point_inf);
771 template<
int codim,
class DataHandle>
772 void CpGridData::scatterCodimData(DataHandle& data, CpGridData* global_data,
773 CpGridData* distributed_data)
775 CpGridData *gather_view, *scatter_view;
776 gather_view=global_data;
777 scatter_view=distributed_data;
779 mover::Mover<DataHandle,codim> mover(data, gather_view, scatter_view);
782 for(
auto index=distributed_data->cellIndexSet().begin(),
783 end = distributed_data->cellIndexSet().end();
786 std::size_t from=index->global();
787 std::size_t to=index->local();
795 template<
int codim,
class T,
class F>
796 void visitInterior(CpGridData& distributed_data, T begin, T endit, F& func)
798 for(T it=begin; it!=endit; ++it)
800 Entity<codim> entity(distributed_data, it-begin,
true);
801 PartitionType pt = entity.partitionType();
802 if(pt==Dune::InteriorEntity)
809 template<
class DataHandle>
810 struct GlobalIndexSizeGatherer
812 GlobalIndexSizeGatherer(DataHandle& data_,
813 std::vector<int>& ownedGlobalIndices_,
814 std::vector<int>& ownedSizes_)
815 : data(data_), ownedGlobalIndices(ownedGlobalIndices_), ownedSizes(ownedSizes_)
818 template<
class T,
class E>
819 void operator()(T& i, E& entity)
821 ownedGlobalIndices.push_back(i);
822 ownedSizes.push_back(data.size(entity));
825 std::vector<int>& ownedGlobalIndices;
826 std::vector<int>& ownedSizes;
829 template<
class DataHandle>
832 DataGatherer(mover::MoveBuffer<typename DataHandle::DataType>& buffer_,
834 : buffer(buffer_), data(data_)
837 template<
class T,
class E>
838 void operator()(T& , E& entity)
840 data.gather(buffer, entity);
842 mover::MoveBuffer<typename DataHandle::DataType>& buffer;
848 template<
class DataHandle>
849 void CpGridData::gatherData(DataHandle& data, CpGridData* global_data,
850 CpGridData* distributed_data)
853 if(data.contains(3,0))
854 gatherCodimData<0>(data, global_data, distributed_data);
855 if(data.contains(3,3))
856 gatherCodimData<3>(data, global_data, distributed_data);
860 template<
int codim,
class DataHandle>
861 void CpGridData::gatherCodimData(DataHandle& data, CpGridData* global_data,
862 CpGridData* distributed_data)
866 const std::vector<int>& mapping =
867 distributed_data->global_id_set_->getMapping<codim>();
871 std::vector<int> owned_global_indices;
872 std::vector<int> owned_sizes;
873 owned_global_indices.reserve(mapping.size());
874 owned_sizes.reserve(mapping.size());
876 GlobalIndexSizeGatherer<DataHandle> gisg(data, owned_global_indices, owned_sizes);
877 visitInterior<codim>(*distributed_data, mapping.begin(), mapping.end(), gisg);
880 int no_indices=owned_sizes.size();
883 if ( owned_global_indices.empty() )
884 owned_global_indices.resize(1);
885 if ( owned_sizes.empty() )
886 owned_sizes.resize(1);
887 std::vector<int> no_indices_to_recv(distributed_data->ccobj_.size());
888 distributed_data->ccobj_.allgather(&no_indices, 1, &(no_indices_to_recv[0]));
892 std::vector<int> displ(distributed_data->ccobj_.size()+1, 0);
893 std::transform(displ.begin(), displ.end()-1, no_indices_to_recv.begin(), displ.begin()+1,
895 int global_size=displ[displ.size()-1];
896 std::vector<int> global_indices(global_size);
897 std::vector<int> global_sizes(global_size);
898 MPI_Allgatherv(&(owned_global_indices[0]), no_indices, MPITraits<int>::getType(),
899 &(global_indices[0]), &(no_indices_to_recv[0]), &(displ[0]),
900 MPITraits<int>::getType(),
901 distributed_data->ccobj_);
902 MPI_Allgatherv(&(owned_sizes[0]), no_indices, MPITraits<int>::getType(),
903 &(global_sizes[0]), &(no_indices_to_recv[0]), &(displ[0]),
904 MPITraits<int>::getType(),
905 distributed_data->ccobj_);
906 std::vector<int>().swap(owned_global_indices);
908 std::vector<int> no_data_send(distributed_data->ccobj_.size());
909 for(
typename std::vector<int>::iterator begin=no_data_send.begin(),
910 i=begin, end=no_data_send.end(); i!=end; ++i)
911 *i = std::accumulate(global_sizes.begin()+displ[i-begin],
912 global_sizes.begin()+displ[i-begin+1], std::size_t());
914 std::vector<int>().swap(owned_sizes);
917 std::transform(displ.begin(), displ.end()-1, no_data_send.begin(), displ.begin()+1,
918 std::plus<std::size_t>());
920 int no_data_recv = displ[displ.size()-1];
923 mover::MoveBuffer<typename DataHandle::DataType> local_data_buffer, global_data_buffer;
924 if ( no_data_send[distributed_data->ccobj_.rank()] )
926 local_data_buffer.resize(no_data_send[distributed_data->ccobj_.rank()]);
930 local_data_buffer.resize(1);
932 global_data_buffer.resize(no_data_recv);
934 DataGatherer<DataHandle> gatherer(local_data_buffer, data);
935 visitInterior<codim>(*distributed_data, mapping.begin(), mapping.end(), gatherer);
936 MPI_Allgatherv(&(local_data_buffer.buffer_[0]), no_data_send[distributed_data->ccobj_.rank()],
937 MPITraits<typename DataHandle::DataType>::getType(),
938 &(global_data_buffer.buffer_[0]), &(no_data_send[0]), &(displ[0]),
939 MPITraits<typename DataHandle::DataType>::getType(),
940 distributed_data->ccobj_);
941 Entity2IndexDataHandle<DataHandle, codim> edata(*global_data, data);
943 for(
int i=0; i< codim; ++i)
944 offset+=global_data->size(i);
946 typename std::vector<int>::const_iterator s=global_sizes.begin();
947 for(
typename std::vector<int>::const_iterator i=global_indices.begin(),
948 end=global_indices.end();
951 edata.scatter(global_data_buffer, *i-offset, *s);
[ provides Dune::Grid ]
Definition: CpGrid.hpp:207
Struct that hods all the data needed to represent a Cpgrid.
Definition: CpGridData.hpp:123
@ MAX_DATA_PER_CELL
The maximum data items allowed per cell (DUNE < 2.5.2)
Definition: CpGridData.hpp:140
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition: CpGridData.hpp:166
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir)
communicate objects for all codims on a given level
Definition: CpGridData.hpp:602
bool uniqueBoundaryIds() const
Is the grid currently using unique boundary ids?
Definition: CpGridData.hpp:252
CpGridData()
Constructor.
Definition: CpGridData.cpp:39
void readSintefLegacyFormat(const std::string &grid_prefix)
Read the Sintef legacy grid format ('topogeom').
Definition: readSintefLegacyFormat.cpp:66
int size(int codim) const
number of leaf entities per codim in this process
Definition: CpGridData.cpp:145
void readEclipseFormat(const std::string &filename, bool periodic_extension, bool turn_normals=false)
Read the Eclipse grid format ('grdecl').
~CpGridData()
Destructor.
Definition: CpGridData.cpp:94
void getIJK(int c, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell.
Definition: CpGridData.hpp:238
const IndexSet & indexSet() const
Get the index set.
Definition: CpGridData.hpp:277
const std::array< int, 3 > & logicalCartesianSize() const
The logical cartesian size of the grid.
Definition: CpGridData.hpp:284
void distributeGlobalGrid(CpGrid &grid, const CpGridData &view_data, const std::vector< int > &cell_part)
Redistribute a global grid.
Definition: CpGridData.cpp:1525
void processEclipseFormat(const grdecl &input_data, Opm::EclipseState *ecl_state, std::array< std::set< std::pair< int, int >>, 2 > &nnc, bool remove_ij_boundary, bool turn_normals, bool pinchActive)
Read the Eclipse grid format ('grdecl').
Definition: processEclipseFormat.cpp:249
AttributeSet
The type of the set of the attributes.
Definition: CpGridData.hpp:367
void setUniqueBoundaryIds(bool uids)
Set whether we want to have unique boundary ids.
Definition: CpGridData.hpp:259
void writeSintefLegacyFormat(const std::string &grid_prefix) const
Write the Sintef legacy grid format ('topogeom').
Definition: writeSintefLegacyFormat.cpp:72
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:270
Definition: DefaultGeometryPolicy.hpp:49
const EntityVariable< cpgrid::Geometry< 3 - codim, 3 >, codim > & geomVector() const
Definition: DefaultGeometryPolicy.hpp:74
Wrapper that turns a data handle suitable for dune-grid into one based on integers instead of entitie...
Definition: Entity2IndexDataHandle.hpp:54
This class encapsulates geometry for both vertices, intersections and cells.
Definition: Geometry.hpp:68
The global id set for Dune.
Definition: Indexsets.hpp:325
Definition: Indexsets.hpp:198
Definition: Indexsets.hpp:54
Definition: Indexsets.hpp:257
Definition: PartitionTypeIndicator.hpp:50
Copyright 2019 Equinor AS.
Definition: CartesianIndexMapper.hpp:10
Holds the implementation of the CpGrid as a pimple.
Definition: CellQuadrature.hpp:29
Low-level corner-point processing routines and supporting data structures.
Definition: CpGridData.hpp:115
Raw corner-point specification of a particular geological model.
Definition: preprocess.h:56