36 template<
class RdeltaTType,
class RhoType,
class SpType,
class SuType>
39 const RdeltaTType& rDeltaT,
61 mesh.
Vsc0()().
field()*rho.oldTime().field()
65 )/(rho.field()*rDeltaT - Sp.field());
71 rho.oldTime().field()*psi0*rDeltaT
74 )/(rho.field()*rDeltaT - Sp.field());
81 template<
class RhoType>
93 template<
class RhoType,
class SpType,
class SuType>
105 if (fv::localEulerDdt::enabled(mesh))
107 const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
118 template<
class RhoType,
class PsiMaxType,
class PsiMinType>
125 const PsiMaxType& psiMax,
126 const PsiMinType& psiMin
159 const PsiMaxType& psiMax,
160 const PsiMinType& psiMin
167 if (fv::localEulerDdt::enabled(mesh))
169 const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
170 limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin,
false);
176 limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin,
false);
194 const RdeltaTType& rDeltaT,
201 const PsiMaxType& psiMax,
202 const PsiMinType& psiMin
212 const label nLimiterIter
217 const scalar smoothLimiter
222 const scalar extremaCoeff
227 const scalar boundaryExtremaCoeff
231 "boundaryExtremaCoeff",
236 const scalar boundaryDeltaExtremaCoeff
238 max(boundaryExtremaCoeff - extremaCoeff, 0)
249 const surfaceScalarField::Boundary& phiBDBf =
253 const surfaceScalarField::Boundary& phiCorrBf =
289 const label own = owner[facei];
290 const label nei = neighb[facei];
292 psiMaxn[own] =
max(psiMaxn[own], psiIf[nei]);
293 psiMinn[own] =
min(psiMinn[own], psiIf[nei]);
295 psiMaxn[nei] =
max(psiMaxn[nei], psiIf[own]);
296 psiMinn[nei] =
min(psiMinn[nei], psiIf[own]);
298 sumPhiBD[own] += phiBDIf[facei];
299 sumPhiBD[nei] -= phiBDIf[facei];
301 const scalar phiCorrf = phiCorrIf[facei];
305 sumPhip[own] += phiCorrf;
306 mSumPhim[nei] += phiCorrf;
310 mSumPhim[own] -= phiCorrf;
311 sumPhip[nei] -= phiCorrf;
329 const label pfCelli = pFaceCells[pFacei];
331 psiMaxn[pfCelli] =
max(psiMaxn[pfCelli], psiPNf[pFacei]);
332 psiMinn[pfCelli] =
min(psiMinn[pfCelli], psiPNf[pFacei]);
339 const label pfCelli = pFaceCells[pFacei];
341 psiMaxn[pfCelli] =
max(psiMaxn[pfCelli], psiPf[pFacei]);
342 psiMinn[pfCelli] =
min(psiMinn[pfCelli], psiPf[pFacei]);
348 if (boundaryDeltaExtremaCoeff > 0)
352 const label pfCelli = pFaceCells[pFacei];
354 const scalar extrema =
355 boundaryDeltaExtremaCoeff
356 *(psiMax[pfCelli] - psiMin[pfCelli]);
358 psiMaxn[pfCelli] += extrema;
359 psiMinn[pfCelli] -= extrema;
366 const label pfCelli = pFaceCells[pFacei];
368 sumPhiBD[pfCelli] += phiBDPf[pFacei];
370 const scalar phiCorrf = phiCorrPf[pFacei];
374 sumPhip[pfCelli] += phiCorrf;
378 mSumPhim[pfCelli] -= phiCorrf;
383 psiMaxn =
min(psiMaxn + extremaCoeff*(psiMax - psiMin), psiMax);
384 psiMinn =
max(psiMinn - extremaCoeff*(psiMax - psiMin), psiMin);
386 if (smoothLimiter > small)
389 min(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMaxn, psiMax);
391 max(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMinn, psiMin);
401 (rho.field()*rDeltaT - Sp.field())*psiMaxn
404 - (V0().field()*rDeltaT)*rho.oldTime().field()*psi0
411 - (rho.field()*rDeltaT - Sp.field())*psiMinn
413 + (V0().field()*rDeltaT)*rho.oldTime().field()*psi0
421 (rho.field()*rDeltaT - Sp.field())*psiMaxn
423 - (rho.oldTime().field()*rDeltaT)*psi0
431 - (rho.field()*rDeltaT - Sp.field())*psiMinn
432 + (rho.oldTime().field()*rDeltaT)*psi0
440 for (
int j=0; j<nLimiterIter; j++)
447 const label own = owner[facei];
448 const label nei = neighb[facei];
450 scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
452 if (lambdaPhiCorrf > 0)
454 sumlPhip[own] += lambdaPhiCorrf;
455 mSumlPhim[nei] += lambdaPhiCorrf;
459 mSumlPhim[own] -= lambdaPhiCorrf;
460 sumlPhip[nei] -= lambdaPhiCorrf;
473 const label pfCelli = pFaceCells[pFacei];
474 const scalar lambdaPhiCorrf =
475 lambdaPf[pFacei]*phiCorrfPf[pFacei];
477 if (lambdaPhiCorrf > 0)
479 sumlPhip[pfCelli] += lambdaPhiCorrf;
483 mSumlPhim[pfCelli] -= lambdaPhiCorrf;
493 (sumlPhip[celli] + psiMaxn[celli])
494 /(mSumPhim[celli] + rootVSmall),
501 (mSumlPhim[celli] + psiMinn[celli])
502 /(sumPhip[celli] + rootVSmall),
512 if (phiCorrIf[facei] > 0)
514 lambdaIf[facei] =
min 517 min(lambdap[owner[facei]], lambdam[neighb[facei]])
522 lambdaIf[facei] =
min 525 min(lambdam[owner[facei]], lambdap[neighb[facei]])
547 const label pfCelli = pFaceCells[pFacei];
549 if (phiCorrfPf[pFacei] > 0)
552 min(lambdaPf[pFacei], lambdap[pfCelli]);
557 min(lambdaPf[pFacei], lambdam[pfCelli]);
579 const RdeltaTType& rDeltaT,
586 const PsiMaxType& psiMax,
587 const PsiMinType& psiMin,
588 const bool returnCorr
596 const surfaceScalarField::Boundary& phiPsiBf = phiPsi.
boundaryField();
604 phiBDPf = phiPsiBf[
patchi];
650 phiPsi = phiBD + lambda*phiCorr;
671 const PsiMaxType& psiMax,
672 const PsiMinType& psiMin,
678 if (fv::localEulerDdt::enabled(mesh))
680 const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
681 limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin, rtnCorr);
686 limit(rDeltaT, rho, psi, phi, phiPsi, Sp, Su, psiMax, psiMin, rtnCorr);
691 template<
template<
class>
class AlphaList,
template<
class>
class PhiList>
694 const AlphaList<const volScalarField>& alphas,
695 PhiList<surfaceScalarField>& phiPsis,
723 const surfaceScalarField::Boundary& phibf = phi.
boundaryField();
727 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.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & name() const
Return name.
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
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...
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
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 ...
void limiter(scalarField &allLambda, 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)
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
scalar deltaTValue() const
Return time step value.
virtual bool coupled() const
Return true if this patch field is coupled.
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-cycl old-time cell volumes.
const Mesh & mesh() const
Return mesh.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
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.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Upwind differencing scheme class.
Mesh data needed to do the Finite Volume discretisation.
void correctBoundaryConditions()
Correct boundary field.
MULES: Multidimensional universal limiter for explicit solution.
const volScalarField & psi
A class for managing temporary objects.
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Specialisation of GeometricField which holds slices of given complete fields in a form that they act ...
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.