24#ifndef DUMUX_IO_PORENETWORKGRID_DATA_HH
25#define DUMUX_IO_PORENETWORKGRID_DATA_HH
31#include <unordered_map>
34#include <dune/common/exceptions.hh>
35#include <dune/grid/common/gridfactory.hh>
36#include <dune/grid/common/mcmgmapper.hh>
37#include <dune/grid/io/file/dgfparser/gridptr.hh>
38#include <dune/grid/io/file/dgfparser/parser.hh>
39#include <dune/grid/utility/persistentcontainer.hh>
40#include <dune/geometry/axisalignedcubegeometry.hh>
56 static constexpr int dim = Grid::dimension;
57 static constexpr int dimWorld = Grid::dimensionworld;
58 using Intersection =
typename Grid::LeafIntersection;
59 using Element =
typename Grid::template Codim<0>::Entity;
60 using Vertex =
typename Grid::template Codim<dim>::Entity;
61 using GridView =
typename Grid::LeafGridView;
62 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
64 using StringVector = std::vector<std::string>;
66 using Scalar =
typename Grid::ctype;
67 using PersistentParameterContainer = Dune::PersistentContainer<Grid, std::vector<typename Grid::ctype>>;
79 setParameterIndices_();
88 numSubregions_ = getParamFromGroup<std::size_t>(paramGroup_,
"Grid.NumSubregions", 0);
89 setParameterIndices_();
90 parametersForGeneratedGrid_ = std::make_unique<ParametersForGeneratedGrid>(gridView_(),
paramGroup);
97 const std::vector<double>&
parameters(
const Vertex& vertex)
const
99 if (isDgfData_ && !useCopiedDgfData_)
100 return dgfGrid_.parameters(vertex);
103 assert(!(*vertexParameters_)[vertex].empty() &&
"No parameters available.");
104 return (*vertexParameters_)[vertex];
111 const std::vector<double>&
parameters(
const Element& element)
const
113 if (isDgfData_ && !useCopiedDgfData_)
115 if (element.hasFather())
117 auto level0Element = element;
118 while(level0Element.hasFather())
119 level0Element = level0Element.father();
121 return dgfGrid_.parameters(level0Element);
125 return dgfGrid_.parameters(element);
130 assert(!(*elementParameters_)[element].empty() &&
"No parameters available.");
131 return (*elementParameters_)[element];
138 template <
class Gr
idImp,
class IntersectionImp>
139 const Dune::DGFBoundaryParameter::type&
parameters(
const Dune::Intersection<GridImp, IntersectionImp>& intersection)
const
142 return dgfGrid_.parameters(intersection);
144 DUNE_THROW(Dune::InvalidStateException,
"The parameters method is only available if the grid was constructed with a DGF file.");
152 std::vector<SmallLocalIndex> coordinationNumbers(gridView_().size(dim), 0);
154 for (
const auto &element : elements(gridView_()))
156 for (SmallLocalIndex vIdxLocal = 0; vIdxLocal < 2; ++vIdxLocal)
158 const auto vIdxGlobal = gridView_().indexSet().subIndex(element, vIdxLocal, dim);
159 coordinationNumbers[vIdxGlobal] += 1;
163 if (std::any_of(coordinationNumbers.begin(), coordinationNumbers.end(), [](
auto i){ return i == 0; }))
164 DUNE_THROW(Dune::InvalidStateException,
"One of the pores is not connected to another pore. SanitizeGrid will not help in this case. Check your grid file");
166 return coordinationNumbers;
175 DUNE_THROW(Dune::InvalidStateException,
"Assigning parameter not possible for dgf gids");
177 const auto numVertexParams = vertexParameterNames_.size();
178 const auto numElementParams = elementParameterNames_.size();
179 vertexParameters_ = makeParamContainer_(*factoryGrid_, numVertexParams, 1);
180 elementParameters_ = makeParamContainer_(*factoryGrid_, numElementParams, 0);
182 auto setParamHelper = [&](
const auto& entity,
const std::string& param,
const Scalar value)
184 setParameter_(entity, param, value);
187 auto getParamHelper = [&](
const auto& entity,
const std::string& param)
192 parametersForGeneratedGrid_->assignParameters(setParamHelper, getParamHelper, numSubregions_);
198 vertexParameters_->resize();
199 elementParameters_->resize();
200 vertexParameters_->shrinkToFit();
201 elementParameters_->shrinkToFit();
207 DUNE_THROW(Dune::InvalidStateException,
"copying dgf data only works when a dgf grid is actually used");
209 useCopiedDgfData_ =
true;
210 const auto someVertex = *(vertices(gridView_()).begin());
211 const auto someElement = *(elements(gridView_()).begin());
212 const auto numVertexParams = dgfGrid_.parameters(someVertex).size();
213 const auto numElementParams = dgfGrid_.parameters(someElement).size();
214 vertexParameters_ = makeParamContainer_(*dgfGrid_, numVertexParams, 1);
215 elementParameters_ = makeParamContainer_(*dgfGrid_, numElementParams, 0);
217 for (
const auto& element : elements(gridView_()))
219 for (
int i = 0; i < numElementParams; ++i)
220 (*elementParameters_)[element][i] = dgfGrid_.parameters(element)[i];
223 for (
const auto& vertex : vertices(gridView_()))
225 for (
int i = 0; i < numVertexParams; ++i)
226 (*vertexParameters_)[vertex][i] = dgfGrid_.parameters(vertex)[i];
237 if (
static_cast<bool>(parameterIndex_.count(paramName)))
238 return parameterIndex_.at(paramName);
241 std::stringstream list;
242 list <<
"List of parameter indices:\n";
243 for (
const auto& entry : parameterIndex_)
244 list << entry.first <<
" " << entry.second <<
"\n";
247 if (paramName.find(
"Throat") != std::string::npos)
248 msg =
"Make sure to include it in the vector of parameter names ElementParameters = " + paramName +
" ... ...";
249 else if (paramName.find(
"Pore") != std::string::npos)
250 msg =
"Make sure to include it in the vector of parameter names VertexParameters = " + paramName +
" ... ...";
260 {
return paramGroup_; }
267 return std::any_of(elementParameterNames_.begin(), elementParameterNames_.end(), [¶m](
const auto& i ){ return (i == param); });
275 return std::any_of(vertexParameterNames_.begin(), vertexParameterNames_.end(), [¶m](
const auto& i ){ return (i == param); });
281 Scalar
getParameter(
const Element& element,
const std::string& param)
const
287 Scalar
getParameter(
const Vertex& vertex,
const std::string& param)
const
295 {
return parametersForGeneratedGrid_->boundaryFaceMarkerAtPos(pos); }
301 {
return vertexParameterNames_; }
307 {
return elementParameterNames_; }
311 void setParameter_(
const Element& element,
const std::string& param,
const Scalar value)
312 { (*elementParameters_)[element][
parameterIndex(param)] = value; }
314 void setParameter_(
const Vertex& vertex,
const std::string& param,
const Scalar value)
317 void setParameterIndices_()
321 vertexParameterNames_ = StringVector{
"PoreInscribedRadius",
"PoreVolume",
"PoreLabel"};
322 elementParameterNames_ = StringVector{
"ThroatInscribedRadius",
"ThroatLength",
"ThroatLabel"};
323 if (numSubregions_ > 0)
325 vertexParameterNames_.push_back(
"PoreRegionId");
326 elementParameterNames_.push_back(
"ThroatRegionId");
332 if (
auto dgfFileVertexParameterNames = dgfFileParameterNames_(
"Vertex"); dgfFileVertexParameterNames.empty())
333 DUNE_THROW(Dune::InvalidStateException,
"No vertex parameter names specified in dgf file. Set '% Vertex parameters: Param1 Param2 ...'");
335 vertexParameterNames_ = std::move(dgfFileVertexParameterNames);
338 if (
auto dgfFileElementParameterNames = dgfFileParameterNames_(
"Element"); dgfFileElementParameterNames.empty())
339 DUNE_THROW(Dune::InvalidStateException,
"No element parameter names specified in dgf file. Set '% Element parameters: Param1 Param2 ...'");
341 elementParameterNames_ = std::move(dgfFileElementParameterNames);
344 if (
const auto& someElement = *(elements(gridView_()).begin()); elementParameterNames_.size() != dgfGrid_.nofParameters(someElement))
345 DUNE_THROW(Dune::InvalidStateException,
"Number of user-specified element parameters (" << elementParameterNames_.size()
346 <<
") does not match number of element paramters in dgf file (" << dgfGrid_.nofParameters(someElement) <<
")");
348 if (
const auto& someVertex = *(vertices(gridView_()).begin()); vertexParameterNames_.size() != dgfGrid_.nofParameters(someVertex))
349 DUNE_THROW(Dune::InvalidStateException,
"Number of user-specified vertex parameters (" << vertexParameterNames_.size()
350 <<
") does not match number of vertex paramters in dgf file (" << dgfGrid_.nofParameters(someVertex) <<
")");
353 for (
int i = 0; i < vertexParameterNames_.size(); ++i)
355 std::cout << vertexParameterNames_[i] <<
" is vertex parameter " << i << std::endl;
356 parameterIndex_[vertexParameterNames_[i]] = i;
359 for (
int i = 0; i < elementParameterNames_.size(); ++i)
361 std::cout << elementParameterNames_[i] <<
" is element parameter " << i << std::endl;
362 parameterIndex_[elementParameterNames_[i]] = i;
373 auto makeParamContainer_(
const Grid& grid,
int numParams,
int codim)
const
375 auto parameters = std::make_unique<PersistentParameterContainer>(grid, codim);
376 (*parameters).resize();
382 StringVector dgfFileParameterNames_(
const std::string& entity)
const
384 std::ifstream gridFile(getParamFromGroup<std::string>(paramGroup_,
"Grid.File"));
386 while (getline(gridFile,
line))
388 if (
line.find(entity +
" parameters:", 0) != std::string::npos)
390 std::string args =
line.substr(
line.find(
":")+1, std::string::npos);
391 StringVector paramNames;
392 std::istringstream iss(args);
394 while (std::getline(iss, item,
' '))
396 *std::back_inserter(paramNames)++ = item;
401 return StringVector();
407 const GridView gridView_()
const
408 {
return isDgfData_ ? dgfGrid_->leafGridView() : factoryGrid_->leafGridView(); }
411 Dune::GridPtr<Grid> dgfGrid_;
413 std::shared_ptr<Grid> factoryGrid_;
415 bool isDgfData_ =
false;
416 bool useCopiedDgfData_ =
false;
417 std::string paramGroup_;
419 std::unique_ptr<ParametersForGeneratedGrid> parametersForGeneratedGrid_;
421 std::vector<std::string> vertexParameterNames_;
422 std::vector<std::string> elementParameterNames_;
424 std::unique_ptr<PersistentParameterContainer> vertexParameters_;
425 std::unique_ptr<PersistentParameterContainer> elementParameters_;
427 std::size_t numSubregions_;
429 std::unordered_map<std::string, int> parameterIndex_;
Defines the index types used for grid and local indices.
Helper class to assign parameters to a generated grid.
This file contains functions related to calculate pore-throat properties.
Definition: discretization/porenetwork/fvelementgeometry.hh:34
constexpr Line line
Definition: couplingmanager1d3d_line.hh:43
Exception thrown if a run-time parameter is not specified correctly.
Definition: exceptions.hh:60
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:41
Class for grid data attached to dgf or gmsh grid files.
Definition: porenetwork/griddata.hh:55
bool gridHasElementParameter(const std::string ¶m) const
Return if a given element parameter is provided by the grid.
Definition: porenetwork/griddata.hh:265
const std::vector< std::string > & vertexParameterNames() const
Returns the names of the vertex parameters.
Definition: porenetwork/griddata.hh:300
void resizeParameterContainers()
Definition: porenetwork/griddata.hh:195
GridData(Dune::GridPtr< Grid > grid, const std::string ¶mGroup)
constructor for dgf grid data
Definition: porenetwork/griddata.hh:73
const std::vector< double > & parameters(const Element &element) const
Call the parameters function of the DGF grid pointer if available for element data.
Definition: porenetwork/griddata.hh:111
const Dune::DGFBoundaryParameter::type & parameters(const Dune::Intersection< GridImp, IntersectionImp > &intersection) const
Call the parameters function of the DGF grid pointer if available.
Definition: porenetwork/griddata.hh:139
std::vector< SmallLocalIndex > getCoordinationNumbers() const
Returns the coordination numbers for all pore bodies.
Definition: porenetwork/griddata.hh:150
Scalar getParameter(const Element &element, const std::string ¶m) const
Returns the value of an element parameter.
Definition: porenetwork/griddata.hh:281
const std::vector< std::string > & elementParameterNames() const
Returns the names of the element parameters.
Definition: porenetwork/griddata.hh:306
void copyDgfData()
Definition: porenetwork/griddata.hh:204
int parameterIndex(const std::string ¶mName) const
Return the index for a given parameter name.
Definition: porenetwork/griddata.hh:233
int poreLabelAtPosForGenericGrid(const GlobalPosition &pos) const
Returns the pore label at a given position for a generic grid. This is needed by the grid creator in ...
Definition: porenetwork/griddata.hh:294
const std::string & paramGroup() const
Return the parameter group.
Definition: porenetwork/griddata.hh:259
GridData(std::shared_ptr< Grid > grid, const std::string ¶mGroup)
constructor for non-dgf grid data
Definition: porenetwork/griddata.hh:83
Scalar getParameter(const Vertex &vertex, const std::string ¶m) const
Returns the value of an vertex parameter.
Definition: porenetwork/griddata.hh:287
void assignParameters()
Assign parameters for generically created grids.
Definition: porenetwork/griddata.hh:172
const std::vector< double > & parameters(const Vertex &vertex) const
Call the parameters function of the DGF grid pointer if available for vertex data.
Definition: porenetwork/griddata.hh:97
bool gridHasVertexParameter(const std::string ¶m) const
Return if a given vertex parameter is provided by the grid.
Definition: porenetwork/griddata.hh:273
Helper class to assign parameters to a generated grid.
Definition: parametersforgeneratedgrid.hh:53