12#ifndef DUMUX_IO_LOADSOLUTION_HH
13#define DUMUX_IO_LOADSOLUTION_HH
18#include <unordered_set>
19#include <unordered_map>
23#include <dune/common/exceptions.hh>
24#include <dune/common/indices.hh>
25#include <dune/grid/common/partitionset.hh>
42template <
class Container,
class EntityMapper,
int codim>
44:
public Dune::CommDataHandleIF< LoadSolutionDataHandle<Container, EntityMapper, codim>,
45 std::decay_t<decltype(std::declval<Container>()[0])> >
47 using FieldType = std::decay_t<decltype(std::declval<Container>()[0])>;
50 const EntityMapper& mapper)
52 , container_(container)
56 {
return cd == codim; }
62 template<
class EntityType>
63 std::size_t
size (
const EntityType &e)
const
66 template<
class MessageBufferImp,
class EntityType>
67 void gather(MessageBufferImp& buff,
const EntityType& e)
const
69 const auto vIdx = mapper_.index(e);
70 buff.write(container_[vIdx]);
73 template<
class MessageBufferImp,
class EntityType>
74 void scatter(MessageBufferImp& buff,
const EntityType& e, std::size_t n)
76 const auto vIdx = mapper_.index(e);
79 container_[vIdx] = tmp;
84 Container& container_;
91template <
class SolutionVector,
class PvNameFunc,
class Gr
idGeometry>
93 const std::string fileName,
94 PvNameFunc&& targetPvNameFunc,
95 const GridGeometry& gridGeometry,
101 using PrimaryVariables =
typename SolutionVector::block_type;
102 using Scalar =
typename PrimaryVariables::field_type;
103 constexpr auto dim = GridGeometry::GridView::dimension;
104 const std::size_t targetSolutionSize = PrimaryVariables::dimension;
106 std::size_t matchingLoadedArrays = 0;
107 for (std::size_t i = 0; i < targetSolutionSize; i++)
108 if (vtu.
hasData(targetPvNameFunc(i,0), dataType))
109 matchingLoadedArrays++;
111 if (matchingLoadedArrays < targetSolutionSize)
112 std::cout <<
"The loaded solution does not provide a data array for each of the primary variables. \n"
113 <<
"The target solution has "<< targetSolutionSize <<
" entries, "
114 <<
"whereas the loaded solution provides only " << matchingLoadedArrays <<
" data array(s). \n"
115 <<
"Make sure that the model concepts are compatible, "
116 <<
"and be sure to provide initial conditions for the missing primary variables. \n";
118 for (std::size_t targetPvIdx = 0; targetPvIdx < targetSolutionSize; ++targetPvIdx)
120 std::vector<Scalar> vec;
121 const auto targetPvName = targetPvNameFunc(targetPvIdx, 0);
123 if (vtu.
hasData(targetPvName, dataType))
124 vec = vtu.
readData<std::vector<Scalar>>(targetPvName, dataType);
127 std::cout <<
"The loaded solution does not have a field named \"" << targetPvName <<
"\". "
128 <<
"Make sure this field is filled using the initial method in the problem definition. \n";
135 for (
const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
137 const auto eIdx = gridGeometry.elementMapper().index(element);
138 sol[eIdx][targetPvIdx] = vec[i++];
144 if (sol.size() != vec.size())
145 DUNE_THROW(Dune::InvalidStateException,
"Solution size (" << sol.size() <<
") does not match input size (" << vec.size() <<
")!");
147 for (std::size_t i = 0; i < sol.size(); ++i)
148 sol[i][targetPvIdx] = vec[i];
153 std::vector<bool> visited(gridGeometry.gridView().size(dim),
false);
154 for (
const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
156 for (
int vIdxLocal = 0; vIdxLocal < element.subEntities(dim); ++vIdxLocal)
158 const auto vIdxGlobal = gridGeometry.vertexMapper().subIndex(element, vIdxLocal, dim);
159 if (!visited[vIdxGlobal])
161 sol[vIdxGlobal][targetPvIdx] = vec[i++];
162 visited[vIdxGlobal] =
true;
174template <
class SolutionVector,
class PvNameFunc,
class Gr
idGeometry>
176 const std::string fileName,
177 PvNameFunc&& targetPvNameFunc,
178 const GridGeometry& gridGeometry,
180->
typename std::enable_if_t<
decltype(
isValid(Detail::hasState())(sol[0]))::value, void>
182 VTKReader vtu(fileName);
184 const auto stateAtDof = vtu.readData<std::vector<int>>(
"phase presence", dataType);
187 std::unordered_set<int> states;
188 for (std::size_t i = 0; i < stateAtDof.size(); ++i)
189 states.insert(stateAtDof[i]);
191 using PrimaryVariables =
typename SolutionVector::block_type;
192 using Scalar =
typename PrimaryVariables::field_type;
193 const std::size_t targetSolutionSize = PrimaryVariables::dimension;
195 std::unordered_set<std::string> matchingNames;
196 for (std::size_t i = 0; i < targetSolutionSize; i++)
197 for (
const auto& state : states)
198 if ( vtu.hasData(targetPvNameFunc(i,state), dataType))
199 matchingNames.insert(targetPvNameFunc(i,state));
201 const std::size_t matchingLoadedArrays = matchingNames.size() - (states.size()-1);
203 if (matchingLoadedArrays < targetSolutionSize)
204 std::cout <<
"The loaded solution does not provide a data array for each of the primary variables. \n"
205 <<
"The target solution has "<< targetSolutionSize <<
" entries, "
206 <<
"whereas the loaded solution provides only " << matchingLoadedArrays <<
" data array(s). \n"
207 <<
"Make sure that the model concepts are compatible, "
208 <<
"and be sure to provide initial conditions for the missing primary variables. \n";
210 for (std::size_t targetPvIdx = 0; targetPvIdx < targetSolutionSize; ++targetPvIdx)
212 std::unordered_map<int, std::vector<Scalar>> data;
213 for (
const auto& state : states)
215 const auto targetPvName = targetPvNameFunc(targetPvIdx, state);
217 if (vtu.hasData(targetPvName, dataType))
218 data[state] = vtu.readData<std::vector<Scalar>>(targetPvName, dataType);
221 std::cout <<
"Loaded Solution does not have a field named \"" << targetPvName <<
"\". "
222 <<
"Make sure this field is filled using the initial method in the problem definition. \n";
230 for (
const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
232 const auto eIdx = gridGeometry.elementMapper().index(element);
233 const auto state = stateAtDof[i];
234 sol[eIdx][targetPvIdx] = data[state][i++];
235 sol[eIdx].setState(state);
241 constexpr int dim = GridGeometry::GridView::dimension;
242 std::vector<bool> visited(gridGeometry.gridView().size(dim),
false);
243 for (
const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
245 for (
int vIdxLocal = 0; vIdxLocal <
element.subEntities(dim); ++vIdxLocal)
247 const auto vIdxGlobal = gridGeometry.vertexMapper().subIndex(element, vIdxLocal, dim);
248 if (!visited[vIdxGlobal])
250 const auto state = stateAtDof[i];
251 sol[vIdxGlobal][targetPvIdx] = data[state][i++];
252 sol[vIdxGlobal].setState(state);
253 visited[vIdxGlobal] =
true;
266template<
class IOFields,
class PrimaryVariables,
class ModelTraits =
void,
class Flu
idSystem =
void,
class Sol
idSystem =
void>
268->
typename std::enable_if_t<
decltype(
isValid(
Detail::hasState())(PrimaryVariables(0)))::value, std::function<std::string(int,int)>>
270 return [paramGroup](
int pvIdx,
int state = 0)
272 static auto numStates = (1 << ModelTraits::numFluidPhases()) - 1;
273 const auto paramNameWithState =
"LoadSolution.PriVarNamesState" + std::to_string(state);
276 DUNE_THROW(Dune::NotImplemented,
"please provide LoadSolution.PriVarNamesState1..." << numStates
277 <<
" or remove LoadSolution.PriVarNames to use the model's default primary variable names");
281 const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup, paramNameWithState);
282 return pvName[pvIdx];
285 return IOFields::template primaryVariableName<ModelTraits, FluidSystem, SolidSystem>(pvIdx, state);
294template<
class IOFields,
class PrimaryVariables,
class ModelTraits =
void,
class Flu
idSystem =
void,
class Sol
idSystem =
void>
296->
typename std::enable_if_t<!
decltype(
isValid(Detail::hasState())(PrimaryVariables(0)))::value, std::function<std::string(
int,
int)>>
300 const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup,
"LoadSolution.PriVarNames");
301 return [n = std::move(pvName)](
int pvIdx,
int state = 0){
return n[pvIdx]; };
304 return [](
int pvIdx,
int state = 0){
return IOFields::template primaryVariableName<ModelTraits, FluidSystem, SolidSystem>(pvIdx, state); };
317template <
class SolutionVector,
class PvNameFunc,
class Gr
idGeometry>
319 const std::string& fileName,
320 PvNameFunc&& targetPvNameFunc,
321 const GridGeometry& gridGeometry)
323 const auto extension = fileName.substr(fileName.find_last_of(
".") + 1);
325 VTKReader::DataType::pointData : VTKReader::DataType::cellData;
327 if (extension ==
"vtu" || extension ==
"vtp")
330 dataType = VTKReader::DataType::pointData;
334 else if (extension ==
"pvtu" || extension ==
"pvtp")
337 DUNE_THROW(Dune::NotImplemented,
"reading staggered solution from a parallel vtk file");
342 DUNE_THROW(Dune::NotImplemented,
"loadSolution for file with extension " << extension);
345 if (gridGeometry.gridView().comm().size() > 1)
347 using GridView =
typename GridGeometry::GridView;
348 if (dataType == VTKReader::DataType::cellData)
351 dataHandle(sol, gridGeometry.elementMapper());
353 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, 0>)
354 gridGeometry.gridView().communicate(dataHandle,
355 Dune::InteriorBorder_All_Interface,
356 Dune::ForwardCommunication);
358 DUNE_THROW(Dune::InvalidStateException,
"Cannot call loadSolution on multiple processes for a grid that cannot communicate codim-" << 0 <<
"-entities.");
363 dataHandle(sol, gridGeometry.vertexMapper());
365 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, GridView::dimension>)
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-" << GridView::dimension <<
"-entities.");
a data handle to communicate the solution on ghosts and overlaps when reading from vtk file in parall...
Definition: loadsolution.hh:46
LoadSolutionDataHandle(Container &container, const EntityMapper &mapper)
Definition: loadsolution.hh:49
bool contains(int dim, int cd) const
Definition: loadsolution.hh:55
std::size_t size(const EntityType &e) const
Definition: loadsolution.hh:63
void scatter(MessageBufferImp &buff, const EntityType &e, std::size_t n)
Definition: loadsolution.hh:74
void gather(MessageBufferImp &buff, const EntityType &e) const
Definition: loadsolution.hh:67
bool fixedSize(int dim, int cd) const
returns true if size per entity of given dim and codim is a constant
Definition: loadsolution.hh:59
A vtk file reader using tinyxml2 as xml backend.
Definition: vtkreader.hh:39
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:72
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:95
DataType
The data array types.
Definition: vtkreader.hh:44
dune-grid capabilities compatibility layer
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:165
constexpr auto isValid(const Expression &t)
A function that creates a test functor to do class member introspection at compile time.
Definition: isvalid.hh:81
A helper function for class member function introspection.
The available discretization methods in Dumux.
constexpr Box box
Definition: method.hh:147
constexpr Staggered staggered
Definition: method.hh:149
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Type traits to be used with matrix types.
helper struct detecting if a PrimaryVariables object has a state() function
Definition: state.hh:21
Type traits to be used with vector types.
A vtk file reader using tinyxml2 as xml backend.