35 template<
class RdeltaTType,
class RhoType,
class SpType,
class SuType>
38 const RdeltaTType& rDeltaT,
60 mesh.
Vsc0()().
field()*rho.oldTime().field()
64 )/(rho.field()*rDeltaT - Sp.field());
70 rho.oldTime().field()*psi0*rDeltaT
73 )/(rho.field()*rDeltaT - Sp.field());
80 template<
class RhoType>
92 template<
class RhoType,
class SpType,
class SuType>
104 if (fv::localEulerDdt::enabled(mesh))
106 const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
117 template<
class RhoType,
class PsiMaxType,
class PsiMinType>
124 const PsiMaxType& psiMax,
125 const PsiMinType& psiMin
158 const PsiMaxType& psiMax,
159 const PsiMinType& psiMin
166 if (fv::localEulerDdt::enabled(mesh))
168 const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
169 limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin,
false);
175 limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin,
false);
193 const RdeltaTType& rDeltaT,
200 const PsiMaxType& psiMax,
201 const PsiMinType& psiMin
211 const label nLimiterIter
216 const scalar smoothLimiter
221 const scalar extremaCoeff
226 const scalar boundaryExtremaCoeff
230 "boundaryExtremaCoeff",
235 const scalar boundaryDeltaExtremaCoeff
237 max(boundaryExtremaCoeff - extremaCoeff, 0)
271 const label own = owner[facei];
272 const label nei = neighb[facei];
274 psiMaxn[own] =
max(psiMaxn[own], psiIf[nei]);
275 psiMinn[own] =
min(psiMinn[own], psiIf[nei]);
277 psiMaxn[nei] =
max(psiMaxn[nei], psiIf[own]);
278 psiMinn[nei] =
min(psiMinn[nei], psiIf[own]);
280 sumPhiBD[own] += phiBDIf[facei];
281 sumPhiBD[nei] -= phiBDIf[facei];
283 const scalar phiCorrf = phiCorrIf[facei];
287 sumPhip[own] += phiCorrf;
288 mSumPhim[nei] += phiCorrf;
292 mSumPhim[own] -= phiCorrf;
293 sumPhip[nei] -= phiCorrf;
311 const label pfCelli = pFaceCells[pFacei];
313 psiMaxn[pfCelli] =
max(psiMaxn[pfCelli], psiPNf[pFacei]);
314 psiMinn[pfCelli] =
min(psiMinn[pfCelli], psiPNf[pFacei]);
321 const label pfCelli = pFaceCells[pFacei];
323 psiMaxn[pfCelli] =
max(psiMaxn[pfCelli], psiPf[pFacei]);
324 psiMinn[pfCelli] =
min(psiMinn[pfCelli], psiPf[pFacei]);
330 if (boundaryDeltaExtremaCoeff > 0)
334 const label pfCelli = pFaceCells[pFacei];
336 const scalar extrema =
337 boundaryDeltaExtremaCoeff
338 *(psiMax[pfCelli] - psiMin[pfCelli]);
340 psiMaxn[pfCelli] += extrema;
341 psiMinn[pfCelli] -= extrema;
348 const label pfCelli = pFaceCells[pFacei];
350 sumPhiBD[pfCelli] += phiBDPf[pFacei];
352 const scalar phiCorrf = phiCorrPf[pFacei];
356 sumPhip[pfCelli] += phiCorrf;
360 mSumPhim[pfCelli] -= phiCorrf;
365 psiMaxn =
min(psiMaxn + extremaCoeff*(psiMax - psiMin), psiMax);
366 psiMinn =
max(psiMinn - extremaCoeff*(psiMax - psiMin), psiMin);
368 if (smoothLimiter > small)
371 min(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMaxn, psiMax);
373 max(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMinn, psiMin);
383 (rho.field()*rDeltaT - Sp.field())*psiMaxn
386 - (V0().field()*rDeltaT)*rho.oldTime().field()*psi0
393 - (rho.field()*rDeltaT - Sp.field())*psiMinn
395 + (V0().field()*rDeltaT)*rho.oldTime().field()*psi0
403 (rho.field()*rDeltaT - Sp.field())*psiMaxn
405 - (rho.oldTime().field()*rDeltaT)*psi0
413 - (rho.field()*rDeltaT - Sp.field())*psiMinn
414 + (rho.oldTime().field()*rDeltaT)*psi0
422 for (
int j=0; j<nLimiterIter; j++)
429 const label own = owner[facei];
430 const label nei = neighb[facei];
432 scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
434 if (lambdaPhiCorrf > 0)
436 sumlPhip[own] += lambdaPhiCorrf;
437 mSumlPhim[nei] += lambdaPhiCorrf;
441 mSumlPhim[own] -= lambdaPhiCorrf;
442 sumlPhip[nei] -= lambdaPhiCorrf;
455 const label pfCelli = pFaceCells[pFacei];
456 const scalar lambdaPhiCorrf =
457 lambdaPf[pFacei]*phiCorrfPf[pFacei];
459 if (lambdaPhiCorrf > 0)
461 sumlPhip[pfCelli] += lambdaPhiCorrf;
465 mSumlPhim[pfCelli] -= lambdaPhiCorrf;
475 (sumlPhip[celli] + psiMaxn[celli])
476 /(mSumPhim[celli] + rootVSmall),
483 (mSumlPhim[celli] + psiMinn[celli])
484 /(sumPhip[celli] + rootVSmall),
494 if (phiCorrIf[facei] > 0)
496 lambdaIf[facei] =
min 499 min(lambdap[owner[facei]], lambdam[neighb[facei]])
504 lambdaIf[facei] =
min 507 min(lambdam[owner[facei]], lambdap[neighb[facei]])
529 const label pfCelli = pFaceCells[pFacei];
531 if (phiCorrfPf[pFacei] > 0)
534 min(lambdaPf[pFacei], lambdap[pfCelli]);
539 min(lambdaPf[pFacei], lambdam[pfCelli]);
548 surfaceScalarField::Internal::null(),
557 lambdaPf =
min(lambdaPf, lambdaNbrPf);
575 const RdeltaTType& rDeltaT,
582 const PsiMaxType& psiMax,
583 const PsiMinType& psiMin,
584 const bool returnCorr
600 phiBDPf = phiPsiBf[
patchi];
642 phiPsi = phiBD + lambda*phiCorr;
663 const PsiMaxType& psiMax,
664 const PsiMinType& psiMin,
670 if (fv::localEulerDdt::enabled(mesh))
672 const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
673 limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin, rtnCorr);
678 limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin, rtnCorr);
683 template<
template<
class>
class AlphaList,
template<
class>
class PhiList>
686 const AlphaList<const volScalarField>& alphas,
687 PhiList<surfaceScalarField>& phiPsis,
719 if (phibf[
patchi].coupled())
void limit(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin, const bool returnCorr)
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
#define forAll(list, i)
Loop across all elements in list.
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const word & name() const
Return name.
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
bool moving() const
Is mesh moving.
virtual bool fixesValue() const
Return true if this patch field fixes a value.
A list of keyword definitions, which are a keyword followed by any number of values (e...
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void size(const label)
Override size to be inconsistent with allocated storage.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar lambda(viscosity->lookup("lambda"))
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
const fvSolution & solution() const
Return the fvSchemes.
const dimensionSet dimless
const Time & time() const
Return the top-level database.
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
const labelUList & neighbour() const
Internal face neighbour.
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
virtual bool coupled() const
Return true if this patch field is coupled.
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
scalar deltaTValue() const
Return time step value.
const volScalarField & psi
virtual bool coupled() const
Return true if this patch field is coupled.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void limiter(surfaceScalarField &lambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phiBD, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const PsiMaxType &psiMax, const PsiMinType &psiMin)
tmp< GeometricBoundaryField > boundaryNeighbourField() const
Return BoundaryField of the values on the other side of couples.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
tmp< DimensionedField< scalar, volMesh > > Vsc0() const
Return sub-cycle old-time cell volumes.
const Mesh & mesh() const
Return mesh.
const labelUList & owner() const
Internal face owner.
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Mesh data needed to do the Finite Volume discretisation.
void correctBoundaryConditions()
Correct boundary field.
MULES: Multidimensional universal limiter for explicit solution.
A class for managing temporary objects.
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
virtual tmp< Field< Type > > patchNeighbourField(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking) const
Return patchField on the opposite patch of a coupled patch.
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.