25#ifndef DUMUX_STAGGERED_LOCAL_RESIDUAL_HH
26#define DUMUX_STAGGERED_LOCAL_RESIDUAL_HH
38template<
class TypeTag>
46 using Element =
typename GridView::template Codim<0>::Entity;
51 using SubControlVolume =
typename FVElementGeometry::SubControlVolume;
52 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
75 const FVElementGeometry& fvGeometry,
76 const ElementVolumeVariables& elemVolVars,
77 const ElementFaceVariables& elemFaceVars,
78 const ElementBoundaryTypes& bcTypes,
79 const ElementFluxVariablesCache& elemFluxVarsCache)
const
84 for (
auto&& scv : scvs(fvGeometry))
85 asImp().evalSourceForCellCenter(residual, this->
problem(), element, fvGeometry, elemVolVars, elemFaceVars, scv);
88 for (
auto&& scvf : scvfs(fvGeometry))
89 asImp().evalFluxForCellCenter(residual, this->
problem(), element, fvGeometry, elemVolVars, elemFaceVars, bcTypes, elemFluxVarsCache, scvf);
97 const Element& element,
98 const FVElementGeometry& fvGeometry,
99 const ElementVolumeVariables& elemVolVars,
100 const ElementFaceVariables& elemFaceVars,
101 const ElementBoundaryTypes& elemBcTypes,
102 const ElementFluxVariablesCache& elemFluxVarsCache,
103 const SubControlVolumeFace& scvf)
const
105 if (!scvf.boundary())
106 residual +=
asImp_().computeFluxForCellCenter(
problem, element, fvGeometry, elemVolVars, elemFaceVars, scvf, elemFluxVarsCache);
108 residual +=
asImp_().computeBoundaryFluxForCellCenter(
problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars, elemBcTypes, elemFluxVarsCache);
114 const Element& element,
115 const FVElementGeometry& fvGeometry,
116 const ElementVolumeVariables& curElemVolVars,
117 const ElementFaceVariables& curElemFaceVars,
118 const SubControlVolume& scv)
const
120 const auto curExtrusionFactor = curElemVolVars[scv].extrusionFactor();
123 auto source =
asImp_().computeSourceForCellCenter(
problem, element, fvGeometry, curElemVolVars, curElemFaceVars, scv);
124 source *= scv.volume()*curExtrusionFactor;
130 const FVElementGeometry& fvGeometry,
131 const ElementVolumeVariables& prevElemVolVars,
132 const ElementVolumeVariables& curElemVolVars)
const
134 assert(timeLoop_ &&
"no time loop set for storage term evaluation");
137 for (
auto&& scv : scvs(fvGeometry))
138 asImp().evalStorageForCellCenter(storage,
problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, scv);
146 const Element &element,
147 const FVElementGeometry& fvGeometry,
148 const ElementVolumeVariables& prevElemVolVars,
149 const ElementVolumeVariables& curElemVolVars,
150 const SubControlVolume& scv)
const
153 const auto& curVolVars = curElemVolVars[scv];
154 const auto& prevVolVars = prevElemVolVars[scv];
162 auto prevCCStorage =
asImp_().computeStorageForCellCenter(
problem, scv, prevVolVars);
163 auto curCCStorage =
asImp_().computeStorageForCellCenter(
problem, scv, curVolVars);
165 prevCCStorage *= prevVolVars.extrusionFactor();
166 curCCStorage *= curVolVars.extrusionFactor();
168 storage = std::move(curCCStorage);
169 storage -= std::move(prevCCStorage);
170 storage *= scv.volume();
177 template<
class... Args>
184 template<
class... Args>
199 const FVElementGeometry& fvGeometry,
200 const ElementVolumeVariables& elemVolVars,
201 const ElementFaceVariables& elemFaceVars,
202 const ElementBoundaryTypes& bcTypes,
203 const ElementFluxVariablesCache& elemFluxVarsCache,
204 const SubControlVolumeFace& scvf)
const
207 asImp().evalSourceForFace(residual, this->
problem(), element, fvGeometry, elemVolVars, elemFaceVars, scvf);
208 asImp().evalFluxForFace(residual, this->
problem(), element, fvGeometry, elemVolVars, elemFaceVars, bcTypes, elemFluxVarsCache, scvf);
216 const Element& element,
217 const FVElementGeometry& fvGeometry,
218 const ElementVolumeVariables& elemVolVars,
219 const ElementFaceVariables& elemFaceVars,
220 const ElementBoundaryTypes& elemBcTypes,
221 const ElementFluxVariablesCache& elemFluxVarsCache,
222 const SubControlVolumeFace& scvf)
const
224 if (!scvf.boundary())
225 residual +=
asImp_().computeFluxForFace(
problem, element, scvf, fvGeometry, elemVolVars, elemFaceVars, elemFluxVarsCache);
227 residual +=
asImp_().computeBoundaryFluxForFace(
problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars, elemBcTypes, elemFluxVarsCache);
233 const Element& element,
234 const FVElementGeometry& fvGeometry,
235 const ElementVolumeVariables& elemVolVars,
236 const ElementFaceVariables& elemFaceVars,
237 const SubControlVolumeFace& scvf)
const
240 auto source =
asImp_().computeSourceForFace(
problem, element, fvGeometry, scvf, elemVolVars, elemFaceVars);
241 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
242 const auto extrusionFactor = elemVolVars[scv].extrusionFactor();
245 source *= 0.5*scv.volume()*extrusionFactor;
251 const FVElementGeometry& fvGeometry,
252 const ElementVolumeVariables& prevElemVolVars,
253 const ElementVolumeVariables& curElemVolVars,
254 const ElementFaceVariables& prevElemFaceVars,
255 const ElementFaceVariables& curElemFaceVars,
256 const SubControlVolumeFace& scvf)
const
258 assert(timeLoop_ &&
"no time loop set for storage term evaluation");
260 asImp().evalStorageForFace(storage,
problem(), element, fvGeometry, prevElemVolVars, curElemVolVars, prevElemFaceVars, curElemFaceVars, scvf);
267 const Element& element,
268 const FVElementGeometry& fvGeometry,
269 const ElementVolumeVariables& prevElemVolVars,
270 const ElementVolumeVariables& curElemVolVars,
271 const ElementFaceVariables& prevElemFaceVars,
272 const ElementFaceVariables& curElemFaceVars,
273 const SubControlVolumeFace& scvf)
const
276 const auto& scv = fvGeometry.scv(scvf.insideScvIdx());
277 auto prevFaceStorage =
asImp_().computeStorageForFace(
problem, scvf, prevElemVolVars[scv], prevElemFaceVars);
278 auto curFaceStorage =
asImp_().computeStorageForFace(
problem, scvf, curElemVolVars[scv], curElemFaceVars);
280 storage = std::move(curFaceStorage);
281 storage -= std::move(prevFaceStorage);
283 const auto extrusionFactor = curElemVolVars[scv].extrusionFactor();
286 storage *= 0.5*scv.volume()*extrusionFactor;
294 {
return !timeLoop_; }
298 {
return *problem_; }
303 {
return *
static_cast<Implementation*
>(
this); }
306 {
return *
static_cast<const Implementation*
>(
this); }
310 {
return *timeLoop_; }
313 {
return *timeLoop_; }
316 {
return *
static_cast<Implementation*
>(
this); }
319 {
return *
static_cast<const Implementation*
>(
this); }
322 const Problem* problem_;
Manages the handling of time dependent problems.
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:149
Calculates the element-wise residual for the staggered FV scheme.
Definition: staggeredlocalresidual.hh:40
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:129
GetPropType< TypeTag, Properties::CellCenterPrimaryVariables > CellCenterResidualValue
Definition: staggeredlocalresidual.hh:62
Implementation & asImp()
Definition: staggeredlocalresidual.hh:315
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:144
Implementation & asImp_()
Definition: staggeredlocalresidual.hh:302
const TimeLoop & timeLoop() const
Definition: staggeredlocalresidual.hh:312
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:198
GetPropType< TypeTag, Properties::FacePrimaryVariables > FaceResidualValue
Definition: staggeredlocalresidual.hh:63
const Problem & problem() const
the problem
Definition: staggeredlocalresidual.hh:297
const Implementation & asImp_() const
Definition: staggeredlocalresidual.hh:305
CellCenterResidualValue evalStorage(Args &&... args) const
for compatibility with FVLocalAssemblerBase
Definition: staggeredlocalresidual.hh:185
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:231
bool isStationary() const
If no solution has been set, we treat the problem as stationary.
Definition: staggeredlocalresidual.hh:293
TimeLoop & timeLoop()
Definition: staggeredlocalresidual.hh:309
StaggeredLocalResidual(const Problem *problem, const TimeLoop *timeLoop=nullptr)
the constructor
Definition: staggeredlocalresidual.hh:67
CellCenterResidualValue evalFluxAndSource(Args &&... args) const
for compatibility with FVLocalAssemblerBase
Definition: staggeredlocalresidual.hh:178
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:265
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:74
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:214
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:112
CellCenterResidualValue ElementResidualVector
Definition: staggeredlocalresidual.hh:64
const Implementation & asImp() const
Definition: staggeredlocalresidual.hh:318
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:95
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:250
Manages the handling of time dependent problems.
Definition: timeloop.hh:72
virtual Scalar timeStepSize() const =0
Returns the suggested time step length .
The default time loop for instationary simulations.
Definition: timeloop.hh:117
Declares all properties used in Dumux.