3.3.0
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
45
46namespace Dumux {
47
53template <class Container, class EntityMapper, int codim>
55: public Dune::CommDataHandleIF< LoadSolutionDataHandle<Container, EntityMapper, codim>,
56 std::decay_t<decltype(std::declval<Container>()[0])> >
57{
58 using FieldType = std::decay_t<decltype(std::declval<Container>()[0])>;
59public:
60 LoadSolutionDataHandle(Container& container,
61 const EntityMapper& mapper)
62 : mapper_(mapper)
63 , container_(container)
64 {}
65
66 bool contains(int dim, int cd) const
67 { return cd == codim; }
68
70 bool fixedSize(int dim, int cd) const
71 { return true; }
72
73 template<class EntityType>
74 std::size_t size (const EntityType &e) const
75 { return 1; }
76
77 template<class MessageBufferImp, class EntityType>
78 void gather(MessageBufferImp& buff, const EntityType& e) const
79 {
80 const auto vIdx = mapper_.index(e);
81 buff.write(container_[vIdx]);
82 }
83
84 template<class MessageBufferImp, class EntityType>
85 void scatter(MessageBufferImp& buff, const EntityType& e, std::size_t n)
86 {
87 const auto vIdx = mapper_.index(e);
88 FieldType tmp;
89 buff.read(tmp);
90 container_[vIdx] = tmp;
91 }
92
93private:
94 EntityMapper mapper_;
95 Container& container_;
96};
97
102template <class SolutionVector, class PvNameFunc, class GridGeometry>
103auto loadSolutionFromVtkFile(SolutionVector& sol,
104 const std::string fileName,
105 PvNameFunc&& targetPvNameFunc,
106 const GridGeometry& gridGeometry,
107 const VTKReader::DataType& dataType)
108-> typename std::enable_if_t<!decltype(isValid(Detail::hasState())(sol[0]))::value, void>
109{
110 VTKReader vtu(fileName);
111
112 using PrimaryVariables = typename SolutionVector::block_type;
113 using Scalar = typename PrimaryVariables::field_type;
114 constexpr auto dim = GridGeometry::GridView::dimension;
115 const std::size_t targetSolutionSize = PrimaryVariables::dimension;
116
117 std::size_t matchingLoadedArrays = 0;
118 for (std::size_t i = 0; i < targetSolutionSize; i++)
119 if (vtu.hasData(targetPvNameFunc(i,0), dataType))
120 matchingLoadedArrays++;
121
122 if (matchingLoadedArrays < targetSolutionSize)
123 std::cout << "The loaded solution does not provide a data array for each of the primary variables. \n"
124 << "The target solution has "<< targetSolutionSize << " entries, "
125 << "whereas the loaded solution provides only " << matchingLoadedArrays << " data array(s). \n"
126 << "Make sure that the model concepts are compatible, "
127 << "and be sure to provide initial conditions for the missing primary variables. \n";
128
129 for (std::size_t targetPvIdx = 0; targetPvIdx < targetSolutionSize; ++targetPvIdx)
130 {
131 std::vector<Scalar> vec;
132 const auto targetPvName = targetPvNameFunc(targetPvIdx, 0);
133
134 if (vtu.hasData(targetPvName, dataType))
135 vec = vtu.readData<std::vector<Scalar>>(targetPvName, dataType);
136 else
137 {
138 std::cout << "The loaded solution does not have a field named \"" << targetPvName << "\". "
139 << "Make sure this field is filled using the initial method in the problem definition. \n";
140 continue;
141 }
142
143 if (dataType == VTKReader::DataType::cellData)
144 {
145 std::size_t i = 0;
146 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
147 {
148 const auto eIdx = gridGeometry.elementMapper().index(element);
149 sol[eIdx][targetPvIdx] = vec[i++];
150 }
151 }
152 // for staggered face data (which is written out as VTK point data) we just read in the vector
153 else if (dataType == VTKReader::DataType::pointData && GridGeometry::discMethod == DiscretizationMethod::staggered)
154 {
155 if (sol.size() != vec.size())
156 DUNE_THROW(Dune::InvalidStateException, "Solution size (" << sol.size() << ") does not match input size (" << vec.size() << ")!");
157
158 for (std::size_t i = 0; i < sol.size(); ++i)
159 sol[i][targetPvIdx] = vec[i];
160 }
161 else
162 {
163 std::size_t i = 0;
164 std::vector<bool> visited(gridGeometry.gridView().size(dim), false);
165 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
166 {
167 for (int vIdxLocal = 0; vIdxLocal < element.subEntities(dim); ++vIdxLocal)
168 {
169 const auto vIdxGlobal = gridGeometry.vertexMapper().subIndex(element, vIdxLocal, dim);
170 if (!visited[vIdxGlobal])
171 {
172 sol[vIdxGlobal][targetPvIdx] = vec[i++];
173 visited[vIdxGlobal] = true;
174 }
175 }
176 }
177 }
178 }
179}
180
185template <class SolutionVector, class PvNameFunc, class GridGeometry>
186auto loadSolutionFromVtkFile(SolutionVector& sol,
187 const std::string fileName,
188 PvNameFunc&& targetPvNameFunc,
189 const GridGeometry& gridGeometry,
190 const VTKReader::DataType& dataType)
191-> typename std::enable_if_t<decltype(isValid(Detail::hasState())(sol[0]))::value, void>
192{
193 VTKReader vtu(fileName);
194 // get states at each dof location
195 const auto stateAtDof = vtu.readData<std::vector<int>>("phase presence", dataType);
196
197 // determine all states that are present
198 std::unordered_set<int> states;
199 for (std::size_t i = 0; i < stateAtDof.size(); ++i)
200 states.insert(stateAtDof[i]);
201
202 using PrimaryVariables = typename SolutionVector::block_type;
203 using Scalar = typename PrimaryVariables::field_type;
204 const std::size_t targetSolutionSize = PrimaryVariables::dimension;
205
206 std::unordered_set<std::string> matchingNames;
207 for (std::size_t i = 0; i < targetSolutionSize; i++)
208 for (const auto& state : states)
209 if ( vtu.hasData(targetPvNameFunc(i,state), dataType))
210 matchingNames.insert(targetPvNameFunc(i,state));
211
212 const std::size_t matchingLoadedArrays = matchingNames.size() - (states.size()-1);
213
214 if (matchingLoadedArrays < targetSolutionSize)
215 std::cout << "The loaded solution does not provide a data array for each of the primary variables. \n"
216 << "The target solution has "<< targetSolutionSize << " entries, "
217 << "whereas the loaded solution provides only " << matchingLoadedArrays << " data array(s). \n"
218 << "Make sure that the model concepts are compatible, "
219 << "and be sure to provide initial conditions for the missing primary variables. \n";
220
221 for (std::size_t targetPvIdx = 0; targetPvIdx < targetSolutionSize; ++targetPvIdx)
222 {
223 std::unordered_map<int, std::vector<Scalar>> data;
224 for (const auto& state : states)
225 {
226 const auto targetPvName = targetPvNameFunc(targetPvIdx, state);
227
228 if (vtu.hasData(targetPvName, dataType))
229 data[state] = vtu.readData<std::vector<Scalar>>(targetPvName, dataType);
230 else
231 {
232 std::cout << "Loaded Solution does not have a field named \"" << targetPvName << "\". "
233 << "Make sure this field is filled using the initial method in the problem definition. \n";
234 continue;
235 }
236 }
237
238 if (dataType == VTKReader::DataType::cellData)
239 {
240 std::size_t i = 0;
241 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
242 {
243 const auto eIdx = gridGeometry.elementMapper().index(element);
244 const auto state = stateAtDof[i];
245 sol[eIdx][targetPvIdx] = data[state][i++];
246 sol[eIdx].setState(state);
247 }
248 }
249 else
250 {
251 std::size_t i = 0;
252 constexpr int dim = GridGeometry::GridView::dimension;
253 std::vector<bool> visited(gridGeometry.gridView().size(dim), false);
254 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
255 {
256 for (int vIdxLocal = 0; vIdxLocal < element.subEntities(dim); ++vIdxLocal)
257 {
258 const auto vIdxGlobal = gridGeometry.vertexMapper().subIndex(element, vIdxLocal, dim);
259 if (!visited[vIdxGlobal])
260 {
261 const auto state = stateAtDof[i];
262 sol[vIdxGlobal][targetPvIdx] = data[state][i++];
263 sol[vIdxGlobal].setState(state);
264 visited[vIdxGlobal] = true;
265 }
266 }
267 }
268 }
269 }
270}
271
277template<class IOFields, class PrimaryVariables, class ModelTraits = void, class FluidSystem = void, class SolidSystem = void>
278auto createPVNameFunction(const std::string& paramGroup = "")
279-> typename std::enable_if_t<decltype(isValid(Detail::hasState())(PrimaryVariables(0)))::value, std::function<std::string(int,int)>>
280{
281 return [paramGroup](int pvIdx, int state = 0)
282 {
283 static auto numStates = (1 << ModelTraits::numFluidPhases()) - 1;
284 const auto paramNameWithState = "LoadSolution.PriVarNamesState" + std::to_string(state);
285
286 if (hasParamInGroup(paramGroup, "LoadSolution.PriVarNames") && !hasParamInGroup(paramGroup, paramNameWithState))
287 DUNE_THROW(Dune::NotImplemented, "please provide LoadSolution.PriVarNamesState1..." << numStates
288 << " or remove LoadSolution.PriVarNames to use the model's default primary variable names");
289
290 else if (hasParamInGroup(paramGroup, paramNameWithState))
291 {
292 const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup, paramNameWithState);
293 return pvName[pvIdx];
294 }
295 else
296 return IOFields::template primaryVariableName<ModelTraits, FluidSystem, SolidSystem>(pvIdx, state);
297 };
298}
299
305template<class IOFields, class PrimaryVariables, class ModelTraits = void, class FluidSystem = void, class SolidSystem = void>
306auto createPVNameFunction(const std::string& paramGroup = "")
307-> typename std::enable_if_t<!decltype(isValid(Detail::hasState())(PrimaryVariables(0)))::value, std::function<std::string(int,int)>>
308{
309 if (hasParamInGroup(paramGroup, "LoadSolution.PriVarNames"))
310 {
311 const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup, "LoadSolution.PriVarNames");
312 return [n = std::move(pvName)](int pvIdx, int state = 0){ return n[pvIdx]; };
313 }
314 else
315 return [](int pvIdx, int state = 0){ return IOFields::template primaryVariableName<ModelTraits, FluidSystem, SolidSystem>(pvIdx, state); };
316}
317
328template <class SolutionVector, class PvNameFunc, class GridGeometry>
329void loadSolution(SolutionVector& sol,
330 const std::string& fileName,
331 PvNameFunc&& targetPvNameFunc,
332 const GridGeometry& gridGeometry)
333{
334 const auto extension = fileName.substr(fileName.find_last_of(".") + 1);
335 auto dataType = GridGeometry::discMethod == DiscretizationMethod::box ?
336 VTKReader::DataType::pointData : VTKReader::DataType::cellData;
337
338 if (extension == "vtu" || extension == "vtp")
339 {
340 if (GridGeometry::discMethod == DiscretizationMethod::staggered && extension == "vtp")
341 dataType = VTKReader::DataType::pointData;
342
343 loadSolutionFromVtkFile(sol, fileName, targetPvNameFunc, gridGeometry, dataType);
344 }
345 else if (extension == "pvtu" || extension == "pvtp")
346 {
347 if (GridGeometry::discMethod == DiscretizationMethod::staggered)
348 DUNE_THROW(Dune::NotImplemented, "reading staggered solution from a parallel vtk file");
349
350 loadSolutionFromVtkFile(sol, fileName, targetPvNameFunc, gridGeometry, dataType);
351 }
352 else
353 DUNE_THROW(Dune::NotImplemented, "loadSolution for file with extension " << extension);
354
355 // communicate solution on ghost and overlap dofs
356 if (gridGeometry.gridView().comm().size() > 1)
357 {
358 using GridView = typename GridGeometry::GridView;
359 if (dataType == VTKReader::DataType::cellData)
360 {
362 dataHandle(sol, gridGeometry.elementMapper());
363 gridGeometry.gridView().communicate(dataHandle,
364 Dune::InteriorBorder_All_Interface,
365 Dune::ForwardCommunication);
366 }
367 else
368 {
370 dataHandle(sol, gridGeometry.vertexMapper());
371 gridGeometry.gridView().communicate(dataHandle,
372 Dune::InteriorBorder_All_Interface,
373 Dune::ForwardCommunication);
374 }
375 }
376}
377} // end namespace Dumux
378
379#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.
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:390
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:103
void loadSolution(SolutionVector &sol, const std::string &fileName, PvNameFunc &&targetPvNameFunc, const GridGeometry &gridGeometry)
load a solution vector from file
Definition: loadsolution.hh:329
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:278
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
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:57
LoadSolutionDataHandle(Container &container, const EntityMapper &mapper)
Definition: loadsolution.hh:60
bool contains(int dim, int cd) const
Definition: loadsolution.hh:66
std::size_t size(const EntityType &e) const
Definition: loadsolution.hh:74
void scatter(MessageBufferImp &buff, const EntityType &e, std::size_t n)
Definition: loadsolution.hh:85
void gather(MessageBufferImp &buff, const EntityType &e) const
Definition: loadsolution.hh:78
bool fixedSize(int dim, int cd) const
returns true if size per entity of given dim and codim is a constant
Definition: loadsolution.hh:70
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