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