version 3.10-dev
porousmediumflow/2p/griddatatransfer.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//
13#ifndef DUMUX_TWOP_GRIDDATA_TRANSFER_HH
14#define DUMUX_TWOP_GRIDDATA_TRANSFER_HH
15
16#include <memory>
17
18#include <dune/grid/common/partitionset.hh>
19#include <dune/grid/utility/persistentcontainer.hh>
20
22
29
30namespace Dumux {
31
36template<class TypeTag>
37class TwoPGridDataTransfer : public GridDataTransfer<GetPropType<TypeTag, Properties::Grid>>
38{
44 using GridView = typename GridGeometry::GridView;
45 using FVElementGeometry = typename GridGeometry::LocalView;
46 using SubControlVolume = typename FVElementGeometry::SubControlVolume;
47 using Extrusion = Extrusion_t<GridGeometry>;
52 using Element = typename Grid::template Codim<0>::Entity;
53 using ElementSolution = std::decay_t<decltype(elementSolution(std::declval<Element>(),
54 std::declval<SolutionVector>(),
55 std::declval<GridGeometry>()))>;
58 using Indices = typename ModelTraits::Indices;
59
60 struct AdaptedValues
61 {
62 AdaptedValues() : associatedMass(0.0) {}
63 ElementSolution u;
64 int count = 0;
65 PrimaryVariables associatedMass;
66 bool wasLeaf = false;
67 };
68
69 using PersistentContainer = Dune::PersistentContainer<Grid, AdaptedValues>;
70
71 static constexpr int dim = Grid::dimension;
72 static constexpr int dimWorld = Grid::dimensionworld;
74
75 // saturation primary variable index
76 enum { saturationIdx = Indices::saturationIdx };
77
78 // phase indices
79 enum
80 {
81 phase0Idx = FluidSystem::phase0Idx,
82 phase1Idx = FluidSystem::phase1Idx,
83 };
84
85 // formulations
86 static constexpr auto p0s1 = TwoPFormulation::p0s1;
87 static constexpr auto p1s0 = TwoPFormulation::p1s0;
88
89 // the formulation that is actually used
90 static constexpr auto formulation = ModelTraits::priVarFormulation();
91
92 // This won't work (mass conservative) for compressible fluids
93 static_assert(!FluidSystem::isCompressible(phase0Idx)
94 && !FluidSystem::isCompressible(phase1Idx),
95 "This adaption helper is only mass conservative for incompressible fluids!");
96
97 // check if the used formulation is implemented here
98 static_assert(formulation == p0s1 || formulation == p1s0, "Chosen formulation not known to the TwoPGridDataTransfer");
99
100public:
109 TwoPGridDataTransfer(std::shared_ptr<const Problem> problem,
110 std::shared_ptr<GridGeometry> gridGeometry,
111 std::shared_ptr<const GridVariables> gridVariables,
112 SolutionVector& sol)
113 : ParentType()
114 , problem_(problem)
115 , gridGeometry_(gridGeometry)
116 , gridVariables_(gridVariables)
117 , sol_(sol)
118 , adaptionMap_(gridGeometry->gridView().grid(), 0)
119 {}
120
129 void store(const Grid& grid) override
130 {
131 adaptionMap_.resize();
132
133 for (auto level = grid.maxLevel(); level >= 0; level--)
134 {
135 auto fvGeometry = localView(*gridGeometry_);
136 for (const auto& element : elements(grid.levelGridView(level)))
137 {
138 // get map entry
139 auto& adaptedValues = adaptionMap_[element];
140
141 // put values in the map for leaf elements
142 if (element.isLeaf())
143 {
144 fvGeometry.bindElement(element);
145
146 // store current element solution
147 adaptedValues.u = ElementSolution(element, sol_, *gridGeometry_);
148
149 // compute mass in the scvs
150 for (const auto& scv : scvs(fvGeometry))
151 {
152 VolumeVariables volVars;
153 volVars.update(adaptedValues.u, *problem_, element, scv);
154
155 const auto poreVolume = Extrusion::volume(fvGeometry, scv)*volVars.porosity();
156 adaptedValues.associatedMass[phase1Idx] += poreVolume * volVars.density(phase1Idx) * volVars.saturation(phase1Idx);
157 adaptedValues.associatedMass[phase0Idx] += poreVolume * volVars.density(phase0Idx) * volVars.saturation(phase0Idx);
158 }
159
160 // leaf elements always start with count = 1
161 adaptedValues.count = 1;
162 adaptedValues.wasLeaf = true;
163 }
164 // Average in father elements
165 if (element.level() > 0)
166 {
167 auto& adaptedValuesFather = adaptionMap_[element.father()];
168 // For some grids the father element is identical to the son element.
169 // In that case averaging is not necessary.
170 if(&adaptedValues != &adaptedValuesFather)
171 storeAdaptionValues(adaptedValues, adaptedValuesFather);
172 }
173
174 // The vertices of the non-leaf elements exist on the leaf as well
175 // This element solution constructor uses the vertex mapper to obtain
176 // the privars at the vertices, thus, this works for non-leaf elements!
177 if(isBox && !element.isLeaf())
178 adaptedValues.u = ElementSolution(element, sol_, *gridGeometry_);
179 }
180 }
181 }
182
194 void reconstruct(const Grid& grid) override
195 {
196 gridGeometry_->update(grid.leafGridView());
197 reconstruct_();
198 }
199
200 private:
201
202 void reconstruct_()
203 {
204 // resize stuff (grid might have changed)
205 adaptionMap_.resize();
206 sol_.resize(gridGeometry_->numDofs());
207
208 // vectors storing the mass associated with each vertex, when using the box method
209 std::vector<Scalar> massCoeff;
210 std::vector<Scalar> associatedMass;
211
212 if(isBox)
213 {
214 massCoeff.resize(gridGeometry_->numDofs(), 0.0);
215 associatedMass.resize(gridGeometry_->numDofs(), 0.0);
216 }
217
218 // iterate over leaf and reconstruct the solution
219 auto fvGeometry = localView(*gridGeometry_);
220 for (const auto& element : elements(gridGeometry_->gridView().grid().leafGridView(), Dune::Partitions::interior))
221 {
222 if (!element.isNew())
223 {
224 const auto& adaptedValues = adaptionMap_[element];
225 fvGeometry.bindElement(element);
226
227 // obtain element solution from map (divide by count!)
228 auto elemSol = adaptedValues.u;
229 if (!isBox)
230 elemSol[0] /= adaptedValues.count;
231
232 const auto elementVolume = volume(element.geometry(), Extrusion{});
233 for (const auto& scv : scvs(fvGeometry))
234 {
235 VolumeVariables volVars;
236 volVars.update(elemSol, *problem_, element, scv);
237
238 // write solution at dof in current solution vector
239 sol_[scv.dofIndex()] = elemSol[scv.localDofIndex()];
240
241 const auto dofIdxGlobal = scv.dofIndex();
242 // For cc schemes, overwrite the saturation by a mass conservative one here
243 if (!isBox)
244 {
245 // only recalculate the saturations if element hasn't been leaf before adaptation
246 if (!adaptedValues.wasLeaf)
247 {
248 if (formulation == p0s1)
249 {
250 sol_[dofIdxGlobal][saturationIdx] = adaptedValues.associatedMass[phase1Idx];
251 sol_[dofIdxGlobal][saturationIdx] /= elementVolume * volVars.density(phase1Idx) * volVars.porosity();
252 }
253 else if (formulation == p1s0)
254 {
255 sol_[dofIdxGlobal][saturationIdx] = adaptedValues.associatedMass[phase0Idx];
256 sol_[dofIdxGlobal][saturationIdx] /= elementVolume * volVars.density(phase0Idx) * volVars.porosity();
257 }
258 }
259 }
260
261 // For the box scheme, add mass & mass coefficient to container (saturations are recalculated at the end)
262 else
263 {
264 const auto scvVolume = Extrusion::volume(fvGeometry, scv);
265 if (formulation == p0s1)
266 {
267 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase1Idx) * volVars.porosity();
268 associatedMass[dofIdxGlobal] += scvVolume / elementVolume * adaptedValues.associatedMass[phase1Idx];
269 }
270 else if (formulation == p1s0)
271 {
272 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase0Idx) * volVars.porosity();
273 associatedMass[dofIdxGlobal] += scvVolume / elementVolume * adaptedValues.associatedMass[phase0Idx];
274 }
275 }
276 }
277 }
278 else
279 {
280 // value is not in map, interpolate from father element
281 assert(element.hasFather() && "new element does not have a father element!");
282
283 // find the ancestor element that existed on the old grid already
284 auto fatherElement = element.father();
285 while(fatherElement.isNew() && fatherElement.level() > 0)
286 fatherElement = fatherElement.father();
287
288 if(!isBox)
289 {
290 const auto& adaptedValuesFather = adaptionMap_[fatherElement];
291
292 // obtain the mass contained in father
293 Scalar massFather = 0.0;
294 if (formulation == p0s1)
295 massFather = adaptedValuesFather.associatedMass[phase1Idx];
296 else if (formulation == p1s0)
297 massFather = adaptedValuesFather.associatedMass[phase0Idx];
298
299 // obtain the element solution through the father
300 auto elemSolSon = adaptedValuesFather.u;
301 elemSolSon[0] /= adaptedValuesFather.count;
302
303 fvGeometry.bindElement(element);
304
305 for (const auto& scv : scvs(fvGeometry))
306 {
307 VolumeVariables volVars;
308 volVars.update(elemSolSon, *problem_, element, scv);
309
310 // store constructed values of son in the current solution
311 sol_[scv.dofIndex()] = elemSolSon[0];
312
313 // overwrite the saturation by a mass conservative one here
314 Scalar massCoeffSon = 0.0;
315 if (formulation == p0s1)
316 massCoeffSon = Extrusion::volume(fvGeometry, scv) * volVars.density(phase1Idx) * volVars.porosity();
317 else if (formulation == p1s0)
318 massCoeffSon = Extrusion::volume(fvGeometry, scv) * volVars.density(phase0Idx) * volVars.porosity();
319 sol_[scv.dofIndex()][saturationIdx] =
320 ( Extrusion::volume(fvGeometry, scv)/volume(fatherElement.geometry(), Extrusion{})*massFather )/massCoeffSon;
321 }
322 }
323 else
324 {
325 auto& adaptedValuesFather = adaptionMap_[fatherElement];
326
327 fvGeometry.bindElement(element);
328
329 // interpolate solution in the father to the vertices of the new son
330 ElementSolution elemSolSon(element, sol_, *gridGeometry_);
331 const auto fatherGeometry = fatherElement.geometry();
332 for (const auto& scv : scvs(fvGeometry))
333 elemSolSon[scv.localDofIndex()] = evalSolution(fatherElement,
334 fatherGeometry,
335 adaptedValuesFather.u,
336 scv.dofPosition());
337
338 // compute mass & mass coefficients for the scvs (saturations are recalculated at the end)
339 const auto fatherElementVolume = volume(fatherGeometry, Extrusion{});
340 for (const auto& scv : scvs(fvGeometry))
341 {
342 VolumeVariables volVars;
343 volVars.update(elemSolSon, *problem_, element, scv);
344
345 const auto dofIdxGlobal = scv.dofIndex();
346 const auto scvVolume = Extrusion::volume(fvGeometry, scv);
347 if (formulation == p0s1)
348 {
349 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase1Idx) * volVars.porosity();
350 associatedMass[dofIdxGlobal] += scvVolume / fatherElementVolume * adaptedValuesFather.associatedMass[phase1Idx];
351 }
352 else if (formulation == p1s0)
353 {
354 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase0Idx) * volVars.porosity();
355 associatedMass[dofIdxGlobal] += scvVolume / fatherElementVolume * adaptedValuesFather.associatedMass[phase0Idx];
356 }
357
358 // store constructed (pressure) values of son in the current solution (saturation comes later)
359 sol_[dofIdxGlobal] = elemSolSon[scv.localDofIndex()];
360 }
361 }
362 }
363 }
364
365 if(isBox)
366 {
367 for(std::size_t dofIdxGlobal = 0; dofIdxGlobal < gridGeometry_->numDofs(); dofIdxGlobal++)
368 sol_[dofIdxGlobal][saturationIdx] = associatedMass[dofIdxGlobal] / massCoeff[dofIdxGlobal];
369 }
370
371 // reset entries in adaptation map
372 adaptionMap_.resize( typename PersistentContainer::Value() );
373 adaptionMap_.shrinkToFit();
374 adaptionMap_.fill( typename PersistentContainer::Value() );
375
377//#if HAVE_MPI
378// // communicate ghost data
379// using SolutionTypes = typename GetProp<TypeTag, SolutionTypes>;
380// using ElementMapper = typename SolutionTypes::ElementMapper;
381// using DataHandle = VectorExchange<ElementMapper, std::vector<CellData> >;
382// DataHandle dataHandle(problem.elementMapper(), this->cellDataGlobal());
383// problem.gridView().template communicate<DataHandle>(dataHandle,
384// Dune::InteriorBorder_All_Interface,
385// Dune::ForwardCommunication);
386//#endif
387 }
388
399 static void storeAdaptionValues(AdaptedValues& adaptedValues,
400 AdaptedValues& adaptedValuesFather)
401 {
402 // Add associated mass of the child to the one of the father
403 adaptedValuesFather.associatedMass += adaptedValues.associatedMass;
404
405 if(!isBox)
406 {
407 // add the child's primary variables to the ones of father
408 // we have to divide the child's ones in case it was composed
409 // of several children as well!
410 auto values = adaptedValues.u[0];
411 values /= adaptedValues.count;
412 adaptedValuesFather.u[0] += values;
413
414 // keep track of the number of children that composed this father
415 adaptedValuesFather.count += 1;
416
417 // A father element is never leaf
418 adaptedValuesFather.wasLeaf = false;
419 }
420 else
421 {
422 // For the box scheme, scaling of primary variables by count is obsolete
423 // Thus, we always want count = 1
424 adaptedValuesFather.count = 1;
425
426 // A father element is never leaf
427 adaptedValuesFather.wasLeaf = false;
428 }
429 }
430
431 std::shared_ptr<const Problem> problem_;
432 std::shared_ptr<GridGeometry> gridGeometry_;
433 std::shared_ptr<const GridVariables> gridVariables_;
434 SolutionVector& sol_;
435 PersistentContainer adaptionMap_;
436};
437
438} // end namespace Dumux
439
440#endif
Interface to be used by classes transferring grid data on adaptive grids.
Interface to be used by classes transferring grid data on adaptive grids.
Definition: adaptive/griddatatransfer.hh:23
Class performing the transfer of data on a grid from before to after adaptation.
Definition: porousmediumflow/2p/griddatatransfer.hh:38
void store(const Grid &grid) override
Stores primary variables and additional data.
Definition: porousmediumflow/2p/griddatatransfer.hh:129
TwoPGridDataTransfer(std::shared_ptr< const Problem > problem, std::shared_ptr< GridGeometry > gridGeometry, std::shared_ptr< const GridVariables > gridVariables, SolutionVector &sol)
Constructor.
Definition: porousmediumflow/2p/griddatatransfer.hh:109
void reconstruct(const Grid &grid) override
Reconstruct missing primary variables (where elements are created/deleted)
Definition: porousmediumflow/2p/griddatatransfer.hh:194
Defines all properties used in Dumux.
Element solution classes and factory functions.
Helper classes to compute the integration elements.
Defines an enumeration for the formulations accepted by the two-phase model.
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
PrimaryVariables evalSolution(const Element &element, const typename Element::Geometry &geometry, const typename FVElementGeometry::GridGeometry &gridGeometry, const CVFEElementSolution< FVElementGeometry, PrimaryVariables > &elemSol, const typename Element::Geometry::GlobalCoordinate &globalPos, bool ignoreState=false)
Interpolates a given box element solution at a given global position. Uses the finite element cache o...
Definition: evalsolution.hh:152
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
@ p1s0
first phase saturation and second phase pressure as primary variables
@ p0s1
first phase pressure and second phase saturation as primary variables
The available discretization methods in Dumux.
constexpr Box box
Definition: method.hh:147
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
Compute the volume of several common geometry types.