44 using GridView =
typename GridGeometry::GridView;
45 using FVElementGeometry =
typename GridGeometry::LocalView;
46 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
52 using Element =
typename Grid::template Codim<0>::Entity;
53 using ElementSolution = std::decay_t<decltype(elementSolution(std::declval<Element>(),
54 std::declval<SolutionVector>(),
55 std::declval<GridGeometry>()))>;
58 using Indices =
typename ModelTraits::Indices;
62 AdaptedValues() : associatedMass(0.0) {}
65 PrimaryVariables associatedMass;
69 using PersistentContainer = Dune::PersistentContainer<Grid, AdaptedValues>;
71 static constexpr int dim = Grid::dimension;
72 static constexpr int dimWorld = Grid::dimensionworld;
76 enum { saturationIdx = Indices::saturationIdx };
81 phase0Idx = FluidSystem::phase0Idx,
82 phase1Idx = FluidSystem::phase1Idx,
90 static constexpr auto formulation = ModelTraits::priVarFormulation();
93 static_assert(!FluidSystem::isCompressible(phase0Idx)
94 && !FluidSystem::isCompressible(phase1Idx),
95 "This adaption helper is only mass conservative for incompressible fluids!");
98 static_assert(formulation == p0s1 || formulation == p1s0,
"Chosen formulation not known to the TwoPGridDataTransfer");
110 std::shared_ptr<GridGeometry> gridGeometry,
111 std::shared_ptr<const GridVariables> gridVariables,
115 , gridGeometry_(gridGeometry)
116 , gridVariables_(gridVariables)
118 , adaptionMap_(gridGeometry->gridView().grid(), 0)
129 void store(
const Grid& grid)
override
131 adaptionMap_.resize();
133 for (
auto level = grid.maxLevel(); level >= 0; level--)
135 auto fvGeometry =
localView(*gridGeometry_);
136 for (
const auto& element : elements(grid.levelGridView(level)))
139 auto& adaptedValues = adaptionMap_[element];
142 if (element.isLeaf())
144 fvGeometry.bindElement(element);
147 adaptedValues.u = ElementSolution(element, sol_, *gridGeometry_);
150 for (
const auto& scv : scvs(fvGeometry))
152 VolumeVariables volVars;
153 volVars.update(adaptedValues.u, *problem_, element, scv);
155 const auto poreVolume = Extrusion::volume(fvGeometry, scv)*volVars.porosity();
156 adaptedValues.associatedMass[phase1Idx] += poreVolume * volVars.density(phase1Idx) * volVars.saturation(phase1Idx);
157 adaptedValues.associatedMass[phase0Idx] += poreVolume * volVars.density(phase0Idx) * volVars.saturation(phase0Idx);
161 adaptedValues.count = 1;
162 adaptedValues.wasLeaf =
true;
165 if (element.level() > 0)
167 auto& adaptedValuesFather = adaptionMap_[element.father()];
170 if(&adaptedValues != &adaptedValuesFather)
171 storeAdaptionValues(adaptedValues, adaptedValuesFather);
177 if(isBox && !element.isLeaf())
178 adaptedValues.u = ElementSolution(element, sol_, *gridGeometry_);
196 gridGeometry_->update(grid.leafGridView());
205 adaptionMap_.resize();
206 sol_.resize(gridGeometry_->numDofs());
209 std::vector<Scalar> massCoeff;
210 std::vector<Scalar> associatedMass;
214 massCoeff.resize(gridGeometry_->numDofs(), 0.0);
215 associatedMass.resize(gridGeometry_->numDofs(), 0.0);
219 auto fvGeometry =
localView(*gridGeometry_);
220 for (
const auto& element : elements(gridGeometry_->gridView().grid().leafGridView(),
Dune::Partitions::interior))
222 if (!element.isNew())
224 const auto& adaptedValues = adaptionMap_[element];
225 fvGeometry.bindElement(element);
228 auto elemSol = adaptedValues.u;
230 elemSol[0] /= adaptedValues.count;
232 const auto elementVolume =
volume(element.geometry(), Extrusion{});
233 for (
const auto& scv : scvs(fvGeometry))
235 VolumeVariables volVars;
236 volVars.update(elemSol, *problem_, element, scv);
239 sol_[scv.dofIndex()] = elemSol[scv.localDofIndex()];
241 const auto dofIdxGlobal = scv.dofIndex();
246 if (!adaptedValues.wasLeaf)
248 if (formulation == p0s1)
250 sol_[dofIdxGlobal][saturationIdx] = adaptedValues.associatedMass[phase1Idx];
251 sol_[dofIdxGlobal][saturationIdx] /= elementVolume * volVars.density(phase1Idx) * volVars.porosity();
253 else if (formulation == p1s0)
255 sol_[dofIdxGlobal][saturationIdx] = adaptedValues.associatedMass[phase0Idx];
256 sol_[dofIdxGlobal][saturationIdx] /= elementVolume * volVars.density(phase0Idx) * volVars.porosity();
264 const auto scvVolume = Extrusion::volume(fvGeometry, scv);
265 if (formulation == p0s1)
267 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase1Idx) * volVars.porosity();
268 associatedMass[dofIdxGlobal] += scvVolume / elementVolume * adaptedValues.associatedMass[phase1Idx];
270 else if (formulation == p1s0)
272 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase0Idx) * volVars.porosity();
273 associatedMass[dofIdxGlobal] += scvVolume / elementVolume * adaptedValues.associatedMass[phase0Idx];
281 assert(
element.hasFather() &&
"new element does not have a father element!");
284 auto fatherElement =
element.father();
285 while(fatherElement.isNew() && fatherElement.level() > 0)
286 fatherElement = fatherElement.father();
290 const auto& adaptedValuesFather = adaptionMap_[fatherElement];
293 Scalar massFather = 0.0;
294 if (formulation == p0s1)
295 massFather = adaptedValuesFather.associatedMass[phase1Idx];
296 else if (formulation == p1s0)
297 massFather = adaptedValuesFather.associatedMass[phase0Idx];
300 auto elemSolSon = adaptedValuesFather.u;
301 elemSolSon[0] /= adaptedValuesFather.count;
303 fvGeometry.bindElement(element);
305 for (
const auto& scv : scvs(fvGeometry))
307 VolumeVariables volVars;
308 volVars.update(elemSolSon, *problem_, element, scv);
311 sol_[scv.dofIndex()] = elemSolSon[0];
314 Scalar massCoeffSon = 0.0;
315 if (formulation == p0s1)
316 massCoeffSon = Extrusion::volume(fvGeometry, scv) * volVars.density(phase1Idx) * volVars.porosity();
317 else if (formulation == p1s0)
318 massCoeffSon = Extrusion::volume(fvGeometry, scv) * volVars.density(phase0Idx) * volVars.porosity();
319 sol_[scv.dofIndex()][saturationIdx] =
320 ( Extrusion::volume(fvGeometry, scv)/
volume(fatherElement.geometry(), Extrusion{})*massFather )/massCoeffSon;
325 auto& adaptedValuesFather = adaptionMap_[fatherElement];
327 fvGeometry.bindElement(element);
330 ElementSolution elemSolSon(element, sol_, *gridGeometry_);
331 const auto fatherGeometry = fatherElement.geometry();
332 for (
const auto& scv : scvs(fvGeometry))
333 elemSolSon[scv.localDofIndex()] =
evalSolution(fatherElement,
335 adaptedValuesFather.u,
339 const auto fatherElementVolume =
volume(fatherGeometry, Extrusion{});
340 for (
const auto& scv : scvs(fvGeometry))
342 VolumeVariables volVars;
343 volVars.update(elemSolSon, *problem_, element, scv);
345 const auto dofIdxGlobal = scv.dofIndex();
346 const auto scvVolume = Extrusion::volume(fvGeometry, scv);
347 if (formulation == p0s1)
349 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase1Idx) * volVars.porosity();
350 associatedMass[dofIdxGlobal] += scvVolume / fatherElementVolume * adaptedValuesFather.associatedMass[phase1Idx];
352 else if (formulation == p1s0)
354 massCoeff[dofIdxGlobal] += scvVolume * volVars.density(phase0Idx) * volVars.porosity();
355 associatedMass[dofIdxGlobal] += scvVolume / fatherElementVolume * adaptedValuesFather.associatedMass[phase0Idx];
359 sol_[dofIdxGlobal] = elemSolSon[scv.localDofIndex()];
367 for(std::size_t dofIdxGlobal = 0; dofIdxGlobal < gridGeometry_->numDofs(); dofIdxGlobal++)
368 sol_[dofIdxGlobal][saturationIdx] = associatedMass[dofIdxGlobal] / massCoeff[dofIdxGlobal];
372 adaptionMap_.resize(
typename PersistentContainer::Value() );
373 adaptionMap_.shrinkToFit();
374 adaptionMap_.fill(
typename PersistentContainer::Value() );
399 static void storeAdaptionValues(AdaptedValues& adaptedValues,
400 AdaptedValues& adaptedValuesFather)
403 adaptedValuesFather.associatedMass += adaptedValues.associatedMass;
410 auto values = adaptedValues.u[0];
411 values /= adaptedValues.count;
412 adaptedValuesFather.u[0] += values;
415 adaptedValuesFather.count += 1;
418 adaptedValuesFather.wasLeaf =
false;
424 adaptedValuesFather.count = 1;
427 adaptedValuesFather.wasLeaf =
false;
431 std::shared_ptr<const Problem> problem_;
432 std::shared_ptr<GridGeometry> gridGeometry_;
433 std::shared_ptr<const GridVariables> gridVariables_;
434 SolutionVector& sol_;
435 PersistentContainer adaptionMap_;