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
83 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
84 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
86
87public:
91 explicit MultiPhaseBaseProblem(Simulator& simulator)
92 : ParentType(simulator)
93 { init_(); }
94
98 static void registerParameters()
99 {
100 ParentType::registerParameters();
101
102 Parameters::Register<Parameters::EnableGravity>
103 ("Use the gravity correction for the pressure gradients.");
104 }
105
115 template <class Context>
116 void intersectionIntrinsicPermeability(DimMatrix& result,
117 const Context& context,
118 unsigned intersectionIdx,
119 unsigned timeIdx) const
120 {
121 const auto& scvf = context.stencil(timeIdx).interiorFace(intersectionIdx);
122
123 const DimMatrix& K1 = asImp_().intrinsicPermeability(context, scvf.interiorIndex(), timeIdx);
124 const DimMatrix& K2 = asImp_().intrinsicPermeability(context, scvf.exteriorIndex(), timeIdx);
125
126 // entry-wise harmonic mean. this is almost certainly wrong if
127 // you have off-main diagonal entries in your permeabilities!
128 for (unsigned i = 0; i < dimWorld; ++i) {
129 for (unsigned j = 0; j < dimWorld; ++j) {
130 result[i][j] = harmonicMean(K1[i][j], K2[i][j]);
131 }
132 }
133 }
134
138 // \{
139
148 template <class Context>
149 const DimMatrix& intrinsicPermeability(const Context&,
150 unsigned,
151 unsigned) const
152 {
153 throw std::logic_error("Not implemented: Problem::intrinsicPermeability()");
154 }
155
165 template <class Context>
166 Scalar porosity(const Context&,
167 unsigned,
168 unsigned) const
169 {
170 throw std::logic_error("Not implemented: Problem::porosity()");
171 }
172
182 template <class Context>
183 const SolidEnergyLawParams&
184 solidEnergyParams(const Context&,
185 unsigned,
186 unsigned) const
187 {
188 throw std::logic_error("Not implemented: Problem::solidEnergyParams()");
189 }
190
200 template <class Context>
201 const ThermalConductionLawParams&
203 unsigned,
204 unsigned) const
205 {
206 throw std::logic_error("Not implemented: Problem::thermalConductionParams()");
207 }
208
217 template <class Context>
218 Scalar tortuosity(const Context&,
219 unsigned,
220 unsigned) const
221 {
222 throw std::logic_error("Not implemented: Problem::tortuosity()");
223 }
224
233 template <class Context>
234 Scalar dispersivity(const Context&,
235 unsigned,
236 unsigned) const
237 {
238 throw std::logic_error("Not implemented: Problem::dispersivity()");
239 }
240
254 template <class Context>
255 const MaterialLawParams &
256 materialLawParams(const Context&,
257 unsigned,
258 unsigned) const
259 {
260 static MaterialLawParams dummy;
261 return dummy;
262 }
263
264 template <class FluidState>
265 void updateRelperms([[maybe_unused]] std::array<Evaluation,numPhases>& mobility,
266 [[maybe_unused]] DirectionalMobilityPtr& dirMob,
267 [[maybe_unused]] FluidState& fluidState,
268 [[maybe_unused]] unsigned globalSpaceIdx) const
269 {}
270
279 template <class Context>
280 Scalar temperature(const Context&,
281 unsigned,
282 unsigned) const
283 { return asImp_().temperature(); }
284
292 Scalar temperature() const
293 { throw std::logic_error("Not implemented:temperature() method not implemented by the actual problem"); }
294
295
304 template <class Context>
305 const DimVector& gravity(const Context&,
306 unsigned,
307 unsigned) const
308 { return asImp_().gravity(); }
309
319 const DimVector& gravity() const
320 { return gravity_; }
321
328 {
329 using Toolbox = MathToolbox<Evaluation>;
330
331 unsigned numMarked = 0;
332 ElementContext elemCtx( this->simulator() );
333 auto gridView = this->simulator().vanguard().gridView();
334 auto& grid = this->simulator().vanguard().grid();
335 for (const auto& element : elements(gridView, Dune::Partitions::interior)) {
336 elemCtx.updateAll(element);
337
338 // HACK: this should better be part of an AdaptionCriterion class
339 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
340 Scalar minSat = 1e100 ;
341 Scalar maxSat = -1e100;
342 size_t nDofs = elemCtx.numDof(/*timeIdx=*/0);
343 for (unsigned dofIdx = 0; dofIdx < nDofs; ++dofIdx) {
344 const auto& intQuant = elemCtx.intensiveQuantities( dofIdx, /*timeIdx=*/0);
345 minSat = std::min(minSat,
346 Toolbox::value(intQuant.fluidState().saturation(phaseIdx)));
347 maxSat = std::max(maxSat,
348 Toolbox::value(intQuant.fluidState().saturation(phaseIdx)));
349 }
350
351 const Scalar indicator =
352 (maxSat - minSat) / (std::max<Scalar>(0.01, maxSat + minSat) / 2);
353 if (indicator > 0.2 && element.level() < 2) {
354 grid.mark(1, element);
355 ++numMarked;
356 }
357 else if (indicator < 0.025) {
358 grid.mark(-1, element);
359 ++numMarked;
360 }
361 else {
362 grid.mark(0, element);
363 }
364 }
365 }
366
367 // get global sum so that every proc is on the same page
368 return this->simulator().vanguard().grid().comm().sum(numMarked);
369 }
370
371 // \}
372
373protected:
384 DimMatrix toDimMatrix_(Scalar val) const
385 {
386 DimMatrix ret(0.0);
387 for (unsigned i = 0; i < DimMatrix::rows; ++i) {
388 ret[i][i] = val;
389 }
390 return ret;
391 }
392
393 DimVector gravity_;
394
395private:
397 Implementation& asImp_()
398 { return *static_cast<Implementation *>(this); }
399
401 const Implementation& asImp_() const
402 { return *static_cast<const Implementation *>(this); }
403
404 void init_()
405 {
406 gravity_ = 0.0;
407 if (Parameters::Get<Parameters::EnableGravity>()) {
408 gravity_[dimWorld-1] = -9.81;
409 }
410 }
411};
412
413} // namespace Opm
414
415#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:692
const GridView & gridView() const
The GridView which used by the problem.
Definition: fvbaseproblem.hh:660
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:116
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:202
DimVector gravity_
Definition: multiphasebaseproblem.hh:393
const MaterialLawParams & materialLawParams(const Context &, unsigned, unsigned) const
Returns the material law parameters within a control volume.
Definition: multiphasebaseproblem.hh:256
DimMatrix toDimMatrix_(Scalar val) const
Converts a Scalar value to an isotropic Tensor.
Definition: multiphasebaseproblem.hh:384
Scalar temperature() const
Returns the temperature for an isothermal problem.
Definition: multiphasebaseproblem.hh:292
Scalar tortuosity(const Context &, unsigned, unsigned) const
Define the tortuosity.
Definition: multiphasebaseproblem.hh:218
Scalar temperature(const Context &, unsigned, unsigned) const
Returns the temperature within a control volume.
Definition: multiphasebaseproblem.hh:280
const DimVector & gravity() const
Returns the acceleration due to gravity .
Definition: multiphasebaseproblem.hh:319
const DimVector & gravity(const Context &, unsigned, unsigned) const
Returns the acceleration due to gravity .
Definition: multiphasebaseproblem.hh:305
Scalar porosity(const Context &, unsigned, unsigned) const
Returns the porosity [] of the porous medium for a given control volume.
Definition: multiphasebaseproblem.hh:166
MultiPhaseBaseProblem(Simulator &simulator)
Definition: multiphasebaseproblem.hh:91
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:184
Scalar dispersivity(const Context &, unsigned, unsigned) const
Define the dispersivity.
Definition: multiphasebaseproblem.hh:234
void updateRelperms(std::array< Evaluation, numPhases > &mobility, DirectionalMobilityPtr &dirMob, FluidState &fluidState, unsigned globalSpaceIdx) const
Definition: multiphasebaseproblem.hh:265
const DimMatrix & intrinsicPermeability(const Context &, unsigned, unsigned) const
Returns the intrinsic permeability tensor at a given position.
Definition: multiphasebaseproblem.hh:149
unsigned markForGridAdaptation()
Mark grid cells for refinement or coarsening.
Definition: multiphasebaseproblem.hh:327
static void registerParameters()
Register all run-time parameters for the problem and the model.
Definition: multiphasebaseproblem.hh:98
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: blackoilbioeffectsmodules.hh:43
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