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>
57class MultiPhaseBaseModel;
58}
59
60namespace Opm::Properties {
61
63// Create new type tags
64namespace TTag {
65
67
68} // end namespace TTag
69
71template<class TypeTag>
72struct Splices<TypeTag, TTag::MultiPhaseBaseModel>
73{
74 using type = std::tuple<GetSplicePropType<TypeTag,
77};
78
82template<class TypeTag>
85
87template<class TypeTag>
88struct NumEq<TypeTag, TTag::MultiPhaseBaseModel>
89{ static constexpr int value = GetPropType<TypeTag, Properties::Indices>::numEq; };
90
92template<class TypeTag>
93struct NumPhases<TypeTag, TTag::MultiPhaseBaseModel>
94{ static constexpr int value = GetPropType<TypeTag, Properties::FluidSystem>::numPhases; };
95
97template<class TypeTag>
100
102template<class TypeTag>
105
107template<class TypeTag>
108struct FluxModule<TypeTag, TTag::MultiPhaseBaseModel>
109{ using type = DarcyFluxModule<TypeTag>; };
110
114template<class TypeTag>
115struct MaterialLaw<TypeTag, TTag::MultiPhaseBaseModel>
116{
117private:
121
122public:
123 using type = NullMaterial<Traits>;
124};
125
130template<class TypeTag>
133
136template<class TypeTag>
139
142template<class TypeTag>
145
147template<class TypeTag>
150
153template<class TypeTag>
156
157} // namespace Opm::Properties
158
159namespace Opm {
160
166template <class TypeTag>
167class MultiPhaseBaseModel : public GetPropType<TypeTag, Properties::Discretization>
168{
170 using Implementation = GetPropType<TypeTag, Properties::Model>;
179
180 using ElementIterator = typename GridView::template Codim<0>::Iterator;
181 using Element = typename GridView::template Codim<0>::Entity;
182
184 enum { numComponents = FluidSystem::numComponents };
185
186public:
187 explicit MultiPhaseBaseModel(Simulator& simulator)
188 : ParentType(simulator)
189 {}
190
194 static void registerParameters()
195 {
196 ParentType::registerParameters();
197
198 // register runtime parameters of the VTK output modules
201 }
202
208 bool phaseIsConsidered(unsigned) const
209 { return true; }
210
218 void globalPhaseStorage(EqVector& storage, unsigned phaseIdx)
219 {
220 assert(phaseIdx < numPhases);
221
222 storage = 0;
223
224 ThreadedEntityIterator<GridView, /*codim=*/0> threadedElemIt(this->gridView());
225 std::mutex mutex;
226#ifdef _OPENMP
227#pragma omp parallel
228#endif
229 {
230 // Attention: the variables below are thread specific and thus cannot be
231 // moved in front of the #pragma!
232 const unsigned threadId = ThreadManager::threadId();
233 ElementContext elemCtx(this->simulator_);
234 ElementIterator elemIt = threadedElemIt.beginParallel();
235 EqVector tmp;
236
237 for (; !threadedElemIt.isFinished(elemIt); elemIt = threadedElemIt.increment()) {
238 const Element& elem = *elemIt;
239 if (elem.partitionType() != Dune::InteriorEntity) {
240 continue; // ignore ghost and overlap elements
241 }
242
243 elemCtx.updateStencil(elem);
244 elemCtx.updateIntensiveQuantities(/*timeIdx=*/0);
245
246 const auto& stencil = elemCtx.stencil(/*timeIdx=*/0);
247
248 for (unsigned dofIdx = 0; dofIdx < elemCtx.numDof(/*timeIdx=*/0); ++dofIdx) {
249 const auto& scv = stencil.subControlVolume(dofIdx);
250 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, /*timeIdx=*/0);
251
252 tmp = 0;
253 this->localResidual(threadId).addPhaseStorage(tmp,
254 elemCtx,
255 dofIdx,
256 /*timeIdx=*/0,
257 phaseIdx);
258 tmp *= scv.volume() * intQuants.extrusionFactor();
259
260 mutex.lock();
261 storage += tmp;
262 mutex.unlock();
263 }
264 }
265 }
266
267 storage = this->gridView_.comm().sum(storage);
268 }
269
270 void registerOutputModules_()
271 {
272 ParentType::registerOutputModules_();
273
274 // add the VTK output modules which make sense for all multi-phase models
275 this->addOutputModule(std::make_unique<VtkMultiPhaseModule<TypeTag>>(this->simulator_));
276 this->addOutputModule(std::make_unique<VtkTemperatureModule<TypeTag>>(this->simulator_));
277 }
278
279private:
280 const Implementation& asImp_() const
281 { return *static_cast<const Implementation*>(this); }
282};
283
284} // namespace Opm
285
286#endif
A base class for fully-implicit multi-phase porous-media flow models which assume multiple fluid phas...
Definition multiphasebasemodel.hh:168
void globalPhaseStorage(EqVector &storage, unsigned phaseIdx)
Compute the total storage inside one phase of all conservation quantities.
Definition multiphasebasemodel.hh:218
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition multiphasebasemodel.hh:194
bool phaseIsConsidered(unsigned) const
Returns true iff a fluid phase is used by the model.
Definition multiphasebasemodel.hh:208
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.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
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
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 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.