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>
75 using type = DenseAd::Evaluation<Scalar, numDerivatives>;
90template<
class TypeTag>
91class FvBaseAdLocalLinearizer
103 using Element =
typename GridView::template Codim<0>::Entity;
107 using ScalarVectorBlock = Dune::FieldVector<Scalar, numEq>;
111 using ScalarLocalBlockVector = Dune::BlockVector<ScalarVectorBlock>;
112 using ScalarLocalBlockMatrix = Dune::Matrix<ScalarMatrixBlock>;
115 FvBaseAdLocalLinearizer() =
default;
119 FvBaseAdLocalLinearizer(
const FvBaseAdLocalLinearizer&) =
delete;
135 void init(Simulator& simulator)
137 simulatorPtr_ = &simulator;
138 internalElemContext_ = std::make_unique<ElementContext>(simulator);
154 linearize(*internalElemContext_, element);
172 void linearize(ElementContext& elemCtx,
const Element& elem)
174 elemCtx.updateStencil(elem);
175 elemCtx.updateAllIntensiveQuantities();
178 model_().updatePVWeights(elemCtx);
184 const unsigned numPrimaryDof = elemCtx.numPrimaryDof(0);
185 for (
unsigned focusDofIdx = 0; focusDofIdx < numPrimaryDof; ++focusDofIdx) {
186 elemCtx.setFocusDofIndex(focusDofIdx);
187 elemCtx.updateAllExtensiveQuantities();
190 localResidual_.eval(elemCtx);
203 {
return localResidual_; }
209 {
return localResidual_; }
219 const ScalarMatrixBlock&
jacobian(
unsigned domainScvIdx,
unsigned rangeScvIdx)
const
220 {
return jacobian_[domainScvIdx][rangeScvIdx]; }
227 const ScalarVectorBlock&
residual(
unsigned dofIdx)
const
228 {
return residual_[dofIdx]; }
231 Implementation& asImp_()
232 {
return *
static_cast<Implementation*
>(
this); }
234 const Implementation& asImp_()
const
235 {
return *
static_cast<const Implementation*
>(
this); }
237 const Simulator& simulator_()
const
238 {
return *simulatorPtr_; }
240 const Problem& problem_()
const
241 {
return simulatorPtr_->problem(); }
243 const Model& model_()
const
244 {
return simulatorPtr_->model(); }
252 const std::size_t numDof = elemCtx.numDof(0);
253 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(0);
255 residual_.resize(numDof);
256 if (jacobian_.N() != numDof || jacobian_.M() != numPrimaryDof) {
257 jacobian_.setSize(numDof, numPrimaryDof);
275 unsigned focusDofIdx)
277 const auto& resid = localResidual_.residual();
279 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
280 residual_[focusDofIdx][eqIdx] = resid[focusDofIdx][eqIdx].value();
283 const std::size_t numDof = elemCtx.numDof(0);
284 for (
unsigned dofIdx = 0; dofIdx < numDof; ++dofIdx) {
285 for (
unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
286 for (
unsigned pvIdx = 0; pvIdx < numEq; ++pvIdx) {
291 jacobian_[dofIdx][focusDofIdx][eqIdx][pvIdx] = resid[dofIdx][eqIdx].derivative(pvIdx);
292 Valgrind::CheckDefined(jacobian_[dofIdx][focusDofIdx][eqIdx][pvIdx]);
300 std::unique_ptr<ElementContext> internalElemContext_{};
302 LocalResidual localResidual_{};
304 ScalarLocalBlockVector residual_{};
305 ScalarLocalBlockMatrix jacobian_{};
Calculates the local residual and its Jacobian for a single element of the grid.
Definition fvbaseadlocallinearizer.hh:92
const ScalarVectorBlock & residual(unsigned dofIdx) const
Returns the local residual of a sub-control volume.
Definition fvbaseadlocallinearizer.hh:227
void linearize(const Element &element)
Compute an element's local Jacobian matrix and evaluate its residual.
Definition fvbaseadlocallinearizer.hh:152
LocalResidual & localResidual()
Return reference to the local residual.
Definition fvbaseadlocallinearizer.hh:202
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:274
void reset_(const ElementContext &)
Reset the all relevant internal attributes to 0.
Definition fvbaseadlocallinearizer.hh:264
void resize_(const ElementContext &elemCtx)
Resize all internal attributes to the size of the element.
Definition fvbaseadlocallinearizer.hh:250
const LocalResidual & localResidual() const
Return reference to the local residual.
Definition fvbaseadlocallinearizer.hh:208
void init(Simulator &simulator)
Initialize the local Jacobian object.
Definition fvbaseadlocallinearizer.hh:135
static void registerParameters()
Register all run-time parameters for the local jacobian.
Definition fvbaseadlocallinearizer.hh:124
void linearize(ElementContext &elemCtx, const Element &elem)
Compute an element's local Jacobian matrix and evaluate its residual.
Definition fvbaseadlocallinearizer.hh:172
const ScalarMatrixBlock & jacobian(unsigned domainScvIdx, unsigned rangeScvIdx) const
Returns the local Jacobian matrix of the residual of a sub-control volume.
Definition fvbaseadlocallinearizer.hh:219
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:82
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:45
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