3.2-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/version.hh>
36#include <dune/common/exceptions.hh>
37#include <dune/common/indices.hh>
38#include <dune/grid/common/partitionset.hh>
39
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
70#if DUNE_VERSION_GTE(DUNE_GRID,2,7)
72 bool fixedSize(int dim, int cd) const
73 { return true; }
74#else
76 bool fixedsize(int dim, int cd) const
77 { return true; }
78#endif
79
80 template<class EntityType>
81 size_t size (const EntityType &e) const
82 { return 1; }
83
84 template<class MessageBufferImp, class EntityType>
85 void gather(MessageBufferImp& buff, const EntityType& e) const
86 {
87 const auto vIdx = mapper_.index(e);
88 buff.write(container_[vIdx]);
89 }
90
91 template<class MessageBufferImp, class EntityType>
92 void scatter(MessageBufferImp& buff, const EntityType& e, size_t n)
93 {
94 const auto vIdx = mapper_.index(e);
95 FieldType tmp;
96 buff.read(tmp);
97 container_[vIdx] = tmp;
98 }
99
100private:
101 EntityMapper mapper_;
102 Container& container_;
103};
104
109template <class SolutionVector, class PvNameFunc, class GridGeometry>
110auto loadSolutionFromVtkFile(SolutionVector& sol,
111 const std::string fileName,
112 PvNameFunc&& pvNameFunc,
113 const GridGeometry& gridGeometry,
114 const VTKReader::DataType& dataType)
115-> typename std::enable_if_t<!decltype(isValid(Detail::hasState())(sol[0]))::value, void>
116{
117 VTKReader vtu(fileName);
118 using PrimaryVariables = typename SolutionVector::block_type;
119 using Scalar = typename PrimaryVariables::field_type;
120 constexpr auto dim = GridGeometry::GridView::dimension;
121
122 for (size_t pvIdx = 0; pvIdx < PrimaryVariables::dimension; ++pvIdx)
123 {
124 const auto pvName = pvNameFunc(pvIdx, 0);
125 auto vec = vtu.readData<std::vector<Scalar>>(pvName, dataType);
126
127 if (dataType == VTKReader::DataType::cellData)
128 {
129 std::size_t i = 0;
130 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
131 {
132 const auto eIdx = gridGeometry.elementMapper().index(element);
133 sol[eIdx][pvIdx] = vec[i++];
134 }
135 }
136 // for staggered face data (which is written out as VTK point data) we just read in the vector
137 else if (dataType == VTKReader::DataType::pointData && GridGeometry::discMethod == DiscretizationMethod::staggered)
138 {
139 if (sol.size() != vec.size())
140 DUNE_THROW(Dune::InvalidStateException, "Solution size (" << sol.size() << ") does not match input size (" << vec.size() << ")!");
141
142 for (std::size_t i = 0; i < sol.size(); ++i)
143 sol[i][pvIdx] = vec[i];
144 }
145 else
146 {
147 std::size_t i = 0;
148 std::vector<bool> visited(gridGeometry.gridView().size(dim), false);
149 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
150 {
151 for (int vIdxLocal = 0; vIdxLocal < element.subEntities(dim); ++vIdxLocal)
152 {
153 const auto vIdxGlobal = gridGeometry.vertexMapper().subIndex(element, vIdxLocal, dim);
154 if (!visited[vIdxGlobal])
155 {
156 sol[vIdxGlobal][pvIdx] = vec[i++];
157 visited[vIdxGlobal] = true;
158 }
159 }
160 }
161 }
162 }
163}
164
169template <class SolutionVector, class PvNameFunc, class GridGeometry>
170auto loadSolutionFromVtkFile(SolutionVector& sol,
171 const std::string fileName,
172 PvNameFunc&& pvNameFunc,
173 const GridGeometry& gridGeometry,
174 const VTKReader::DataType& dataType)
175-> typename std::enable_if_t<decltype(isValid(Detail::hasState())(sol[0]))::value, void>
176{
177 VTKReader vtu(fileName);
178
179 // get states at each dof location
180 const auto stateAtDof = vtu.readData<std::vector<int>>("phase presence", dataType);
181
182 // determine all states that are present
183 std::unordered_set<int> states;
184 for (size_t i = 0; i < stateAtDof.size(); ++i)
185 states.insert(stateAtDof[i]);
186
187 using PrimaryVariables = typename SolutionVector::block_type;
188 using Scalar = typename PrimaryVariables::field_type;
189 for (size_t pvIdx = 0; pvIdx < PrimaryVariables::dimension; ++pvIdx)
190 {
191 std::unordered_map<int, std::vector<Scalar>> data;
192 for (const auto& state : states)
193 data[state] = vtu.readData<std::vector<Scalar>>(pvNameFunc(pvIdx, state), dataType);
194
195 if (dataType == VTKReader::DataType::cellData)
196 {
197 std::size_t i = 0;
198 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
199 {
200 const auto eIdx = gridGeometry.elementMapper().index(element);
201 const auto state = stateAtDof[i];
202 sol[eIdx][pvIdx] = data[state][i++];
203 sol[eIdx].setState(state);
204 }
205 }
206 else
207 {
208 std::size_t i = 0;
209 constexpr int dim = GridGeometry::GridView::dimension;
210 std::vector<bool> visited(gridGeometry.gridView().size(dim), false);
211 for (const auto& element : elements(gridGeometry.gridView(), Dune::Partitions::interior))
212 {
213 for (int vIdxLocal = 0; vIdxLocal < element.subEntities(dim); ++vIdxLocal)
214 {
215 const auto vIdxGlobal = gridGeometry.vertexMapper().subIndex(element, vIdxLocal, dim);
216 if (!visited[vIdxGlobal])
217 {
218 const auto state = stateAtDof[i];
219 sol[vIdxGlobal][pvIdx] = data[state][i++];
220 sol[vIdxGlobal].setState(state);
221 visited[vIdxGlobal] = true;
222 }
223 }
224 }
225 }
226 }
227}
228
234template<class IOFields, class PrimaryVariables, class ModelTraits = void, class FluidSystem = void, class SolidSystem = void>
235auto createPVNameFunction(const std::string& paramGroup = "")
236-> typename std::enable_if_t<decltype(isValid(Detail::hasState())(PrimaryVariables(0)))::value, std::function<std::string(int,int)>>
237{
238 return [paramGroup](int pvIdx, int state = 0)
239 {
240 static auto numStates = (1 << ModelTraits::numFluidPhases()) - 1;
241 const auto paramNameWithState = "LoadSolution.PriVarNamesState" + std::to_string(state);
242
243 if (hasParamInGroup(paramGroup, "LoadSolution.PriVarNames") && !hasParamInGroup(paramGroup, paramNameWithState))
244 DUNE_THROW(Dune::NotImplemented, "please provide LoadSolution.PriVarNamesState1..." << numStates
245 << " or remove LoadSolution.PriVarNames to use the model's default primary variable names");
246
247 else if (hasParamInGroup(paramGroup, paramNameWithState))
248 {
249 const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup, paramNameWithState);
250 return pvName[pvIdx];
251 }
252 else
253 return IOFields::template primaryVariableName<ModelTraits, FluidSystem, SolidSystem>(pvIdx, state);
254 };
255}
256
262template<class IOFields, class PrimaryVariables, class ModelTraits = void, class FluidSystem = void, class SolidSystem = void>
263auto createPVNameFunction(const std::string& paramGroup = "")
264-> typename std::enable_if_t<!decltype(isValid(Detail::hasState())(PrimaryVariables(0)))::value, std::function<std::string(int,int)>>
265{
266 if (hasParamInGroup(paramGroup, "LoadSolution.PriVarNames"))
267 {
268 const auto pvName = getParamFromGroup<std::vector<std::string>>(paramGroup, "LoadSolution.PriVarNames");
269 return [n = std::move(pvName)](int pvIdx, int state = 0){ return n[pvIdx]; };
270 }
271 else
272 return [](int pvIdx, int state = 0){ return IOFields::template primaryVariableName<ModelTraits, FluidSystem, SolidSystem>(pvIdx, state); };
273}
274
285template <class SolutionVector, class PvNameFunc, class GridGeometry>
286void loadSolution(SolutionVector& sol,
287 const std::string& fileName,
288 PvNameFunc&& pvNameFunc,
289 const GridGeometry& gridGeometry)
290{
291 const auto extension = fileName.substr(fileName.find_last_of(".") + 1);
292 auto dataType = GridGeometry::discMethod == DiscretizationMethod::box ?
293 VTKReader::DataType::pointData : VTKReader::DataType::cellData;
294
295 if (extension == "vtu" || extension == "vtp")
296 {
297 if (GridGeometry::discMethod == DiscretizationMethod::staggered && extension == "vtp")
298 dataType = VTKReader::DataType::pointData;
299
300 loadSolutionFromVtkFile(sol, fileName, pvNameFunc, gridGeometry, dataType);
301 }
302 else if (extension == "pvtu" || extension == "pvtp")
303 {
304 if (GridGeometry::discMethod == DiscretizationMethod::staggered)
305 DUNE_THROW(Dune::NotImplemented, "reading staggered solution from a parallel vtk file");
306
307 loadSolutionFromVtkFile(sol, fileName, pvNameFunc, gridGeometry, dataType);
308 }
309 else
310 DUNE_THROW(Dune::NotImplemented, "loadSolution for file with extension " << extension);
311
312 // communicate solution on ghost and overlap dofs
313 if (gridGeometry.gridView().comm().size() > 1)
314 {
315 using GridView = typename GridGeometry::GridView;
316 if (dataType == VTKReader::DataType::cellData)
317 {
319 dataHandle(sol, gridGeometry.elementMapper());
320 gridGeometry.gridView().communicate(dataHandle,
321 Dune::InteriorBorder_All_Interface,
322 Dune::ForwardCommunication);
323 }
324 else
325 {
327 dataHandle(sol, gridGeometry.vertexMapper());
328 gridGeometry.gridView().communicate(dataHandle,
329 Dune::InteriorBorder_All_Interface,
330 Dune::ForwardCommunication);
331 }
332 }
333}
334
335} // end namespace Dumux
336
337#endif
A helper function for class member function introspection.
Type traits to be used with matrix types.
Type traits to be used with vector types.
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:391
auto loadSolutionFromVtkFile(SolutionVector &sol, const std::string fileName, PvNameFunc &&pvNameFunc, 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:110
void loadSolution(SolutionVector &sol, const std::string &fileName, PvNameFunc &&pvNameFunc, const GridGeometry &gridGeometry)
load a solution vector from file
Definition: loadsolution.hh:286
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:235
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:58
LoadSolutionDataHandle(Container &container, const EntityMapper &mapper)
Definition: loadsolution.hh:61
bool contains(int dim, int cd) const
Definition: loadsolution.hh:67
size_t size(const EntityType &e) const
Definition: loadsolution.hh:81
void scatter(MessageBufferImp &buff, const EntityType &e, size_t n)
Definition: loadsolution.hh:92
void gather(MessageBufferImp &buff, const EntityType &e) const
Definition: loadsolution.hh:85
bool fixedsize(int dim, int cd) const
returns true if size per entity of given dim and codim is a constant
Definition: loadsolution.hh:76
A vtk file reader using tinyxml2 as xml backend.
Definition: vtkreader.hh:49
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:79
DataType
The data array types.
Definition: vtkreader.hh:54