34 template<
class RdeltaTType,
class RhoType,
class SpType,
class SuType>
37 const RdeltaTType& rDeltaT,
54 psi.primitiveFieldRef() =
56 rho.primitiveField()*
psi.primitiveField()*rDeltaT
59 )/(
rho.primitiveField()*rDeltaT -
Sp.primitiveField());
63 psi.primitiveFieldRef() =
65 rho.primitiveField()*
psi.primitiveField()*rDeltaT
68 )/(
rho.primitiveField()*rDeltaT -
Sp.primitiveField());
71 psi.correctBoundaryConditions();
75 template<
class RhoType>
87 template<
class RhoType,
class SpType,
class SuType>
99 if (fv::localEulerDdt::enabled(mesh))
101 const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
112 template<
class RhoType,
class PsiMaxType,
class PsiMinType>
119 const PsiMaxType& psiMax,
120 const PsiMinType& psiMin
143 const PsiMaxType& psiMax,
144 const PsiMinType& psiMin
149 if (fv::localEulerDdt::enabled(mesh))
151 const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
202 const RdeltaTType& rDeltaT,
209 const PsiMaxType& psiMax,
210 const PsiMinType& psiMin
220 const label nLimiterIter
225 const scalar smoothLimiter
230 const scalar extremaCoeff
235 const scalar boundaryExtremaCoeff
239 "boundaryExtremaCoeff",
244 const scalar boundaryDeltaExtremaCoeff
246 max(boundaryExtremaCoeff - extremaCoeff, 0)
263 lambda.boundaryFieldRef();
276 const label own = owner[facei];
277 const label nei = neighb[facei];
279 psiMaxn[own] =
max(psiMaxn[own], psiIf[nei]);
280 psiMinn[own] =
min(psiMinn[own], psiIf[nei]);
282 psiMaxn[nei] =
max(psiMaxn[nei], psiIf[own]);
283 psiMinn[nei] =
min(psiMinn[nei], psiIf[own]);
285 const scalar phiCorrf = phiCorrIf[facei];
289 sumPhip[own] += phiCorrf;
290 mSumPhim[nei] += phiCorrf;
294 mSumPhim[own] -= phiCorrf;
295 sumPhip[nei] -= phiCorrf;
312 label pfCelli = pFaceCells[pFacei];
314 psiMaxn[pfCelli] =
max(psiMaxn[pfCelli], psiPNf[pFacei]);
315 psiMinn[pfCelli] =
min(psiMinn[pfCelli], psiPNf[pFacei]);
322 const label pfCelli = pFaceCells[pFacei];
324 psiMaxn[pfCelli] =
max(psiMaxn[pfCelli], psiPf[pFacei]);
325 psiMinn[pfCelli] =
min(psiMinn[pfCelli], psiPf[pFacei]);
331 if (boundaryDeltaExtremaCoeff > 0)
335 const label pfCelli = pFaceCells[pFacei];
337 const scalar extrema =
338 boundaryDeltaExtremaCoeff
339 *(psiMax[pfCelli] - psiMin[pfCelli]);
341 psiMaxn[pfCelli] += extrema;
342 psiMinn[pfCelli] -= extrema;
349 const label pfCelli = pFaceCells[pFacei];
351 const scalar phiCorrf = phiCorrPf[pFacei];
355 sumPhip[pfCelli] += phiCorrf;
359 mSumPhim[pfCelli] -= phiCorrf;
364 if (extremaCoeff > 0)
366 psiMaxn =
min(psiMaxn + extremaCoeff*(psiMax - psiMin), psiMax);
367 psiMinn =
max(psiMinn - extremaCoeff*(psiMax - psiMin), psiMin);
371 psiMaxn =
min(psiMaxn, psiMax);
372 psiMinn =
max(psiMinn, psiMin);
375 if (smoothLimiter > small)
378 min(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMaxn, psiMax);
380 max(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMinn, psiMin);
386 (
rho.primitiveField()*rDeltaT -
Sp.primitiveField())*psiMaxn
387 -
Su.primitiveField()
388 -
rho.primitiveField()*
psi.primitiveField()*rDeltaT
395 - (
rho.primitiveField()*rDeltaT -
Sp.primitiveField())*psiMinn
396 +
rho.primitiveField()*
psi.primitiveField()*rDeltaT
402 for (
int j=0; j<nLimiterIter; j++)
409 const label own = owner[facei];
410 const label nei = neighb[facei];
412 const scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
414 if (lambdaPhiCorrf > 0)
416 sumlPhip[own] += lambdaPhiCorrf;
417 mSumlPhim[nei] += lambdaPhiCorrf;
421 mSumlPhim[own] -= lambdaPhiCorrf;
422 sumlPhip[nei] -= lambdaPhiCorrf;
435 label pfCelli = pFaceCells[pFacei];
437 scalar lambdaPhiCorrf = lambdaPf[pFacei]*phiCorrfPf[pFacei];
439 if (lambdaPhiCorrf > 0)
441 sumlPhip[pfCelli] += lambdaPhiCorrf;
445 mSumlPhim[pfCelli] -= lambdaPhiCorrf;
455 (sumlPhip[celli] + psiMaxn[celli])
456 /(mSumPhim[celli] + rootVSmall),
463 (mSumlPhim[celli] + psiMinn[celli])
464 /(sumPhip[celli] + rootVSmall),
474 if (phiCorrIf[facei] > 0)
476 lambdaIf[facei] =
min
479 min(lambdap[owner[facei]], lambdam[neighb[facei]])
484 lambdaIf[facei] =
min
487 min(lambdam[owner[facei]], lambdap[neighb[facei]])
510 const label pfCelli = pFaceCells[pFacei];
512 if (phiCorrfPf[pFacei] > 0)
515 min(lambdaPf[pFacei], lambdap[pfCelli]);
520 min(lambdaPf[pFacei], lambdam[pfCelli]);
533 if ((phiPf[pFacei] + phiCorrfPf[pFacei]) > small*small)
535 const label pfCelli = pFaceCells[pFacei];
537 if (phiCorrfPf[pFacei] > 0)
540 min(lambdaPf[pFacei], lambdap[pfCelli]);
545 min(lambdaPf[pFacei], lambdam[pfCelli]);
555 surfaceScalarField::Internal::null(),
564 lambdaPf =
min(lambdaPf, lambdaNbrPf);
582 const RdeltaTType& rDeltaT,
589 const PsiMaxType& psiMax,
590 const PsiMinType& psiMin
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
#define forAll(list, i)
Loop across all elements in list.
Generic GeometricBoundaryField class.
tmp< GeometricBoundaryField > boundaryNeighbourField() const
Return BoundaryField of the values on the other side of couples.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
void size(const label)
Override size to be inconsistent with allocated storage.
scalar deltaTValue() const
Return time step value.
A list of keyword definitions, which are a keyword followed by any number of values (e....
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T, if not found return the given default.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
const word & name() const
Return const reference to name.
Mesh data needed to do the Finite Volume discretisation.
const Time & time() const
Return the top-level database.
const labelUList & owner() const
Internal face owner.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
const fvSolution & solution() const
Return the fvSolution.
const labelUList & neighbour() const
Internal face neighbour.
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual bool fixesValue() const
Return true if this patch field fixes a value.
virtual bool coupled() const
Return true if this patch field is coupled.
virtual tmp< Field< Type > > patchNeighbourField(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking) const
Return patchField on the opposite patch of a coupled patch.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
virtual bool coupled() const
Return true if this patch field is coupled.
bool moving() const
Is mesh moving.
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
A class for managing temporary objects.
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
const volScalarField & psi
dimensionedScalar lambda(viscosity->lookup("lambda"))
void limiterCorr(surfaceScalarField &allLambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin)
void limitCorr(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin)
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
void surfaceIntegrate(Field< Type > &ivf, const SurfaceField< Type > &ssf)
tmp< VolField< Type > > Su(const VolField< Type > &su, const VolField< Type > &vf)
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const dimensionSet dimless
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)