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/fvector.hh>
32#include <dune/common/fmatrix.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
49namespace Opm {
56template<class TypeTag>
58 : public FvBaseProblem<TypeTag>
59 , public GetPropType<TypeTag, Properties::FluxModule>::FluxBaseProblem
60{
62 using ParentType = FvBaseProblem<TypeTag>;
63
64 using Implementation = GetPropType<TypeTag, Properties::Problem>;
71 using ThermalConductionLawParams = GetPropType<TypeTag, Properties::ThermalConductionLawParams>;
72 using MaterialLawParams = typename GetPropType<TypeTag, Properties::MaterialLaw>::Params;
73 using DirectionalMobilityPtr = Opm::Utility::CopyablePtr<DirectionalMobility<TypeTag, Evaluation>>;
74
75 enum { dimWorld = GridView::dimensionworld };
76 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
77 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
78 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
80
81public:
86 : ParentType(simulator)
87 { init_(); }
88
92 static void registerParameters()
93 {
94 ParentType::registerParameters();
95
96 Parameters::Register<Parameters::EnableGravity>
97 ("Use the gravity correction for the pressure gradients.");
98 }
99
109 template <class Context>
110 void intersectionIntrinsicPermeability(DimMatrix& result,
111 const Context& context,
112 unsigned intersectionIdx,
113 unsigned timeIdx) const
114 {
115 const auto& scvf = context.stencil(timeIdx).interiorFace(intersectionIdx);
116
117 const DimMatrix& K1 = asImp_().intrinsicPermeability(context, scvf.interiorIndex(), timeIdx);
118 const DimMatrix& K2 = asImp_().intrinsicPermeability(context, scvf.exteriorIndex(), timeIdx);
119
120 // entry-wise harmonic mean. this is almost certainly wrong if
121 // you have off-main diagonal entries in your permeabilities!
122 for (unsigned i = 0; i < dimWorld; ++i)
123 for (unsigned j = 0; j < dimWorld; ++j)
124 result[i][j] = harmonicMean(K1[i][j], K2[i][j]);
125 }
126
130 // \{
131
140 template <class Context>
141 const DimMatrix& intrinsicPermeability(const Context&,
142 unsigned,
143 unsigned) const
144 {
145 throw std::logic_error("Not implemented: Problem::intrinsicPermeability()");
146 }
147
157 template <class Context>
158 Scalar porosity(const Context&,
159 unsigned,
160 unsigned) const
161 {
162 throw std::logic_error("Not implemented: Problem::porosity()");
163 }
164
174 template <class Context>
175 const SolidEnergyLawParams&
176 solidEnergyParams(const Context&,
177 unsigned,
178 unsigned) const
179 {
180 throw std::logic_error("Not implemented: Problem::solidEnergyParams()");
181 }
182
192 template <class Context>
193 const ThermalConductionLawParams&
195 unsigned,
196 unsigned) const
197 {
198 throw std::logic_error("Not implemented: Problem::thermalConductionParams()");
199 }
200
209 template <class Context>
210 Scalar tortuosity(const Context&,
211 unsigned,
212 unsigned) const
213 {
214 throw std::logic_error("Not implemented: Problem::tortuosity()");
215 }
216
225 template <class Context>
226 Scalar dispersivity(const Context&,
227 unsigned,
228 unsigned) const
229 {
230 throw std::logic_error("Not implemented: Problem::dispersivity()");
231 }
232
246 template <class Context>
247 const MaterialLawParams &
248 materialLawParams(const Context&,
249 unsigned,
250 unsigned) const
251 {
252 static MaterialLawParams dummy;
253 return dummy;
254 }
255
256 template <class FluidState>
257 void updateRelperms([[maybe_unused]] std::array<Evaluation,numPhases>& mobility,
258 [[maybe_unused]] DirectionalMobilityPtr& dirMob,
259 [[maybe_unused]] FluidState& fluidState,
260 [[maybe_unused]] unsigned globalSpaceIdx) const
261 {}
262
271 template <class Context>
272 Scalar temperature(const Context&,
273 unsigned,
274 unsigned) const
275 { return asImp_().temperature(); }
276
284 Scalar temperature() const
285 { throw std::logic_error("Not implemented:temperature() method not implemented by the actual problem"); }
286
287
296 template <class Context>
297 const DimVector& gravity(const Context&,
298 unsigned,
299 unsigned) const
300 { return asImp_().gravity(); }
301
311 const DimVector& gravity() const
312 { return gravity_; }
313
320 {
321 using Toolbox = MathToolbox<Evaluation>;
322
323 unsigned numMarked = 0;
324 ElementContext elemCtx( this->simulator() );
325 auto gridView = this->simulator().vanguard().gridView();
326 auto& grid = this->simulator().vanguard().grid();
327 for (const auto& element : elements(gridView, Dune::Partitions::interior)) {
328 elemCtx.updateAll(element);
329
330 // HACK: this should better be part of an AdaptionCriterion class
331 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
332 Scalar minSat = 1e100 ;
333 Scalar maxSat = -1e100;
334 size_t nDofs = elemCtx.numDof(/*timeIdx=*/0);
335 for (unsigned dofIdx = 0; dofIdx < nDofs; ++dofIdx)
336 {
337 const auto& intQuant = elemCtx.intensiveQuantities( dofIdx, /*timeIdx=*/0 );
338 minSat = std::min(minSat,
339 Toolbox::value(intQuant.fluidState().saturation(phaseIdx)));
340 maxSat = std::max(maxSat,
341 Toolbox::value(intQuant.fluidState().saturation(phaseIdx)));
342 }
343
344 const Scalar indicator =
345 (maxSat - minSat)/(std::max<Scalar>(0.01, maxSat+minSat)/2);
346 if( indicator > 0.2 && element.level() < 2 ) {
347 grid.mark( 1, element );
348 ++ numMarked;
349 }
350 else if ( indicator < 0.025 ) {
351 grid.mark( -1, element );
352 ++ numMarked;
353 }
354 else
355 {
356 grid.mark( 0, element );
357 }
358 }
359 }
360
361 // get global sum so that every proc is on the same page
362 numMarked = this->simulator().vanguard().grid().comm().sum( numMarked );
363
364 return numMarked;
365 }
366
367 // \}
368
369protected:
380 DimMatrix toDimMatrix_(Scalar val) const
381 {
382 DimMatrix ret(0.0);
383 for (unsigned i = 0; i < DimMatrix::rows; ++i)
384 ret[i][i] = val;
385 return ret;
386 }
387
388 DimVector gravity_;
389
390private:
392 Implementation& asImp_()
393 { return *static_cast<Implementation *>(this); }
395 const Implementation& asImp_() const
396 { return *static_cast<const Implementation *>(this); }
397
398 void init_()
399 {
400 gravity_ = 0.0;
401 if (Parameters::Get<Parameters::EnableGravity>()) {
402 gravity_[dimWorld-1] = -9.81;
403 }
404 }
405};
406
407} // namespace Opm
408
409#endif
Base class for all problems which use a finite volume spatial discretization.
Definition: fvbaseproblem.hh:66
Simulator & simulator()
Returns Simulator object used by the simulation.
Definition: fvbaseproblem.hh:686
const GridView & gridView() const
The GridView which used by the problem.
Definition: fvbaseproblem.hh:654
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
Definition: multiphasebaseproblem.hh:60
void intersectionIntrinsicPermeability(DimMatrix &result, const Context &context, unsigned intersectionIdx, unsigned timeIdx) const
Returns the intrinsic permeability of an intersection.
Definition: multiphasebaseproblem.hh:110
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:194
DimVector gravity_
Definition: multiphasebaseproblem.hh:388
const MaterialLawParams & materialLawParams(const Context &, unsigned, unsigned) const
Returns the material law parameters within a control volume.
Definition: multiphasebaseproblem.hh:248
DimMatrix toDimMatrix_(Scalar val) const
Converts a Scalar value to an isotropic Tensor.
Definition: multiphasebaseproblem.hh:380
Scalar temperature() const
Returns the temperature for an isothermal problem.
Definition: multiphasebaseproblem.hh:284
Scalar tortuosity(const Context &, unsigned, unsigned) const
Define the tortuosity.
Definition: multiphasebaseproblem.hh:210
Scalar temperature(const Context &, unsigned, unsigned) const
Returns the temperature within a control volume.
Definition: multiphasebaseproblem.hh:272
const DimVector & gravity() const
Returns the acceleration due to gravity .
Definition: multiphasebaseproblem.hh:311
const DimVector & gravity(const Context &, unsigned, unsigned) const
Returns the acceleration due to gravity .
Definition: multiphasebaseproblem.hh:297
Scalar porosity(const Context &, unsigned, unsigned) const
Returns the porosity [] of the porous medium for a given control volume.
Definition: multiphasebaseproblem.hh:158
MultiPhaseBaseProblem(Simulator &simulator)
Definition: multiphasebaseproblem.hh:85
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:176
Scalar dispersivity(const Context &, unsigned, unsigned) const
Define the dispersivity.
Definition: multiphasebaseproblem.hh:226
void updateRelperms(std::array< Evaluation, numPhases > &mobility, DirectionalMobilityPtr &dirMob, FluidState &fluidState, unsigned globalSpaceIdx) const
Definition: multiphasebaseproblem.hh:257
const DimMatrix & intrinsicPermeability(const Context &, unsigned, unsigned) const
Returns the intrinsic permeability tensor at a given position.
Definition: multiphasebaseproblem.hh:141
unsigned markForGridAdaptation()
Mark grid cells for refinement or coarsening.
Definition: multiphasebaseproblem.hh:319
static void registerParameters()
Register all run-time parameters for the problem and the model.
Definition: multiphasebaseproblem.hh:92
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:37
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:235