47#ifndef OPM_CPGRIDDATA_HEADER
48#define OPM_CPGRIDDATA_HEADER
51#include <dune/common/parallel/mpihelper.hh>
53#include <dune/istl/owneroverlapcopy.hh>
56#include <dune/common/parallel/communication.hh>
57#include <dune/common/parallel/variablesizecommunicator.hh>
58#include <dune/grid/common/gridenums.hh>
61#include <opm/input/eclipse/EclipseState/Grid/NNC.hpp>
66#include "Entity2IndexDataHandle.hpp"
67#include "CpGridDataTraits.hpp"
70#include "Geometry.hpp"
85class LevelGlobalIdSet;
86class PartitionTypeIndicator;
87template<
int,
int>
class Geometry;
88template<
int>
class Entity;
89template<
int>
class EntityRep;
94 const std::array<int, 3>&,
96void refinePatch_and_check(
const std::array<int,3>&,
97 const std::array<int,3>&,
98 const std::array<int,3>&);
101 const std::array<int,3>&,
102 const std::array<int,3>&,
103 const std::array<int,3>&);
114template<
class T,
int i>
struct Mover;
131 const std::array<int, 3>&,
134 void::refinePatch_and_check(
const std::array<int,3>&,
135 const std::array<int,3>&,
136 const std::array<int,3>&);
140 const std::array<int,3>&,
141 const std::array<int,3>&,
142 const std::array<int,3>&);
153#ifndef MAX_DATA_COMMUNICATED_PER_ENTITY
173 explicit CpGridData(MPIHelper::MPICommunicator comm, std::vector<std::shared_ptr<CpGridData>>& data);
178 CpGridData(std::vector<std::shared_ptr<CpGridData>>& data);
186 int size(
int codim)
const;
189 int size (GeometryType type)
const
192 return size(3 - type.dim());
212 void readEclipseFormat(
const std::string& filename,
bool periodic_extension,
bool turn_normals =
false);
223 void processEclipseFormat(
const Opm::Deck& deck,
bool periodic_extension,
bool turn_normals =
false,
bool clip_z =
false,
224 const std::vector<double>& poreVolume = std::vector<double>());
237 std::vector<std::size_t>
processEclipseFormat(
const Opm::EclipseGrid* ecl_grid, Opm::EclipseState* ecl_state,
238 bool periodic_extension,
bool turn_normals =
false,
bool clip_z =
false,
239 bool pinchActive =
true);
254 Opm::EclipseState* ecl_state,
256 std::array<std::set<std::pair<int, int>>, 2>& nnc,
257 bool remove_ij_boundary,
bool turn_normals,
bool pinchActive);
266 void getIJK(
int c, std::array<int,3>& ijk)
const
268 int gc = global_cell_[c];
269 ijk[0] = gc % logical_cartesian_size_[0]; gc /= logical_cartesian_size_[0];
270 ijk[1] = gc % logical_cartesian_size_[1];
271 ijk[2] = gc / logical_cartesian_size_[1];
282 const std::array<int,3> getPatchDim(
const std::array<int,3>& startIJK,
const std::array<int,3>& endIJK)
const;
291 const std::array<std::vector<int>,3> getPatchGeomIndices(
const std::array<int,3>& startIJK,
const std::array<int,3>& endIJK)
const;
307 const Geometry<3,3> cellifyPatch(
const std::array<int,3>& startIJK,
const std::array<int,3>& endIJK,
309 std::array<int,8>& cellifiedPatch_to_point,
310 std::array<int,8>& allcorners_cellifiedPatch)
const;
334 const std::tuple< const std::shared_ptr<CpGridData>,
335 const std::vector<std::array<int,2>>,
336 const std::vector<std::tuple<int,std::vector<int>>>,
337 const std::tuple<int, std::vector<int>>,
338 const std::vector<std::array<int,2>>,
339 const std::vector<std::array<int,2>>>
340 refineSingleCell(
const std::array<int,3>& cells_per_dim,
const int& parent_idx)
const;
354 const std::tuple< std::shared_ptr<CpGridData>,
355 const std::vector<std::array<int,2>>,
356 const std::vector<std::tuple<int,std::vector<int>>>,
357 const std::vector<std::tuple<int,std::vector<int>>>,
358 const std::vector<std::tuple<int,std::vector<int>>>,
359 const std::vector<std::array<int,2>>,
360 const std::vector<std::array<int,2>>>
361 refinePatch(
const std::array<int,3>& cells_per_dim,
const std::array<int,3>& startIJK,
const std::array<int,3>& endIJK)
const;
364 void computeUniqueBoundaryIds();
371 return use_unique_boundary_ids_;
378 use_unique_boundary_ids_ = uids;
379 if (use_unique_boundary_ids_ && unique_boundary_ids_.empty()) {
380 computeUniqueBoundaryIds();
403 return logical_cartesian_size_;
411 const std::vector<int>& cell_part);
418 template<
class DataHandle>
419 void communicate(DataHandle& data, InterfaceType iftype, CommunicationDirection dir);
425 using CollectiveCommunication = CpGridDataTraits::CollectiveCommunication;
431 using Communicator = CpGridDataTraits::Communicator;
434 using InterfaceMap = CpGridDataTraits::InterfaceMap;
437 using CommunicationType = CpGridDataTraits::CommunicationType;
440 using ParallelIndexSet = CpGridDataTraits::ParallelIndexSet;
443 using RemoteIndices = CpGridDataTraits::RemoteIndices;
448 CommunicationType& cellCommunication()
456 const CommunicationType& cellCommunication()
const
461 ParallelIndexSet& cellIndexSet()
463 return cellCommunication().indexSet();
466 const ParallelIndexSet& cellIndexSet()
const
468 return cellCommunication().indexSet();
471 RemoteIndices& cellRemoteIndices()
473 return cellCommunication().remoteIndices();
476 const RemoteIndices& cellRemoteIndices()
const
478 return cellCommunication().remoteIndices();
485 return aquifer_cells_;
491 void populateGlobalCellIndexSet();
500 template<
class DataHandle>
501 void gatherData(DataHandle& data,
CpGridData* global_view,
511 template<
int codim,
class DataHandle>
512 void gatherCodimData(DataHandle& data,
CpGridData* global_data,
521 template<
class DataHandle>
522 void scatterData(DataHandle& data,
CpGridData* global_data,
523 CpGridData* distributed_data,
const InterfaceMap& cell_inf,
524 const InterfaceMap& point_inf);
533 template<
int codim,
class DataHandle>
534 void scatterCodimData(DataHandle& data,
CpGridData* global_data,
545 template<
int codim,
class DataHandle>
547 const Interface& interface);
557 template<
int codim,
class DataHandle>
559 const InterfaceMap& interface);
563 void computeGeometry(
CpGrid& grid,
565 const std::vector<int>& globalAquiferCells,
568 std::vector<int>& aquiferCells,
570 const std::vector< std::array<int,8> >& cell2Points);
586 std::vector< std::array<int,8> > cell_to_point_;
593 std::array<int, 3> logical_cartesian_size_;
600 std::vector<int> global_cell_;
606 typedef FieldVector<double, 3> PointType;
612 std::unique_ptr<cpgrid::IndexSet> index_set_;
614 std::shared_ptr<const cpgrid::IdSet> local_id_set_;
616 std::shared_ptr<LevelGlobalIdSet> global_id_set_;
618 std::shared_ptr<PartitionTypeIndicator> partition_type_indicator_;
622 std::vector<std::shared_ptr<CpGridData>>* level_data_ptr_;
625 std::map<int,int> level_to_leaf_cells_;
627 std::vector<std::tuple<int,std::vector<int>>> parent_to_children_cells_;
630 std::vector<std::array<int,2>> leaf_to_level_cells_;
633 std::vector<std::array<int,2>> child_to_parent_cells_;
639 bool use_unique_boundary_ids_;
646 std::vector<double> zcorn;
649 std::vector<int> aquifer_cells_;
654 CommunicationType cell_comm_;
657 std::tuple<Interface,Interface,Interface,Interface,Interface> cell_interfaces_;
666 std::tuple<InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap,InterfaceMap>
679 template<
int>
friend class Entity;
680 template<
int>
friend class EntityRep;
681 friend class Intersection;
682 friend class PartitionTypeIndicator;
696T& getInterface(InterfaceType iftype,
697 std::tuple<T,T,T,T,T>& interfaces)
702 return std::get<0>(interfaces);
704 return std::get<1>(interfaces);
706 return std::get<2>(interfaces);
708 return std::get<3>(interfaces);
710 return std::get<4>(interfaces);
712 OPM_THROW(std::runtime_error,
"Invalid Interface type was used during communication");
717template<
int codim,
class DataHandle>
718void CpGridData::communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data, CommunicationDirection dir,
719 const Interface& interface)
721 this->
template communicateCodim<codim>(data, dir, interface.interfaces());
724template<
int codim,
class DataHandle>
725void CpGridData::communicateCodim(Entity2IndexDataHandle<DataHandle, codim>& data_wrapper, CommunicationDirection dir,
726 const InterfaceMap& interface)
728 Communicator comm(ccobj_, interface);
730 if(dir==ForwardCommunication)
731 comm.forward(data_wrapper);
733 comm.backward(data_wrapper);
737template<
class DataHandle>
739 CommunicationDirection dir)
742 if(data.contains(3,0))
745 communicateCodim<0>(data_wrapper, dir, getInterface(iftype, cell_interfaces_));
747 if(data.contains(3,3))
750 communicateCodim<3>(data_wrapper, dir, getInterface(iftype, point_interfaces_));
762#include <opm/grid/cpgrid/Entity.hpp>
763#include <opm/grid/cpgrid/Indexsets.hpp>
777 data=buffer_[index_++];
779 void write(
const T& data)
781 buffer_[index_++]=data;
787 void resize(std::size_t size)
789 buffer_.resize(size);
793 std::vector<T> buffer_;
794 typename std::vector<T>::size_type index_;
796template<
class DataHandle,
int codim>
801template<
class DataHandle>
804 BaseMover(DataHandle& data)
808 void moveData(
const E& from,
const E& to)
810 std::size_t size=data_.size(from);
812 data_.gather(buffer, from);
814 data_.scatter(buffer, to, size);
817 MoveBuffer<typename DataHandle::DataType> buffer;
821template<
class DataHandle>
822struct Mover<DataHandle,0> :
public BaseMover<DataHandle>
824 Mover<DataHandle,0>(DataHandle& data, CpGridData* gatherView,
825 CpGridData* scatterView)
826 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
829 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
831 Entity<0> from_entity=Entity<0>(*gatherView_, from_cell_index,
true);
832 Entity<0> to_entity=Entity<0>(*scatterView_, to_cell_index,
true);
833 this->moveData(from_entity, to_entity);
835 CpGridData* gatherView_;
836 CpGridData* scatterView_;
839template<
class DataHandle>
840struct Mover<DataHandle,1> :
public BaseMover<DataHandle>
842 Mover<DataHandle,1>(DataHandle& data, CpGridData* gatherView,
843 CpGridData* scatterView)
844 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
847 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
849 typedef typename OrientedEntityTable<0,1>::row_type row_type;
850 EntityRep<0> from_cell=EntityRep<0>(from_cell_index,
true);
851 EntityRep<0> to_cell=EntityRep<0>(to_cell_index,
true);
852 OrientedEntityTable<0,1>& table = gatherView_->cell_to_face_;
853 row_type from_faces=table.operator[](from_cell);
854 row_type to_faces=scatterView_->cell_to_face_[to_cell];
856 for(
int i=0; i<from_faces.size(); ++i)
857 this->moveData(from_faces[i], to_faces[i]);
859 CpGridData *gatherView_;
860 CpGridData *scatterView_;
863template<
class DataHandle>
864struct Mover<DataHandle,3> :
public BaseMover<DataHandle>
866 Mover<DataHandle,3>(DataHandle& data, CpGridData* gatherView,
867 CpGridData* scatterView)
868 : BaseMover<DataHandle>(data), gatherView_(gatherView), scatterView_(scatterView)
870 void operator()(std::size_t from_cell_index,std::size_t to_cell_index)
872 const std::array<int,8>& from_cell_points=
873 gatherView_->cell_to_point_[from_cell_index];
874 const std::array<int,8>& to_cell_points=
875 scatterView_->cell_to_point_[to_cell_index];
876 for(std::size_t i=0; i<8; ++i)
878 this->moveData(Entity<3>(*gatherView_, from_cell_points[i],
true),
879 Entity<3>(*scatterView_, to_cell_points[i],
true));
882 CpGridData* gatherView_;
883 CpGridData* scatterView_;
888template<
class DataHandle>
889void CpGridData::scatterData(DataHandle& data, CpGridData* global_data,
890 CpGridData* distributed_data,
const InterfaceMap& cell_inf,
891 const InterfaceMap& point_inf)
894 if(data.contains(3,0))
896 Entity2IndexDataHandle<DataHandle, 0> data_wrapper(*global_data, *distributed_data, data);
897 communicateCodim<0>(data_wrapper, ForwardCommunication, cell_inf);
899 if(data.contains(3,3))
901 Entity2IndexDataHandle<DataHandle, 3> data_wrapper(*global_data, *distributed_data, data);
902 communicateCodim<3>(data_wrapper, ForwardCommunication, point_inf);
907template<
int codim,
class DataHandle>
908void CpGridData::scatterCodimData(DataHandle& data, CpGridData* global_data,
909 CpGridData* distributed_data)
911 CpGridData *gather_view, *scatter_view;
912 gather_view=global_data;
913 scatter_view=distributed_data;
915 mover::Mover<DataHandle,codim> mover(data, gather_view, scatter_view);
918 for(
auto index=distributed_data->cellIndexSet().begin(),
919 end = distributed_data->cellIndexSet().end();
922 std::size_t from=index->global();
923 std::size_t to=index->local();
931template<
int codim,
class T,
class F>
932void visitInterior(CpGridData& distributed_data, T begin, T endit, F& func)
934 for(T it=begin; it!=endit; ++it)
936 Entity<codim> entity(distributed_data, it-begin,
true);
937 PartitionType pt = entity.partitionType();
938 if(pt==Dune::InteriorEntity)
945template<
class DataHandle>
946struct GlobalIndexSizeGatherer
948 GlobalIndexSizeGatherer(DataHandle& data_,
949 std::vector<int>& ownedGlobalIndices_,
950 std::vector<int>& ownedSizes_)
951 : data(data_), ownedGlobalIndices(ownedGlobalIndices_), ownedSizes(ownedSizes_)
954 template<
class T,
class E>
955 void operator()(T& i, E& entity)
957 ownedGlobalIndices.push_back(i);
958 ownedSizes.push_back(data.size(entity));
961 std::vector<int>& ownedGlobalIndices;
962 std::vector<int>& ownedSizes;
965template<
class DataHandle>
968 DataGatherer(mover::MoveBuffer<typename DataHandle::DataType>& buffer_,
970 : buffer(buffer_), data(data_)
973 template<
class T,
class E>
974 void operator()(T& , E& entity)
976 data.gather(buffer, entity);
978 mover::MoveBuffer<typename DataHandle::DataType>& buffer;
984template<
class DataHandle>
985void CpGridData::gatherData(DataHandle& data, CpGridData* global_data,
986 CpGridData* distributed_data)
989 if(data.contains(3,0))
990 gatherCodimData<0>(data, global_data, distributed_data);
991 if(data.contains(3,3))
992 gatherCodimData<3>(data, global_data, distributed_data);
996template<
int codim,
class DataHandle>
997void CpGridData::gatherCodimData(DataHandle& data, CpGridData* global_data,
998 CpGridData* distributed_data)
1002 const std::vector<int>& mapping =
1003 distributed_data->global_id_set_->getMapping<codim>();
1007 std::vector<int> owned_global_indices;
1008 std::vector<int> owned_sizes;
1009 owned_global_indices.reserve(mapping.size());
1010 owned_sizes.reserve(mapping.size());
1012 GlobalIndexSizeGatherer<DataHandle> gisg(data, owned_global_indices, owned_sizes);
1013 visitInterior<codim>(*distributed_data, mapping.begin(), mapping.end(), gisg);
1016 int no_indices=owned_sizes.size();
1019 if ( owned_global_indices.empty() )
1020 owned_global_indices.resize(1);
1021 if ( owned_sizes.empty() )
1022 owned_sizes.resize(1);
1023 std::vector<int> no_indices_to_recv(distributed_data->ccobj_.size());
1024 distributed_data->ccobj_.allgather(&no_indices, 1, &(no_indices_to_recv[0]));
1028 std::vector<int> displ(distributed_data->ccobj_.size()+1, 0);
1029 std::transform(displ.begin(), displ.end()-1, no_indices_to_recv.begin(), displ.begin()+1,
1031 int global_size=displ[displ.size()-1];
1032 std::vector<int> global_indices(global_size);
1033 std::vector<int> global_sizes(global_size);
1034 MPI_Allgatherv(&(owned_global_indices[0]), no_indices, MPITraits<int>::getType(),
1035 &(global_indices[0]), &(no_indices_to_recv[0]), &(displ[0]),
1036 MPITraits<int>::getType(),
1037 distributed_data->ccobj_);
1038 MPI_Allgatherv(&(owned_sizes[0]), no_indices, MPITraits<int>::getType(),
1039 &(global_sizes[0]), &(no_indices_to_recv[0]), &(displ[0]),
1040 MPITraits<int>::getType(),
1041 distributed_data->ccobj_);
1042 std::vector<int>().swap(owned_global_indices);
1044 std::vector<int> no_data_send(distributed_data->ccobj_.size());
1045 for(
typename std::vector<int>::iterator begin=no_data_send.begin(),
1046 i=begin, end=no_data_send.end(); i!=end; ++i)
1047 *i = std::accumulate(global_sizes.begin()+displ[i-begin],
1048 global_sizes.begin()+displ[i-begin+1], std::size_t());
1050 std::vector<int>().swap(owned_sizes);
1053 std::transform(displ.begin(), displ.end()-1, no_data_send.begin(), displ.begin()+1,
1054 std::plus<std::size_t>());
1056 int no_data_recv = displ[displ.size()-1];
1059 mover::MoveBuffer<typename DataHandle::DataType> local_data_buffer, global_data_buffer;
1060 if ( no_data_send[distributed_data->ccobj_.rank()] )
1062 local_data_buffer.resize(no_data_send[distributed_data->ccobj_.rank()]);
1066 local_data_buffer.resize(1);
1068 global_data_buffer.resize(no_data_recv);
1070 DataGatherer<DataHandle> gatherer(local_data_buffer, data);
1071 visitInterior<codim>(*distributed_data, mapping.begin(), mapping.end(), gatherer);
1072 MPI_Allgatherv(&(local_data_buffer.buffer_[0]), no_data_send[distributed_data->ccobj_.rank()],
1073 MPITraits<typename DataHandle::DataType>::getType(),
1074 &(global_data_buffer.buffer_[0]), &(no_data_send[0]), &(displ[0]),
1075 MPITraits<typename DataHandle::DataType>::getType(),
1076 distributed_data->ccobj_);
1077 Entity2IndexDataHandle<DataHandle, codim> edata(*global_data, data);
1079 for(
int i=0; i< codim; ++i)
1080 offset+=global_data->size(i);
1082 typename std::vector<int>::const_iterator s=global_sizes.begin();
1083 for(
typename std::vector<int>::const_iterator i=global_indices.begin(),
1084 end=global_indices.end();
1087 edata.scatter(global_data_buffer, *i-offset, *s);
[ provides Dune::Grid ]
Definition: CpGrid.hpp:214
Struct that hods all the data needed to represent a Cpgrid.
Definition: CpGridData.hpp:122
@ MAX_DATA_PER_CELL
The maximum data items allowed per cell (DUNE < 2.5.2)
Definition: CpGridData.hpp:160
int size(GeometryType type) const
number of leaf entities per geometry type in this process
Definition: CpGridData.hpp:189
const std::array< int, 3 > & logicalCartesianSize() const
The logical cartesian size of the grid.
Definition: CpGridData.hpp:401
void communicate(DataHandle &data, InterfaceType iftype, CommunicationDirection dir)
communicate objects for all codims on a given level
Definition: CpGridData.hpp:738
const 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:1820
bool uniqueBoundaryIds() const
Is the grid currently using unique boundary ids?
Definition: CpGridData.hpp:369
const 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:1949
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)
Read the Eclipse grid format ('grdecl').
Definition: processEclipseFormat.cpp:252
void readSintefLegacyFormat(const std::string &grid_prefix)
Read the Sintef legacy grid format ('topogeom').
Definition: readSintefLegacyFormat.cpp:67
int size(int codim) const
number of leaf entities per codim in this process
Definition: CpGridData.cpp:144
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:424
~CpGridData()
Destructor.
Definition: CpGridData.cpp:97
void getIJK(int c, std::array< int, 3 > &ijk) const
Extract Cartesian index triplet (i,j,k) of an active cell.
Definition: CpGridData.hpp:266
const IndexSet & indexSet() const
Get the index set.
Definition: CpGridData.hpp:394
CpGridDataTraits::MPICommunicator MPICommunicator
The type of the mpi communicator.
Definition: CpGridData.hpp:422
void distributeGlobalGrid(CpGrid &grid, const CpGridData &view_data, const std::vector< int > &cell_part)
Redistribute a global grid.
Definition: CpGridData.cpp:1513
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:387
void setUniqueBoundaryIds(bool uids)
Set whether we want to have unique boundary ids.
Definition: CpGridData.hpp:376
void writeSintefLegacyFormat(const std::string &grid_prefix) const
Write the Sintef legacy grid format ('topogeom').
Definition: writeSintefLegacyFormat.cpp:73
const std::vector< int > & sortedNumAquiferCells() const
Get sorted active cell indices of numerical aquifer.
Definition: CpGridData.hpp:483
Definition: DefaultGeometryPolicy.hpp:53
const EntityVariable< cpgrid::Geometry< 3 - codim, 3 >, codim > & geomVector() const
Definition: DefaultGeometryPolicy.hpp:81
Wrapper that turns a data handle suitable for dune-grid into one based on integers instead of entitie...
Definition: Entity2IndexDataHandle.hpp:56
This class encapsulates geometry for vertices, intersections, and cells.
Definition: Geometry.hpp:74
The global id set for Dune.
Definition: Indexsets.hpp:304
Only needs to provide interface for doing nothing.
Definition: Iterators.hpp:108
Definition: Indexsets.hpp:55
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.
MPIHelper::MPICommunicator MPICommunicator
The type of the collective communication.
Definition: CpGridDataTraits.hpp:55
AttributeSet
The type of the set of the attributes.
Definition: CpGridDataTraits.hpp:65
Definition: CpGridData.hpp:114
Raw corner-point specification of a particular geological model.
Definition: preprocess.h:56