24#ifndef DUMUX_IO_LOADSOLUTION_HH
25#define DUMUX_IO_LOADSOLUTION_HH
30#include <unordered_set>
31#include <unordered_map>
35#include <dune/common/exceptions.hh>
36#include <dune/common/indices.hh>
37#include <dune/grid/common/partitionset.hh>
54template <
class Container,
class EntityMapper,
int codim>
56:
public Dune::CommDataHandleIF< LoadSolutionDataHandle<Container, EntityMapper, codim>,
57 std::decay_t<decltype(std::declval<Container>()[0])> >
59 using FieldType = std::decay_t<decltype(std::declval<Container>()[0])>;
62 const EntityMapper& mapper)
64 , container_(container)
68 {
return cd == codim; }
74 template<
class EntityType>
75 std::size_t
size (
const EntityType &e)
const
78 template<
class MessageBufferImp,
class EntityType>
79 void gather(MessageBufferImp& buff,
const EntityType& e)
const
81 const auto vIdx = mapper_.index(e);
82 buff.write(container_[vIdx]);
85 template<
class MessageBufferImp,
class EntityType>
86 void scatter(MessageBufferImp& buff,
const EntityType& e, std::size_t n)
88 const auto vIdx = mapper_.index(e);
91 container_[vIdx] = tmp;
96 Container& container_;
103template <
class SolutionVector,
class PvNameFunc,
class Gr
idGeometry>
105 const std::string fileName,
106 PvNameFunc&& targetPvNameFunc,
107 const GridGeometry& gridGeometry,
113 using PrimaryVariables =
typename SolutionVector::block_type;
114 using Scalar =
typename PrimaryVariables::field_type;
115 constexpr auto dim = GridGeometry::GridView::dimension;
116 const std::size_t targetSolutionSize = PrimaryVariables::dimension;
118 std::size_t matchingLoadedArrays = 0;
119 for (std::size_t i = 0; i < targetSolutionSize; i++)
120 if (vtu.
hasData(targetPvNameFunc(i,0), dataType))
121 matchingLoadedArrays++;
123 if (matchingLoadedArrays < targetSolutionSize)
124 std::cout <<
"The loaded solution does not provide a data array for each of the primary variables. \n"
125 <<
"The target solution has "<< targetSolutionSize <<
" entries, "
126 <<
"whereas the loaded solution provides only " << matchingLoadedArrays <<
" data array(s). \n"
127 <<
"Make sure that the model concepts are compatible, "
128 <<
"and be sure to provide initial conditions for the missing primary variables. \n";
130 for (std::size_t targetPvIdx = 0; targetPvIdx < targetSolutionSize; ++targetPvIdx)
132 std::vector<Scalar> vec;
133 const auto targetPvName = targetPvNameFunc(targetPvIdx, 0);
135 if (vtu.
hasData(targetPvName, dataType))
136 vec = vtu.
readData<std::vector<Scalar>>(targetPvName, dataType);
139 std::cout <<
"The loaded solution does not have a field named \"" << targetPvName <<
"\". "
140 <<
"Make sure this field is filled using the initial method in the problem definition. \n";
147 for (
const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
149 const auto eIdx = gridGeometry.elementMapper().index(element);
150 sol[eIdx][targetPvIdx] = vec[i++];
156 if (sol.size() != vec.size())
157 DUNE_THROW(Dune::InvalidStateException,
"Solution size (" << sol.size() <<
") does not match input size (" << vec.size() <<
")!");
159 for (std::size_t i = 0; i < sol.size(); ++i)
160 sol[i][targetPvIdx] = vec[i];
165 std::vector<bool> visited(gridGeometry.gridView().size(dim),
false);
166 for (
const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
168 for (
int vIdxLocal = 0; vIdxLocal < element.subEntities(dim); ++vIdxLocal)
170 const auto vIdxGlobal = gridGeometry.vertexMapper().subIndex(element, vIdxLocal, dim);
171 if (!visited[vIdxGlobal])
173 sol[vIdxGlobal][targetPvIdx] = vec[i++];
174 visited[vIdxGlobal] =
true;
186template <
class SolutionVector,
class PvNameFunc,
class Gr
idGeometry>
188 const std::string fileName,
189 PvNameFunc&& targetPvNameFunc,
190 const GridGeometry& gridGeometry,
192->
typename std::enable_if_t<
decltype(
isValid(Detail::hasState())(sol[0]))::value, void>
194 VTKReader vtu(fileName);
196 const auto stateAtDof = vtu.readData<std::vector<int>>(
"phase presence", dataType);
199 std::unordered_set<int> states;
200 for (std::size_t i = 0; i < stateAtDof.size(); ++i)
201 states.insert(stateAtDof[i]);
203 using PrimaryVariables =
typename SolutionVector::block_type;
204 using Scalar =
typename PrimaryVariables::field_type;
205 const std::size_t targetSolutionSize = PrimaryVariables::dimension;
207 std::unordered_set<std::string> matchingNames;
208 for (std::size_t i = 0; i < targetSolutionSize; i++)
209 for (
const auto& state : states)
210 if ( vtu.hasData(targetPvNameFunc(i,state), dataType))
211 matchingNames.insert(targetPvNameFunc(i,state));
213 const std::size_t matchingLoadedArrays = matchingNames.size() - (states.size()-1);
215 if (matchingLoadedArrays < targetSolutionSize)
216 std::cout <<
"The loaded solution does not provide a data array for each of the primary variables. \n"
217 <<
"The target solution has "<< targetSolutionSize <<
" entries, "
218 <<
"whereas the loaded solution provides only " << matchingLoadedArrays <<
" data array(s). \n"
219 <<
"Make sure that the model concepts are compatible, "
220 <<
"and be sure to provide initial conditions for the missing primary variables. \n";
222 for (std::size_t targetPvIdx = 0; targetPvIdx < targetSolutionSize; ++targetPvIdx)
224 std::unordered_map<int, std::vector<Scalar>> data;
225 for (
const auto& state : states)
227 const auto targetPvName = targetPvNameFunc(targetPvIdx, state);
229 if (vtu.hasData(targetPvName, dataType))
230 data[state] = vtu.readData<std::vector<Scalar>>(targetPvName, dataType);
233 std::cout <<
"Loaded Solution does not have a field named \"" << targetPvName <<
"\". "
234 <<
"Make sure this field is filled using the initial method in the problem definition. \n";
242 for (
const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
244 const auto eIdx = gridGeometry.elementMapper().index(element);
245 const auto state = stateAtDof[i];
246 sol[eIdx][targetPvIdx] = data[state][i++];
247 sol[eIdx].setState(state);
253 constexpr int dim = GridGeometry::GridView::dimension;
254 std::vector<bool> visited(gridGeometry.gridView().size(dim),
false);
255 for (
const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
257 for (
int vIdxLocal = 0; vIdxLocal <
element.subEntities(dim); ++vIdxLocal)
259 const auto vIdxGlobal = gridGeometry.vertexMapper().subIndex(element, vIdxLocal, dim);
260 if (!visited[vIdxGlobal])
262 const auto state = stateAtDof[i];
263 sol[vIdxGlobal][targetPvIdx] = data[state][i++];
264 sol[vIdxGlobal].setState(state);
265 visited[vIdxGlobal] =
true;
278template<
class IOFields,
class PrimaryVariables,
class ModelTraits =
void,
class Flu
idSystem =
void,
class Sol
idSystem =
void>
280->
typename std::enable_if_t<
decltype(
isValid(
Detail::hasState())(PrimaryVariables(0)))::value, std::function<std::string(int,int)>>
282 return [paramGroup](
int pvIdx,
int state = 0)
284 static auto numStates = (1 << ModelTraits::numFluidPhases()) - 1;
285 const auto paramNameWithState =
"LoadSolution.PriVarNamesState" + std::to_string(state);
288 DUNE_THROW(Dune::NotImplemented,
"please provide LoadSolution.PriVarNamesState1..." << numStates
289 <<
" or remove LoadSolution.PriVarNames to use the model's default primary variable names");
293 const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup, paramNameWithState);
294 return pvName[pvIdx];
297 return IOFields::template primaryVariableName<ModelTraits, FluidSystem, SolidSystem>(pvIdx, state);
306template<
class IOFields,
class PrimaryVariables,
class ModelTraits =
void,
class Flu
idSystem =
void,
class Sol
idSystem =
void>
308->
typename std::enable_if_t<!
decltype(
isValid(Detail::hasState())(PrimaryVariables(0)))::value, std::function<std::string(
int,
int)>>
312 const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup,
"LoadSolution.PriVarNames");
313 return [n = std::move(pvName)](
int pvIdx,
int state = 0){
return n[pvIdx]; };
316 return [](
int pvIdx,
int state = 0){
return IOFields::template primaryVariableName<ModelTraits, FluidSystem, SolidSystem>(pvIdx, state); };
329template <
class SolutionVector,
class PvNameFunc,
class Gr
idGeometry>
331 const std::string& fileName,
332 PvNameFunc&& targetPvNameFunc,
333 const GridGeometry& gridGeometry)
335 const auto extension = fileName.substr(fileName.find_last_of(
".") + 1);
337 VTKReader::DataType::pointData : VTKReader::DataType::cellData;
339 if (extension ==
"vtu" || extension ==
"vtp")
342 dataType = VTKReader::DataType::pointData;
346 else if (extension ==
"pvtu" || extension ==
"pvtp")
349 DUNE_THROW(Dune::NotImplemented,
"reading staggered solution from a parallel vtk file");
354 DUNE_THROW(Dune::NotImplemented,
"loadSolution for file with extension " << extension);
357 if (gridGeometry.gridView().comm().size() > 1)
359 using GridView =
typename GridGeometry::GridView;
360 if (dataType == VTKReader::DataType::cellData)
363 dataHandle(sol, gridGeometry.elementMapper());
365 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, 0>)
366 gridGeometry.gridView().communicate(dataHandle,
367 Dune::InteriorBorder_All_Interface,
368 Dune::ForwardCommunication);
370 DUNE_THROW(Dune::InvalidStateException,
"Cannot call loadSolution on multiple processes for a grid that cannot communicate codim-" << 0 <<
"-entities.");
375 dataHandle(sol, gridGeometry.vertexMapper());
377 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, GridView::dimension>)
378 gridGeometry.gridView().communicate(dataHandle,
379 Dune::InteriorBorder_All_Interface,
380 Dune::ForwardCommunication);
382 DUNE_THROW(Dune::InvalidStateException,
"Cannot call loadSolution on multiple processes for a grid that cannot communicate codim-" << GridView::dimension <<
"-entities.");
Type traits to be used with matrix types.
Type traits to be used with vector types.
A helper function for class member function introspection.
dune-grid capabilities compatibility layer
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
A vtk file reader using tinyxml2 as xml backend.
bool hasParamInGroup(const std::string ¶mGroup, const std::string ¶m)
Check whether a key exists in the parameter tree with a model group prefix.
Definition: parameters.hh:177
auto loadSolutionFromVtkFile(SolutionVector &sol, const std::string fileName, PvNameFunc &&targetPvNameFunc, const GridGeometry &gridGeometry, const VTKReader::DataType &dataType) -> typename std::enable_if_t<!decltype(isValid(Detail::hasState())(sol[0]))::value, void >
read from a vtk file into a solution vector with primary variables without state
Definition: loadsolution.hh:104
void loadSolution(SolutionVector &sol, const std::string &fileName, PvNameFunc &&targetPvNameFunc, const GridGeometry &gridGeometry)
load a solution vector from file
Definition: loadsolution.hh:330
auto createPVNameFunction(const std::string ¶mGroup="") -> typename std::enable_if_t< decltype(isValid(Detail::hasState())(PrimaryVariables(0)))::value, std::function< std::string(int, int)> >
helper function to determine the primary variable names of a model with privar state
Definition: loadsolution.hh:279
constexpr auto isValid(const Expression &t)
A function that creates a test functor to do class member introspection at compile time.
Definition: isvalid.hh:93
constexpr Box box
Definition: method.hh:139
constexpr Staggered staggered
Definition: method.hh:140
helper struct detecting if a PrimaryVariables object has a state() function
Definition: state.hh:32
a data handle to communicate the solution on ghosts and overlaps when reading from vtk file in parall...
Definition: loadsolution.hh:58
LoadSolutionDataHandle(Container &container, const EntityMapper &mapper)
Definition: loadsolution.hh:61
bool contains(int dim, int cd) const
Definition: loadsolution.hh:67
std::size_t size(const EntityType &e) const
Definition: loadsolution.hh:75
void scatter(MessageBufferImp &buff, const EntityType &e, std::size_t n)
Definition: loadsolution.hh:86
void gather(MessageBufferImp &buff, const EntityType &e) const
Definition: loadsolution.hh:79
bool fixedSize(int dim, int cd) const
returns true if size per entity of given dim and codim is a constant
Definition: loadsolution.hh:71
A vtk file reader using tinyxml2 as xml backend.
Definition: vtkreader.hh:51
bool hasData(const std::string &name, const DataType &type) const
Reviews data from the vtk file to check if there is a data array with a specified name.
Definition: vtkreader.hh:84
Container readData(const std::string &name, const DataType &type) const
read data from the vtk file to a container, e.g. std::vector<double>
Definition: vtkreader.hh:107
DataType
The data array types.
Definition: vtkreader.hh:56