discretefracturemodel.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_DISCRETE_FRACTURE_MODEL_HH
29#define EWOMS_DISCRETE_FRACTURE_MODEL_HH
30
31#include <opm/material/densead/Math.hpp>
32
41
42#include <memory>
43#include <stdexcept>
44#include <string>
45#include <tuple>
46
47namespace Opm {
48template <class TypeTag>
49class DiscreteFractureModel;
50}
51
52namespace Opm::Properties {
53
54// Create new type tags
55namespace TTag {
56
59{ using InheritsFrom = std::tuple<ImmiscibleTwoPhaseModel>; };
60
61} // end namespace TTag
62
64template<class TypeTag>
65struct Model<TypeTag, TTag::DiscreteFractureModel>
67
69template<class TypeTag>
70struct BaseProblem<TypeTag, TTag::DiscreteFractureModel>
72
74template<class TypeTag>
77
78// The type of the base base class for actual problems.
79// TODO!?
80// template<class TypeTag>
81// struct BaseProblem<TypeTag, TTag::DiscreteFractureModel>
82// { using type = DiscreteFractureBaseProblem<TypeTag>; };
83
85template<class TypeTag>
88
90template<class TypeTag>
93
95template<class TypeTag>
98
101template<class TypeTag>
103{ static constexpr bool value = true; };
104
105} // namespace Opm::Properties
106
107namespace Opm {
108
126template <class TypeTag>
128{
129 using ParentType = ImmiscibleModel<TypeTag>;
131
132public:
133 explicit DiscreteFractureModel(Simulator& simulator)
134 : ParentType(simulator)
135 {
136 if (Parameters::Get<Parameters::EnableIntensiveQuantityCache>()) {
137 throw std::runtime_error("The discrete fracture model does not work in conjunction "
138 "with intensive quantities caching");
139 }
140 }
141
145 static void registerParameters()
146 {
148
149 // register runtime parameters of the VTK output modules
151
152 // The intensive quantity cache cannot be used by the discrete fracture model, because
153 // the intensive quantities of a control degree of freedom are not identical to the
154 // intensive quantities of the other intensive quantities of the same of the same degree
155 // of freedom. This is because the fracture properties (volume, permeability, etc) are
156 // specific for each...
157 Parameters::SetDefault<Parameters::EnableIntensiveQuantityCache>(false);
158 }
159
163 static std::string name()
164 { return "discretefracture"; }
165
167 {
169
170 this->addOutputModule(std::make_unique<VtkDiscreteFractureModule<TypeTag>>(this->simulator_));
171 }
172};
173
174} // namespace Opm
175
176#endif
This class expresses all intensive quantities of the discrete fracture model.
Definition: discretefractureextensivequantities.hh:49
Contains the quantities which are are constant within a finite volume in the discret fracture immisci...
Definition: discretefractureintensivequantities.hh:50
Calculates the local residual of the discrete fracture immiscible multi-phase model.
Definition: discretefracturelocalresidual.hh:43
A fully-implicit multi-phase flow model which assumes immiscibility of the phases and is able to incl...
Definition: discretefracturemodel.hh:128
static std::string name()
Definition: discretefracturemodel.hh:163
DiscreteFractureModel(Simulator &simulator)
Definition: discretefracturemodel.hh:133
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: discretefracturemodel.hh:145
void registerOutputModules_()
Definition: discretefracturemodel.hh:166
Represents the primary variables used by the discrete fracture multi-phase model.
Definition: discretefractureprimaryvariables.hh:48
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
Definition: discretefractureproblem.hh:50
A fully-implicit multi-phase flow model which assumes immiscibility of the phases.
Definition: immisciblemodel.hh:211
static void registerParameters()
Register all run-time parameters for the immiscible model.
Definition: immisciblemodel.hh:236
void registerOutputModules_()
Definition: immisciblemodel.hh:354
VTK output module for quantities which make sense for all models which deal with discrete fractures i...
Definition: vtkdiscretefracturemodule.hpp:66
static void registerParameters()
Register all run-time parameters for the multi-phase VTK output module.
Definition: vtkdiscretefracturemodule.hpp:101
Defines the properties required for the immiscible multi-phase model which considers discrete fractur...
Definition: blackoilmodel.hh:79
Definition: blackoilboundaryratevector.hh:39
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
The type of the base class for all problems which use this model.
Definition: fvbaseproperties.hh:84
Data required to calculate a flux over a face.
Definition: fvbaseproperties.hh:149
The secondary variables within a sub-control volume.
Definition: fvbaseproperties.hh:133
The type of the local residual function.
Definition: fvbaseproperties.hh:94
The type of the model.
Definition: basicproperties.hh:88
A vector of primary variables within a sub-control volume.
Definition: fvbaseproperties.hh:130
The generic type tag for problems using the immiscible multi-phase model.
Definition: discretefracturemodel.hh:59
std::tuple< ImmiscibleTwoPhaseModel > InheritsFrom
Definition: discretefracturemodel.hh:59
Definition: discretefractureproperties.hh:38