version 3.8
loadsolution.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_IO_LOADSOLUTION_HH
13#define DUMUX_IO_LOADSOLUTION_HH
14
15#include <string>
16#include <iostream>
17#include <vector>
18#include <unordered_set>
19#include <unordered_map>
20#include <type_traits>
21#include <functional>
22
23#include <dune/common/exceptions.hh>
24#include <dune/common/indices.hh>
25#include <dune/grid/common/partitionset.hh>
26
34
35namespace Dumux {
36
42template <class Container, class EntityMapper, int codim>
44: public Dune::CommDataHandleIF< LoadSolutionDataHandle<Container, EntityMapper, codim>,
45 std::decay_t<decltype(std::declval<Container>()[0])> >
46{
47 using FieldType = std::decay_t<decltype(std::declval<Container>()[0])>;
48public:
49 LoadSolutionDataHandle(Container& container,
50 const EntityMapper& mapper)
51 : mapper_(mapper)
52 , container_(container)
53 {}
54
55 bool contains(int dim, int cd) const
56 { return cd == codim; }
57
59 bool fixedSize(int dim, int cd) const
60 { return true; }
61
62 template<class EntityType>
63 std::size_t size (const EntityType &e) const
64 { return 1; }
65
66 template<class MessageBufferImp, class EntityType>
67 void gather(MessageBufferImp& buff, const EntityType& e) const
68 {
69 const auto vIdx = mapper_.index(e);
70 buff.write(container_[vIdx]);
71 }
72
73 template<class MessageBufferImp, class EntityType>
74 void scatter(MessageBufferImp& buff, const EntityType& e, std::size_t n)
75 {
76 const auto vIdx = mapper_.index(e);
77 FieldType tmp;
78 buff.read(tmp);
79 container_[vIdx] = tmp;
80 }
81
82private:
83 EntityMapper mapper_;
84 Container& container_;
85};
86
91template <class SolutionVector, class PvNameFunc, class GridGeometry>
92auto loadSolutionFromVtkFile(SolutionVector& sol,
93 const std::string fileName,
94 PvNameFunc&& targetPvNameFunc,
95 const GridGeometry& gridGeometry,
96 const VTKReader::DataType& dataType)
97-> typename std::enable_if_t<!decltype(isValid(Detail::hasState())(sol[0]))::value, void>
98{
99 VTKReader vtu(fileName);
100
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;
105
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++;
110
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";
117
118 for (std::size_t targetPvIdx = 0; targetPvIdx < targetSolutionSize; ++targetPvIdx)
119 {
120 std::vector<Scalar> vec;
121 const auto targetPvName = targetPvNameFunc(targetPvIdx, 0);
122
123 if (vtu.hasData(targetPvName, dataType))
124 vec = vtu.readData<std::vector<Scalar>>(targetPvName, dataType);
125 else
126 {
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";
129 continue;
130 }
131
132 if (dataType == VTKReader::DataType::cellData)
133 {
134 std::size_t i = 0;
135 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
136 {
137 const auto eIdx = gridGeometry.elementMapper().index(element);
138 sol[eIdx][targetPvIdx] = vec[i++];
139 }
140 }
141 // for staggered face data (which is written out as VTK point data) we just read in the vector
142 else if (dataType == VTKReader::DataType::pointData && GridGeometry::discMethod == DiscretizationMethods::staggered)
143 {
144 if (sol.size() != vec.size())
145 DUNE_THROW(Dune::InvalidStateException, "Solution size (" << sol.size() << ") does not match input size (" << vec.size() << ")!");
146
147 for (std::size_t i = 0; i < sol.size(); ++i)
148 sol[i][targetPvIdx] = vec[i];
149 }
150 else
151 {
152 std::size_t i = 0;
153 std::vector<bool> visited(gridGeometry.gridView().size(dim), false);
154 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
155 {
156 for (int vIdxLocal = 0; vIdxLocal < element.subEntities(dim); ++vIdxLocal)
157 {
158 const auto vIdxGlobal = gridGeometry.vertexMapper().subIndex(element, vIdxLocal, dim);
159 if (!visited[vIdxGlobal])
160 {
161 sol[vIdxGlobal][targetPvIdx] = vec[i++];
162 visited[vIdxGlobal] = true;
163 }
164 }
165 }
166 }
167 }
168}
169
174template <class SolutionVector, class PvNameFunc, class GridGeometry>
175auto loadSolutionFromVtkFile(SolutionVector& sol,
176 const std::string fileName,
177 PvNameFunc&& targetPvNameFunc,
178 const GridGeometry& gridGeometry,
179 const VTKReader::DataType& dataType)
180-> typename std::enable_if_t<decltype(isValid(Detail::hasState())(sol[0]))::value, void>
181{
182 VTKReader vtu(fileName);
183 // get states at each dof location
184 const auto stateAtDof = vtu.readData<std::vector<int>>("phase presence", dataType);
185
186 // determine all states that are present
187 std::unordered_set<int> states;
188 for (std::size_t i = 0; i < stateAtDof.size(); ++i)
189 states.insert(stateAtDof[i]);
190
191 using PrimaryVariables = typename SolutionVector::block_type;
192 using Scalar = typename PrimaryVariables::field_type;
193 const std::size_t targetSolutionSize = PrimaryVariables::dimension;
194
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));
200
201 const std::size_t matchingLoadedArrays = matchingNames.size() - (states.size()-1);
202
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";
209
210 for (std::size_t targetPvIdx = 0; targetPvIdx < targetSolutionSize; ++targetPvIdx)
211 {
212 std::unordered_map<int, std::vector<Scalar>> data;
213 for (const auto& state : states)
214 {
215 const auto targetPvName = targetPvNameFunc(targetPvIdx, state);
216
217 if (vtu.hasData(targetPvName, dataType))
218 data[state] = vtu.readData<std::vector<Scalar>>(targetPvName, dataType);
219 else
220 {
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";
223 continue;
224 }
225 }
226
227 if (dataType == VTKReader::DataType::cellData)
228 {
229 std::size_t i = 0;
230 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
231 {
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);
236 }
237 }
238 else
239 {
240 std::size_t i = 0;
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))
244 {
245 for (int vIdxLocal = 0; vIdxLocal < element.subEntities(dim); ++vIdxLocal)
246 {
247 const auto vIdxGlobal = gridGeometry.vertexMapper().subIndex(element, vIdxLocal, dim);
248 if (!visited[vIdxGlobal])
249 {
250 const auto state = stateAtDof[i];
251 sol[vIdxGlobal][targetPvIdx] = data[state][i++];
252 sol[vIdxGlobal].setState(state);
253 visited[vIdxGlobal] = true;
254 }
255 }
256 }
257 }
258 }
259}
260
266template<class IOFields, class PrimaryVariables, class ModelTraits = void, class FluidSystem = void, class SolidSystem = void>
267auto createPVNameFunction(const std::string& paramGroup = "")
268-> typename std::enable_if_t<decltype(isValid(Detail::hasState())(PrimaryVariables(0)))::value, std::function<std::string(int,int)>>
269{
270 return [paramGroup](int pvIdx, int state = 0)
271 {
272 static auto numStates = (1 << ModelTraits::numFluidPhases()) - 1;
273 const auto paramNameWithState = "LoadSolution.PriVarNamesState" + std::to_string(state);
274
275 if (hasParamInGroup(paramGroup, "LoadSolution.PriVarNames") && !hasParamInGroup(paramGroup, paramNameWithState))
276 DUNE_THROW(Dune::NotImplemented, "please provide LoadSolution.PriVarNamesState1..." << numStates
277 << " or remove LoadSolution.PriVarNames to use the model's default primary variable names");
278
279 else if (hasParamInGroup(paramGroup, paramNameWithState))
280 {
281 const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup, paramNameWithState);
282 return pvName[pvIdx];
283 }
284 else
285 return IOFields::template primaryVariableName<ModelTraits, FluidSystem, SolidSystem>(pvIdx, state);
286 };
287}
288
294template<class IOFields, class PrimaryVariables, class ModelTraits = void, class FluidSystem = void, class SolidSystem = void>
295auto createPVNameFunction(const std::string& paramGroup = "")
296-> typename std::enable_if_t<!decltype(isValid(Detail::hasState())(PrimaryVariables(0)))::value, std::function<std::string(int,int)>>
297{
298 if (hasParamInGroup(paramGroup, "LoadSolution.PriVarNames"))
299 {
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]; };
302 }
303 else
304 return [](int pvIdx, int state = 0){ return IOFields::template primaryVariableName<ModelTraits, FluidSystem, SolidSystem>(pvIdx, state); };
305}
306
317template <class SolutionVector, class PvNameFunc, class GridGeometry>
318void loadSolution(SolutionVector& sol,
319 const std::string& fileName,
320 PvNameFunc&& targetPvNameFunc,
321 const GridGeometry& gridGeometry)
322{
323 const auto extension = fileName.substr(fileName.find_last_of(".") + 1);
324 auto dataType = GridGeometry::discMethod == DiscretizationMethods::box ?
325 VTKReader::DataType::pointData : VTKReader::DataType::cellData;
326
327 if (extension == "vtu" || extension == "vtp")
328 {
329 if (GridGeometry::discMethod == DiscretizationMethods::staggered && extension == "vtp")
330 dataType = VTKReader::DataType::pointData;
331
332 loadSolutionFromVtkFile(sol, fileName, targetPvNameFunc, gridGeometry, dataType);
333 }
334 else if (extension == "pvtu" || extension == "pvtp")
335 {
336 if (GridGeometry::discMethod == DiscretizationMethods::staggered)
337 DUNE_THROW(Dune::NotImplemented, "reading staggered solution from a parallel vtk file");
338
339 loadSolutionFromVtkFile(sol, fileName, targetPvNameFunc, gridGeometry, dataType);
340 }
341 else
342 DUNE_THROW(Dune::NotImplemented, "loadSolution for file with extension " << extension);
343
344 // communicate solution on ghost and overlap dofs
345 if (gridGeometry.gridView().comm().size() > 1)
346 {
347 using GridView = typename GridGeometry::GridView;
348 if (dataType == VTKReader::DataType::cellData)
349 {
351 dataHandle(sol, gridGeometry.elementMapper());
352
353 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, 0>)
354 gridGeometry.gridView().communicate(dataHandle,
355 Dune::InteriorBorder_All_Interface,
356 Dune::ForwardCommunication);
357 else
358 DUNE_THROW(Dune::InvalidStateException, "Cannot call loadSolution on multiple processes for a grid that cannot communicate codim-" << 0 << "-entities.");
359 }
360 else
361 {
363 dataHandle(sol, gridGeometry.vertexMapper());
364
365 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, GridView::dimension>)
366 gridGeometry.gridView().communicate(dataHandle,
367 Dune::InteriorBorder_All_Interface,
368 Dune::ForwardCommunication);
369 else
370 DUNE_THROW(Dune::InvalidStateException, "Cannot call loadSolution on multiple processes for a grid that cannot communicate codim-" << GridView::dimension << "-entities.");
371 }
372 }
373}
374} // end namespace Dumux
375
376#endif
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
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:92
void loadSolution(SolutionVector &sol, const std::string &fileName, PvNameFunc &&targetPvNameFunc, const GridGeometry &gridGeometry)
load a solution vector from file
Definition: loadsolution.hh:318
auto createPVNameFunction(const std::string &paramGroup="") -> 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:267
bool hasParamInGroup(const std::string &paramGroup, const std::string &param)
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
Definition: adapt.hh:17
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.