fvbasediscretizationfemadapt.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_FV_BASE_DISCRETIZATION_FEMADAPT_HH
29#define EWOMS_FV_BASE_DISCRETIZATION_FEMADAPT_HH
30
31#include <dune/fem/function/blockvectorfunction.hh>
32#include <dune/fem/misc/capabilities.hh>
33#include <dune/fem/space/common/adaptationmanager.hh>
34#include <dune/fem/space/common/restrictprolongtuple.hh>
35
37
38#include <memory>
39#include <stdexcept>
40
41namespace Opm {
42
43template<class TypeTag>
44class FvBaseDiscretizationFemAdapt;
45
46namespace Properties {
47
48template<class TypeTag>
49struct BaseDiscretizationType<TypeTag, TTag::FvBaseDiscretization>
51
52template<class TypeTag>
53struct DiscreteFunction<TypeTag, TTag::FvBaseDiscretization>
54{
57 using type = Dune::Fem::ISTLBlockVectorDiscreteFunction<DiscreteFunctionSpace, PrimaryVariables>;
58};
59
60} // namespace Properties
61
68template <class TypeTag>
70{
76
77 static constexpr unsigned historySize = getPropValue<TypeTag, Properties::TimeDiscHistorySize>();
78
80 // discrete function storing solution data
81 using DiscreteFunction = Dune::Fem::ISTLBlockVectorDiscreteFunction<DiscreteFunctionSpace, PrimaryVariables>;
82
83 // problem restriction and prolongation operator for adaptation
84 using ProblemRestrictProlongOperator = typename Problem::RestrictProlongOperator;
85
86 // discrete function restriction and prolongation operator for adaptation
87 using DiscreteFunctionRestrictProlong = Dune::Fem::RestrictProlongDefault<DiscreteFunction>;
88 using RestrictProlong
89 = Dune::Fem::RestrictProlongTuple<DiscreteFunctionRestrictProlong, ProblemRestrictProlongOperator>;
90
91 // adaptation classes
92 using AdaptationManager = Dune::Fem::AdaptationManager<Grid, RestrictProlong>;
93
94public:
95 template<class Serializer>
97 {
98 template<class SolutionType>
99 static void serializeOp(Serializer& serializer,
100 SolutionType& solution)
101 {
102 for (auto& sol : solution) {
103 serializer(sol->blockVector());
104 }
105 }
106 };
107
108 explicit FvBaseDiscretizationFemAdapt(Simulator& simulator)
109 : ParentType(simulator)
110 , space_(simulator.vanguard().gridPart())
111 {
112 if (this->enableGridAdaptation_ && !Dune::Fem::Capabilities::isLocallyAdaptive<Grid>::v) {
113 throw std::invalid_argument("Grid adaptation enabled, but chosen Grid is not capable"
114 " of adaptivity");
115 }
116
117 for (unsigned timeIdx = 0; timeIdx < historySize; ++timeIdx) {
118 this->solution_[timeIdx] = std::make_unique<DiscreteFunction>("solution", space_);
119 }
120 }
121
123 {
124 // adapt the grid if enabled and if all dependencies are available
125 // adaptation is only done if markForGridAdaptation returns true
126 if (this->enableGridAdaptation_) {
127 // check if problem allows for adaptation and cells were marked
128 if (this->simulator_.problem().markForGridAdaptation()) {
129 // adapt the grid and load balance if necessary
130 adaptationManager().adapt();
131
132 // if the grid has potentially changed, we need to re-create the
133 // supporting data structures.
134 this->elementMapper_.update(this->gridView_);
135 this->vertexMapper_.update(this->gridView_);
136 this->resetLinearizer();
137 // this is a bit hacky because it supposes that Problem::finishInit()
138 // works fine multiple times in a row.
139 //
140 // TODO: move this to Problem::gridChanged()
141 this->finishInit();
142
143 // notify the problem that the grid has changed
144 //
145 // TODO: come up with a mechanism to access the unadapted data structures
146 // outside of the problem (i.e., grid, mappers, solutions)
147 this->simulator_.problem().gridChanged();
148
149 // notify the modules for visualization output
150 for (auto& module : this->outputModules_) {
151 module->allocBuffers();
152 }
153 }
154 }
155 }
156
157 AdaptationManager& adaptationManager()
158 {
159 if (!adaptationManager_) {
160 // create adaptation objects here, because when doing so in constructor
161 // problem is not yet intialized, aka seg fault
162 restrictProlong_ =
163 std::make_unique<RestrictProlong>(DiscreteFunctionRestrictProlong(*(this->solution_[/*timeIdx=*/0])),
164 this->simulator_.problem().restrictProlongOperator());
165 adaptationManager_ =
166 std::make_unique<AdaptationManager>(this->simulator_.vanguard().grid(), *restrictProlong_);
167 }
168 return *adaptationManager_;
169 }
170
171private:
172 DiscreteFunctionSpace space_;
173 std::unique_ptr<RestrictProlong> restrictProlong_;
174 std::unique_ptr<AdaptationManager> adaptationManager_;
175};
176
177} // namespace Opm
178
179#endif
The base class for the finite volume discretization schemes.
Definition: fvbasediscretizationfemadapt.hh:70
AdaptationManager & adaptationManager()
Definition: fvbasediscretizationfemadapt.hh:157
FvBaseDiscretizationFemAdapt(Simulator &simulator)
Definition: fvbasediscretizationfemadapt.hh:108
void adaptGrid()
Definition: fvbasediscretizationfemadapt.hh:122
The base class for the finite volume discretization schemes without adaptation.
Definition: fvbasediscretization.hh:1952
The base class for the finite volume discretization schemes.
Definition: fvbasediscretization.hh:298
void finishInit()
Apply the initial conditions to the model.
Definition: fvbasediscretization.hh:468
bool enableGridAdaptation_
Definition: fvbasediscretization.hh:1939
std::array< std::unique_ptr< DiscreteFunction >, historySize > solution_
Definition: fvbasediscretization.hh:1929
std::list< std::unique_ptr< BaseOutputModule< TypeTag > > > outputModules_
Definition: fvbasediscretization.hh:1931
ElementMapper elementMapper_
Definition: fvbasediscretization.hh:1903
GridView gridView_
Definition: fvbasediscretization.hh:1900
const SolutionVector & solution(unsigned timeIdx) const
Reference to the solution at a given history index as a block vector.
Definition: fvbasediscretization.hh:1144
Simulator & simulator_
Definition: fvbasediscretization.hh:1897
VertexMapper vertexMapper_
Definition: fvbasediscretization.hh:1904
void resetLinearizer()
Resets the Jacobian matrix linearizer, so that the boundary types can be altered.
Definition: fvbasediscretization.hh:1537
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
Definition: fvbasediscretizationfemadapt.hh:97
static void serializeOp(Serializer &serializer, SolutionType &solution)
Definition: fvbasediscretizationfemadapt.hh:99
GetPropType< TypeTag, Properties::DiscreteFunctionSpace > DiscreteFunctionSpace
Definition: fvbasediscretizationfemadapt.hh:55
GetPropType< TypeTag, Properties::PrimaryVariables > PrimaryVariables
Definition: fvbasediscretizationfemadapt.hh:56
typename BaseDiscretization::BlockVectorWrapper type
Definition: fvbasediscretization.hh:283
Definition: fvbaseproperties.hh:77