25#ifndef DUMUX_TWOP_GRIDDATA_TRANSFER_HH
26#define DUMUX_TWOP_GRIDDATA_TRANSFER_HH
30#include <dune/grid/common/partitionset.hh>
31#include <dune/grid/utility/persistentcontainer.hh>
48template<
class TypeTag>
56 using GridView =
typename GridGeometry::GridView;
57 using FVElementGeometry =
typename GridGeometry::LocalView;
58 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
64 using Element =
typename Grid::template Codim<0>::Entity;
65 using ElementSolution = std::decay_t<decltype(elementSolution(std::declval<Element>(),
66 std::declval<SolutionVector>(),
67 std::declval<GridGeometry>()))>;
70 using Indices =
typename ModelTraits::Indices;
74 AdaptedValues() : associatedMass(0.0) {}
77 PrimaryVariables associatedMass;
81 using PersistentContainer = Dune::PersistentContainer<Grid, AdaptedValues>;
83 static constexpr int dim = Grid::dimension;
84 static constexpr int dimWorld = Grid::dimensionworld;
88 enum { saturationIdx = Indices::saturationIdx };
93 phase0Idx = FluidSystem::phase0Idx,
94 phase1Idx = FluidSystem::phase1Idx,
102 static constexpr auto formulation = ModelTraits::priVarFormulation();
105 static_assert(!FluidSystem::isCompressible(phase0Idx)
106 && !FluidSystem::isCompressible(phase1Idx),
107 "This adaption helper is only mass conservative for incompressible fluids!");
110 static_assert(formulation == p0s1 || formulation == p1s0,
"Chosen formulation not known to the TwoPGridDataTransfer");
122 std::shared_ptr<GridGeometry> gridGeometry,
123 std::shared_ptr<const GridVariables> gridVariables,
127 , gridGeometry_(gridGeometry)
128 , gridVariables_(gridVariables)
130 , adaptionMap_(gridGeometry->gridView().grid(), 0)
141 void store(
const Grid& grid)
override
143 adaptionMap_.resize();
145 for (
auto level = grid.maxLevel(); level >= 0; level--)
147 auto fvGeometry =
localView(*gridGeometry_);
148 for (
const auto& element : elements(grid.levelGridView(level)))
151 auto& adaptedValues = adaptionMap_[element];
154 if (element.isLeaf())
156 fvGeometry.bindElement(element);
159 adaptedValues.u = ElementSolution(element, sol_, *gridGeometry_);
162 for (
const auto& scv : scvs(fvGeometry))
164 VolumeVariables volVars;
165 volVars.update(adaptedValues.u, *problem_, element, scv);
168 adaptedValues.associatedMass[phase1Idx] += poreVolume * volVars.density(phase1Idx) * volVars.saturation(phase1Idx);
169 adaptedValues.associatedMass[phase0Idx] += poreVolume * volVars.density(phase0Idx) * volVars.saturation(phase0Idx);
173 adaptedValues.count = 1;
174 adaptedValues.wasLeaf =
true;
177 if (element.level() > 0)
179 auto& adaptedValuesFather = adaptionMap_[element.father()];
182 if(&adaptedValues != &adaptedValuesFather)
183 storeAdaptionValues(adaptedValues, adaptedValuesFather);
189 if(isBox && !element.isLeaf())
190 adaptedValues.u = ElementSolution(element, sol_, *gridGeometry_);
208 gridGeometry_->update(grid.leafGridView());
217 adaptionMap_.resize();
218 sol_.resize(gridGeometry_->numDofs());
221 std::vector<Scalar> massCoeff;
222 std::vector<Scalar> associatedMass;
226 massCoeff.resize(gridGeometry_->numDofs(), 0.0);
227 associatedMass.resize(gridGeometry_->numDofs(), 0.0);
231 auto fvGeometry =
localView(*gridGeometry_);
232 for (
const auto& element : elements(gridGeometry_->gridView().grid().leafGridView(), Dune::Partitions::interior))
234 if (!element.isNew())
236 const auto& adaptedValues = adaptionMap_[element];
237 fvGeometry.bindElement(element);
240 auto elemSol = adaptedValues.u;
242 elemSol[0] /= adaptedValues.count;
244 const auto elementVolume =
volume(element.geometry(), Extrusion{});
245 for (
const auto& scv : scvs(fvGeometry))
247 VolumeVariables volVars;
248 volVars.update(elemSol, *problem_, element, scv);
251 sol_[scv.dofIndex()] = elemSol[scv.localDofIndex()];
253 const auto dofIdxGlobal = scv.dofIndex();
258 if (!adaptedValues.wasLeaf)
260 if (formulation == p0s1)
262 sol_[dofIdxGlobal][saturationIdx] = adaptedValues.associatedMass[phase1Idx];
263 sol_[dofIdxGlobal][saturationIdx] /= elementVolume * volVars.density(phase1Idx) * volVars.porosity();
265 else if (formulation == p1s0)
267 sol_[dofIdxGlobal][saturationIdx] = adaptedValues.associatedMass[phase0Idx];
268 sol_[dofIdxGlobal][saturationIdx] /= elementVolume * volVars.density(phase0Idx) * volVars.porosity();
277 if (formulation == p0s1)
279 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase1Idx) * volVars.porosity();
280 associatedMass[dofIdxGlobal] += scvVolume / elementVolume * adaptedValues.associatedMass[phase1Idx];
282 else if (formulation == p1s0)
284 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase0Idx) * volVars.porosity();
285 associatedMass[dofIdxGlobal] += scvVolume / elementVolume * adaptedValues.associatedMass[phase0Idx];
293 assert(
element.hasFather() &&
"new element does not have a father element!");
296 auto fatherElement =
element.father();
297 while(fatherElement.isNew() && fatherElement.level() > 0)
298 fatherElement = fatherElement.father();
302 const auto& adaptedValuesFather = adaptionMap_[fatherElement];
305 Scalar massFather = 0.0;
306 if (formulation == p0s1)
307 massFather = adaptedValuesFather.associatedMass[phase1Idx];
308 else if (formulation == p1s0)
309 massFather = adaptedValuesFather.associatedMass[phase0Idx];
312 auto elemSolSon = adaptedValuesFather.u;
313 elemSolSon[0] /= adaptedValuesFather.count;
315 fvGeometry.bindElement(element);
317 for (
const auto& scv : scvs(fvGeometry))
319 VolumeVariables volVars;
320 volVars.update(elemSolSon, *problem_, element, scv);
323 sol_[scv.dofIndex()] = elemSolSon[0];
326 Scalar massCoeffSon = 0.0;
327 if (formulation == p0s1)
328 massCoeffSon =
Extrusion::volume(fvGeometry, scv) * volVars.density(phase1Idx) * volVars.porosity();
329 else if (formulation == p1s0)
330 massCoeffSon =
Extrusion::volume(fvGeometry, scv) * volVars.density(phase0Idx) * volVars.porosity();
331 sol_[scv.dofIndex()][saturationIdx] =
337 auto& adaptedValuesFather = adaptionMap_[fatherElement];
339 fvGeometry.bindElement(element);
342 ElementSolution elemSolSon(element, sol_, *gridGeometry_);
343 const auto fatherGeometry = fatherElement.geometry();
344 for (
const auto& scv : scvs(fvGeometry))
345 elemSolSon[scv.localDofIndex()] =
evalSolution(fatherElement,
347 adaptedValuesFather.u,
351 const auto fatherElementVolume =
volume(fatherGeometry, Extrusion{});
352 for (
const auto& scv : scvs(fvGeometry))
354 VolumeVariables volVars;
355 volVars.update(elemSolSon, *problem_, element, scv);
357 const auto dofIdxGlobal = scv.dofIndex();
359 if (formulation == p0s1)
361 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase1Idx) * volVars.porosity();
362 associatedMass[dofIdxGlobal] += scvVolume / fatherElementVolume * adaptedValuesFather.associatedMass[phase1Idx];
364 else if (formulation == p1s0)
366 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase0Idx) * volVars.porosity();
367 associatedMass[dofIdxGlobal] += scvVolume / fatherElementVolume * adaptedValuesFather.associatedMass[phase0Idx];
371 sol_[dofIdxGlobal] = elemSolSon[scv.localDofIndex()];
379 for(std::size_t dofIdxGlobal = 0; dofIdxGlobal < gridGeometry_->numDofs(); dofIdxGlobal++)
380 sol_[dofIdxGlobal][saturationIdx] = associatedMass[dofIdxGlobal] / massCoeff[dofIdxGlobal];
384 adaptionMap_.resize(
typename PersistentContainer::Value() );
385 adaptionMap_.shrinkToFit();
386 adaptionMap_.fill(
typename PersistentContainer::Value() );
411 static void storeAdaptionValues(AdaptedValues& adaptedValues,
412 AdaptedValues& adaptedValuesFather)
415 adaptedValuesFather.associatedMass += adaptedValues.associatedMass;
422 auto values = adaptedValues.u[0];
423 values /= adaptedValues.count;
424 adaptedValuesFather.u[0] += values;
427 adaptedValuesFather.count += 1;
430 adaptedValuesFather.wasLeaf =
false;
436 adaptedValuesFather.count = 1;
439 adaptedValuesFather.wasLeaf =
false;
443 std::shared_ptr<const Problem> problem_;
444 std::shared_ptr<GridGeometry> gridGeometry_;
445 std::shared_ptr<const GridVariables> gridVariables_;
446 SolutionVector& sol_;
447 PersistentContainer adaptionMap_;
Element solution classes and factory functions.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
Defines an enumeration for the formulations accepted by the two-phase model.
Compute the volume of several common geometry types.
@ p1s0
first phase saturation and second phase pressure as primary variables
@ p0s1
first phase pressure and second phase saturation as primary variables
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
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:165
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:251
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
constexpr Box box
Definition: method.hh:136
Interface to be used by classes transferring grid data on adaptive grids.
Definition: adaptive/griddatatransfer.hh:35
Class performing the transfer of data on a grid from before to after adaptation.
Definition: porousmediumflow/2p/griddatatransfer.hh:50
void store(const Grid &grid) override
Stores primary variables and additional data.
Definition: porousmediumflow/2p/griddatatransfer.hh:141
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:121
void reconstruct(const Grid &grid) override
Reconstruct missing primary variables (where elements are created/deleted)
Definition: porousmediumflow/2p/griddatatransfer.hh:206
Interface to be used by classes transferring grid data on adaptive grids.
Declares all properties used in Dumux.