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
32
36
37#include <opm/material/fluidmatrixinteractions/NullMaterial.hpp>
38#include <opm/material/common/Means.hpp>
39#include <opm/material/densead/Evaluation.hpp>
40
41#include <opm/utility/CopyablePtr.hpp>
42
43#include <dune/common/fvector.hh>
44#include <dune/common/fmatrix.hh>
45
46namespace Opm {
53template<class TypeTag>
55 : public FvBaseProblem<TypeTag>
56 , public GetPropType<TypeTag, Properties::FluxModule>::FluxBaseProblem
57{
59 using ParentType = FvBaseProblem<TypeTag>;
60
61 using Implementation = GetPropType<TypeTag, Properties::Problem>;
68 using ThermalConductionLawParams = GetPropType<TypeTag, Properties::ThermalConductionLawParams>;
69 using MaterialLawParams = typename GetPropType<TypeTag, Properties::MaterialLaw>::Params;
70 using DirectionalMobilityPtr = Opm::Utility::CopyablePtr<DirectionalMobility<TypeTag, Evaluation>>;
71
72 enum { dimWorld = GridView::dimensionworld };
73 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
74 using DimVector = Dune::FieldVector<Scalar, dimWorld>;
75 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
77
78public:
83 : ParentType(simulator)
84 { init_(); }
85
89 static void registerParameters()
90 {
91 ParentType::registerParameters();
92
93 Parameters::registerParam<TypeTag, Properties::EnableGravity>
94 ("Use the gravity correction for the pressure gradients.");
95 }
96
106 template <class Context>
107 void intersectionIntrinsicPermeability(DimMatrix& result,
108 const Context& context,
109 unsigned intersectionIdx,
110 unsigned timeIdx) const
111 {
112 const auto& scvf = context.stencil(timeIdx).interiorFace(intersectionIdx);
113
114 const DimMatrix& K1 = asImp_().intrinsicPermeability(context, scvf.interiorIndex(), timeIdx);
115 const DimMatrix& K2 = asImp_().intrinsicPermeability(context, scvf.exteriorIndex(), timeIdx);
116
117 // entry-wise harmonic mean. this is almost certainly wrong if
118 // you have off-main diagonal entries in your permeabilities!
119 for (unsigned i = 0; i < dimWorld; ++i)
120 for (unsigned j = 0; j < dimWorld; ++j)
121 result[i][j] = harmonicMean(K1[i][j], K2[i][j]);
122 }
123
127 // \{
128
137 template <class Context>
138 const DimMatrix& intrinsicPermeability(const Context&,
139 unsigned,
140 unsigned) const
141 {
142 throw std::logic_error("Not implemented: Problem::intrinsicPermeability()");
143 }
144
154 template <class Context>
155 Scalar porosity(const Context&,
156 unsigned,
157 unsigned) const
158 {
159 throw std::logic_error("Not implemented: Problem::porosity()");
160 }
161
171 template <class Context>
172 const SolidEnergyLawParams&
173 solidEnergyParams(const Context&,
174 unsigned,
175 unsigned) const
176 {
177 throw std::logic_error("Not implemented: Problem::solidEnergyParams()");
178 }
179
189 template <class Context>
190 const ThermalConductionLawParams&
192 unsigned,
193 unsigned) const
194 {
195 throw std::logic_error("Not implemented: Problem::thermalConductionParams()");
196 }
197
206 template <class Context>
207 Scalar tortuosity(const Context&,
208 unsigned,
209 unsigned) const
210 {
211 throw std::logic_error("Not implemented: Problem::tortuosity()");
212 }
213
222 template <class Context>
223 Scalar dispersivity(const Context&,
224 unsigned,
225 unsigned) const
226 {
227 throw std::logic_error("Not implemented: Problem::dispersivity()");
228 }
229
243 template <class Context>
244 const MaterialLawParams &
245 materialLawParams(const Context&,
246 unsigned,
247 unsigned) const
248 {
249 static MaterialLawParams dummy;
250 return dummy;
251 }
252
253 template <class FluidState>
254 void updateRelperms([[maybe_unused]] std::array<Evaluation,numPhases>& mobility,
255 [[maybe_unused]] DirectionalMobilityPtr& dirMob,
256 [[maybe_unused]] FluidState& fluidState,
257 [[maybe_unused]] unsigned globalSpaceIdx) const
258 {}
259
268 template <class Context>
269 Scalar temperature(const Context&,
270 unsigned,
271 unsigned) const
272 { return asImp_().temperature(); }
273
281 Scalar temperature() const
282 { throw std::logic_error("Not implemented:temperature() method not implemented by the actual problem"); }
283
284
293 template <class Context>
294 const DimVector& gravity(const Context&,
295 unsigned,
296 unsigned) const
297 { return asImp_().gravity(); }
298
308 const DimVector& gravity() const
309 { return gravity_; }
310
317 {
318 using Toolbox = MathToolbox<Evaluation>;
319
320 unsigned numMarked = 0;
321 ElementContext elemCtx( this->simulator() );
322 auto gridView = this->simulator().vanguard().gridView();
323 auto& grid = this->simulator().vanguard().grid();
324 for (const auto& element : elements(gridView, Dune::Partitions::interior)) {
325 elemCtx.updateAll(element);
326
327 // HACK: this should better be part of an AdaptionCriterion class
328 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
329 Scalar minSat = 1e100 ;
330 Scalar maxSat = -1e100;
331 size_t nDofs = elemCtx.numDof(/*timeIdx=*/0);
332 for (unsigned dofIdx = 0; dofIdx < nDofs; ++dofIdx)
333 {
334 const auto& intQuant = elemCtx.intensiveQuantities( dofIdx, /*timeIdx=*/0 );
335 minSat = std::min(minSat,
336 Toolbox::value(intQuant.fluidState().saturation(phaseIdx)));
337 maxSat = std::max(maxSat,
338 Toolbox::value(intQuant.fluidState().saturation(phaseIdx)));
339 }
340
341 const Scalar indicator =
342 (maxSat - minSat)/(std::max<Scalar>(0.01, maxSat+minSat)/2);
343 if( indicator > 0.2 && element.level() < 2 ) {
344 grid.mark( 1, element );
345 ++ numMarked;
346 }
347 else if ( indicator < 0.025 ) {
348 grid.mark( -1, element );
349 ++ numMarked;
350 }
351 else
352 {
353 grid.mark( 0, element );
354 }
355 }
356 }
357
358 // get global sum so that every proc is on the same page
359 numMarked = this->simulator().vanguard().grid().comm().sum( numMarked );
360
361 return numMarked;
362 }
363
364 // \}
365
366protected:
377 DimMatrix toDimMatrix_(Scalar val) const
378 {
379 DimMatrix ret(0.0);
380 for (unsigned i = 0; i < DimMatrix::rows; ++i)
381 ret[i][i] = val;
382 return ret;
383 }
384
385 DimVector gravity_;
386
387private:
389 Implementation& asImp_()
390 { return *static_cast<Implementation *>(this); }
392 const Implementation& asImp_() const
393 { return *static_cast<const Implementation *>(this); }
394
395 void init_()
396 {
397 gravity_ = 0.0;
398 if (Parameters::get<TypeTag, Properties::EnableGravity>())
399 gravity_[dimWorld-1] = -9.81;
400 }
401};
402
403} // namespace Opm
404
405#endif
Base class for all problems which use a finite volume spatial discretization.
Definition: fvbaseproblem.hh:65
Simulator & simulator()
Returns Simulator object used by the simulation.
Definition: fvbaseproblem.hh:685
const GridView & gridView() const
The GridView which used by the problem.
Definition: fvbaseproblem.hh:653
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
Definition: multiphasebaseproblem.hh:57
void intersectionIntrinsicPermeability(DimMatrix &result, const Context &context, unsigned intersectionIdx, unsigned timeIdx) const
Returns the intrinsic permeability of an intersection.
Definition: multiphasebaseproblem.hh:107
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:191
DimVector gravity_
Definition: multiphasebaseproblem.hh:385
const MaterialLawParams & materialLawParams(const Context &, unsigned, unsigned) const
Returns the material law parameters within a control volume.
Definition: multiphasebaseproblem.hh:245
DimMatrix toDimMatrix_(Scalar val) const
Converts a Scalar value to an isotropic Tensor.
Definition: multiphasebaseproblem.hh:377
Scalar temperature() const
Returns the temperature for an isothermal problem.
Definition: multiphasebaseproblem.hh:281
Scalar tortuosity(const Context &, unsigned, unsigned) const
Define the tortuosity.
Definition: multiphasebaseproblem.hh:207
Scalar temperature(const Context &, unsigned, unsigned) const
Returns the temperature within a control volume.
Definition: multiphasebaseproblem.hh:269
const DimVector & gravity() const
Returns the acceleration due to gravity .
Definition: multiphasebaseproblem.hh:308
const DimVector & gravity(const Context &, unsigned, unsigned) const
Returns the acceleration due to gravity .
Definition: multiphasebaseproblem.hh:294
Scalar porosity(const Context &, unsigned, unsigned) const
Returns the porosity [] of the porous medium for a given control volume.
Definition: multiphasebaseproblem.hh:155
MultiPhaseBaseProblem(Simulator &simulator)
Definition: multiphasebaseproblem.hh:82
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:173
Scalar dispersivity(const Context &, unsigned, unsigned) const
Define the dispersivity.
Definition: multiphasebaseproblem.hh:223
void updateRelperms(std::array< Evaluation, numPhases > &mobility, DirectionalMobilityPtr &dirMob, FluidState &fluidState, unsigned globalSpaceIdx) const
Definition: multiphasebaseproblem.hh:254
const DimMatrix & intrinsicPermeability(const Context &, unsigned, unsigned) const
Returns the intrinsic permeability tensor at a given position.
Definition: multiphasebaseproblem.hh:138
unsigned markForGridAdaptation()
Mark grid cells for refinement or coarsening.
Definition: multiphasebaseproblem.hh:316
static void registerParameters()
Register all run-time parameters for the problem and the model.
Definition: multiphasebaseproblem.hh:89
This file contains definitions related to directional mobilities.
Declare the properties used by the infrastructure code of the finite volume discretizations.
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:242