opm-simulators
Loading...
Searching...
No Matches
fvbaseadlocallinearizer.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_AD_LOCAL_LINEARIZER_HH
29#define EWOMS_FV_BASE_AD_LOCAL_LINEARIZER_HH
30
31#include <dune/common/fmatrix.hh>
32#include <dune/common/fvector.hh>
33
34#include <dune/istl/bvector.hh>
35#include <dune/istl/matrix.hh>
36
37#include <opm/material/common/Valgrind.hpp>
38#include <opm/material/densead/Math.hpp>
39
41
43
44#include <cstddef>
45#include <memory>
46
47namespace Opm {
48// forward declaration
49template<class TypeTag>
51}
52
53namespace Opm::Properties {
54
55// declare the property tags required for the finite differences local linearizer
56
57namespace TTag {
59} // namespace TTag
60
61// set the properties to be spliced in
62template<class TypeTag>
63struct LocalLinearizer<TypeTag, TTag::AutoDiffLocalLinearizer>
64{ using type = FvBaseAdLocalLinearizer<TypeTag>; };
65
67template<class TypeTag>
68struct Evaluation<TypeTag, TTag::AutoDiffLocalLinearizer>
69{
70private:
71 static constexpr unsigned numEq = getPropValue<TypeTag, Properties::NumEq>();
72
74
75public:
76 using type = DenseAd::Evaluation<Scalar, numEq>;
77};
78
79} // namespace Opm::Properties
80
81namespace Opm {
82
91template<class TypeTag>
92class FvBaseAdLocalLinearizer
93{
94private:
104 using Element = typename GridView::template Codim<0>::Entity;
105
107
108 using ScalarVectorBlock = Dune::FieldVector<Scalar, numEq>;
109 // extract local matrices from jacobian matrix for consistency
111
112 using ScalarLocalBlockVector = Dune::BlockVector<ScalarVectorBlock>;
113 using ScalarLocalBlockMatrix = Dune::Matrix<ScalarMatrixBlock>;
114
115public:
116 FvBaseAdLocalLinearizer() = default;
117
118 // copying local linearizer objects around is a very bad idea, so we explicitly
119 // prevent it...
120 FvBaseAdLocalLinearizer(const FvBaseAdLocalLinearizer&) = delete;
121
125 static void registerParameters()
126 {}
127
136 void init(Simulator& simulator)
137 {
138 simulatorPtr_ = &simulator;
139 internalElemContext_ = std::make_unique<ElementContext>(simulator);
140 }
141
153 void linearize(const Element& element)
154 {
155 linearize(*internalElemContext_, element);
156 }
157
173 void linearize(ElementContext& elemCtx, const Element& elem)
174 {
175 elemCtx.updateStencil(elem);
176 elemCtx.updateAllIntensiveQuantities();
177
178 // update the weights of the primary variables for the context
179 model_().updatePVWeights(elemCtx);
180
181 // resize the internal arrays of the linearizer
182 resize_(elemCtx);
183
184 // compute the local residual and its Jacobian
185 const unsigned numPrimaryDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
186 for (unsigned focusDofIdx = 0; focusDofIdx < numPrimaryDof; ++focusDofIdx) {
187 elemCtx.setFocusDofIndex(focusDofIdx);
188 elemCtx.updateAllExtensiveQuantities();
189
190 // calculate the local residual
191 localResidual_.eval(elemCtx);
192
193 // convert the local Jacobian matrix and the right hand side from the data
194 // structures used by the automatic differentiation code to the conventional
195 // ones used by the linear solver.
196 updateLocalLinearization_(elemCtx, focusDofIdx);
197 }
198 }
199
203 LocalResidual& localResidual()
204 { return localResidual_; }
205
209 const LocalResidual& localResidual() const
210 { return localResidual_; }
211
220 const ScalarMatrixBlock& jacobian(unsigned domainScvIdx, unsigned rangeScvIdx) const
221 { return jacobian_[domainScvIdx][rangeScvIdx]; }
222
228 const ScalarVectorBlock& residual(unsigned dofIdx) const
229 { return residual_[dofIdx]; }
230
231protected:
232 Implementation& asImp_()
233 { return *static_cast<Implementation*>(this); }
234
235 const Implementation& asImp_() const
236 { return *static_cast<const Implementation*>(this); }
237
238 const Simulator& simulator_() const
239 { return *simulatorPtr_; }
240
241 const Problem& problem_() const
242 { return simulatorPtr_->problem(); }
243
244 const Model& model_() const
245 { return simulatorPtr_->model(); }
246
251 void resize_(const ElementContext& elemCtx)
252 {
253 const std::size_t numDof = elemCtx.numDof(/*timeIdx=*/0);
254 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
255
256 residual_.resize(numDof);
257 if (jacobian_.N() != numDof || jacobian_.M() != numPrimaryDof) {
258 jacobian_.setSize(numDof, numPrimaryDof);
259 }
260 }
261
265 void reset_(const ElementContext&)
266 {
267 residual_ = 0.0;
268 jacobian_ = 0.0;
269 }
270
275 void updateLocalLinearization_(const ElementContext& elemCtx,
276 unsigned focusDofIdx)
277 {
278 const auto& resid = localResidual_.residual();
279
280 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
281 residual_[focusDofIdx][eqIdx] = resid[focusDofIdx][eqIdx].value();
282 }
283
284 const std::size_t numDof = elemCtx.numDof(/*timeIdx=*/0);
285 for (unsigned dofIdx = 0; dofIdx < numDof; ++dofIdx) {
286 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
287 for (unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
288 // A[dofIdx][focusDofIdx][eqIdx][pvIdx] is the partial derivative of
289 // the residual function 'eqIdx' for the degree of freedom 'dofIdx'
290 // with regard to the focus variable 'pvIdx' of the degree of freedom
291 // 'focusDofIdx'
292 jacobian_[dofIdx][focusDofIdx][eqIdx][pvIdx] = resid[dofIdx][eqIdx].derivative(pvIdx);
293 Valgrind::CheckDefined(jacobian_[dofIdx][focusDofIdx][eqIdx][pvIdx]);
294 }
295 }
296 }
297 }
298
299 Simulator* simulatorPtr_{};
300
301 std::unique_ptr<ElementContext> internalElemContext_{};
302
303 LocalResidual localResidual_{};
304
305 ScalarLocalBlockVector residual_{};
306 ScalarLocalBlockMatrix jacobian_{};
307};
308
309} // namespace Opm
310
311#endif
Calculates the local residual and its Jacobian for a single element of the grid.
Definition fvbaseadlocallinearizer.hh:93
const ScalarVectorBlock & residual(unsigned dofIdx) const
Returns the local residual of a sub-control volume.
Definition fvbaseadlocallinearizer.hh:228
void linearize(const Element &element)
Compute an element's local Jacobian matrix and evaluate its residual.
Definition fvbaseadlocallinearizer.hh:153
LocalResidual & localResidual()
Return reference to the local residual.
Definition fvbaseadlocallinearizer.hh:203
void updateLocalLinearization_(const ElementContext &elemCtx, unsigned focusDofIdx)
Updates the current local Jacobian matrix with the partial derivatives of all equations for the degre...
Definition fvbaseadlocallinearizer.hh:275
void reset_(const ElementContext &)
Reset the all relevant internal attributes to 0.
Definition fvbaseadlocallinearizer.hh:265
void resize_(const ElementContext &elemCtx)
Resize all internal attributes to the size of the element.
Definition fvbaseadlocallinearizer.hh:251
const LocalResidual & localResidual() const
Return reference to the local residual.
Definition fvbaseadlocallinearizer.hh:209
void init(Simulator &simulator)
Initialize the local Jacobian object.
Definition fvbaseadlocallinearizer.hh:136
static void registerParameters()
Register all run-time parameters for the local jacobian.
Definition fvbaseadlocallinearizer.hh:125
void linearize(ElementContext &elemCtx, const Element &elem)
Compute an element's local Jacobian matrix and evaluate its residual.
Definition fvbaseadlocallinearizer.hh:173
const ScalarMatrixBlock & jacobian(unsigned domainScvIdx, unsigned rangeScvIdx) const
Returns the local Jacobian matrix of the residual of a sub-control volume.
Definition fvbaseadlocallinearizer.hh:220
Manages the initializing and running of time dependent problems.
Definition simulator.hh:84
Declare the properties used by the infrastructure code of the finite volume discretizations.
Declares the properties required by the black oil model.
The generic type tag for problems using the immiscible multi-phase model.
Definition blackoilmodel.hh:81
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
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
Representation of a function evaluation and all necessary derivatives with regard to the intensive qu...
Definition fvbaseproperties.hh:66
The type of the local linearizer.
Definition fvbaseproperties.hh:97
Definition fvbaseadlocallinearizer.hh:58