28#ifndef EWOMS_FV_BASE_FD_LOCAL_LINEARIZER_HH
29#define EWOMS_FV_BASE_FD_LOCAL_LINEARIZER_HH
31#include <dune/common/fmatrix.hh>
32#include <dune/common/fvector.hh>
34#include <dune/istl/bvector.hh>
35#include <dune/istl/matrix.hh>
37#include <opm/material/common/MathToolbox.hpp>
38#include <opm/material/common/Valgrind.hpp>
55template<
class TypeTag>
60namespace Opm::Properties {
68template<
class TypeTag,
class MyTypeTag>
72template<
class TypeTag>
76template<
class TypeTag>
81template<
class TypeTag>
85 static constexpr type value = std::max<type>(0.9123e-10,
86 std::numeric_limits<type>::epsilon() * 1.23e3);
91namespace Opm::Parameters {
142template<
class TypeTag>
143class FvBaseFdLocalLinearizer
155 using Element =
typename GridView::template Codim<0>::Entity;
161 using ScalarVectorBlock = Dune::FieldVector<Scalar, numEq>;
163 using ScalarLocalBlockVector = Dune::BlockVector<ScalarVectorBlock>;
164 using ScalarLocalBlockMatrix = Dune::Matrix<ScalarMatrixBlock>;
166 using LocalEvalBlockVector =
typename LocalResidual::LocalEvalBlockVector;
169 FvBaseFdLocalLinearizer() =
default;
172 FvBaseFdLocalLinearizer(
const FvBaseFdLocalLinearizer&) =
delete;
180 (
"The method used for numeric differentiation (-1: backward "
181 "differences, 0: central differences, 1: forward differences)");
192 void init(Simulator& simulator)
194 simulatorPtr_ = &simulator;
195 internalElemContext_ = std::make_unique<ElementContext>(simulator);
211 linearize(*internalElemContext_, element);
229 void linearize(ElementContext& elemCtx,
const Element& elem)
231 elemCtx.updateAll(elem);
234 model_().updatePVWeights(elemCtx);
240 localResidual_.eval(residual_, elemCtx);
243 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(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);
274 unsigned pvIdx)
const
276 const unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, 0);
277 const Scalar pvWeight = elemCtx.model().primaryVarWeight(globalIdx, pvIdx);
278 assert(pvWeight > 0 && std::isfinite(pvWeight));
279 Valgrind::CheckDefined(pvWeight);
288 {
return localResidual_; }
294 {
return localResidual_; }
304 const ScalarMatrixBlock&
jacobian(
unsigned domainScvIdx,
unsigned rangeScvIdx)
const
305 {
return jacobian_[domainScvIdx][rangeScvIdx]; }
312 const ScalarVectorBlock&
residual(
unsigned dofIdx)
const
313 {
return residual_[dofIdx]; }
316 Implementation& asImp_()
317 {
return *
static_cast<Implementation*
>(
this); }
319 const Implementation& asImp_()
const
320 {
return *
static_cast<const Implementation*
>(
this); }
322 const Simulator& simulator_()
const
323 {
return *simulatorPtr_; }
325 const Problem& problem_()
const
326 {
return simulatorPtr_->problem(); }
328 const Model& model_()
const
329 {
return simulatorPtr_->model(); }
346 const std::size_t numDof = elemCtx.numDof(0);
347 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(0);
349 residual_.resize(numDof);
350 jacobian_.setSize(numDof, numPrimaryDof);
352 derivResidual_.resize(numDof);
358 void reset_(
const ElementContext& elemCtx)
360 const std::size_t numDof = elemCtx.numDof(0);
361 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(0);
362 for (
unsigned primaryDofIdx = 0; primaryDofIdx < numPrimaryDof; ++primaryDofIdx) {
363 for (
unsigned dof2Idx = 0; dof2Idx < numDof; ++dof2Idx) {
364 jacobian_[dof2Idx][primaryDofIdx] = 0.0;
368 for (
unsigned primaryDofIdx = 0; primaryDofIdx < numDof; ++primaryDofIdx) {
369 residual_[primaryDofIdx] = 0.0;
421 elemCtx.stashIntensiveQuantities(dofIdx);
423 PrimaryVariables priVars(elemCtx.primaryVars(dofIdx, 0));
424 const Scalar eps = asImp_().numericEpsilon(elemCtx, dofIdx, pvIdx);
432 priVars[pvIdx] += eps;
436 elemCtx.updateIntensiveQuantities(priVars, dofIdx, 0);
437 elemCtx.updateAllExtensiveQuantities();
438 localResidual_.eval(derivResidual_, elemCtx);
444 derivResidual_ = residual_;
452 priVars[pvIdx] -= delta + eps;
457 elemCtx.updateIntensiveQuantities(priVars, dofIdx, 0);
458 elemCtx.updateAllExtensiveQuantities();
459 localResidual_.eval(elemCtx);
461 derivResidual_ -= localResidual_.residual();
467 derivResidual_ -= residual_;
474 derivResidual_ /= delta;
478 elemCtx.restoreIntensiveQuantities(dofIdx);
481 for (
unsigned i = 0; i < derivResidual_.size(); ++i) {
482 Valgrind::CheckDefined(derivResidual_[i]);
493 unsigned focusDofIdx,
496 const std::size_t numDof = elemCtx.numDof(0);
497 for (
unsigned dofIdx = 0; dofIdx < numDof; ++dofIdx) {
498 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
503 jacobian_[dofIdx][focusDofIdx][eqIdx][pvIdx] = derivResidual_[dofIdx][eqIdx];
504 Valgrind::CheckDefined(jacobian_[dofIdx][focusDofIdx][eqIdx][pvIdx]);
511 std::unique_ptr<ElementContext> internalElemContext_{};
513 LocalEvalBlockVector residual_{};
514 LocalEvalBlockVector derivResidual_{};
515 ScalarLocalBlockMatrix jacobian_{};
517 LocalResidual localResidual_{};
Calculates the Jacobian of the local residual for finite volume spatial discretizations using a finit...
Definition fvbasefdlocallinearizer.hh:144
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
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
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 LocalResidual & localResidual() const
Return reference to the local residual.
Definition fvbasefdlocallinearizer.hh:293
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
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
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
This file provides the infrastructure to retrieve run-time parameters.
void Register(const char *usageString)
Register a run-time parameter.
Definition parametersystem.hpp:292
auto Get(bool errorIfNotRegistered=true)
Retrieve a runtime parameter.
Definition parametersystem.hpp:187
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
Definition fvbasefdlocallinearizer.hh:69
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