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>
54 using FVElementGeometry =
typename GridGeometry::LocalView;
55 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
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;
71 AdaptedValues() : associatedMass(0.0) {}
74 PrimaryVariables associatedMass;
78 using PersistentContainer = Dune::PersistentContainer<Grid, AdaptedValues>;
80 static constexpr int dim = Grid::dimension;
81 static constexpr int dimWorld = Grid::dimensionworld;
85 enum { saturationIdx = Indices::saturationIdx };
90 phase0Idx = FluidSystem::phase0Idx,
91 phase1Idx = FluidSystem::phase1Idx,
99 static constexpr auto formulation = ModelTraits::priVarFormulation();
102 static_assert(!FluidSystem::isCompressible(phase0Idx)
103 && !FluidSystem::isCompressible(phase1Idx),
104 "This adaption helper is only mass conservative for incompressible fluids!");
107 static_assert(formulation == p0s1 || formulation == p1s0,
"Chosen formulation not known to the TwoPGridDataTransfer");
119 std::shared_ptr<GridGeometry> gridGeometry,
120 std::shared_ptr<const GridVariables> gridVariables,
124 , gridGeometry_(gridGeometry)
125 , gridVariables_(gridVariables)
127 , adaptionMap_(gridGeometry->gridView().grid(), 0)
140 adaptionMap_.resize();
142 const auto& grid = gridGeometry_->gridView().grid();
143 for (
auto level = grid.maxLevel(); level >= 0; level--)
145 for (
const auto& element : elements(grid.levelGridView(level)))
148 auto& adaptedValues = adaptionMap_[element];
151 if (element.isLeaf())
153 auto fvGeometry =
localView(*gridGeometry_);
154 fvGeometry.bindElement(element);
157 adaptedValues.u = ElementSolution(element, sol_, *gridGeometry_);
160 for (
const auto& scv : scvs(fvGeometry))
162 VolumeVariables volVars;
163 volVars.update(adaptedValues.u, *problem_, element, scv);
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);
171 adaptedValues.count = 1;
172 adaptedValues.wasLeaf =
true;
175 if (element.level() > 0)
177 auto& adaptedValuesFather = adaptionMap_[element.father()];
180 if(&adaptedValues != &adaptedValuesFather)
181 storeAdaptionValues(adaptedValues, adaptedValuesFather);
187 if(isBox && !element.isLeaf())
188 adaptedValues.u = ElementSolution(element, sol_, *gridGeometry_);
207 adaptionMap_.resize();
208 gridGeometry_->update();
209 sol_.resize(gridGeometry_->numDofs());
212 std::vector<Scalar> massCoeff;
213 std::vector<Scalar> associatedMass;
217 massCoeff.resize(gridGeometry_->numDofs(), 0.0);
218 associatedMass.resize(gridGeometry_->numDofs(), 0.0);
222 for (
const auto& element : elements(gridGeometry_->gridView().grid().leafGridView(), Dune::Partitions::interior))
224 if (!element.isNew())
226 const auto& adaptedValues = adaptionMap_[element];
228 auto fvGeometry =
localView(*gridGeometry_);
229 fvGeometry.bindElement(element);
232 auto elemSol = adaptedValues.u;
234 elemSol[0] /= adaptedValues.count;
236 const auto elementVolume = Extrusion::volume(element.geometry());
237 for (
const auto& scv : scvs(fvGeometry))
239 VolumeVariables volVars;
240 volVars.update(elemSol, *problem_, element, scv);
243 sol_[scv.dofIndex()] = elemSol[scv.localDofIndex()];
245 const auto dofIdxGlobal = scv.dofIndex();
250 if (!adaptedValues.wasLeaf)
252 if (formulation == p0s1)
254 sol_[dofIdxGlobal][saturationIdx] = adaptedValues.associatedMass[phase1Idx];
255 sol_[dofIdxGlobal][saturationIdx] /= elementVolume * volVars.density(phase1Idx) * volVars.porosity();
257 else if (formulation == p1s0)
259 sol_[dofIdxGlobal][saturationIdx] = adaptedValues.associatedMass[phase0Idx];
260 sol_[dofIdxGlobal][saturationIdx] /= elementVolume * volVars.density(phase0Idx) * volVars.porosity();
268 const auto scvVolume = Extrusion::volume(scv);
269 if (formulation == p0s1)
271 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase1Idx) * volVars.porosity();
272 associatedMass[dofIdxGlobal] += scvVolume / elementVolume * adaptedValues.associatedMass[phase1Idx];
274 else if (formulation == p1s0)
276 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase0Idx) * volVars.porosity();
277 associatedMass[dofIdxGlobal] += scvVolume / elementVolume * adaptedValues.associatedMass[phase0Idx];
285 assert(element.hasFather() &&
"new element does not have a father element!");
288 auto fatherElement = element.father();
289 while(fatherElement.isNew() && fatherElement.level() > 0)
290 fatherElement = fatherElement.father();
294 const auto& adaptedValuesFather = adaptionMap_[fatherElement];
297 Scalar massFather = 0.0;
298 if (formulation == p0s1)
299 massFather = adaptedValuesFather.associatedMass[phase1Idx];
300 else if (formulation == p1s0)
301 massFather = adaptedValuesFather.associatedMass[phase0Idx];
304 auto elemSolSon = adaptedValuesFather.u;
305 elemSolSon[0] /= adaptedValuesFather.count;
307 auto fvGeometry =
localView(*gridGeometry_);
308 fvGeometry.bindElement(element);
310 for (
const auto& scv : scvs(fvGeometry))
312 VolumeVariables volVars;
313 volVars.update(elemSolSon, *problem_, element, scv);
316 sol_[scv.dofIndex()] = elemSolSon[0];
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;
330 auto& adaptedValuesFather = adaptionMap_[fatherElement];
332 auto fvGeometry =
localView(*gridGeometry_);
333 fvGeometry.bindElement(element);
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,
341 adaptedValuesFather.u,
345 const auto fatherElementVolume = Extrusion::volume(fatherGeometry);
346 for (
const auto& scv : scvs(fvGeometry))
348 VolumeVariables volVars;
349 volVars.update(elemSolSon, *problem_, element, scv);
351 const auto dofIdxGlobal = scv.dofIndex();
352 const auto scvVolume = Extrusion::volume(scv);
353 if (formulation == p0s1)
355 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase1Idx) * volVars.porosity();
356 associatedMass[dofIdxGlobal] += scvVolume / fatherElementVolume * adaptedValuesFather.associatedMass[phase1Idx];
358 else if (formulation == p1s0)
360 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase0Idx) * volVars.porosity();
361 associatedMass[dofIdxGlobal] += scvVolume / fatherElementVolume * adaptedValuesFather.associatedMass[phase0Idx];
365 sol_[dofIdxGlobal] = elemSolSon[scv.localDofIndex()];
373 for(std::size_t dofIdxGlobal = 0; dofIdxGlobal < gridGeometry_->numDofs(); dofIdxGlobal++)
374 sol_[dofIdxGlobal][saturationIdx] = associatedMass[dofIdxGlobal] / massCoeff[dofIdxGlobal];
378 adaptionMap_.resize(
typename PersistentContainer::Value() );
379 adaptionMap_.shrinkToFit();
380 adaptionMap_.fill(
typename PersistentContainer::Value() );
407 static void storeAdaptionValues(AdaptedValues& adaptedValues,
408 AdaptedValues& adaptedValuesFather)
411 adaptedValuesFather.associatedMass += adaptedValues.associatedMass;
418 auto values = adaptedValues.u[0];
419 values /= adaptedValues.count;
420 adaptedValuesFather.u[0] += values;
423 adaptedValuesFather.count += 1;
426 adaptedValuesFather.wasLeaf =
false;
432 adaptedValuesFather.count = 1;
435 adaptedValuesFather.wasLeaf =
false;
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_;
Element solution classes and factory functions.
The available discretization methods in Dumux.
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
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.