3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_IO_LOADSOLUTION_HH
25#define DUMUX_IO_LOADSOLUTION_HH
26
27#include <string>
28#include <iostream>
29#include <vector>
30#include <unordered_set>
31#include <unordered_map>
32#include <type_traits>
33#include <functional>
34
35#include <dune/common/exceptions.hh>
36#include <dune/common/indices.hh>
37#include <dune/grid/common/partitionset.hh>
38
46
47namespace Dumux {
48
54template <class Container, class EntityMapper, int codim>
56: public Dune::CommDataHandleIF< LoadSolutionDataHandle<Container, EntityMapper, codim>,
57 std::decay_t<decltype(std::declval<Container>()[0])> >
58{
59 using FieldType = std::decay_t<decltype(std::declval<Container>()[0])>;
60public:
61 LoadSolutionDataHandle(Container& container,
62 const EntityMapper& mapper)
63 : mapper_(mapper)
64 , container_(container)
65 {}
66
67 bool contains(int dim, int cd) const
68 { return cd == codim; }
69
71 bool fixedSize(int dim, int cd) const
72 { return true; }
73
74 template<class EntityType>
75 std::size_t size (const EntityType &e) const
76 { return 1; }
77
78 template<class MessageBufferImp, class EntityType>
79 void gather(MessageBufferImp& buff, const EntityType& e) const
80 {
81 const auto vIdx = mapper_.index(e);
82 buff.write(container_[vIdx]);
83 }
84
85 template<class MessageBufferImp, class EntityType>
86 void scatter(MessageBufferImp& buff, const EntityType& e, std::size_t n)
87 {
88 const auto vIdx = mapper_.index(e);
89 FieldType tmp;
90 buff.read(tmp);
91 container_[vIdx] = tmp;
92 }
93
94private:
95 EntityMapper mapper_;
96 Container& container_;
97};
98
103template <class SolutionVector, class PvNameFunc, class GridGeometry>
104auto loadSolutionFromVtkFile(SolutionVector& sol,
105 const std::string fileName,
106 PvNameFunc&& targetPvNameFunc,
107 const GridGeometry& gridGeometry,
108 const VTKReader::DataType& dataType)
109-> typename std::enable_if_t<!decltype(isValid(Detail::hasState())(sol[0]))::value, void>
110{
111 VTKReader vtu(fileName);
112
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;
117
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++;
122
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";
129
130 for (std::size_t targetPvIdx = 0; targetPvIdx < targetSolutionSize; ++targetPvIdx)
131 {
132 std::vector<Scalar> vec;
133 const auto targetPvName = targetPvNameFunc(targetPvIdx, 0);
134
135 if (vtu.hasData(targetPvName, dataType))
136 vec = vtu.readData<std::vector<Scalar>>(targetPvName, dataType);
137 else
138 {
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";
141 continue;
142 }
143
144 if (dataType == VTKReader::DataType::cellData)
145 {
146 std::size_t i = 0;
147 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
148 {
149 const auto eIdx = gridGeometry.elementMapper().index(element);
150 sol[eIdx][targetPvIdx] = vec[i++];
151 }
152 }
153 // for staggered face data (which is written out as VTK point data) we just read in the vector
154 else if (dataType == VTKReader::DataType::pointData && GridGeometry::discMethod == DiscretizationMethods::staggered)
155 {
156 if (sol.size() != vec.size())
157 DUNE_THROW(Dune::InvalidStateException, "Solution size (" << sol.size() << ") does not match input size (" << vec.size() << ")!");
158
159 for (std::size_t i = 0; i < sol.size(); ++i)
160 sol[i][targetPvIdx] = vec[i];
161 }
162 else
163 {
164 std::size_t i = 0;
165 std::vector<bool> visited(gridGeometry.gridView().size(dim), false);
166 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
167 {
168 for (int vIdxLocal = 0; vIdxLocal < element.subEntities(dim); ++vIdxLocal)
169 {
170 const auto vIdxGlobal = gridGeometry.vertexMapper().subIndex(element, vIdxLocal, dim);
171 if (!visited[vIdxGlobal])
172 {
173 sol[vIdxGlobal][targetPvIdx] = vec[i++];
174 visited[vIdxGlobal] = true;
175 }
176 }
177 }
178 }
179 }
180}
181
186template <class SolutionVector, class PvNameFunc, class GridGeometry>
187auto loadSolutionFromVtkFile(SolutionVector& sol,
188 const std::string fileName,
189 PvNameFunc&& targetPvNameFunc,
190 const GridGeometry& gridGeometry,
191 const VTKReader::DataType& dataType)
192-> typename std::enable_if_t<decltype(isValid(Detail::hasState())(sol[0]))::value, void>
193{
194 VTKReader vtu(fileName);
195 // get states at each dof location
196 const auto stateAtDof = vtu.readData<std::vector<int>>("phase presence", dataType);
197
198 // determine all states that are present
199 std::unordered_set<int> states;
200 for (std::size_t i = 0; i < stateAtDof.size(); ++i)
201 states.insert(stateAtDof[i]);
202
203 using PrimaryVariables = typename SolutionVector::block_type;
204 using Scalar = typename PrimaryVariables::field_type;
205 const std::size_t targetSolutionSize = PrimaryVariables::dimension;
206
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));
212
213 const std::size_t matchingLoadedArrays = matchingNames.size() - (states.size()-1);
214
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";
221
222 for (std::size_t targetPvIdx = 0; targetPvIdx < targetSolutionSize; ++targetPvIdx)
223 {
224 std::unordered_map<int, std::vector<Scalar>> data;
225 for (const auto& state : states)
226 {
227 const auto targetPvName = targetPvNameFunc(targetPvIdx, state);
228
229 if (vtu.hasData(targetPvName, dataType))
230 data[state] = vtu.readData<std::vector<Scalar>>(targetPvName, dataType);
231 else
232 {
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";
235 continue;
236 }
237 }
238
239 if (dataType == VTKReader::DataType::cellData)
240 {
241 std::size_t i = 0;
242 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
243 {
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);
248 }
249 }
250 else
251 {
252 std::size_t i = 0;
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))
256 {
257 for (int vIdxLocal = 0; vIdxLocal < element.subEntities(dim); ++vIdxLocal)
258 {
259 const auto vIdxGlobal = gridGeometry.vertexMapper().subIndex(element, vIdxLocal, dim);
260 if (!visited[vIdxGlobal])
261 {
262 const auto state = stateAtDof[i];
263 sol[vIdxGlobal][targetPvIdx] = data[state][i++];
264 sol[vIdxGlobal].setState(state);
265 visited[vIdxGlobal] = true;
266 }
267 }
268 }
269 }
270 }
271}
272
278template<class IOFields, class PrimaryVariables, class ModelTraits = void, class FluidSystem = void, class SolidSystem = void>
279auto createPVNameFunction(const std::string& paramGroup = "")
280-> typename std::enable_if_t<decltype(isValid(Detail::hasState())(PrimaryVariables(0)))::value, std::function<std::string(int,int)>>
281{
282 return [paramGroup](int pvIdx, int state = 0)
283 {
284 static auto numStates = (1 << ModelTraits::numFluidPhases()) - 1;
285 const auto paramNameWithState = "LoadSolution.PriVarNamesState" + std::to_string(state);
286
287 if (hasParamInGroup(paramGroup, "LoadSolution.PriVarNames") && !hasParamInGroup(paramGroup, paramNameWithState))
288 DUNE_THROW(Dune::NotImplemented, "please provide LoadSolution.PriVarNamesState1..." << numStates
289 << " or remove LoadSolution.PriVarNames to use the model's default primary variable names");
290
291 else if (hasParamInGroup(paramGroup, paramNameWithState))
292 {
293 const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup, paramNameWithState);
294 return pvName[pvIdx];
295 }
296 else
297 return IOFields::template primaryVariableName<ModelTraits, FluidSystem, SolidSystem>(pvIdx, state);
298 };
299}
300
306template<class IOFields, class PrimaryVariables, class ModelTraits = void, class FluidSystem = void, class SolidSystem = void>
307auto createPVNameFunction(const std::string& paramGroup = "")
308-> typename std::enable_if_t<!decltype(isValid(Detail::hasState())(PrimaryVariables(0)))::value, std::function<std::string(int,int)>>
309{
310 if (hasParamInGroup(paramGroup, "LoadSolution.PriVarNames"))
311 {
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]; };
314 }
315 else
316 return [](int pvIdx, int state = 0){ return IOFields::template primaryVariableName<ModelTraits, FluidSystem, SolidSystem>(pvIdx, state); };
317}
318
329template <class SolutionVector, class PvNameFunc, class GridGeometry>
330void loadSolution(SolutionVector& sol,
331 const std::string& fileName,
332 PvNameFunc&& targetPvNameFunc,
333 const GridGeometry& gridGeometry)
334{
335 const auto extension = fileName.substr(fileName.find_last_of(".") + 1);
336 auto dataType = GridGeometry::discMethod == DiscretizationMethods::box ?
337 VTKReader::DataType::pointData : VTKReader::DataType::cellData;
338
339 if (extension == "vtu" || extension == "vtp")
340 {
341 if (GridGeometry::discMethod == DiscretizationMethods::staggered && extension == "vtp")
342 dataType = VTKReader::DataType::pointData;
343
344 loadSolutionFromVtkFile(sol, fileName, targetPvNameFunc, gridGeometry, dataType);
345 }
346 else if (extension == "pvtu" || extension == "pvtp")
347 {
348 if (GridGeometry::discMethod == DiscretizationMethods::staggered)
349 DUNE_THROW(Dune::NotImplemented, "reading staggered solution from a parallel vtk file");
350
351 loadSolutionFromVtkFile(sol, fileName, targetPvNameFunc, gridGeometry, dataType);
352 }
353 else
354 DUNE_THROW(Dune::NotImplemented, "loadSolution for file with extension " << extension);
355
356 // communicate solution on ghost and overlap dofs
357 if (gridGeometry.gridView().comm().size() > 1)
358 {
359 using GridView = typename GridGeometry::GridView;
360 if (dataType == VTKReader::DataType::cellData)
361 {
363 dataHandle(sol, gridGeometry.elementMapper());
364
365 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, 0>)
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-" << 0 << "-entities.");
371 }
372 else
373 {
375 dataHandle(sol, gridGeometry.vertexMapper());
376
377 if constexpr (Detail::canCommunicate<typename GridView::Traits::Grid, GridView::dimension>)
378 gridGeometry.gridView().communicate(dataHandle,
379 Dune::InteriorBorder_All_Interface,
380 Dune::ForwardCommunication);
381 else
382 DUNE_THROW(Dune::InvalidStateException, "Cannot call loadSolution on multiple processes for a grid that cannot communicate codim-" << GridView::dimension << "-entities.");
383 }
384 }
385}
386} // end namespace Dumux
387
388#endif
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 &paramGroup, const std::string &param)
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 &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:279
Definition: adapt.hh:29
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