28#ifndef EWOMS_FV_BASE_AD_LOCAL_LINEARIZER_HH
29#define EWOMS_FV_BASE_AD_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/Valgrind.hpp>
38#include <opm/material/densead/Math.hpp>
49template<
class TypeTag>
53namespace Opm::Properties {
62template<
class TypeTag>
67template<
class TypeTag>
76 using type = DenseAd::Evaluation<Scalar, numEq>;
91template<
class TypeTag>
92class FvBaseAdLocalLinearizer
104 using Element =
typename GridView::template Codim<0>::Entity;
108 using ScalarVectorBlock = Dune::FieldVector<Scalar, numEq>;
112 using ScalarLocalBlockVector = Dune::BlockVector<ScalarVectorBlock>;
113 using ScalarLocalBlockMatrix = Dune::Matrix<ScalarMatrixBlock>;
116 FvBaseAdLocalLinearizer() =
default;
120 FvBaseAdLocalLinearizer(
const FvBaseAdLocalLinearizer&) =
delete;
136 void init(Simulator& simulator)
138 simulatorPtr_ = &simulator;
139 internalElemContext_ = std::make_unique<ElementContext>(simulator);
155 linearize(*internalElemContext_, element);
173 void linearize(ElementContext& elemCtx,
const Element& elem)
175 elemCtx.updateStencil(elem);
176 elemCtx.updateAllIntensiveQuantities();
179 model_().updatePVWeights(elemCtx);
185 const unsigned numPrimaryDof = elemCtx.numPrimaryDof(0);
186 for (
unsigned focusDofIdx = 0; focusDofIdx < numPrimaryDof; ++focusDofIdx) {
187 elemCtx.setFocusDofIndex(focusDofIdx);
188 elemCtx.updateAllExtensiveQuantities();
191 localResidual_.eval(elemCtx);
204 {
return localResidual_; }
210 {
return localResidual_; }
220 const ScalarMatrixBlock&
jacobian(
unsigned domainScvIdx,
unsigned rangeScvIdx)
const
221 {
return jacobian_[domainScvIdx][rangeScvIdx]; }
228 const ScalarVectorBlock&
residual(
unsigned dofIdx)
const
229 {
return residual_[dofIdx]; }
232 Implementation& asImp_()
233 {
return *
static_cast<Implementation*
>(
this); }
235 const Implementation& asImp_()
const
236 {
return *
static_cast<const Implementation*
>(
this); }
238 const Simulator& simulator_()
const
239 {
return *simulatorPtr_; }
241 const Problem& problem_()
const
242 {
return simulatorPtr_->problem(); }
244 const Model& model_()
const
245 {
return simulatorPtr_->model(); }
253 const std::size_t numDof = elemCtx.numDof(0);
254 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(0);
256 residual_.resize(numDof);
257 if (jacobian_.N() != numDof || jacobian_.M() != numPrimaryDof) {
258 jacobian_.setSize(numDof, numPrimaryDof);
276 unsigned focusDofIdx)
278 const auto& resid = localResidual_.residual();
280 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
281 residual_[focusDofIdx][eqIdx] = resid[focusDofIdx][eqIdx].value();
284 const std::size_t numDof = elemCtx.numDof(0);
285 for (
unsigned dofIdx = 0; dofIdx < numDof; ++dofIdx) {
286 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
287 for (
unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
292 jacobian_[dofIdx][focusDofIdx][eqIdx][pvIdx] = resid[dofIdx][eqIdx].derivative(pvIdx);
293 Valgrind::CheckDefined(jacobian_[dofIdx][focusDofIdx][eqIdx][pvIdx]);
301 std::unique_ptr<ElementContext> internalElemContext_{};
303 LocalResidual localResidual_{};
305 ScalarLocalBlockVector residual_{};
306 ScalarLocalBlockMatrix jacobian_{};
Calculates the local residual and its Jacobian for a single element of the grid.
Definition fvbaseadlocallinearizer.hh:93
const ScalarVectorBlock & residual(unsigned dofIdx) const
Returns the local residual of a sub-control volume.
Definition fvbaseadlocallinearizer.hh:228
void linearize(const Element &element)
Compute an element's local Jacobian matrix and evaluate its residual.
Definition fvbaseadlocallinearizer.hh:153
LocalResidual & localResidual()
Return reference to the local residual.
Definition fvbaseadlocallinearizer.hh:203
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:275
void reset_(const ElementContext &)
Reset the all relevant internal attributes to 0.
Definition fvbaseadlocallinearizer.hh:265
void resize_(const ElementContext &elemCtx)
Resize all internal attributes to the size of the element.
Definition fvbaseadlocallinearizer.hh:251
const LocalResidual & localResidual() const
Return reference to the local residual.
Definition fvbaseadlocallinearizer.hh:209
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:173
const ScalarMatrixBlock & jacobian(unsigned domainScvIdx, unsigned rangeScvIdx) const
Returns the local Jacobian matrix of the residual of a sub-control volume.
Definition fvbaseadlocallinearizer.hh:220
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
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