opm-simulators
Loading...
Searching...
No Matches
multiphasebasemodel.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_MODEL_HH
29#define EWOMS_MULTI_PHASE_BASE_MODEL_HH
30
31#include <opm/material/densead/Math.hpp>
32
33#include <opm/material/fluidmatrixinteractions/MaterialTraits.hpp>
34#include <opm/material/fluidmatrixinteractions/NullMaterial.hpp>
35
36#include <opm/material/thermal/NullSolidEnergyLaw.hpp>
37#include <opm/material/thermal/NullThermalConductionLaw.hpp>
38
44
46
49
50#include <cassert>
51#include <memory>
52#include <mutex>
53#include <tuple>
54
55namespace Opm {
56template <class TypeTag>
58}
59
60namespace Opm::Properties {
61
63// Create new type tags
64namespace TTag {
65
67
68} // end namespace TTag
69
71template<class TypeTag>
73{
74 using type = std::tuple<GetSplicePropType<TypeTag,
77};
78
82template<class TypeTag>
85
87template<class TypeTag>
89{ static constexpr int value = GetPropType<TypeTag, Properties::Indices>::numEq; };
90
92// Default is a fully implicit approach where numDerivaties = numEq;
93template<class TypeTag>
95{ static constexpr int value = GetPropType<TypeTag, Properties::Indices>::numEq; };
96
97
99template<class TypeTag>
102
104template<class TypeTag>
107
109template<class TypeTag>
112
114template<class TypeTag>
116{ using type = DarcyFluxModule<TypeTag>; };
117
121template<class TypeTag>
123{
124private:
127 using Traits = NullMaterialTraits<Scalar, FluidSystem::numPhases>;
128
129public:
130 using type = NullMaterial<Traits>;
131};
132
137template<class TypeTag>
140
143template<class TypeTag>
145{ using type = NullSolidEnergyLaw<GetPropType<TypeTag, Properties::Scalar>>; };
146
149template<class TypeTag>
152
154template<class TypeTag>
156{ using type = NullThermalConductionLaw<GetPropType<TypeTag, Properties::Scalar>>; };
157
160template<class TypeTag>
163
164} // namespace Opm::Properties
165
166namespace Opm {
167
173template <class TypeTag>
174class MultiPhaseBaseModel : public GetPropType<TypeTag, Properties::Discretization>
175{
177 using Implementation = GetPropType<TypeTag, Properties::Problem>;
183
184 using ElementIterator = typename GridView::template Codim<0>::Iterator;
185 using Element = typename GridView::template Codim<0>::Entity;
186
188
189public:
190 explicit MultiPhaseBaseModel(Simulator& simulator)
191 : ParentType(simulator)
192 {}
193
197 static void registerParameters()
198 {
199 ParentType::registerParameters();
200
201 // register runtime parameters of the VTK output modules
204 }
205
211 bool phaseIsConsidered(unsigned) const
212 { return true; }
213
221 void globalPhaseStorage(EqVector& storage, unsigned phaseIdx)
222 {
223 assert(phaseIdx < numPhases);
224
225 storage = 0;
226
227 ThreadedEntityIterator<GridView, /*codim=*/0> threadedElemIt(this->gridView());
228 std::mutex mutex;
229#ifdef _OPENMP
230#pragma omp parallel
231#endif
232 {
233 // Attention: the variables below are thread specific and thus cannot be
234 // moved in front of the #pragma!
235 const unsigned threadId = ThreadManager::threadId();
236 ElementContext elemCtx(this->simulator_);
237 ElementIterator elemIt = threadedElemIt.beginParallel();
238 EqVector tmp;
239
240 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
241 const Element& elem = *elemIt;
242 if (elem.partitionType() != Dune::InteriorEntity) {
243 continue; // ignore ghost and overlap elements
244 }
245
246 elemCtx.updateStencil(elem);
247 elemCtx.updateIntensiveQuantities(/*timeIdx=*/0);
248
249 const auto& stencil = elemCtx.stencil(/*timeIdx=*/0);
250
251 for (unsigned dofIdx = 0; dofIdx < elemCtx.numDof(/*timeIdx=*/0); ++dofIdx) {
252 const auto& scv = stencil.subControlVolume(dofIdx);
253 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
254
255 tmp = 0;
256 this->localResidual(threadId).addPhaseStorage(tmp,
257 elemCtx,
258 dofIdx,
259 /*timeIdx=*/0,
260 phaseIdx);
261 tmp *= scv.volume() * intQuants.extrusionFactor();
262
263 mutex.lock();
264 storage += tmp;
265 mutex.unlock();
266 }
267 }
268 }
269
270 storage = this->gridView_.comm().sum(storage);
271 }
272
273 void registerOutputModules_()
274 {
275 ParentType::registerOutputModules_();
276
277 // add the VTK output modules which make sense for all multi-phase models
278 this->addOutputModule(std::make_unique<VtkMultiPhaseModule<TypeTag>>(this->simulator_));
279 this->addOutputModule(std::make_unique<VtkTemperatureModule<TypeTag>>(this->simulator_));
280 }
281
282private:
283 const Implementation& asImp_() const
284 { return *static_cast<const Implementation*>(this); }
285};
286
287} // namespace Opm
288
289#endif
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
Definition multiphasebasemodel.hh:175
void globalPhaseStorage(EqVector &storage, unsigned phaseIdx)
Compute the total storage inside one phase of all conservation quantities.
Definition multiphasebasemodel.hh:221
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition multiphasebasemodel.hh:197
bool phaseIsConsidered(unsigned) const
Returns true iff a fluid phase is used by the model.
Definition multiphasebasemodel.hh:211
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
Definition multiphasebaseproblem.hh:65
static unsigned threadId()
Return the index of the current OpenMP thread.
Definition threadmanager.cpp:84
Provides an STL-iterator like interface to iterate over the enties of a GridView in OpenMP threaded a...
Definition threadedentityiterator.hh:42
VTK output module for quantities which make sense for all models which deal with multiple fluid phase...
Definition vtkmultiphasemodule.hpp:73
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition vtkmultiphasemodule.hpp:110
VTK output module for the temperature in which assume thermal equilibrium.
Definition vtktemperaturemodule.hpp:51
static void registerParameters()
Register all run-time parameters for the Vtk output module.
Definition vtktemperaturemodule.hpp:75
This file contains the necessary classes to calculate the velocity out of a pressure potential gradie...
This class calculates the pressure potential gradients and the filter velocities for multi-phase flow...
Defines the common parameters for the porous medium multi-phase models.
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
Defines the common properties required by the porous medium multi-phase models.
The generic type tag for problems using the immiscible multi-phase model.
Definition blackoilmodel.hh:82
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
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
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
Specifies a flux module which uses the Darcy relation.
Definition darcyfluxmodule.hh:67
The type of the base class for all problems which use this model.
Definition fvbaseproperties.hh:84
Specifies the relation used for velocity.
Definition multiphasebaseproperties.hh:83
The context material law (extracted from the spatial parameters).
Definition multiphasebaseproperties.hh:59
The material law which ought to be used (extracted from the spatial parameters).
Definition multiphasebaseproperties.hh:55
Number of chemical species in the system.
Definition multiphasebaseproperties.hh:47
Number of derivatives in the system of PDEs.
Definition basicproperties.hh:84
Number of equations in the system of PDEs.
Definition basicproperties.hh:80
Number of fluid phases in the system.
Definition multiphasebaseproperties.hh:43
The parameters of the material law for energy storage of the solid.
Definition multiphasebaseproperties.hh:67
The material law for the energy stored in the solid matrix.
Definition multiphasebaseproperties.hh:63
The splice to be used for the spatial discretization.
Definition multiphasebaseproperties.hh:39
Definition propertysystem.hh:42
Definition multiphasebasemodel.hh:66
Definition vcfvproperties.hh:41
The parameters of the material law for thermal conduction.
Definition multiphasebaseproperties.hh:75
The material law for thermal conduction.
Definition multiphasebaseproperties.hh:71
The base class for the vertex centered finite volume discretization scheme.
VTK output module for quantities which make sense for all models which deal with multiple fluid phase...
VTK output module for the temperature in which assume thermal equilibrium.