25#ifndef DUMUX_STAGGERED_LOCAL_RESIDUAL_HH
26#define DUMUX_STAGGERED_LOCAL_RESIDUAL_HH
39template<
class TypeTag>
47 using Element =
typename GridView::template Codim<0>::Entity;
52 using FVElementGeometry =
typename GridGeometry::LocalView;
53 using SubControlVolume =
typename GridGeometry::SubControlVolume;
54 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
55 using FaceSubControlVolume =
typename GridGeometry::Traits::FaceSubControlVolume;
79 const FVElementGeometry& fvGeometry,
80 const ElementVolumeVariables& elemVolVars,
81 const ElementFaceVariables& elemFaceVars,
82 const ElementBoundaryTypes& bcTypes,
83 const ElementFluxVariablesCache& elemFluxVarsCache)
const
88 for (
auto&& scv : scvs(fvGeometry))
89 asImp().evalSourceForCellCenter(residual, this->
problem(), element, fvGeometry, elemVolVars, elemFaceVars, scv);
92 for (
auto&& scvf : scvfs(fvGeometry))
93 asImp().evalFluxForCellCenter(residual, this->
problem(), element, fvGeometry, elemVolVars, elemFaceVars, bcTypes, elemFluxVarsCache, scvf);
101 const Element& element,
102 const FVElementGeometry& fvGeometry,
103 const ElementVolumeVariables& elemVolVars,
104 const ElementFaceVariables& elemFaceVars,
105 const ElementBoundaryTypes& elemBcTypes,
106 const ElementFluxVariablesCache& elemFluxVarsCache,
107 const SubControlVolumeFace& scvf)
const
109 if (!scvf.boundary())
110 residual +=
asImp_().computeFluxForCellCenter(
problem, element, fvGeometry, elemVolVars, elemFaceVars, scvf, elemFluxVarsCache);
112 residual +=
asImp_().computeBoundaryFluxForCellCenter(
problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars, elemBcTypes, elemFluxVarsCache);
118 const Element& element,
119 const FVElementGeometry& fvGeometry,
120 const ElementVolumeVariables& curElemVolVars,
121 const ElementFaceVariables& curElemFaceVars,
122 const SubControlVolume& scv)
const
124 const auto curExtrusionFactor = curElemVolVars[scv].extrusionFactor();
127 auto source =
asImp_().computeSourceForCellCenter(
problem, element, fvGeometry, curElemVolVars, curElemFaceVars, scv);
134 const FVElementGeometry& fvGeometry,
135 const ElementVolumeVariables& prevElemVolVars,
136 const ElementVolumeVariables& curElemVolVars)
const
138 assert(timeLoop_ &&
"no time loop set for storage term evaluation");
141 for (
auto&& scv : scvs(fvGeometry))
142 asImp().evalStorageForCellCenter(storage,
problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, scv);
150 const Element &element,
151 const FVElementGeometry& fvGeometry,
152 const ElementVolumeVariables& prevElemVolVars,
153 const ElementVolumeVariables& curElemVolVars,
154 const SubControlVolume& scv)
const
157 const auto& curVolVars = curElemVolVars[scv];
158 const auto& prevVolVars = prevElemVolVars[scv];
166 auto prevCCStorage =
asImp_().computeStorageForCellCenter(
problem, scv, prevVolVars);
167 auto curCCStorage =
asImp_().computeStorageForCellCenter(
problem, scv, curVolVars);
169 prevCCStorage *= prevVolVars.extrusionFactor();
170 curCCStorage *= curVolVars.extrusionFactor();
172 storage = std::move(curCCStorage);
173 storage -= std::move(prevCCStorage);
181 template<
class... Args>
188 template<
class... Args>
203 const FVElementGeometry& fvGeometry,
204 const ElementVolumeVariables& elemVolVars,
205 const ElementFaceVariables& elemFaceVars,
206 const ElementBoundaryTypes& bcTypes,
207 const ElementFluxVariablesCache& elemFluxVarsCache,
208 const SubControlVolumeFace& scvf)
const
211 asImp().evalSourceForFace(residual, this->
problem(), element, fvGeometry, elemVolVars, elemFaceVars, scvf);
212 asImp().evalFluxForFace(residual, this->
problem(), element, fvGeometry, elemVolVars, elemFaceVars, bcTypes, elemFluxVarsCache, scvf);
220 const Element& element,
221 const FVElementGeometry& fvGeometry,
222 const ElementVolumeVariables& elemVolVars,
223 const ElementFaceVariables& elemFaceVars,
224 const ElementBoundaryTypes& elemBcTypes,
225 const ElementFluxVariablesCache& elemFluxVarsCache,
226 const SubControlVolumeFace& scvf)
const
228 if (!scvf.boundary())
229 residual +=
asImp_().computeFluxForFace(
problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, elemFluxVarsCache);
231 residual +=
asImp_().computeBoundaryFluxForFace(
problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars, elemBcTypes, elemFluxVarsCache);
237 const Element& element,
238 const FVElementGeometry& fvGeometry,
239 const ElementVolumeVariables& elemVolVars,
240 const ElementFaceVariables& elemFaceVars,
241 const SubControlVolumeFace& scvf)
const
244 auto source =
asImp_().computeSourceForFace(
problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars);
245 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
246 const auto extrusionFactor = elemVolVars[scv].extrusionFactor();
249 auto faceScvCenter = scvf.center() + scv.center();
250 faceScvCenter *= 0.5;
251 FaceSubControlVolume faceScv(faceScvCenter, 0.5*scv.volume());
259 const FVElementGeometry& fvGeometry,
260 const ElementVolumeVariables& prevElemVolVars,
261 const ElementVolumeVariables& curElemVolVars,
262 const ElementFaceVariables& prevElemFaceVars,
263 const ElementFaceVariables& curElemFaceVars,
264 const SubControlVolumeFace& scvf)
const
266 assert(timeLoop_ &&
"no time loop set for storage term evaluation");
268 asImp().evalStorageForFace(storage,
problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, prevElemFaceVars, curElemFaceVars, scvf);
275 const Element& element,
276 const FVElementGeometry& fvGeometry,
277 const ElementVolumeVariables& prevElemVolVars,
278 const ElementVolumeVariables& curElemVolVars,
279 const ElementFaceVariables& prevElemFaceVars,
280 const ElementFaceVariables& curElemFaceVars,
281 const SubControlVolumeFace& scvf)
const
283 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
285 auto storage =
asImp_().computeStorageForFace(
problem, scvf, curElemVolVars[scv], curElemFaceVars);
286 storage -=
asImp_().computeStorageForFace(
problem, scvf, prevElemVolVars[scv], prevElemFaceVars);
288 const auto extrusionFactor = curElemVolVars[scv].extrusionFactor();
291 auto faceScvCenter = scvf.center() + scv.center();
292 faceScvCenter *= 0.5;
293 FaceSubControlVolume faceScv(faceScvCenter, 0.5*scv.volume());
303 {
return !timeLoop_; }
307 {
return *problem_; }
312 {
return *
static_cast<Implementation*
>(
this); }
315 {
return *
static_cast<const Implementation*
>(
this); }
319 {
return *timeLoop_; }
322 {
return *timeLoop_; }
325 {
return *
static_cast<Implementation*
>(
this); }
328 {
return *
static_cast<const Implementation*
>(
this); }
331 const Problem* problem_;
Helper classes to compute the integration elements.
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Scalar volume(Shape shape, Scalar inscribedRadius)
Returns the volume of a given geometry based on the inscribed radius.
Definition: poreproperties.hh:73
Calculates the element-wise residual for the staggered FV scheme.
Definition: staggeredlocalresidual.hh:41
CellCenterResidualValue evalStorageForCellCenter(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &prevElemVolVars, const ElementVolumeVariables &curElemVolVars) const
Evaluate the storage terms for a cell center residual.
Definition: staggeredlocalresidual.hh:133
GetPropType< TypeTag, Properties::CellCenterPrimaryVariables > CellCenterResidualValue
Definition: staggeredlocalresidual.hh:66
Implementation & asImp()
Definition: staggeredlocalresidual.hh:324
void evalStorageForCellCenter(CellCenterResidualValue &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &prevElemVolVars, const ElementVolumeVariables &curElemVolVars, const SubControlVolume &scv) const
Evaluate the storage terms for a cell center residual.
Definition: staggeredlocalresidual.hh:148
Implementation & asImp_()
Definition: staggeredlocalresidual.hh:311
const TimeLoop & timeLoop() const
Definition: staggeredlocalresidual.hh:321
FaceResidualValue evalFluxAndSourceForFace(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const ElementBoundaryTypes &bcTypes, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Convenience function to evaluate the flux and source terms for the face residual.
Definition: staggeredlocalresidual.hh:202
GetPropType< TypeTag, Properties::FacePrimaryVariables > FaceResidualValue
Definition: staggeredlocalresidual.hh:67
const Problem & problem() const
the problem
Definition: staggeredlocalresidual.hh:306
const Implementation & asImp_() const
Definition: staggeredlocalresidual.hh:314
CellCenterResidualValue evalStorage(Args &&... args) const
for compatibility with FVLocalAssemblerBase
Definition: staggeredlocalresidual.hh:189
void evalSourceForFace(FaceResidualValue &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const SubControlVolumeFace &scvf) const
Evaluate the source terms for a face residual.
Definition: staggeredlocalresidual.hh:235
bool isStationary() const
If no solution has been set, we treat the problem as stationary.
Definition: staggeredlocalresidual.hh:302
TimeLoop & timeLoop()
Definition: staggeredlocalresidual.hh:318
StaggeredLocalResidual(const Problem *problem, const TimeLoop *timeLoop=nullptr)
the constructor
Definition: staggeredlocalresidual.hh:71
CellCenterResidualValue evalFluxAndSource(Args &&... args) const
for compatibility with FVLocalAssemblerBase
Definition: staggeredlocalresidual.hh:182
void evalStorageForFace(FaceResidualValue &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &prevElemVolVars, const ElementVolumeVariables &curElemVolVars, const ElementFaceVariables &prevElemFaceVars, const ElementFaceVariables &curElemFaceVars, const SubControlVolumeFace &scvf) const
Evaluate the storage terms for a face residual.
Definition: staggeredlocalresidual.hh:273
CellCenterResidualValue evalFluxAndSourceForCellCenter(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const ElementBoundaryTypes &bcTypes, const ElementFluxVariablesCache &elemFluxVarsCache) const
Convenience function to evaluate the flux and source terms for the cell center residual.
Definition: staggeredlocalresidual.hh:78
void evalFluxForFace(FaceResidualValue &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Evaluate the flux terms for a face residual.
Definition: staggeredlocalresidual.hh:218
void evalSourceForCellCenter(CellCenterResidualValue &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &curElemVolVars, const ElementFaceVariables &curElemFaceVars, const SubControlVolume &scv) const
Evaluate the source terms for a cell center residual.
Definition: staggeredlocalresidual.hh:116
CellCenterResidualValue ElementResidualVector
Definition: staggeredlocalresidual.hh:68
const Implementation & asImp() const
Definition: staggeredlocalresidual.hh:327
void evalFluxForCellCenter(CellCenterResidualValue &residual, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFaceVariables &elemFaceVars, const ElementBoundaryTypes &elemBcTypes, const ElementFluxVariablesCache &elemFluxVarsCache, const SubControlVolumeFace &scvf) const
Evaluate the flux terms for a cell center residual.
Definition: staggeredlocalresidual.hh:99
FaceResidualValue evalStorageForFace(const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &prevElemVolVars, const ElementVolumeVariables &curElemVolVars, const ElementFaceVariables &prevElemFaceVars, const ElementFaceVariables &curElemFaceVars, const SubControlVolumeFace &scvf) const
Evaluate the storage terms for a face residual.
Definition: staggeredlocalresidual.hh:258
Manages the handling of time dependent problems.
Definition: common/timeloop.hh:69
virtual Scalar timeStepSize() const =0
Returns the suggested time step length .
The default time loop for instationary simulations.
Definition: common/timeloop.hh:114
Declares all properties used in Dumux.
Manages the handling of time dependent problems.