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
34
35#include <opm/material/common/MathToolbox.hpp>
36#include <opm/material/common/Valgrind.hpp>
37
38#include <dune/istl/bvector.hh>
39#include <dune/istl/matrix.hh>
40
41#include <dune/common/fvector.hh>
42#include <dune/common/fmatrix.hh>
43
44#include <limits>
45
46namespace Opm {
47// forward declaration
48template<class TypeTag>
49class FvBaseFdLocalLinearizer;
50
51} // namespace Opm
52
53namespace Opm::Properties {
54
55// declare the property tags required for the finite differences local linearizer
56
57namespace TTag {
59} // namespace TTag
60
61template<class TypeTag, class MyTypeTag>
63template<class TypeTag, class MyTypeTag>
65
66// set the properties to be spliced in
67template<class TypeTag>
68struct LocalLinearizer<TypeTag, TTag::FiniteDifferenceLocalLinearizer>
70
71template<class TypeTag>
72struct Evaluation<TypeTag, TTag::FiniteDifferenceLocalLinearizer>
74
82template<class TypeTag>
83struct NumericDifferenceMethod<TypeTag, TTag::FiniteDifferenceLocalLinearizer> { static constexpr int value = +1; };
84
86template<class TypeTag>
87struct BaseEpsilon<TypeTag, TTag::FiniteDifferenceLocalLinearizer>
88{
90 static constexpr type value = std::max<type>(0.9123e-10, std::numeric_limits<type>::epsilon()*1.23e3);
91};
92
93} // namespace Opm::Properties
94
95namespace Opm {
96
133template<class TypeTag>
135{
136private:
146 using Element = typename GridView::template Codim<0>::Entity;
147
148 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
149
150 // extract local matrices from jacobian matrix for consistency
152 using ScalarVectorBlock = Dune::FieldVector<Scalar, numEq>;
153
154 using ScalarLocalBlockVector = Dune::BlockVector<ScalarVectorBlock>;
155 using ScalarLocalBlockMatrix = Dune::Matrix<ScalarMatrixBlock>;
156
157 using LocalEvalBlockVector = typename LocalResidual::LocalEvalBlockVector;
158
159#if __GNUC__ == 4 && __GNUC_MINOR__ <= 6
160public:
161 // make older GCCs happy by providing a public copy constructor (this is necessary
162 // for their implementation of std::vector, although the method is never called...)
165 {}
166
167#else
168 // copying local residual objects around is a very bad idea, so we explicitly prevent
169 // it...
171#endif
172public:
175 { }
176
178 { delete internalElemContext_; }
179
183 static void registerParameters()
184 {
185 Parameters::registerParam<TypeTag, Properties::NumericDifferenceMethod>
186 ("The method used for numeric differentiation (-1: backward "
187 "differences, 0: central differences, 1: forward differences)");
188 }
189
198 void init(Simulator& simulator)
199 {
200 simulatorPtr_ = &simulator;
202 internalElemContext_ = new ElementContext(simulator);
203 }
204
216 void linearize(const Element& element)
217 {
219 }
220
235 void linearize(ElementContext& elemCtx, const Element& elem)
236 {
237 elemCtx.updateAll(elem);
238
239 // update the weights of the primary variables for the context
240 model_().updatePVWeights(elemCtx);
241
242 resize_(elemCtx);
243 reset_(elemCtx);
244
245 // calculate the local residual
246 localResidual_.eval(residual_, elemCtx);
247
248 // calculate the local jacobian matrix
249 size_t numPrimaryDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
250 for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; dofIdx++) {
251 for (unsigned pvIdx = 0; pvIdx < numEq; pvIdx++) {
252 asImp_().evalPartialDerivative_(elemCtx, dofIdx, pvIdx);
253
254 // incorporate the partial derivatives into the local Jacobian matrix
255 updateLocalJacobian_(elemCtx, dofIdx, pvIdx);
256 }
257 }
258 }
259
264 static Scalar baseEpsilon()
265 { return getPropValue<TypeTag, Properties::BaseEpsilon>(); }
266
278 Scalar numericEpsilon(const ElementContext& elemCtx,
279 unsigned dofIdx,
280 unsigned pvIdx) const
281 {
282 unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
283 Scalar pvWeight = elemCtx.model().primaryVarWeight(globalIdx, pvIdx);
284 assert(pvWeight > 0 && std::isfinite(pvWeight));
285 Valgrind::CheckDefined(pvWeight);
286
287 return baseEpsilon()/pvWeight;
288 }
289
293 LocalResidual& localResidual()
294 { return localResidual_; }
295
299 const LocalResidual& localResidual() const
300 { return localResidual_; }
301
310 const ScalarMatrixBlock& jacobian(unsigned domainScvIdx, unsigned rangeScvIdx) const
311 { return jacobian_[domainScvIdx][rangeScvIdx]; }
312
318 const ScalarVectorBlock& residual(unsigned dofIdx) const
319 { return residual_[dofIdx]; }
320
321protected:
322 Implementation& asImp_()
323 { return *static_cast<Implementation*>(this); }
324 const Implementation& asImp_() const
325 { return *static_cast<const Implementation*>(this); }
326
327 const Simulator& simulator_() const
328 { return *simulatorPtr_; }
329 const Problem& problem_() const
330 { return simulatorPtr_->problem(); }
331 const Model& model_() const
332 { return simulatorPtr_->model(); }
333
338 { return Parameters::get<TypeTag, Properties::NumericDifferenceMethod>(); }
339
344 void resize_(const ElementContext& elemCtx)
345 {
346 size_t numDof = elemCtx.numDof(/*timeIdx=*/0);
347 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 size_t numDof = elemCtx.numDof(/*timeIdx=*/0);
361 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 for (unsigned primaryDofIdx = 0; primaryDofIdx < numDof; ++ primaryDofIdx)
367 residual_[primaryDofIdx] = 0.0;
368 }
369
412 void evalPartialDerivative_(ElementContext& elemCtx,
413 unsigned dofIdx,
414 unsigned pvIdx)
415 {
416 // save all quantities which depend on the specified primary
417 // variable at the given sub control volume
418 elemCtx.stashIntensiveQuantities(dofIdx);
419
420 PrimaryVariables priVars(elemCtx.primaryVars(dofIdx, /*timeIdx=*/0));
421 Scalar eps = asImp_().numericEpsilon(elemCtx, dofIdx, pvIdx);
422 Scalar delta = 0.0;
423
424 if (numericDifferenceMethod_() >= 0) {
425 // we are not using backward differences, i.e. we need to
426 // calculate f(x + \epsilon)
427
428 // deflect primary variables
429 priVars[pvIdx] += eps;
430 delta += eps;
431
432 // calculate the deflected residual
433 elemCtx.updateIntensiveQuantities(priVars, dofIdx, /*timeIdx=*/0);
434 elemCtx.updateAllExtensiveQuantities();
435 localResidual_.eval(derivResidual_, elemCtx);
436 }
437 else {
438 // we are using backward differences, i.e. we don't need
439 // to calculate f(x + \epsilon) and we can recycle the
440 // (already calculated) residual f(x)
442 }
443
444 if (numericDifferenceMethod_() <= 0) {
445 // we are not using forward differences, i.e. we don't
446 // need to calculate f(x - \epsilon)
447
448 // deflect the primary variables
449 priVars[pvIdx] -= delta + eps;
450 delta += eps;
451
452 // calculate the deflected residual again, this time we use the local
453 // residual's internal storage.
454 elemCtx.updateIntensiveQuantities(priVars, dofIdx, /*timeIdx=*/0);
455 elemCtx.updateAllExtensiveQuantities();
456 localResidual_.eval(elemCtx);
457
458 derivResidual_ -= localResidual_.residual();
459 }
460 else {
461 // we are using forward differences, i.e. we don't need to
462 // calculate f(x - \epsilon) and we can recycle the
463 // (already calculated) residual f(x)
465 }
466
467 assert(delta > 0);
468
469 // divide difference in residuals by the magnitude of the
470 // deflections between the two function evaluation
471 derivResidual_ /= delta;
472
473 // restore the original state of the element's volume
474 // variables
475 elemCtx.restoreIntensiveQuantities(dofIdx);
476
477#ifndef NDEBUG
478 for (unsigned i = 0; i < derivResidual_.size(); ++i)
479 Valgrind::CheckDefined(derivResidual_[i]);
480#endif
481 }
482
488 void updateLocalJacobian_(const ElementContext& elemCtx,
489 unsigned focusDofIdx,
490 unsigned pvIdx)
491 {
492 size_t numDof = elemCtx.numDof(/*timeIdx=*/0);
493 for (unsigned dofIdx = 0; dofIdx < numDof; dofIdx++) {
494 for (unsigned eqIdx = 0; eqIdx < numEq; eqIdx++) {
495 // A[dofIdx][focusDofIdx][eqIdx][pvIdx] is the partial derivative of the
496 // residual function 'eqIdx' for the degree of freedom 'dofIdx' with
497 // regard to the primary variable 'pvIdx' of the degree of freedom
498 // 'focusDofIdx'
499 jacobian_[dofIdx][focusDofIdx][eqIdx][pvIdx] = derivResidual_[dofIdx][eqIdx];
500 Valgrind::CheckDefined(jacobian_[dofIdx][focusDofIdx][eqIdx][pvIdx]);
501 }
502 }
503 }
504
505 Simulator *simulatorPtr_;
506 Model *modelPtr_;
507
508 ElementContext *internalElemContext_;
509
510 LocalEvalBlockVector residual_;
511 LocalEvalBlockVector derivResidual_;
512 ScalarLocalBlockMatrix jacobian_;
513
514 LocalResidual localResidual_;
515};
516
517} // namespace Opm
518
519#endif
Calculates the Jacobian of the local residual for finite volume spatial discretizations using a finit...
Definition: fvbasefdlocallinearizer.hh:135
const Implementation & asImp_() const
Definition: fvbasefdlocallinearizer.hh:324
LocalEvalBlockVector derivResidual_
Definition: fvbasefdlocallinearizer.hh:511
void evalPartialDerivative_(ElementContext &elemCtx, unsigned dofIdx, unsigned pvIdx)
Compute the partial derivatives of a context's residual functions.
Definition: fvbasefdlocallinearizer.hh:412
void linearize(const Element &element)
Compute an element's local Jacobian matrix and evaluate its residual.
Definition: fvbasefdlocallinearizer.hh:216
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:278
LocalResidual localResidual_
Definition: fvbasefdlocallinearizer.hh:514
static Scalar baseEpsilon()
Returns the unweighted epsilon value used to calculate the local derivatives.
Definition: fvbasefdlocallinearizer.hh:264
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:488
const Model & model_() const
Definition: fvbasefdlocallinearizer.hh:331
ScalarLocalBlockMatrix jacobian_
Definition: fvbasefdlocallinearizer.hh:512
static int numericDifferenceMethod_()
Returns the numeric difference method which is applied.
Definition: fvbasefdlocallinearizer.hh:337
ElementContext * internalElemContext_
Definition: fvbasefdlocallinearizer.hh:508
const ScalarMatrixBlock & jacobian(unsigned domainScvIdx, unsigned rangeScvIdx) const
Returns the local Jacobian matrix of the residual of a sub-control volume.
Definition: fvbasefdlocallinearizer.hh:310
~FvBaseFdLocalLinearizer()
Definition: fvbasefdlocallinearizer.hh:177
static void registerParameters()
Register all run-time parameters for the local jacobian.
Definition: fvbasefdlocallinearizer.hh:183
Model * modelPtr_
Definition: fvbasefdlocallinearizer.hh:506
const Simulator & simulator_() const
Definition: fvbasefdlocallinearizer.hh:327
const LocalResidual & localResidual() const
Return reference to the local residual.
Definition: fvbasefdlocallinearizer.hh:299
Simulator * simulatorPtr_
Definition: fvbasefdlocallinearizer.hh:505
void init(Simulator &simulator)
Initialize the local Jacobian object.
Definition: fvbasefdlocallinearizer.hh:198
const ScalarVectorBlock & residual(unsigned dofIdx) const
Returns the local residual of a sub-control volume.
Definition: fvbasefdlocallinearizer.hh:318
const Problem & problem_() const
Definition: fvbasefdlocallinearizer.hh:329
LocalEvalBlockVector residual_
Definition: fvbasefdlocallinearizer.hh:510
void linearize(ElementContext &elemCtx, const Element &elem)
Compute an element's local Jacobian matrix and evaluate its residual.
Definition: fvbasefdlocallinearizer.hh:235
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:293
void resize_(const ElementContext &elemCtx)
Resize all internal attributes to the size of the element.
Definition: fvbasefdlocallinearizer.hh:344
Implementation & asImp_()
Definition: fvbasefdlocallinearizer.hh:322
FvBaseFdLocalLinearizer()
Definition: fvbasefdlocallinearizer.hh:173
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilmodel.hh:72
Definition: blackoilboundaryratevector.hh:37
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:242
This file provides the infrastructure to retrieve run-time parameters.
The Opm property system, traits with inheritance.
GetPropType< TypeTag, Properties::Scalar > type
Definition: fvbasefdlocallinearizer.hh:89
Definition: fvbasefdlocallinearizer.hh:64
GetPropType< TypeTag, Properties::Scalar > type
Definition: fvbasefdlocallinearizer.hh:73
Representation of a function evaluation and all necessary derivatives with regard to the intensive qu...
Definition: fvbaseproperties.hh:83
The type of the local linearizer.
Definition: fvbaseproperties.hh:114
Definition: fvbasefdlocallinearizer.hh:62
Definition: fvbasefdlocallinearizer.hh:58
a tag to mark properties as undefined
Definition: propertysystem.hh:40