28#ifndef EWOMS_FV_BASE_LOCAL_RESIDUAL_HH
29#define EWOMS_FV_BASE_LOCAL_RESIDUAL_HH
31#include <dune/common/classname.hh>
32#include <dune/common/fvector.hh>
34#include <dune/grid/common/geometry.hh>
36#include <dune/istl/bvector.hh>
38#include <opm/material/common/MathToolbox.hpp>
39#include <opm/material/common/Valgrind.hpp>
61template<
class TypeTag>
68 using Element =
typename GridView::template Codim<0>::Entity;
86 using EvalVector = Dune::FieldVector<Evaluation, numEq>;
89 using LocalEvalBlockVector = Dune::BlockVector<EvalVector,
aligned_allocator<EvalVector,
alignof(EvalVector)>>;
107 {
return internalResidual_; }
116 {
return internalResidual_[dofIdx]; }
128 void eval(
const Problem& problem,
const Element& element)
130 ElementContext elemCtx(problem);
131 elemCtx.updateAll(element);
144 void eval(ElementContext& elemCtx)
146 std::size_t numDof = elemCtx.numDof(0);
147 internalResidual_.resize(numDof);
148 asImp_().eval(internalResidual_, elemCtx);
159 ElementContext& elemCtx)
const
166 asImp_().evalFluxes(
residual, elemCtx, 0);
169 asImp_().evalVolumeTerms_(
residual, elemCtx);
172 asImp_().evalBoundary_(
residual, elemCtx, 0);
174 if (useVolumetricResidual) {
177 const std::size_t numDof = elemCtx.numDof(0);
178 for (
unsigned dofIdx = 0; dofIdx < numDof; ++dofIdx) {
179 if (elemCtx.dofTotalVolume(dofIdx, 0) > 0.0) {
181 const Scalar dofVolume = elemCtx.dofTotalVolume(dofIdx, 0);
183 assert(std::isfinite(dofVolume));
184 Valgrind::CheckDefined(dofVolume);
206 const ElementContext& elemCtx,
215 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(0);
216 for (
unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
221 elemCtx.stencil(
timeIdx).subControlVolume(dofIdx).volume() *
222 elemCtx.intensiveQuantities(dofIdx,
timeIdx).extrusionFactor();
228 if (dofIdx == elemCtx.focusDofIndex()) {
229 asImp_().computeStorage(
storage[dofIdx],
239 Dune::FieldVector<Scalar, numEq> tmp;
240 asImp_().computeStorage(tmp,
254 if (elemCtx.enableStorageCache()) {
255 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(
timeIdx);
256 for (
unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
257 const unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx,
timeIdx);
258 const auto& cachedStorage = elemCtx.model().cachedStorage(globalDofIdx,
timeIdx);
267 Dune::FieldVector<Scalar, numEq> tmp;
268 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(0);
269 for (
unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
271 asImp_().computeStorage(tmp,
275 tmp *= elemCtx.stencil(
timeIdx).subControlVolume(dofIdx).volume() *
276 elemCtx.intensiveQuantities(dofIdx,
timeIdx).extrusionFactor();
285 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(0);
286 for (
unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
303 const ElementContext& elemCtx,
308 const auto& stencil = elemCtx.stencil(
timeIdx);
310 const std::size_t numInteriorFaces = elemCtx.numInteriorFaces(
timeIdx);
313 const unsigned i =
face.interiorIndex();
314 const unsigned j =
face.exteriorIndex();
316 Valgrind::SetUndefined(
flux);
318 Valgrind::CheckDefined(
flux);
325 const Scalar alpha = elemCtx.extensiveQuantities(
scvfIdx,
timeIdx).extrusionFactor() *
face.area();
326 Valgrind::CheckDefined(alpha);
356 const std::size_t numDof = elemCtx.numDof(
timeIdx);
357 for (
unsigned i = 0; i < numDof; ++i) {
358 for (
unsigned j = 0; j < numEq; ++j) {
360 Valgrind::CheckDefined(
residual[i][j]);
379 const ElementContext&,
383 throw std::logic_error(
"Not implemented: The local residual " + Dune::className<Implementation>() +
384 " does not implement the required method 'computeStorage()'");
395 const ElementContext&,
399 throw std::logic_error(
"Not implemented: The local residual " + Dune::className<Implementation>() +
400 " does not implement the required method 'computeFlux()'");
410 const ElementContext&,
414 throw std::logic_error(
"Not implemented: The local residual " + Dune::className<Implementation>() +
415 " does not implement the required method 'computeSource()'");
423 const ElementContext& elemCtx,
426 if (!elemCtx.onBoundary()) {
437 const std::size_t numBoundaryFaces =
boundaryCtx.numBoundaryFaces(0);
438 for (
unsigned faceIdx = 0; faceIdx < numBoundaryFaces; ++faceIdx,
boundaryCtx.increment()) {
449 const std::size_t numDof = elemCtx.numDof(0);
450 for (
unsigned i = 0; i < numDof; ++i) {
451 for (
unsigned j = 0; j < numEq; ++j) {
453 Valgrind::CheckDefined(
residual[i][j]);
468 BoundaryRateVector values;
470 Valgrind::SetUndefined(values);
472 Valgrind::CheckDefined(values);
475 const unsigned dofIdx = stencil.boundaryFace(
boundaryFaceIdx).interiorIndex();
481 Valgrind::CheckDefined(values[
eqIdx]);
496 ElementContext& elemCtx)
const
506 const std::size_t numPrimaryDof = elemCtx.numPrimaryDof(0);
507 for (
unsigned dofIdx = 0; dofIdx < numPrimaryDof; ++dofIdx) {
508 const Scalar extrusionFactor =
509 elemCtx.intensiveQuantities(dofIdx, 0).extrusionFactor();
510 Valgrind::CheckDefined(extrusionFactor);
512 assert(extrusionFactor > 0.0);
514 elemCtx.stencil(0).subControlVolume(dofIdx).volume() * extrusionFactor;
522 if (!extensiveStorageTerm &&
523 !std::is_same_v<Scalar, Evaluation> &&
524 dofIdx != elemCtx.focusDofIndex())
526 asImp_().computeStorage(
tmp2, elemCtx, dofIdx, 0);
532 asImp_().computeStorage(tmp, elemCtx, dofIdx, 0);
536 Valgrind::CheckDefined(tmp);
542 if (elemCtx.enableStorageCache()) {
543 const auto& model = elemCtx.model();
544 const unsigned globalDofIdx = elemCtx.globalSpaceIndex(dofIdx, 0);
545 if (model.newtonMethod().numIterations() == 0 &&
546 !elemCtx.haveStashedIntensiveQuantities())
548 if (!elemCtx.problem().recycleFirstIterationStorage()) {
553 elemCtx.updatePrimaryIntensiveQuantities(1);
554 asImp_().computeStorage(
tmp2, elemCtx, dofIdx, 1);
569 Valgrind::CheckDefined(
tmp2);
571 model.updateCachedStorage(globalDofIdx, 1,
tmp2);
577 tmp2 = model.cachedStorage(globalDofIdx, 1);
578 Valgrind::CheckDefined(
tmp2);
585 asImp_().computeStorage(
tmp2, elemCtx, dofIdx, 1);
586 Valgrind::CheckDefined(
tmp2);
591 const double dt = elemCtx.simulator().timeStepSize();
599 Valgrind::CheckDefined(
residual[dofIdx]);
602 asImp_().computeSource(
sourceRate, elemCtx, dofIdx, 0);
607 if (!extensiveStorageTerm &&
608 !std::is_same_v<Scalar, Evaluation> &&
609 dofIdx != elemCtx.focusDofIndex())
622 Valgrind::CheckDefined(
residual[dofIdx]);
627 std::size_t numDof = elemCtx.numDof(0);
628 for (
unsigned i = 0; i < numDof; ++i) {
629 for (
unsigned j = 0; j < numEq; ++j) {
631 Valgrind::CheckDefined(
residual[i][j]);
638 Implementation& asImp_()
639 {
return *
static_cast<Implementation*
>(
this); }
641 const Implementation& asImp_()
const
642 {
return *
static_cast<const Implementation*
>(
this); }
644 LocalEvalBlockVector internalResidual_;
This is a stand-alone version of boost::alignment::aligned_allocator from Boost 1....
Element-wise caculation of the residual matrix for models based on a finite volume spatial discretiza...
Definition fvbaselocalresidual.hh:63
void evalVolumeTerms_(LocalEvalBlockVector &residual, ElementContext &elemCtx) const
Add the change in the storage terms and the source term to the local residual of all sub-control volu...
Definition fvbaselocalresidual.hh:495
static void registerParameters()
Register all run-time parameters for the local residual.
Definition fvbaselocalresidual.hh:99
void computeSource(RateVector &, const ElementContext &, unsigned, unsigned) const
Calculate the source term of the equation.
Definition fvbaselocalresidual.hh:409
void evalFluxes(LocalEvalBlockVector &residual, const ElementContext &elemCtx, unsigned timeIdx) const
Add the flux term to a local residual.
Definition fvbaselocalresidual.hh:302
void computeStorage(EqVector &, const ElementContext &, unsigned, unsigned) const
Evaluate the amount all conservation quantities (e.g.
Definition fvbaselocalresidual.hh:378
void evalStorage(LocalEvalBlockVector &storage, const ElementContext &elemCtx, unsigned timeIdx) const
Calculate the amount of all conservation quantities stored in all element's sub-control volumes for a...
Definition fvbaselocalresidual.hh:205
void eval(const Problem &problem, const Element &element)
Compute the local residual, i.e.
Definition fvbaselocalresidual.hh:128
void evalBoundarySegment_(LocalEvalBlockVector &residual, const BoundaryContext &boundaryCtx, unsigned boundaryFaceIdx, unsigned timeIdx) const
Evaluate all boundary conditions for a single sub-control volume face to the local residual.
Definition fvbaselocalresidual.hh:463
void eval(LocalEvalBlockVector &residual, ElementContext &elemCtx) const
Compute the local residual, i.e.
Definition fvbaselocalresidual.hh:158
void computeFlux(RateVector &, const ElementContext &, unsigned, unsigned) const
Evaluates the total mass flux of all conservation quantities over a face of a sub-control volume.
Definition fvbaselocalresidual.hh:394
void evalBoundary_(LocalEvalBlockVector &residual, const ElementContext &elemCtx, unsigned timeIdx) const
Evaluate the boundary conditions of an element.
Definition fvbaselocalresidual.hh:422
const LocalEvalBlockVector & residual() const
Return the result of the eval() call using internal storage.
Definition fvbaselocalresidual.hh:106
void eval(ElementContext &elemCtx)
Compute the local residual, i.e.
Definition fvbaselocalresidual.hh:144
const EvalVector & residual(unsigned dofIdx) const
Return the result of the eval() call using internal storage.
Definition fvbaselocalresidual.hh:115
Definition alignedallocator.hh:97
Declare the properties used by the infrastructure code of the finite volume discretizations.
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilbioeffectsmodules.hh:43
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:240
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
This file provides the infrastructure to retrieve run-time parameters.