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