version 3.11-dev
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-FileCopyrightText: 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
142 // for staggered face data (which is written out as VTK point data) we just read in the vector
143 else if (dataType == VTKReader::DataType::pointData && GridGeometry::discMethod == DiscretizationMethods::staggered)
144 {
145 if (sol.size() != vec.size())
146 DUNE_THROW(Dune::InvalidStateException, "Solution size (" << sol.size() << ") does not match input size (" << vec.size() << ")!");
147
148 for (std::size_t i = 0; i < sol.size(); ++i)
149 sol[i][targetPvIdx] = vec[i];
150 }
151
152 // structured grid vertex data
153 else if (const auto extension = fileName.substr(fileName.find_last_of(".") + 1); extension == "vti" || extension == "pvti")
154 {
155 std::size_t i = 0;
156 for (const auto& vertex : vertices(gridGeometry.gridView(), Dune::Partitions::interior + Dune::Partitions::border))
157 {
158 const auto vIdx = gridGeometry.vertexMapper().index(vertex);
159 sol[vIdx][targetPvIdx] = vec[i++];
160 }
161 }
162
163 // unstructured grid vertex data
164 else
165 {
166 std::size_t i = 0;
167 std::vector<bool> visited(gridGeometry.gridView().size(dim), false);
168 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
169 {
170 for (int vIdxLocal = 0; vIdxLocal < element.subEntities(dim); ++vIdxLocal)
171 {
172 const auto vIdxGlobal = gridGeometry.vertexMapper().subIndex(element, vIdxLocal, dim);
173 if (!visited[vIdxGlobal])
174 {
175 sol[vIdxGlobal][targetPvIdx] = vec[i++];
176 visited[vIdxGlobal] = true;
177 }
178 }
179 }
180 }
181 }
182}
183
188template <class SolutionVector, class PvNameFunc, class GridGeometry>
189auto loadSolutionFromVtkFile(SolutionVector& sol,
190 const std::string fileName,
191 PvNameFunc&& targetPvNameFunc,
192 const GridGeometry& gridGeometry,
193 const VTKReader::DataType& dataType)
194-> typename std::enable_if_t<decltype(isValid(Detail::hasState())(sol[0]))::value, void>
195{
196 VTKReader vtu(fileName);
197 // get states at each dof location
198 const auto stateAtDof = vtu.readData<std::vector<int>>("phase presence", dataType);
199
200 // determine all states that are present
201 std::unordered_set<int> states;
202 for (std::size_t i = 0; i < stateAtDof.size(); ++i)
203 states.insert(stateAtDof[i]);
204
205 using PrimaryVariables = typename SolutionVector::block_type;
206 using Scalar = typename PrimaryVariables::field_type;
207 const std::size_t targetSolutionSize = PrimaryVariables::dimension;
208
209 std::unordered_set<std::string> matchingNames;
210 for (std::size_t i = 0; i < targetSolutionSize; i++)
211 for (const auto& state : states)
212 if ( vtu.hasData(targetPvNameFunc(i,state), dataType))
213 matchingNames.insert(targetPvNameFunc(i,state));
214
215 const std::size_t matchingLoadedArrays = matchingNames.size() - (states.size()-1);
216
217 if (matchingLoadedArrays < targetSolutionSize)
218 std::cout << "The loaded solution does not provide a data array for each of the primary variables. \n"
219 << "The target solution has "<< targetSolutionSize << " entries, "
220 << "whereas the loaded solution provides only " << matchingLoadedArrays << " data array(s). \n"
221 << "Make sure that the model concepts are compatible, "
222 << "and be sure to provide initial conditions for the missing primary variables. \n";
223
224 for (std::size_t targetPvIdx = 0; targetPvIdx < targetSolutionSize; ++targetPvIdx)
225 {
226 std::unordered_map<int, std::vector<Scalar>> data;
227 for (const auto& state : states)
228 {
229 const auto targetPvName = targetPvNameFunc(targetPvIdx, state);
230
231 if (vtu.hasData(targetPvName, dataType))
232 data[state] = vtu.readData<std::vector<Scalar>>(targetPvName, dataType);
233 else
234 {
235 std::cout << "Loaded Solution does not have a field named \"" << targetPvName << "\". "
236 << "Make sure this field is filled using the initial method in the problem definition. \n";
237 continue;
238 }
239 }
240
241 if (dataType == VTKReader::DataType::cellData)
242 {
243 std::size_t i = 0;
244 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
245 {
246 const auto eIdx = gridGeometry.elementMapper().index(element);
247 const auto state = stateAtDof[i];
248 sol[eIdx][targetPvIdx] = data[state][i++];
249 sol[eIdx].setState(state);
250 }
251 }
252 else
253 {
254 std::size_t i = 0;
255 constexpr int dim = GridGeometry::GridView::dimension;
256 std::vector<bool> visited(gridGeometry.gridView().size(dim), false);
257 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
258 {
259 for (int vIdxLocal = 0; vIdxLocal < element.subEntities(dim); ++vIdxLocal)
260 {
261 const auto vIdxGlobal = gridGeometry.vertexMapper().subIndex(element, vIdxLocal, dim);
262 if (!visited[vIdxGlobal])
263 {
264 const auto state = stateAtDof[i];
265 sol[vIdxGlobal][targetPvIdx] = data[state][i++];
266 sol[vIdxGlobal].setState(state);
267 visited[vIdxGlobal] = true;
268 }
269 }
270 }
271 }
272 }
273}
274
280template<class IOFields, class PrimaryVariables, class ModelTraits = void, class FluidSystem = void, class SolidSystem = void>
281auto createPVNameFunction(const std::string& paramGroup = "")
282-> typename std::enable_if_t<decltype(isValid(Detail::hasState())(PrimaryVariables(0)))::value, std::function<std::string(int,int)>>
283{
284 return [paramGroup](int pvIdx, int state = 0)
285 {
286 static auto numStates = (1 << ModelTraits::numFluidPhases()) - 1;
287 const auto paramNameWithState = "LoadSolution.PriVarNamesState" + std::to_string(state);
288
289 if (hasParamInGroup(paramGroup, "LoadSolution.PriVarNames") && !hasParamInGroup(paramGroup, paramNameWithState))
290 DUNE_THROW(Dune::NotImplemented, "please provide LoadSolution.PriVarNamesState1..." << numStates
291 << " or remove LoadSolution.PriVarNames to use the model's default primary variable names");
292
293 else if (hasParamInGroup(paramGroup, paramNameWithState))
294 {
295 const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup, paramNameWithState);
296 return pvName[pvIdx];
297 }
298 else
299 return IOFields::template primaryVariableName<ModelTraits, FluidSystem, SolidSystem>(pvIdx, state);
300 };
301}
302
308template<class IOFields, class PrimaryVariables, class ModelTraits = void, class FluidSystem = void, class SolidSystem = void>
309auto createPVNameFunction(const std::string& paramGroup = "")
310-> typename std::enable_if_t<!decltype(isValid(Detail::hasState())(PrimaryVariables(0)))::value, std::function<std::string(int,int)>>
311{
312 if (hasParamInGroup(paramGroup, "LoadSolution.PriVarNames"))
313 {
314 const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup, "LoadSolution.PriVarNames");
315 return [n = std::move(pvName)](int pvIdx, int state = 0){ return n[pvIdx]; };
316 }
317 else
318 return [](int pvIdx, int state = 0){ return IOFields::template primaryVariableName<ModelTraits, FluidSystem, SolidSystem>(pvIdx, state); };
319}
320
331template <class SolutionVector, class PvNameFunc, class GridGeometry>
332void loadSolution(SolutionVector& sol,
333 const std::string& fileName,
334 PvNameFunc&& targetPvNameFunc,
335 const GridGeometry& gridGeometry)
336{
337 const auto extension = fileName.substr(fileName.find_last_of(".") + 1);
338 auto dataType = GridGeometry::discMethod == DiscretizationMethods::box ?
339 VTKReader::DataType::pointData : VTKReader::DataType::cellData;
340
341 if (extension == "vtu" || extension == "vtp" || extension == "vti")
342 {
343 if (GridGeometry::discMethod == DiscretizationMethods::staggered && extension == "vtp")
344 dataType = VTKReader::DataType::pointData;
345
346 loadSolutionFromVtkFile(sol, fileName, targetPvNameFunc, gridGeometry, dataType);
347 }
348 else if (extension == "pvtu" || extension == "pvtp" || extension == "pvti")
349 {
350 if (GridGeometry::discMethod == DiscretizationMethods::staggered)
351 DUNE_THROW(Dune::NotImplemented, "reading staggered solution from a parallel vtk file");
352
353 loadSolutionFromVtkFile(sol, fileName, targetPvNameFunc, gridGeometry, dataType);
354 }
355 else
356 DUNE_THROW(Dune::NotImplemented, "loadSolution for file with extension " << extension);
357
358 // communicate solution on ghost and overlap dofs
359 if (gridGeometry.gridView().comm().size() > 1)
360 {
361 using GridView = typename GridGeometry::GridView;
362 if (dataType == VTKReader::DataType::cellData)
363 {
365 dataHandle(sol, gridGeometry.elementMapper());
366
367 if constexpr (Dune::Capabilities::canCommunicate<typename GridView::Traits::Grid, 0>::v)
368 gridGeometry.gridView().communicate(dataHandle,
369 Dune::InteriorBorder_All_Interface,
370 Dune::ForwardCommunication);
371 else
372 DUNE_THROW(Dune::InvalidStateException, "Cannot call loadSolution on multiple processes for a grid that cannot communicate codim-" << 0 << "-entities.");
373 }
374 else
375 {
377 dataHandle(sol, gridGeometry.vertexMapper());
378
379 if constexpr (Dune::Capabilities::canCommunicate<typename GridView::Traits::Grid, GridView::dimension>::v)
380 gridGeometry.gridView().communicate(dataHandle,
381 Dune::InteriorBorder_All_Interface,
382 Dune::ForwardCommunication);
383 else
384 DUNE_THROW(Dune::InvalidStateException, "Cannot call loadSolution on multiple processes for a grid that cannot communicate codim-" << GridView::dimension << "-entities.");
385 }
386 }
387}
388} // end namespace Dumux
389
390#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:342
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:379
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:402
DataType
The data array types.
Definition: vtkreader.hh:347
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:332
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:281
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.