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>
63
64// set the properties to be spliced in
65template<class TypeTag>
66struct LocalLinearizer<TypeTag, TTag::FiniteDifferenceLocalLinearizer>
68
69template<class TypeTag>
70struct Evaluation<TypeTag, TTag::FiniteDifferenceLocalLinearizer>
72
74template<class TypeTag>
75struct BaseEpsilon<TypeTag, TTag::FiniteDifferenceLocalLinearizer>
76{
78 static constexpr type value = std::max<type>(0.9123e-10, std::numeric_limits<type>::epsilon()*1.23e3);
79};
80
81} // namespace Opm::Properties
82
83namespace Opm::Parameters {
84
92struct NumericDifferenceMethod { static constexpr int value = +1; };
93
94} // namespace Opm::Parameters
95
96namespace Opm {
97
134template<class TypeTag>
136{
137private:
147 using Element = typename GridView::template Codim<0>::Entity;
148
149 enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
150
151 // extract local matrices from jacobian matrix for consistency
153 using ScalarVectorBlock = Dune::FieldVector<Scalar, numEq>;
154
155 using ScalarLocalBlockVector = Dune::BlockVector<ScalarVectorBlock>;
156 using ScalarLocalBlockMatrix = Dune::Matrix<ScalarMatrixBlock>;
157
158 using LocalEvalBlockVector = typename LocalResidual::LocalEvalBlockVector;
159
160#if __GNUC__ == 4 && __GNUC_MINOR__ <= 6
161public:
162 // make older GCCs happy by providing a public copy constructor (this is necessary
163 // for their implementation of std::vector, although the method is never called...)
166 {}
167
168#else
169 // copying local residual objects around is a very bad idea, so we explicitly prevent
170 // it...
172#endif
173public:
176 { }
177
179 { delete internalElemContext_; }
180
184 static void registerParameters()
185 {
186 Parameters::Register<Parameters::NumericDifferenceMethod>
187 ("The method used for numeric differentiation (-1: backward "
188 "differences, 0: central differences, 1: forward differences)");
189 }
190
199 void init(Simulator& simulator)
200 {
201 simulatorPtr_ = &simulator;
203 internalElemContext_ = new ElementContext(simulator);
204 }
205
217 void linearize(const Element& element)
218 {
220 }
221
236 void linearize(ElementContext& elemCtx, const Element& elem)
237 {
238 elemCtx.updateAll(elem);
239
240 // update the weights of the primary variables for the context
241 model_().updatePVWeights(elemCtx);
242
243 resize_(elemCtx);
244 reset_(elemCtx);
245
246 // calculate the local residual
247 localResidual_.eval(residual_, elemCtx);
248
249 // calculate the local jacobian matrix
250 size_t numPrimaryDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
251 for (unsigned dofIdx = 0; dofIdx < numPrimaryDof; dofIdx++) {
252 for (unsigned pvIdx = 0; pvIdx < numEq; pvIdx++) {
253 asImp_().evalPartialDerivative_(elemCtx, dofIdx, pvIdx);
254
255 // incorporate the partial derivatives into the local Jacobian matrix
256 updateLocalJacobian_(elemCtx, dofIdx, pvIdx);
257 }
258 }
259 }
260
265 static Scalar baseEpsilon()
266 { return getPropValue<TypeTag, Properties::BaseEpsilon>(); }
267
279 Scalar numericEpsilon(const ElementContext& elemCtx,
280 unsigned dofIdx,
281 unsigned pvIdx) const
282 {
283 unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, /*timeIdx=*/0);
284 Scalar pvWeight = elemCtx.model().primaryVarWeight(globalIdx, pvIdx);
285 assert(pvWeight > 0 && std::isfinite(pvWeight));
286 Valgrind::CheckDefined(pvWeight);
287
288 return baseEpsilon()/pvWeight;
289 }
290
294 LocalResidual& localResidual()
295 { return localResidual_; }
296
300 const LocalResidual& localResidual() const
301 { return localResidual_; }
302
311 const ScalarMatrixBlock& jacobian(unsigned domainScvIdx, unsigned rangeScvIdx) const
312 { return jacobian_[domainScvIdx][rangeScvIdx]; }
313
319 const ScalarVectorBlock& residual(unsigned dofIdx) const
320 { return residual_[dofIdx]; }
321
322protected:
323 Implementation& asImp_()
324 { return *static_cast<Implementation*>(this); }
325 const Implementation& asImp_() const
326 { return *static_cast<const Implementation*>(this); }
327
328 const Simulator& simulator_() const
329 { return *simulatorPtr_; }
330 const Problem& problem_() const
331 { return simulatorPtr_->problem(); }
332 const Model& model_() const
333 { return simulatorPtr_->model(); }
334
339 {
340 static int diff = Parameters::Get<Parameters::NumericDifferenceMethod>();
341 return diff;
342 }
343
348 void resize_(const ElementContext& elemCtx)
349 {
350 size_t numDof = elemCtx.numDof(/*timeIdx=*/0);
351 size_t numPrimaryDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
352
353 residual_.resize(numDof);
354 jacobian_.setSize(numDof, numPrimaryDof);
355
356 derivResidual_.resize(numDof);
357 }
358
362 void reset_(const ElementContext& elemCtx)
363 {
364 size_t numDof = elemCtx.numDof(/*timeIdx=*/0);
365 size_t numPrimaryDof = elemCtx.numPrimaryDof(/*timeIdx=*/0);
366 for (unsigned primaryDofIdx = 0; primaryDofIdx < numPrimaryDof; ++ primaryDofIdx)
367 for (unsigned dof2Idx = 0; dof2Idx < numDof; ++ dof2Idx)
368 jacobian_[dof2Idx][primaryDofIdx] = 0.0;
369
370 for (unsigned primaryDofIdx = 0; primaryDofIdx < numDof; ++ primaryDofIdx)
371 residual_[primaryDofIdx] = 0.0;
372 }
373
416 void evalPartialDerivative_(ElementContext& elemCtx,
417 unsigned dofIdx,
418 unsigned pvIdx)
419 {
420 // save all quantities which depend on the specified primary
421 // variable at the given sub control volume
422 elemCtx.stashIntensiveQuantities(dofIdx);
423
424 PrimaryVariables priVars(elemCtx.primaryVars(dofIdx, /*timeIdx=*/0));
425 Scalar eps = asImp_().numericEpsilon(elemCtx, dofIdx, pvIdx);
426 Scalar delta = 0.0;
427
428 if (numericDifferenceMethod_() >= 0) {
429 // we are not using backward differences, i.e. we need to
430 // calculate f(x + \epsilon)
431
432 // deflect primary variables
433 priVars[pvIdx] += eps;
434 delta += eps;
435
436 // calculate the deflected residual
437 elemCtx.updateIntensiveQuantities(priVars, dofIdx, /*timeIdx=*/0);
438 elemCtx.updateAllExtensiveQuantities();
439 localResidual_.eval(derivResidual_, elemCtx);
440 }
441 else {
442 // we are using backward differences, i.e. we don't need
443 // to calculate f(x + \epsilon) and we can recycle the
444 // (already calculated) residual f(x)
446 }
447
448 if (numericDifferenceMethod_() <= 0) {
449 // we are not using forward differences, i.e. we don't
450 // need to calculate f(x - \epsilon)
451
452 // deflect the primary variables
453 priVars[pvIdx] -= delta + eps;
454 delta += eps;
455
456 // calculate the deflected residual again, this time we use the local
457 // residual's internal storage.
458 elemCtx.updateIntensiveQuantities(priVars, dofIdx, /*timeIdx=*/0);
459 elemCtx.updateAllExtensiveQuantities();
460 localResidual_.eval(elemCtx);
461
462 derivResidual_ -= localResidual_.residual();
463 }
464 else {
465 // we are using forward differences, i.e. we don't need to
466 // calculate f(x - \epsilon) and we can recycle the
467 // (already calculated) residual f(x)
469 }
470
471 assert(delta > 0);
472
473 // divide difference in residuals by the magnitude of the
474 // deflections between the two function evaluation
475 derivResidual_ /= delta;
476
477 // restore the original state of the element's volume
478 // variables
479 elemCtx.restoreIntensiveQuantities(dofIdx);
480
481#ifndef NDEBUG
482 for (unsigned i = 0; i < derivResidual_.size(); ++i)
483 Valgrind::CheckDefined(derivResidual_[i]);
484#endif
485 }
486
492 void updateLocalJacobian_(const ElementContext& elemCtx,
493 unsigned focusDofIdx,
494 unsigned pvIdx)
495 {
496 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 Model *modelPtr_;
511
512 ElementContext *internalElemContext_;
513
514 LocalEvalBlockVector residual_;
515 LocalEvalBlockVector derivResidual_;
516 ScalarLocalBlockMatrix jacobian_;
517
518 LocalResidual localResidual_;
519};
520
521} // namespace Opm
522
523#endif
Calculates the Jacobian of the local residual for finite volume spatial discretizations using a finit...
Definition: fvbasefdlocallinearizer.hh:136
const Implementation & asImp_() const
Definition: fvbasefdlocallinearizer.hh:325
LocalEvalBlockVector derivResidual_
Definition: fvbasefdlocallinearizer.hh:515
void evalPartialDerivative_(ElementContext &elemCtx, unsigned dofIdx, unsigned pvIdx)
Compute the partial derivatives of a context's residual functions.
Definition: fvbasefdlocallinearizer.hh:416
void linearize(const Element &element)
Compute an element's local Jacobian matrix and evaluate its residual.
Definition: fvbasefdlocallinearizer.hh:217
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:279
LocalResidual localResidual_
Definition: fvbasefdlocallinearizer.hh:518
static Scalar baseEpsilon()
Returns the unweighted epsilon value used to calculate the local derivatives.
Definition: fvbasefdlocallinearizer.hh:265
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:332
ScalarLocalBlockMatrix jacobian_
Definition: fvbasefdlocallinearizer.hh:516
static int numericDifferenceMethod_()
Returns the numeric difference method which is applied.
Definition: fvbasefdlocallinearizer.hh:338
ElementContext * internalElemContext_
Definition: fvbasefdlocallinearizer.hh:512
const ScalarMatrixBlock & jacobian(unsigned domainScvIdx, unsigned rangeScvIdx) const
Returns the local Jacobian matrix of the residual of a sub-control volume.
Definition: fvbasefdlocallinearizer.hh:311
~FvBaseFdLocalLinearizer()
Definition: fvbasefdlocallinearizer.hh:178
static void registerParameters()
Register all run-time parameters for the local jacobian.
Definition: fvbasefdlocallinearizer.hh:184
Model * modelPtr_
Definition: fvbasefdlocallinearizer.hh:510
const Simulator & simulator_() const
Definition: fvbasefdlocallinearizer.hh:328
const LocalResidual & localResidual() const
Return reference to the local residual.
Definition: fvbasefdlocallinearizer.hh:300
Simulator * simulatorPtr_
Definition: fvbasefdlocallinearizer.hh:509
void init(Simulator &simulator)
Initialize the local Jacobian object.
Definition: fvbasefdlocallinearizer.hh:199
const ScalarVectorBlock & residual(unsigned dofIdx) const
Returns the local residual of a sub-control volume.
Definition: fvbasefdlocallinearizer.hh:319
const Problem & problem_() const
Definition: fvbasefdlocallinearizer.hh:330
LocalEvalBlockVector residual_
Definition: fvbasefdlocallinearizer.hh:514
void linearize(ElementContext &elemCtx, const Element &elem)
Compute an element's local Jacobian matrix and evaluate its residual.
Definition: fvbasefdlocallinearizer.hh:236
void reset_(const ElementContext &elemCtx)
Reset the all relevant internal attributes to 0.
Definition: fvbasefdlocallinearizer.hh:362
LocalResidual & localResidual()
Return reference to the local residual.
Definition: fvbasefdlocallinearizer.hh:294
void resize_(const ElementContext &elemCtx)
Resize all internal attributes to the size of the element.
Definition: fvbasefdlocallinearizer.hh:348
Implementation & asImp_()
Definition: fvbasefdlocallinearizer.hh:323
FvBaseFdLocalLinearizer()
Definition: fvbasefdlocallinearizer.hh:174
Declare the properties used by the infrastructure code of the finite volume discretizations.
Definition: blackoilnewtonmethodparameters.hh:31
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:235
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:92
static constexpr int value
Definition: fvbasefdlocallinearizer.hh:92
GetPropType< TypeTag, Properties::Scalar > type
Definition: fvbasefdlocallinearizer.hh:77
Definition: fvbasefdlocallinearizer.hh:62
GetPropType< TypeTag, Properties::Scalar > type
Definition: fvbasefdlocallinearizer.hh:71
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:58
a tag to mark properties as undefined
Definition: propertysystem.hh:40