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>
46template<
class TypeTag>
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>()))>;
69 AdaptedValues() : associatedMass(0.0) {}
72 PrimaryVariables associatedMass;
76 using PersistentContainer = Dune::PersistentContainer<Grid, AdaptedValues>;
78 static constexpr int dim = Grid::dimension;
79 static constexpr int dimWorld = Grid::dimensionworld;
83 enum { saturationIdx = Indices::saturationIdx };
88 phase0Idx = FluidSystem::phase0Idx,
89 phase1Idx = FluidSystem::phase1Idx,
97 static constexpr auto formulation = ModelTraits::priVarFormulation();
100 static_assert(!FluidSystem::isCompressible(phase0Idx)
101 && !FluidSystem::isCompressible(phase1Idx),
102 "This adaption helper is only mass conservative for incompressible fluids!");
105 static_assert(formulation == p0s1 || formulation == p1s0,
"Chosen formulation not known to the TwoPGridDataTransfer");
117 std::shared_ptr<GridGeometry> gridGeometry,
118 std::shared_ptr<const GridVariables> gridVariables,
122 , gridGeometry_(gridGeometry)
123 , gridVariables_(gridVariables)
125 , adaptionMap_(gridGeometry->gridView().grid(), 0)
138 adaptionMap_.resize();
140 const auto& grid = gridGeometry_->gridView().grid();
141 for (
auto level = grid.maxLevel(); level >= 0; level--)
143 for (
const auto& element : elements(grid.levelGridView(level)))
146 auto& adaptedValues = adaptionMap_[element];
149 if (element.isLeaf())
151 auto fvGeometry =
localView(*gridGeometry_);
152 fvGeometry.bindElement(element);
155 adaptedValues.u = ElementSolution(element, sol_, *gridGeometry_);
158 for (
const auto& scv : scvs(fvGeometry))
160 VolumeVariables volVars;
161 volVars.update(adaptedValues.u, *problem_, element, scv);
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);
169 adaptedValues.count = 1;
170 adaptedValues.wasLeaf =
true;
173 if (element.level() > 0)
175 auto& adaptedValuesFather = adaptionMap_[element.father()];
178 if(&adaptedValues != &adaptedValuesFather)
179 storeAdaptionValues(adaptedValues, adaptedValuesFather);
185 if(isBox && !element.isLeaf())
186 adaptedValues.u = ElementSolution(element, sol_, *gridGeometry_);
205 adaptionMap_.resize();
206 gridGeometry_->update();
207 sol_.resize(gridGeometry_->numDofs());
210 std::vector<Scalar> massCoeff;
211 std::vector<Scalar> associatedMass;
215 massCoeff.resize(gridGeometry_->numDofs(), 0.0);
216 associatedMass.resize(gridGeometry_->numDofs(), 0.0);
220 for (
const auto& element : elements(gridGeometry_->gridView().grid().leafGridView(), Dune::Partitions::interior))
222 if (!element.isNew())
224 const auto& adaptedValues = adaptionMap_[element];
226 auto fvGeometry =
localView(*gridGeometry_);
227 fvGeometry.bindElement(element);
230 auto elemSol = adaptedValues.u;
232 elemSol[0] /= adaptedValues.count;
234 const auto elementVolume = element.geometry().volume();
235 for (
const auto& scv : scvs(fvGeometry))
237 VolumeVariables volVars;
238 volVars.update(elemSol, *problem_, element, scv);
241 sol_[scv.dofIndex()] = elemSol[scv.localDofIndex()];
243 const auto dofIdxGlobal = scv.dofIndex();
248 if (!adaptedValues.wasLeaf)
250 if (formulation == p0s1)
252 sol_[dofIdxGlobal][saturationIdx] = adaptedValues.associatedMass[phase1Idx];
253 sol_[dofIdxGlobal][saturationIdx] /= elementVolume * volVars.density(phase1Idx) * volVars.porosity();
255 else if (formulation == p1s0)
257 sol_[dofIdxGlobal][saturationIdx] = adaptedValues.associatedMass[phase0Idx];
258 sol_[dofIdxGlobal][saturationIdx] /= elementVolume * volVars.density(phase0Idx) * volVars.porosity();
266 const auto scvVolume = scv.volume();
267 if (formulation == p0s1)
269 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase1Idx) * volVars.porosity();
270 associatedMass[dofIdxGlobal] += scvVolume / elementVolume * adaptedValues.associatedMass[phase1Idx];
272 else if (formulation == p1s0)
274 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase0Idx) * volVars.porosity();
275 associatedMass[dofIdxGlobal] += scvVolume / elementVolume * adaptedValues.associatedMass[phase0Idx];
283 assert(element.hasFather() &&
"new element does not have a father element!");
286 auto fatherElement = element.father();
287 while(fatherElement.isNew() && fatherElement.level() > 0)
288 fatherElement = fatherElement.father();
292 const auto& adaptedValuesFather = adaptionMap_[fatherElement];
295 Scalar massFather = 0.0;
296 if (formulation == p0s1)
297 massFather = adaptedValuesFather.associatedMass[phase1Idx];
298 else if (formulation == p1s0)
299 massFather = adaptedValuesFather.associatedMass[phase0Idx];
302 auto elemSolSon = adaptedValuesFather.u;
303 elemSolSon[0] /= adaptedValuesFather.count;
305 auto fvGeometry =
localView(*gridGeometry_);
306 fvGeometry.bindElement(element);
308 for (
const auto& scv : scvs(fvGeometry))
310 VolumeVariables volVars;
311 volVars.update(elemSolSon, *problem_, element, scv);
314 sol_[scv.dofIndex()] = elemSolSon[0];
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;
327 auto& adaptedValuesFather = adaptionMap_[fatherElement];
329 auto fvGeometry =
localView(*gridGeometry_);
330 fvGeometry.bindElement(element);
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,
338 adaptedValuesFather.u,
342 const auto fatherElementVolume = fatherGeometry.volume();
343 for (
const auto& scv : scvs(fvGeometry))
345 VolumeVariables volVars;
346 volVars.update(elemSolSon, *problem_, element, scv);
348 const auto dofIdxGlobal = scv.dofIndex();
349 const auto scvVolume = scv.volume();
350 if (formulation == p0s1)
352 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase1Idx) * volVars.porosity();
353 associatedMass[dofIdxGlobal] += scvVolume / fatherElementVolume * adaptedValuesFather.associatedMass[phase1Idx];
355 else if (formulation == p1s0)
357 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase0Idx) * volVars.porosity();
358 associatedMass[dofIdxGlobal] += scvVolume / fatherElementVolume * adaptedValuesFather.associatedMass[phase0Idx];
362 sol_[dofIdxGlobal] = elemSolSon[scv.localDofIndex()];
370 for(std::size_t dofIdxGlobal = 0; dofIdxGlobal < gridGeometry_->numDofs(); dofIdxGlobal++)
371 sol_[dofIdxGlobal][saturationIdx] = associatedMass[dofIdxGlobal] / massCoeff[dofIdxGlobal];
375 adaptionMap_.resize(
typename PersistentContainer::Value() );
376 adaptionMap_.shrinkToFit();
377 adaptionMap_.fill(
typename PersistentContainer::Value() );
404 static void storeAdaptionValues(AdaptedValues& adaptedValues,
405 AdaptedValues& adaptedValuesFather)
408 adaptedValuesFather.associatedMass += adaptedValues.associatedMass;
415 auto values = adaptedValues.u[0];
416 values /= adaptedValues.count;
417 adaptedValuesFather.u[0] += values;
420 adaptedValuesFather.count += 1;
423 adaptedValuesFather.wasLeaf =
false;
429 adaptedValuesFather.count = 1;
432 adaptedValuesFather.wasLeaf =
false;
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_;
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:95
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
make the local view function available whenever we use the grid geometry
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
Property tag Indices
Definition: porousmediumflow/sequential/properties.hh:59
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.