12#ifndef DUMUX_IO_PORENETWORKGRID_DATA_HH
13#define DUMUX_IO_PORENETWORKGRID_DATA_HH
19#include <unordered_map>
23#include <dune/common/exceptions.hh>
24#include <dune/grid/common/gridfactory.hh>
25#include <dune/grid/common/mcmgmapper.hh>
26#include <dune/grid/io/file/dgfparser/gridptr.hh>
27#include <dune/grid/io/file/dgfparser/parser.hh>
28#include <dune/grid/utility/persistentcontainer.hh>
29#include <dune/geometry/axisalignedcubegeometry.hh>
45 static constexpr int dim = Grid::dimension;
46 static constexpr int dimWorld = Grid::dimensionworld;
47 using Intersection =
typename Grid::LeafIntersection;
48 using Element =
typename Grid::template Codim<0>::Entity;
49 using Vertex =
typename Grid::template Codim<dim>::Entity;
50 using GridView =
typename Grid::LeafGridView;
51 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
53 using StringVector = std::vector<std::string>;
55 using Scalar =
typename Grid::ctype;
56 using PersistentParameterContainer = Dune::PersistentContainer<Grid, std::vector<typename Grid::ctype>>;
68 setParameterIndices_();
77 numSubregions_ = getParamFromGroup<std::size_t>(paramGroup_,
"Grid.NumSubregions", 0);
78 setParameterIndices_();
79 parametersForGeneratedGrid_ = std::make_unique<ParametersForGeneratedGrid>(gridView_(),
paramGroup);
86 const std::vector<double>&
parameters(
const Vertex& vertex)
const
88 if (isDgfData_ && !useCopiedDgfData_)
89 return dgfGrid_.parameters(vertex);
92 assert(!(*vertexParameters_)[vertex].empty() &&
"No parameters available.");
93 return (*vertexParameters_)[vertex];
100 const std::vector<double>&
parameters(
const Element& element)
const
102 if (isDgfData_ && !useCopiedDgfData_)
104 if (element.hasFather())
106 auto level0Element = element;
107 while(level0Element.hasFather())
108 level0Element = level0Element.father();
110 return dgfGrid_.parameters(level0Element);
114 return dgfGrid_.parameters(element);
119 assert(!(*elementParameters_)[element].empty() &&
"No parameters available.");
120 return (*elementParameters_)[element];
127 template <
class Gr
idImp,
class IntersectionImp>
128 const Dune::DGFBoundaryParameter::type&
parameters(
const Dune::Intersection<GridImp, IntersectionImp>& intersection)
const
131 return dgfGrid_.parameters(intersection);
133 DUNE_THROW(Dune::InvalidStateException,
"The parameters method is only available if the grid was constructed with a DGF file.");
141 std::vector<SmallLocalIndex> coordinationNumbers(gridView_().size(dim), 0);
143 for (
const auto &element : elements(gridView_()))
145 for (SmallLocalIndex vIdxLocal = 0; vIdxLocal < 2; ++vIdxLocal)
147 const auto vIdxGlobal = gridView_().indexSet().subIndex(element, vIdxLocal, dim);
148 coordinationNumbers[vIdxGlobal] += 1;
152 if (std::any_of(coordinationNumbers.begin(), coordinationNumbers.end(), [](
auto i){ return i == 0; }))
153 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");
155 return coordinationNumbers;
164 DUNE_THROW(Dune::InvalidStateException,
"Assigning parameter not possible for dgf gids");
166 const auto numVertexParams = vertexParameterNames_.size();
167 const auto numElementParams = elementParameterNames_.size();
168 vertexParameters_ = makeParamContainer_(*factoryGrid_, numVertexParams, 1);
169 elementParameters_ = makeParamContainer_(*factoryGrid_, numElementParams, 0);
171 auto setParamHelper = [&](
const auto& entity,
const std::string& param,
const Scalar value)
173 setParameter_(entity, param, value);
176 auto getParamHelper = [&](
const auto& entity,
const std::string& param)
181 parametersForGeneratedGrid_->assignParameters(setParamHelper, getParamHelper, numSubregions_);
187 vertexParameters_->resize();
188 elementParameters_->resize();
189 vertexParameters_->shrinkToFit();
190 elementParameters_->shrinkToFit();
196 DUNE_THROW(Dune::InvalidStateException,
"copying dgf data only works when a dgf grid is actually used");
198 useCopiedDgfData_ =
true;
199 const auto someVertex = *(vertices(gridView_()).begin());
200 const auto someElement = *(elements(gridView_()).begin());
201 const auto numVertexParams = dgfGrid_.parameters(someVertex).size();
202 const auto numElementParams = dgfGrid_.parameters(someElement).size();
203 vertexParameters_ = makeParamContainer_(*dgfGrid_, numVertexParams, 1);
204 elementParameters_ = makeParamContainer_(*dgfGrid_, numElementParams, 0);
206 for (
const auto& element : elements(gridView_()))
208 for (
int i = 0; i < numElementParams; ++i)
209 (*elementParameters_)[element][i] = dgfGrid_.parameters(element)[i];
212 for (
const auto& vertex : vertices(gridView_()))
214 for (
int i = 0; i < numVertexParams; ++i)
215 (*vertexParameters_)[vertex][i] = dgfGrid_.parameters(vertex)[i];
226 if (
static_cast<bool>(parameterIndex_.count(paramName)))
227 return parameterIndex_.at(paramName);
230 std::stringstream list;
231 list <<
"List of parameter indices:\n";
232 for (
const auto& entry : parameterIndex_)
233 list << entry.first <<
" " << entry.second <<
"\n";
236 if (paramName.find(
"Throat") != std::string::npos)
237 msg =
"Make sure to include it in the vector of parameter names as a comment in the DGF file: Element Parameters = ... " + paramName +
" ...";
238 else if (paramName.find(
"Pore") != std::string::npos)
239 msg =
"Make sure to include it in the vector of parameter names as a comment in the DGF file: Vertex Parameters = ... " + paramName +
" ...";
249 {
return paramGroup_; }
256 return std::any_of(elementParameterNames_.begin(), elementParameterNames_.end(), [¶m](
const auto& i ){ return (i == param); });
264 return std::any_of(vertexParameterNames_.begin(), vertexParameterNames_.end(), [¶m](
const auto& i ){ return (i == param); });
270 Scalar
getParameter(
const Element& element,
const std::string& param)
const
276 Scalar
getParameter(
const Vertex& vertex,
const std::string& param)
const
284 {
return parametersForGeneratedGrid_->boundaryFaceMarkerAtPos(pos); }
290 {
return vertexParameterNames_; }
296 {
return elementParameterNames_; }
300 void setParameter_(
const Element& element,
const std::string& param,
const Scalar value)
301 { (*elementParameters_)[element][
parameterIndex(param)] = value; }
303 void setParameter_(
const Vertex& vertex,
const std::string& param,
const Scalar value)
306 void setParameterIndices_()
310 vertexParameterNames_ = StringVector{
"PoreInscribedRadius",
"PoreVolume",
"PoreLabel"};
311 elementParameterNames_ = StringVector{
"ThroatInscribedRadius",
"ThroatLength",
"ThroatLabel"};
312 if (numSubregions_ > 0)
314 vertexParameterNames_.push_back(
"PoreRegionId");
315 elementParameterNames_.push_back(
"ThroatRegionId");
321 vertexParameterNames_ = dgfFileParameterNames_(
"Vertex");
322 if (vertexParameterNames_.empty())
323 DUNE_THROW(Dune::InvalidStateException,
"No vertex parameter names specified in dgf file. Set '% Vertex parameters: Param1 Param2 ...'");
326 elementParameterNames_ = dgfFileParameterNames_(
"Element");
327 if (elementParameterNames_.empty())
328 DUNE_THROW(Dune::InvalidStateException,
"No element parameter names specified in dgf file. Set '% Element parameters: Param1 Param2 ...'");
331 if (
const auto& someElement = *(elements(gridView_()).begin()); elementParameterNames_.size() != dgfGrid_.nofParameters(someElement))
332 DUNE_THROW(Dune::InvalidStateException,
"Number of user-specified element parameters (" << elementParameterNames_.size()
333 <<
") does not match number of element parameters in dgf file (" << dgfGrid_.nofParameters(someElement) <<
")");
335 if (
const auto& someVertex = *(vertices(gridView_()).begin()); vertexParameterNames_.size() != dgfGrid_.nofParameters(someVertex))
336 DUNE_THROW(Dune::InvalidStateException,
"Number of user-specified vertex parameters (" << vertexParameterNames_.size()
337 <<
") does not match number of vertex parameters in dgf file (" << dgfGrid_.nofParameters(someVertex) <<
")");
340 for (
int i = 0; i < vertexParameterNames_.size(); ++i)
342 std::cout << vertexParameterNames_[i] <<
" is vertex parameter " << i << std::endl;
343 parameterIndex_[vertexParameterNames_[i]] = i;
346 for (
int i = 0; i < elementParameterNames_.size(); ++i)
348 std::cout << elementParameterNames_[i] <<
" is element parameter " << i << std::endl;
349 parameterIndex_[elementParameterNames_[i]] = i;
360 auto makeParamContainer_(
const Grid& grid,
int numParams,
int codim)
const
362 auto parameters = std::make_unique<PersistentParameterContainer>(grid, codim);
363 (*parameters).resize();
369 StringVector dgfFileParameterNames_(
const std::string& entity)
const
371 StringVector paramNames;
373 std::ifstream gridFile(getParamFromGroup<std::string>(paramGroup_,
"Grid.File"));
374 for (std::string
line; std::getline(gridFile,
line); )
380 if (
line.find(entity +
" parameters:", 0) != std::string::npos)
382 std::istringstream iss(
line.substr(
line.find(
":")+1, std::string::npos));
383 for (std::string item; std::getline(iss, item,
' '); )
385 paramNames.push_back(item);
395 const GridView gridView_()
const
396 {
return isDgfData_ ? dgfGrid_->leafGridView() : factoryGrid_->leafGridView(); }
399 Dune::GridPtr<Grid> dgfGrid_;
401 std::shared_ptr<Grid> factoryGrid_;
403 bool isDgfData_ =
false;
404 bool useCopiedDgfData_ =
false;
405 std::string paramGroup_;
407 std::unique_ptr<ParametersForGeneratedGrid> parametersForGeneratedGrid_;
409 std::vector<std::string> vertexParameterNames_;
410 std::vector<std::string> elementParameterNames_;
412 std::unique_ptr<PersistentParameterContainer> vertexParameters_;
413 std::unique_ptr<PersistentParameterContainer> elementParameters_;
415 std::size_t numSubregions_;
417 std::unordered_map<std::string, int> parameterIndex_;
Exception thrown if a run-time parameter is not specified correctly.
Definition: exceptions.hh:48
Class for grid data attached to dgf or gmsh grid files.
Definition: porenetwork/griddata.hh:44
bool gridHasElementParameter(const std::string ¶m) const
Return if a given element parameter is provided by the grid.
Definition: porenetwork/griddata.hh:254
const std::vector< std::string > & vertexParameterNames() const
Returns the names of the vertex parameters.
Definition: porenetwork/griddata.hh:289
void resizeParameterContainers()
Definition: porenetwork/griddata.hh:184
GridData(Dune::GridPtr< Grid > grid, const std::string ¶mGroup)
constructor for dgf grid data
Definition: porenetwork/griddata.hh:62
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:100
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:128
std::vector< SmallLocalIndex > getCoordinationNumbers() const
Returns the coordination numbers for all pore bodies.
Definition: porenetwork/griddata.hh:139
Scalar getParameter(const Element &element, const std::string ¶m) const
Returns the value of an element parameter.
Definition: porenetwork/griddata.hh:270
const std::vector< std::string > & elementParameterNames() const
Returns the names of the element parameters.
Definition: porenetwork/griddata.hh:295
void copyDgfData()
Definition: porenetwork/griddata.hh:193
int parameterIndex(const std::string ¶mName) const
Return the index for a given parameter name.
Definition: porenetwork/griddata.hh:222
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:283
const std::string & paramGroup() const
Return the parameter group.
Definition: porenetwork/griddata.hh:248
GridData(std::shared_ptr< Grid > grid, const std::string ¶mGroup)
constructor for non-dgf grid data
Definition: porenetwork/griddata.hh:72
Scalar getParameter(const Vertex &vertex, const std::string ¶m) const
Returns the value of an vertex parameter.
Definition: porenetwork/griddata.hh:276
void assignParameters()
Assign parameters for generically created grids.
Definition: porenetwork/griddata.hh:161
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:86
bool gridHasVertexParameter(const std::string ¶m) const
Return if a given vertex parameter is provided by the grid.
Definition: porenetwork/griddata.hh:262
Helper class to assign parameters to a generated grid.
Definition: parametersforgeneratedgrid.hh:41
constexpr Line line
Definition: couplingmanager1d3d_line.hh:31
Definition: discretization/porenetwork/fvelementgeometry.hh:24
Helper class to assign parameters to a generated grid.
std::uint_least8_t SmallLocalIndex
Definition: indextraits.hh:29
This file contains functions related to calculate pore-throat properties.