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>
47template<
class TypeTag>
55 using GridView =
typename GridGeometry::GridView;
56 using FVElementGeometry =
typename GridGeometry::LocalView;
57 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
63 using Element =
typename Grid::template Codim<0>::Entity;
64 using ElementSolution = std::decay_t<decltype(elementSolution(std::declval<Element>(),
65 std::declval<SolutionVector>(),
66 std::declval<GridGeometry>()))>;
69 using Indices =
typename ModelTraits::Indices;
73 AdaptedValues() : associatedMass(0.0) {}
76 PrimaryVariables associatedMass;
80 using PersistentContainer = Dune::PersistentContainer<Grid, AdaptedValues>;
82 static constexpr int dim = Grid::dimension;
83 static constexpr int dimWorld = Grid::dimensionworld;
87 enum { saturationIdx = Indices::saturationIdx };
92 phase0Idx = FluidSystem::phase0Idx,
93 phase1Idx = FluidSystem::phase1Idx,
101 static constexpr auto formulation = ModelTraits::priVarFormulation();
104 static_assert(!FluidSystem::isCompressible(phase0Idx)
105 && !FluidSystem::isCompressible(phase1Idx),
106 "This adaption helper is only mass conservative for incompressible fluids!");
109 static_assert(formulation == p0s1 || formulation == p1s0,
"Chosen formulation not known to the TwoPGridDataTransfer");
121 std::shared_ptr<GridGeometry> gridGeometry,
122 std::shared_ptr<const GridVariables> gridVariables,
126 , gridGeometry_(gridGeometry)
127 , gridVariables_(gridVariables)
129 , adaptionMap_(gridGeometry->gridView().grid(), 0)
140 void store(
const Grid& grid)
override
142 adaptionMap_.resize();
144 for (
auto level = grid.maxLevel(); level >= 0; level--)
146 auto fvGeometry =
localView(*gridGeometry_);
147 for (
const auto& element : elements(grid.levelGridView(level)))
150 auto& adaptedValues = adaptionMap_[element];
153 if (element.isLeaf())
155 fvGeometry.bindElement(element);
158 adaptedValues.u = ElementSolution(element, sol_, *gridGeometry_);
161 for (
const auto& scv : scvs(fvGeometry))
163 VolumeVariables volVars;
164 volVars.update(adaptedValues.u, *problem_, element, scv);
167 adaptedValues.associatedMass[phase1Idx] += poreVolume * volVars.density(phase1Idx) * volVars.saturation(phase1Idx);
168 adaptedValues.associatedMass[phase0Idx] += poreVolume * volVars.density(phase0Idx) * volVars.saturation(phase0Idx);
172 adaptedValues.count = 1;
173 adaptedValues.wasLeaf =
true;
176 if (element.level() > 0)
178 auto& adaptedValuesFather = adaptionMap_[element.father()];
181 if(&adaptedValues != &adaptedValuesFather)
182 storeAdaptionValues(adaptedValues, adaptedValuesFather);
188 if(isBox && !element.isLeaf())
189 adaptedValues.u = ElementSolution(element, sol_, *gridGeometry_);
207 gridGeometry_->update(grid.leafGridView());
216 adaptionMap_.resize();
217 sol_.resize(gridGeometry_->numDofs());
220 std::vector<Scalar> massCoeff;
221 std::vector<Scalar> associatedMass;
225 massCoeff.resize(gridGeometry_->numDofs(), 0.0);
226 associatedMass.resize(gridGeometry_->numDofs(), 0.0);
230 auto fvGeometry =
localView(*gridGeometry_);
231 for (
const auto& element : elements(gridGeometry_->gridView().grid().leafGridView(), Dune::Partitions::interior))
233 if (!element.isNew())
235 const auto& adaptedValues = adaptionMap_[element];
236 fvGeometry.bindElement(element);
239 auto elemSol = adaptedValues.u;
241 elemSol[0] /= adaptedValues.count;
244 for (
const auto& scv : scvs(fvGeometry))
246 VolumeVariables volVars;
247 volVars.update(elemSol, *problem_, element, scv);
250 sol_[scv.dofIndex()] = elemSol[scv.localDofIndex()];
252 const auto dofIdxGlobal = scv.dofIndex();
257 if (!adaptedValues.wasLeaf)
259 if (formulation == p0s1)
261 sol_[dofIdxGlobal][saturationIdx] = adaptedValues.associatedMass[phase1Idx];
262 sol_[dofIdxGlobal][saturationIdx] /= elementVolume * volVars.density(phase1Idx) * volVars.porosity();
264 else if (formulation == p1s0)
266 sol_[dofIdxGlobal][saturationIdx] = adaptedValues.associatedMass[phase0Idx];
267 sol_[dofIdxGlobal][saturationIdx] /= elementVolume * volVars.density(phase0Idx) * volVars.porosity();
276 if (formulation == p0s1)
278 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase1Idx) * volVars.porosity();
279 associatedMass[dofIdxGlobal] += scvVolume / elementVolume * adaptedValues.associatedMass[phase1Idx];
281 else if (formulation == p1s0)
283 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase0Idx) * volVars.porosity();
284 associatedMass[dofIdxGlobal] += scvVolume / elementVolume * adaptedValues.associatedMass[phase0Idx];
292 assert(
element.hasFather() &&
"new element does not have a father element!");
295 auto fatherElement =
element.father();
296 while(fatherElement.isNew() && fatherElement.level() > 0)
297 fatherElement = fatherElement.father();
301 const auto& adaptedValuesFather = adaptionMap_[fatherElement];
304 Scalar massFather = 0.0;
305 if (formulation == p0s1)
306 massFather = adaptedValuesFather.associatedMass[phase1Idx];
307 else if (formulation == p1s0)
308 massFather = adaptedValuesFather.associatedMass[phase0Idx];
311 auto elemSolSon = adaptedValuesFather.u;
312 elemSolSon[0] /= adaptedValuesFather.count;
314 fvGeometry.bindElement(element);
316 for (
const auto& scv : scvs(fvGeometry))
318 VolumeVariables volVars;
319 volVars.update(elemSolSon, *problem_, element, scv);
322 sol_[scv.dofIndex()] = elemSolSon[0];
325 Scalar massCoeffSon = 0.0;
326 if (formulation == p0s1)
327 massCoeffSon =
Extrusion::volume(scv) * volVars.density(phase1Idx) * volVars.porosity();
328 else if (formulation == p1s0)
329 massCoeffSon =
Extrusion::volume(scv) * volVars.density(phase0Idx) * volVars.porosity();
330 sol_[scv.dofIndex()][saturationIdx] =
336 auto& adaptedValuesFather = adaptionMap_[fatherElement];
338 fvGeometry.bindElement(element);
341 ElementSolution elemSolSon(element, sol_, *gridGeometry_);
342 const auto fatherGeometry = fatherElement.geometry();
343 for (
const auto& scv : scvs(fvGeometry))
344 elemSolSon[scv.localDofIndex()] =
evalSolution(fatherElement,
346 adaptedValuesFather.u,
351 for (
const auto& scv : scvs(fvGeometry))
353 VolumeVariables volVars;
354 volVars.update(elemSolSon, *problem_, element, scv);
356 const auto dofIdxGlobal = scv.dofIndex();
358 if (formulation == p0s1)
360 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase1Idx) * volVars.porosity();
361 associatedMass[dofIdxGlobal] += scvVolume / fatherElementVolume * adaptedValuesFather.associatedMass[phase1Idx];
363 else if (formulation == p1s0)
365 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase0Idx) * volVars.porosity();
366 associatedMass[dofIdxGlobal] += scvVolume / fatherElementVolume * adaptedValuesFather.associatedMass[phase0Idx];
370 sol_[dofIdxGlobal] = elemSolSon[scv.localDofIndex()];
378 for(std::size_t dofIdxGlobal = 0; dofIdxGlobal < gridGeometry_->numDofs(); dofIdxGlobal++)
379 sol_[dofIdxGlobal][saturationIdx] = associatedMass[dofIdxGlobal] / massCoeff[dofIdxGlobal];
383 adaptionMap_.resize(
typename PersistentContainer::Value() );
384 adaptionMap_.shrinkToFit();
385 adaptionMap_.fill(
typename PersistentContainer::Value() );
410 static void storeAdaptionValues(AdaptedValues& adaptedValues,
411 AdaptedValues& adaptedValuesFather)
414 adaptedValuesFather.associatedMass += adaptedValues.associatedMass;
421 auto values = adaptedValues.u[0];
422 values /= adaptedValues.count;
423 adaptedValuesFather.u[0] += values;
426 adaptedValuesFather.count += 1;
429 adaptedValuesFather.wasLeaf =
false;
435 adaptedValuesFather.count = 1;
438 adaptedValuesFather.wasLeaf =
false;
442 std::shared_ptr<const Problem> problem_;
443 std::shared_ptr<GridGeometry> gridGeometry_;
444 std::shared_ptr<const GridVariables> gridVariables_;
445 SolutionVector& sol_;
446 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.
@ 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:95
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
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
Definition: propertysystem.hh:150
constexpr Box box
Definition: method.hh:139
Scalar volume(Shape shape, Scalar inscribedRadius)
Returns the volume of a given geometry based on the inscribed radius.
Definition: poreproperties.hh:73
Interface to be used by classes transferring grid data on adpative 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:49
void store(const Grid &grid) override
Stores primary variables and additional data.
Definition: porousmediumflow/2p/griddatatransfer.hh:140
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:120
void reconstruct(const Grid &grid) override
Reconstruct missing primary variables (where elements are created/deleted)
Definition: porousmediumflow/2p/griddatatransfer.hh:205
Interface to be used by classes transferring grid data on adpative grids.
Declares all properties used in Dumux.