fvbasefdlocallinearizer.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_FD_LOCAL_LINEARIZER_HH
29#define EWOMS_FV_BASE_FD_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/MathToolbox.hpp>
38#include <opm/material/common/Valgrind.hpp>
39
43
45
46#include <algorithm>
47#include <cmath>
48#include <cstddef>
49#include <limits>
50#include <memory>
51
52namespace Opm {
53
54// forward declaration
55template<class TypeTag>
56class FvBaseFdLocalLinearizer;
57
58} // namespace Opm
59
60namespace Opm::Properties {
61
62// declare the property tags required for the finite differences local linearizer
63
64namespace TTag {
66} // namespace TTag
67
68template<class TypeTag, class MyTypeTag>
70
71// set the properties to be spliced in
72template<class TypeTag>
73struct LocalLinearizer<TypeTag, TTag::FiniteDifferenceLocalLinearizer>
75
76template<class TypeTag>
77struct Evaluation<TypeTag, TTag::FiniteDifferenceLocalLinearizer>
79
81template<class TypeTag>
82struct BaseEpsilon<TypeTag, TTag::FiniteDifferenceLocalLinearizer>
83{
85 static constexpr type value = std::max<type>(0.9123e-10,
86 std::numeric_limits<type>::epsilon() * 1.23e3);
87};
88
89} // namespace Opm::Properties
90
91namespace Opm::Parameters {
92
100struct NumericDifferenceMethod { static constexpr int value = +1; };
101
102} // namespace Opm::Parameters
103
104namespace Opm {
105
142template<class TypeTag>
144{
145private:
155 using Element = typename GridView::template Codim<0>::Entity;
156
157 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
158
159 // extract local matrices from jacobian matrix for consistency
161 using ScalarVectorBlock = Dune::FieldVector<Scalar, numEq>;
162
163 using ScalarLocalBlockVector = Dune::BlockVector<ScalarVectorBlock>;
164 using ScalarLocalBlockMatrix = Dune::Matrix<ScalarMatrixBlock>;
165
166 using LocalEvalBlockVector = typename LocalResidual::LocalEvalBlockVector;
167
168public:
170
171 // Disable copying for efficiency reasons
173
177 static void registerParameters()
178 {
179 Parameters::Register<Parameters::NumericDifferenceMethod>
180 ("The method used for numeric differentiation (-1: backward "
181 "differences, 0: central differences, 1: forward differences)");
182 }
183
192 void init(Simulator& simulator)
193 {
194 simulatorPtr_ = &simulator;
195 internalElemContext_ = std::make_unique<ElementContext>(simulator);
196 }
197
209 void linearize(const Element& element)
210 {
212 }
213
229 void linearize(ElementContext& elemCtx, const Element& elem)
230 {
231 elemCtx.updateAll(elem);
232
233 // update the weights of the primary variables for the context
234 model_().updatePVWeights(elemCtx);
235
236 resize_(elemCtx);
237 reset_(elemCtx);
238
239 // calculate the local residual
240 localResidual_.eval(residual_, elemCtx);
241
242 // calculate the local jacobian matrix
243 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
244 for (auto dofIdx = 0 * numPrimaryDof; dofIdx < numPrimaryDof; ++dofIdx) {
245 for (auto pvIdx = 0 * numEq; pvIdx < numEq; ++pvIdx) {
246 asImp_().evalPartialDerivative_(elemCtx, dofIdx, pvIdx);
247
248 // incorporate the partial derivatives into the local Jacobian matrix
249 updateLocalJacobian_(elemCtx, dofIdx, pvIdx);
250 }
251 }
252 }
253
258 static Scalar baseEpsilon()
259 { return getPropValue<TypeTag, Properties::BaseEpsilon>(); }
260
272 Scalar numericEpsilon(const ElementContext& elemCtx,
273 unsigned dofIdx,
274 unsigned pvIdx) const
275 {
276 const unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
277 const Scalar pvWeight = elemCtx.model().primaryVarWeight(globalIdx, pvIdx);
278 assert(pvWeight > 0 && std::isfinite(pvWeight));
279 Valgrind::CheckDefined(pvWeight);
280
281 return baseEpsilon() / pvWeight;
282 }
283
287 LocalResidual& localResidual()
288 { return localResidual_; }
289
293 const LocalResidual& localResidual() const
294 { return localResidual_; }
295
304 const ScalarMatrixBlock& jacobian(unsigned domainScvIdx, unsigned rangeScvIdx) const
305 { return jacobian_[domainScvIdx][rangeScvIdx]; }
306
312 const ScalarVectorBlock& residual(unsigned dofIdx) const
313 { return residual_[dofIdx]; }
314
315protected:
316 Implementation& asImp_()
317 { return *static_cast<Implementation*>(this); }
318
319 const Implementation& asImp_() const
320 { return *static_cast<const Implementation*>(this); }
321
322 const Simulator& simulator_() const
323 { return *simulatorPtr_; }
324
325 const Problem& problem_() const
326 { return simulatorPtr_->problem(); }
327
328 const Model& model_() const
329 { return simulatorPtr_->model(); }
330
335 {
336 static int diff = Parameters::Get<Parameters::NumericDifferenceMethod>();
337 return diff;
338 }
339
344 void resize_(const ElementContext& elemCtx)
345 {
346 const std::size_t numDof = elemCtx.numDof(/*timeIdx=*/0);
347 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
348
349 residual_.resize(numDof);
350 jacobian_.setSize(numDof, numPrimaryDof);
351
352 derivResidual_.resize(numDof);
353 }
354
358 void reset_(const ElementContext& elemCtx)
359 {
360 const std::size_t numDof = elemCtx.numDof(/*timeIdx=*/0);
361 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
362 for (unsigned primaryDofIdx = 0; primaryDofIdx < numPrimaryDof; ++primaryDofIdx) {
363 for (unsigned dof2Idx = 0; dof2Idx < numDof; ++dof2Idx) {
364 jacobian_[dof2Idx][primaryDofIdx] = 0.0;
365 }
366 }
367
368 for (unsigned primaryDofIdx = 0; primaryDofIdx < numDof; ++primaryDofIdx) {
369 residual_[primaryDofIdx] = 0.0;
370 }
371 }
372
415 void evalPartialDerivative_(ElementContext& elemCtx,
416 unsigned dofIdx,
417 unsigned pvIdx)
418 {
419 // save all quantities which depend on the specified primary
420 // variable at the given sub control volume
421 elemCtx.stashIntensiveQuantities(dofIdx);
422
423 PrimaryVariables priVars(elemCtx.primaryVars(dofIdx, /*timeIdx=*/0));
424 const Scalar eps = asImp_().numericEpsilon(elemCtx, dofIdx, pvIdx);
425 Scalar delta = 0.0;
426
427 if (numericDifferenceMethod_() >= 0) {
428 // we are not using backward differences, i.e. we need to
429 // calculate f(x + \epsilon)
430
431 // deflect primary variables
432 priVars[pvIdx] += eps;
433 delta += eps;
434
435 // calculate the deflected residual
436 elemCtx.updateIntensiveQuantities(priVars, dofIdx, /*timeIdx=*/0);
437 elemCtx.updateAllExtensiveQuantities();
438 localResidual_.eval(derivResidual_, elemCtx);
439 }
440 else {
441 // we are using backward differences, i.e. we don't need
442 // to calculate f(x + \epsilon) and we can recycle the
443 // (already calculated) residual f(x)
445 }
446
447 if (numericDifferenceMethod_() <= 0) {
448 // we are not using forward differences, i.e. we don't
449 // need to calculate f(x - \epsilon)
450
451 // deflect the primary variables
452 priVars[pvIdx] -= delta + eps;
453 delta += eps;
454
455 // calculate the deflected residual again, this time we use the local
456 // residual's internal storage.
457 elemCtx.updateIntensiveQuantities(priVars, dofIdx, /*timeIdx=*/0);
458 elemCtx.updateAllExtensiveQuantities();
459 localResidual_.eval(elemCtx);
460
461 derivResidual_ -= localResidual_.residual();
462 }
463 else {
464 // we are using forward differences, i.e. we don't need to
465 // calculate f(x - \epsilon) and we can recycle the
466 // (already calculated) residual f(x)
468 }
469
470 assert(delta > 0);
471
472 // divide difference in residuals by the magnitude of the
473 // deflections between the two function evaluation
474 derivResidual_ /= delta;
475
476 // restore the original state of the element's volume
477 // variables
478 elemCtx.restoreIntensiveQuantities(dofIdx);
479
480#ifndef NDEBUG
481 for (unsigned i = 0; i < derivResidual_.size(); ++i) {
482 Valgrind::CheckDefined(derivResidual_[i]);
483 }
484#endif
485 }
486
492 void updateLocalJacobian_(const ElementContext& elemCtx,
493 unsigned focusDofIdx,
494 unsigned pvIdx)
495 {
496 const std::size_t numDof = elemCtx.numDof(/*timeIdx=*/0);
497 for (unsigned dofIdx = 0; dofIdx < numDof; ++dofIdx) {
498 for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
499 // A[dofIdx][focusDofIdx][eqIdx][pvIdx] is the partial derivative of the
500 // residual function 'eqIdx' for the degree of freedom 'dofIdx' with
501 // regard to the primary variable 'pvIdx' of the degree of freedom
502 // 'focusDofIdx'
503 jacobian_[dofIdx][focusDofIdx][eqIdx][pvIdx] = derivResidual_[dofIdx][eqIdx];
504 Valgrind::CheckDefined(jacobian_[dofIdx][focusDofIdx][eqIdx][pvIdx]);
505 }
506 }
507 }
508
509 Simulator* simulatorPtr_{};
510
511 std::unique_ptr<ElementContext> internalElemContext_{};
512
513 LocalEvalBlockVector residual_{};
514 LocalEvalBlockVector derivResidual_{};
515 ScalarLocalBlockMatrix jacobian_{};
516
517 LocalResidual localResidual_{};
518};
519
520} // namespace Opm
521
522#endif
Calculates the Jacobian of the local residual for finite volume spatial discretizations using a finit...
Definition: fvbasefdlocallinearizer.hh:144
const Implementation & asImp_() const
Definition: fvbasefdlocallinearizer.hh:319
std::unique_ptr< ElementContext > internalElemContext_
Definition: fvbasefdlocallinearizer.hh:511
LocalEvalBlockVector derivResidual_
Definition: fvbasefdlocallinearizer.hh:514
void evalPartialDerivative_(ElementContext &elemCtx, unsigned dofIdx, unsigned pvIdx)
Compute the partial derivatives of a context's residual functions.
Definition: fvbasefdlocallinearizer.hh:415
void linearize(const Element &element)
Compute an element's local Jacobian matrix and evaluate its residual.
Definition: fvbasefdlocallinearizer.hh:209
Scalar numericEpsilon(const ElementContext &elemCtx, unsigned dofIdx, unsigned pvIdx) const
Returns the epsilon value which is added and removed from the current solution.
Definition: fvbasefdlocallinearizer.hh:272
LocalResidual localResidual_
Definition: fvbasefdlocallinearizer.hh:517
static Scalar baseEpsilon()
Returns the unweighted epsilon value used to calculate the local derivatives.
Definition: fvbasefdlocallinearizer.hh:258
void updateLocalJacobian_(const ElementContext &elemCtx, unsigned focusDofIdx, unsigned pvIdx)
Updates the current local Jacobian matrix with the partial derivatives of all equations for primary v...
Definition: fvbasefdlocallinearizer.hh:492
const Model & model_() const
Definition: fvbasefdlocallinearizer.hh:328
FvBaseFdLocalLinearizer(const FvBaseFdLocalLinearizer &)=delete
ScalarLocalBlockMatrix jacobian_
Definition: fvbasefdlocallinearizer.hh:515
static int numericDifferenceMethod_()
Returns the numeric difference method which is applied.
Definition: fvbasefdlocallinearizer.hh:334
const ScalarMatrixBlock & jacobian(unsigned domainScvIdx, unsigned rangeScvIdx) const
Returns the local Jacobian matrix of the residual of a sub-control volume.
Definition: fvbasefdlocallinearizer.hh:304
static void registerParameters()
Register all run-time parameters for the local jacobian.
Definition: fvbasefdlocallinearizer.hh:177
const Simulator & simulator_() const
Definition: fvbasefdlocallinearizer.hh:322
const LocalResidual & localResidual() const
Return reference to the local residual.
Definition: fvbasefdlocallinearizer.hh:293
Simulator * simulatorPtr_
Definition: fvbasefdlocallinearizer.hh:509
void init(Simulator &simulator)
Initialize the local Jacobian object.
Definition: fvbasefdlocallinearizer.hh:192
const ScalarVectorBlock & residual(unsigned dofIdx) const
Returns the local residual of a sub-control volume.
Definition: fvbasefdlocallinearizer.hh:312
const Problem & problem_() const
Definition: fvbasefdlocallinearizer.hh:325
LocalEvalBlockVector residual_
Definition: fvbasefdlocallinearizer.hh:513
void linearize(ElementContext &elemCtx, const Element &elem)
Compute an element's local Jacobian matrix and evaluate its residual.
Definition: fvbasefdlocallinearizer.hh:229
void reset_(const ElementContext &elemCtx)
Reset the all relevant internal attributes to 0.
Definition: fvbasefdlocallinearizer.hh:358
LocalResidual & localResidual()
Return reference to the local residual.
Definition: fvbasefdlocallinearizer.hh:287
void resize_(const ElementContext &elemCtx)
Resize all internal attributes to the size of the element.
Definition: fvbasefdlocallinearizer.hh:344
Implementation & asImp_()
Definition: fvbasefdlocallinearizer.hh:316
Declare the properties used by the infrastructure code of the finite volume discretizations.
Declares the properties required by the black oil model.
Definition: blackoilnewtonmethodparams.hpp:31
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
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
Specify which kind of method should be used to numerically calculate the partial derivatives of the r...
Definition: fvbasefdlocallinearizer.hh:100
static constexpr int value
Definition: fvbasefdlocallinearizer.hh:100
GetPropType< TypeTag, Properties::Scalar > type
Definition: fvbasefdlocallinearizer.hh:84
Definition: fvbasefdlocallinearizer.hh:69
GetPropType< TypeTag, Properties::Scalar > type
Definition: fvbasefdlocallinearizer.hh:78
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: fvbasefdlocallinearizer.hh:65
a tag to mark properties as undefined
Definition: propertysystem.hh:38