multiphasebaseproblem.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_MULTI_PHASE_BASE_PROBLEM_HH
29#define EWOMS_MULTI_PHASE_BASE_PROBLEM_HH
30
31#include <dune/common/fmatrix.hh>
32#include <dune/common/fvector.hh>
33
34#include <dune/grid/common/partitionset.hh>
35
36#include <opm/material/fluidmatrixinteractions/NullMaterial.hpp>
37#include <opm/material/common/Means.hpp>
38#include <opm/material/densead/Evaluation.hpp>
39
43
46
47#include <opm/utility/CopyablePtr.hpp>
48
49#include <algorithm>
50#include <array>
51#include <stdexcept>
52
53namespace Opm {
54
61template<class TypeTag>
63 : public FvBaseProblem<TypeTag>
64 , public GetPropType<TypeTag, Properties::FluxModule>::FluxBaseProblem
65{
67 using ParentType = FvBaseProblem<TypeTag>;
68
69 using Implementation = GetPropType<TypeTag, Properties::Problem>;
76 using ThermalConductionLawParams = GetPropType<TypeTag, Properties::ThermalConductionLawParams>;
77 using MaterialLawParams = typename GetPropType<TypeTag, Properties::MaterialLaw>::Params;
78 using DirectionalMobilityPtr = Utility::CopyablePtr<DirectionalMobility<TypeTag>>;
79
80 enum { dimWorld = GridView::dimensionworld };
81 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
82 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
83 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
85
86public:
90 explicit MultiPhaseBaseProblem(Simulator& simulator)
91 : ParentType(simulator)
92 { init_(); }
93
97 static void registerParameters()
98 {
99 ParentType::registerParameters();
100
101 Parameters::Register<Parameters::EnableGravity>
102 ("Use the gravity correction for the pressure gradients.");
103 }
104
114 template <class Context>
115 void intersectionIntrinsicPermeability(DimMatrix& result,
116 const Context& context,
117 unsigned intersectionIdx,
118 unsigned timeIdx) const
119 {
120 const auto& scvf = context.stencil(timeIdx).interiorFace(intersectionIdx);
121
122 const DimMatrix& K1 = asImp_().intrinsicPermeability(context, scvf.interiorIndex(), timeIdx);
123 const DimMatrix& K2 = asImp_().intrinsicPermeability(context, scvf.exteriorIndex(), timeIdx);
124
125 // entry-wise harmonic mean. this is almost certainly wrong if
126 // you have off-main diagonal entries in your permeabilities!
127 for (unsigned i = 0; i < dimWorld; ++i) {
128 for (unsigned j = 0; j < dimWorld; ++j) {
129 result[i][j] = harmonicMean(K1[i][j], K2[i][j]);
130 }
131 }
132 }
133
137 // \{
138
147 template <class Context>
148 const DimMatrix& intrinsicPermeability(const Context&,
149 unsigned,
150 unsigned) const
151 {
152 throw std::logic_error("Not implemented: Problem::intrinsicPermeability()");
153 }
154
164 template <class Context>
165 Scalar porosity(const Context&,
166 unsigned,
167 unsigned) const
168 {
169 throw std::logic_error("Not implemented: Problem::porosity()");
170 }
171
181 template <class Context>
182 const SolidEnergyLawParams&
183 solidEnergyParams(const Context&,
184 unsigned,
185 unsigned) const
186 {
187 throw std::logic_error("Not implemented: Problem::solidEnergyParams()");
188 }
189
199 template <class Context>
200 const ThermalConductionLawParams&
202 unsigned,
203 unsigned) const
204 {
205 throw std::logic_error("Not implemented: Problem::thermalConductionParams()");
206 }
207
216 template <class Context>
217 Scalar tortuosity(const Context&,
218 unsigned,
219 unsigned) const
220 {
221 throw std::logic_error("Not implemented: Problem::tortuosity()");
222 }
223
232 template <class Context>
233 Scalar dispersivity(const Context&,
234 unsigned,
235 unsigned) const
236 {
237 throw std::logic_error("Not implemented: Problem::dispersivity()");
238 }
239
253 template <class Context>
254 const MaterialLawParams &
255 materialLawParams(const Context&,
256 unsigned,
257 unsigned) const
258 {
259 static MaterialLawParams dummy;
260 return dummy;
261 }
262
263 template <class FluidState>
264 void updateRelperms([[maybe_unused]] std::array<Evaluation,numPhases>& mobility,
265 [[maybe_unused]] DirectionalMobilityPtr& dirMob,
266 [[maybe_unused]] FluidState& fluidState,
267 [[maybe_unused]] unsigned globalSpaceIdx) const
268 {}
269
278 template <class Context>
279 Scalar temperature(const Context&,
280 unsigned,
281 unsigned) const
282 { return asImp_().temperature(); }
283
291 Scalar temperature() const
292 { throw std::logic_error("Not implemented:temperature() method not implemented by the actual problem"); }
293
294
303 template <class Context>
304 const DimVector& gravity(const Context&,
305 unsigned,
306 unsigned) const
307 { return asImp_().gravity(); }
308
318 const DimVector& gravity() const
319 { return gravity_; }
320
327 {
328 using Toolbox = MathToolbox<Evaluation>;
329
330 unsigned numMarked = 0;
331 ElementContext elemCtx( this->simulator() );
332 auto gridView = this->simulator().vanguard().gridView();
333 auto& grid = this->simulator().vanguard().grid();
334 for (const auto& element : elements(gridView, Dune::Partitions::interior)) {
335 elemCtx.updateAll(element);
336
337 // HACK: this should better be part of an AdaptionCriterion class
338 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
339 Scalar minSat = 1e100 ;
340 Scalar maxSat = -1e100;
341 size_t nDofs = elemCtx.numDof(/*timeIdx=*/0);
342 for (unsigned dofIdx = 0; dofIdx < nDofs; ++dofIdx) {
343 const auto& intQuant = elemCtx.intensiveQuantities( dofIdx, /*timeIdx=*/0);
344 minSat = std::min(minSat,
345 Toolbox::value(intQuant.fluidState().saturation(phaseIdx)));
346 maxSat = std::max(maxSat,
347 Toolbox::value(intQuant.fluidState().saturation(phaseIdx)));
348 }
349
350 const Scalar indicator =
351 (maxSat - minSat) / (std::max<Scalar>(0.01, maxSat + minSat) / 2);
352 if (indicator > 0.2 && element.level() < 2) {
353 grid.mark(1, element);
354 ++numMarked;
355 }
356 else if (indicator < 0.025) {
357 grid.mark(-1, element);
358 ++numMarked;
359 }
360 else {
361 grid.mark(0, element);
362 }
363 }
364 }
365
366 // get global sum so that every proc is on the same page
367 return this->simulator().vanguard().grid().comm().sum(numMarked);
368 }
369
370 // \}
371
372protected:
383 DimMatrix toDimMatrix_(Scalar val) const
384 {
385 DimMatrix ret(0.0);
386 for (unsigned i = 0; i < DimMatrix::rows; ++i) {
387 ret[i][i] = val;
388 }
389 return ret;
390 }
391
392 DimVector gravity_;
393
394private:
396 Implementation& asImp_()
397 { return *static_cast<Implementation *>(this); }
398
400 const Implementation& asImp_() const
401 { return *static_cast<const Implementation *>(this); }
402
403 void init_()
404 {
405 gravity_ = 0.0;
406 if (Parameters::Get<Parameters::EnableGravity>()) {
407 gravity_[dimWorld-1] = -9.81;
408 }
409 }
410};
411
412} // namespace Opm
413
414#endif
Base class for all problems which use a finite volume spatial discretization.
Definition: fvbaseproblem.hh:69
Simulator & simulator()
Returns Simulator object used by the simulation.
Definition: fvbaseproblem.hh:677
const GridView & gridView() const
The GridView which used by the problem.
Definition: fvbaseproblem.hh:645
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
Definition: multiphasebaseproblem.hh:65
void intersectionIntrinsicPermeability(DimMatrix &result, const Context &context, unsigned intersectionIdx, unsigned timeIdx) const
Returns the intrinsic permeability of an intersection.
Definition: multiphasebaseproblem.hh:115
const ThermalConductionLawParams & thermalConductionParams(const Context &, unsigned, unsigned) const
Returns the parameter object for the thermal conductivity law in a sub-control volume.
Definition: multiphasebaseproblem.hh:201
DimVector gravity_
Definition: multiphasebaseproblem.hh:392
const MaterialLawParams & materialLawParams(const Context &, unsigned, unsigned) const
Returns the material law parameters within a control volume.
Definition: multiphasebaseproblem.hh:255
DimMatrix toDimMatrix_(Scalar val) const
Converts a Scalar value to an isotropic Tensor.
Definition: multiphasebaseproblem.hh:383
Scalar temperature() const
Returns the temperature for an isothermal problem.
Definition: multiphasebaseproblem.hh:291
Scalar tortuosity(const Context &, unsigned, unsigned) const
Define the tortuosity.
Definition: multiphasebaseproblem.hh:217
Scalar temperature(const Context &, unsigned, unsigned) const
Returns the temperature within a control volume.
Definition: multiphasebaseproblem.hh:279
const DimVector & gravity() const
Returns the acceleration due to gravity .
Definition: multiphasebaseproblem.hh:318
const DimVector & gravity(const Context &, unsigned, unsigned) const
Returns the acceleration due to gravity .
Definition: multiphasebaseproblem.hh:304
Scalar porosity(const Context &, unsigned, unsigned) const
Returns the porosity [] of the porous medium for a given control volume.
Definition: multiphasebaseproblem.hh:165
MultiPhaseBaseProblem(Simulator &simulator)
Definition: multiphasebaseproblem.hh:90
const SolidEnergyLawParams & solidEnergyParams(const Context &, unsigned, unsigned) const
Returns the parameter object for the energy storage law of the solid in a sub-control volume.
Definition: multiphasebaseproblem.hh:183
Scalar dispersivity(const Context &, unsigned, unsigned) const
Define the dispersivity.
Definition: multiphasebaseproblem.hh:233
void updateRelperms(std::array< Evaluation, numPhases > &mobility, DirectionalMobilityPtr &dirMob, FluidState &fluidState, unsigned globalSpaceIdx) const
Definition: multiphasebaseproblem.hh:264
const DimMatrix & intrinsicPermeability(const Context &, unsigned, unsigned) const
Returns the intrinsic permeability tensor at a given position.
Definition: multiphasebaseproblem.hh:148
unsigned markForGridAdaptation()
Mark grid cells for refinement or coarsening.
Definition: multiphasebaseproblem.hh:326
static void registerParameters()
Register all run-time parameters for the problem and the model.
Definition: multiphasebaseproblem.hh:97
This file contains definitions related to directional mobilities.
Declare the properties used by the infrastructure code of the finite volume discretizations.
Defines the common parameters for the porous medium multi-phase models.
Defines the common properties required by the porous medium multi-phase models.
Definition: blackoilboundaryratevector.hh:39
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:233