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>
50class FvBaseAdLocalLinearizer;
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>
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>
93{
94private:
104 using Element = typename GridView::template Codim<0>::Entity;
105
106 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
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:
117
118 // copying local linearizer objects around is a very bad idea, so we explicitly
119 // prevent it...
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 {
156 }
157
172 void linearize(ElementContext& elemCtx, const Element& elem)
173 {
174 elemCtx.updateStencil(elem);
175 elemCtx.updateAllIntensiveQuantities();
176
177 // update the weights of the primary variables for the context
178 model_().updatePVWeights(elemCtx);
179
180 // resize the internal arrays of the linearizer
181 resize_(elemCtx);
182
183 // compute the local residual and its Jacobian
184 const unsigned numPrimaryDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
185 for (unsigned focusDofIdx = 0; focusDofIdx < numPrimaryDof; ++focusDofIdx) {
186 elemCtx.setFocusDofIndex(focusDofIdx);
187 elemCtx.updateAllExtensiveQuantities();
188
189 // calculate the local residual
190 localResidual_.eval(elemCtx);
191
192 // convert the local Jacobian matrix and the right hand side from the data
193 // structures used by the automatic differentiation code to the conventional
194 // ones used by the linear solver.
195 updateLocalLinearization_(elemCtx, focusDofIdx);
196 }
197 }
198
202 LocalResidual& localResidual()
203 { return localResidual_; }
204
208 const LocalResidual& localResidual() const
209 { return localResidual_; }
210
219 const ScalarMatrixBlock& jacobian(unsigned domainScvIdx, unsigned rangeScvIdx) const
220 { return jacobian_[domainScvIdx][rangeScvIdx]; }
221
227 const ScalarVectorBlock& residual(unsigned dofIdx) const
228 { return residual_[dofIdx]; }
229
230protected:
231 Implementation& asImp_()
232 { return *static_cast<Implementation*>(this); }
233
234 const Implementation& asImp_() const
235 { return *static_cast<const Implementation*>(this); }
236
237 const Simulator& simulator_() const
238 { return *simulatorPtr_; }
239
240 const Problem& problem_() const
241 { return simulatorPtr_->problem(); }
242
243 const Model& model_() const
244 { return simulatorPtr_->model(); }
245
250 void resize_(const ElementContext& elemCtx)
251 {
252 const std::size_t numDof = elemCtx.numDof(/*timeIdx=*/0);
253 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
254
255 residual_.resize(numDof);
256 if (jacobian_.N() != numDof || jacobian_.M() != numPrimaryDof) {
257 jacobian_.setSize(numDof, numPrimaryDof);
258 }
259 }
260
264 void reset_(const ElementContext&)
265 {
266 residual_ = 0.0;
267 jacobian_ = 0.0;
268 }
269
274 void updateLocalLinearization_(const ElementContext& elemCtx,
275 unsigned focusDofIdx)
276 {
277 const auto& resid = localResidual_.residual();
278
279 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
280 residual_[focusDofIdx][eqIdx] = resid[focusDofIdx][eqIdx].value();
281 }
282
283 const std::size_t numDof = elemCtx.numDof(/*timeIdx=*/0);
284 for (unsigned dofIdx = 0; dofIdx < numDof; ++dofIdx) {
285 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
286 for (unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
287 // A[dofIdx][focusDofIdx][eqIdx][pvIdx] is the partial derivative of
288 // the residual function 'eqIdx' for the degree of freedom 'dofIdx'
289 // with regard to the focus variable 'pvIdx' of the degree of freedom
290 // 'focusDofIdx'
291 jacobian_[dofIdx][focusDofIdx][eqIdx][pvIdx] = resid[dofIdx][eqIdx].derivative(pvIdx);
292 Valgrind::CheckDefined(jacobian_[dofIdx][focusDofIdx][eqIdx][pvIdx]);
293 }
294 }
295 }
296 }
297
298 Simulator* simulatorPtr_{};
299
300 std::unique_ptr<ElementContext> internalElemContext_{};
301
302 LocalResidual localResidual_{};
303
304 ScalarLocalBlockVector residual_{};
305 ScalarLocalBlockMatrix jacobian_{};
306};
307
308} // namespace Opm
309
310#endif
Calculates the local residual and its Jacobian for a single element of the grid.
Definition: fvbaseadlocallinearizer.hh:93
Simulator * simulatorPtr_
Definition: fvbaseadlocallinearizer.hh:298
const ScalarVectorBlock & residual(unsigned dofIdx) const
Returns the local residual of a sub-control volume.
Definition: fvbaseadlocallinearizer.hh:227
void linearize(const Element &element)
Compute an element's local Jacobian matrix and evaluate its residual.
Definition: fvbaseadlocallinearizer.hh:153
Implementation & asImp_()
Definition: fvbaseadlocallinearizer.hh:231
const Problem & problem_() const
Definition: fvbaseadlocallinearizer.hh:240
ScalarLocalBlockVector residual_
Definition: fvbaseadlocallinearizer.hh:304
LocalResidual & localResidual()
Return reference to the local residual.
Definition: fvbaseadlocallinearizer.hh:202
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:274
void reset_(const ElementContext &)
Reset the all relevant internal attributes to 0.
Definition: fvbaseadlocallinearizer.hh:264
const Implementation & asImp_() const
Definition: fvbaseadlocallinearizer.hh:234
std::unique_ptr< ElementContext > internalElemContext_
Definition: fvbaseadlocallinearizer.hh:300
FvBaseAdLocalLinearizer(const FvBaseAdLocalLinearizer &)=delete
ScalarLocalBlockMatrix jacobian_
Definition: fvbaseadlocallinearizer.hh:305
void resize_(const ElementContext &elemCtx)
Resize all internal attributes to the size of the element.
Definition: fvbaseadlocallinearizer.hh:250
const Simulator & simulator_() const
Definition: fvbaseadlocallinearizer.hh:237
const Model & model_() const
Definition: fvbaseadlocallinearizer.hh:243
const LocalResidual & localResidual() const
Return reference to the local residual.
Definition: fvbaseadlocallinearizer.hh:208
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:172
const ScalarMatrixBlock & jacobian(unsigned domainScvIdx, unsigned rangeScvIdx) const
Returns the local Jacobian matrix of the residual of a sub-control volume.
Definition: fvbaseadlocallinearizer.hh:219
LocalResidual localResidual_
Definition: fvbaseadlocallinearizer.hh:302
Declare the properties used by the infrastructure code of the finite volume discretizations.
Declares the properties required by the black oil model.
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
DenseAd::Evaluation< Scalar, numEq > type
Definition: fvbaseadlocallinearizer.hh:76
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