39template<
class RdeltaTType,
class RhoType,
class SpType,
class SuType>
42 const RdeltaTType& rDeltaT,
65 *psi0*rDeltaT/
mesh.Vsc()().field()
68 )/(
rho.field()*rDeltaT -
Sp.field());
77 )/(
rho.field()*rDeltaT -
Sp.field());
80 psi.correctBoundaryConditions();
84template<
class RhoType>
89 const surfaceScalarField& phiPsi
96template<
class RhoType,
class SpType,
class SuType>
111 explicitSolve(rDeltaT,
rho,
psi, phiPsi,
Sp,
Su);
115 const scalar rDeltaT = 1.0/
mesh.time().deltaTValue();
121template<
class RhoType,
class PsiMaxType,
class PsiMinType>
126 const surfaceScalarField& phiBD,
127 surfaceScalarField& phiPsi,
128 const PsiMaxType& psiMax,
129 const PsiMinType& psiMin
162 const PsiMaxType& psiMax,
163 const PsiMinType& psiMin
168 psi.correctBoundaryConditions();
178 const scalar rDeltaT = 1.0/
mesh.time().deltaTValue();
197 const RdeltaTType& rDeltaT,
204 const PsiMaxType& psiMax,
205 const PsiMinType& psiMin
215 const label nLimiterIter
220 const scalar smoothLimiter
225 const scalar extremaCoeff
230 const scalar boundaryExtremaCoeff
234 "boundaryExtremaCoeff",
239 const scalar boundaryDeltaExtremaCoeff
241 max(boundaryExtremaCoeff - extremaCoeff, 0)
264 mesh.time().timeName(),
292 const label own = owner[facei];
293 const label nei = neighb[facei];
295 psiMaxn[own] =
max(psiMaxn[own], psiIf[nei]);
296 psiMinn[own] =
min(psiMinn[own], psiIf[nei]);
298 psiMaxn[nei] =
max(psiMaxn[nei], psiIf[own]);
299 psiMinn[nei] =
min(psiMinn[nei], psiIf[own]);
301 sumPhiBD[own] += phiBDIf[facei];
302 sumPhiBD[nei] -= phiBDIf[facei];
304 const scalar phiCorrf = phiCorrIf[facei];
308 sumPhip[own] += phiCorrf;
309 mSumPhim[nei] += phiCorrf;
313 mSumPhim[own] -= phiCorrf;
314 sumPhip[nei] -= phiCorrf;
324 const labelUList& pFaceCells =
mesh.boundary()[patchi].faceCells();
332 const label pfCelli = pFaceCells[pFacei];
334 psiMaxn[pfCelli] =
max(psiMaxn[pfCelli], psiPNf[pFacei]);
335 psiMinn[pfCelli] =
min(psiMinn[pfCelli], psiPNf[pFacei]);
342 const label pfCelli = pFaceCells[pFacei];
344 psiMaxn[pfCelli] =
max(psiMaxn[pfCelli], psiPf[pFacei]);
345 psiMinn[pfCelli] =
min(psiMinn[pfCelli], psiPf[pFacei]);
351 if (boundaryDeltaExtremaCoeff > 0)
355 const label pfCelli = pFaceCells[pFacei];
357 const scalar extrema =
358 boundaryDeltaExtremaCoeff
359 *(psiMax[pfCelli] - psiMin[pfCelli]);
361 psiMaxn[pfCelli] += extrema;
362 psiMinn[pfCelli] -= extrema;
369 const label pfCelli = pFaceCells[pFacei];
371 sumPhiBD[pfCelli] += phiBDPf[pFacei];
373 const scalar phiCorrf = phiCorrPf[pFacei];
377 sumPhip[pfCelli] += phiCorrf;
381 mSumPhim[pfCelli] -= phiCorrf;
386 psiMaxn =
min(psiMaxn + extremaCoeff*(psiMax - psiMin), psiMax);
387 psiMinn =
max(psiMinn - extremaCoeff*(psiMax - psiMin), psiMin);
389 if (smoothLimiter > SMALL)
392 min(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMaxn, psiMax);
394 max(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMinn, psiMin);
404 (
rho.field()*rDeltaT -
Sp.field())*psiMaxn
407 - (V0().field()*rDeltaT)*
rho.oldTime().field()*psi0
414 - (
rho.field()*rDeltaT -
Sp.field())*psiMinn
416 + (V0().field()*rDeltaT)*
rho.oldTime().field()*psi0
424 (
rho.field()*rDeltaT -
Sp.field())*psiMaxn
426 - (
rho.oldTime().field()*rDeltaT)*psi0
434 - (
rho.field()*rDeltaT -
Sp.field())*psiMinn
435 + (
rho.oldTime().field()*rDeltaT)*psi0
443 for (
int j=0; j<nLimiterIter; j++)
450 const label own = owner[facei];
451 const label nei = neighb[facei];
453 scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
455 if (lambdaPhiCorrf > 0)
457 sumlPhip[own] += lambdaPhiCorrf;
458 mSumlPhim[nei] += lambdaPhiCorrf;
462 mSumlPhim[own] -= lambdaPhiCorrf;
463 sumlPhip[nei] -= lambdaPhiCorrf;
472 const labelUList& pFaceCells =
mesh.boundary()[patchi].faceCells();
476 const label pfCelli = pFaceCells[pFacei];
477 const scalar lambdaPhiCorrf =
478 lambdaPf[pFacei]*phiCorrfPf[pFacei];
480 if (lambdaPhiCorrf > 0)
482 sumlPhip[pfCelli] += lambdaPhiCorrf;
486 mSumlPhim[pfCelli] -= lambdaPhiCorrf;
496 (sumlPhip[celli] + psiMaxn[celli])
497 /(mSumPhim[celli] + ROOTVSMALL),
504 (mSumlPhim[celli] + psiMinn[celli])
505 /(sumPhip[celli] + ROOTVSMALL),
515 if (phiCorrIf[facei] > 0)
517 lambdaIf[facei] =
min
520 min(lambdap[owner[facei]], lambdam[neighb[facei]])
525 lambdaIf[facei] =
min
528 min(lambdam[owner[facei]], lambdap[neighb[facei]])
546 mesh.boundary()[patchi].faceCells();
550 const label pfCelli = pFaceCells[pFacei];
552 if (phiCorrfPf[pFacei] > 0)
555 min(lambdaPf[pFacei], lambdap[pfCelli]);
560 min(lambdaPf[pFacei], lambdam[pfCelli]);
582 const RdeltaTType& rDeltaT,
589 const PsiMaxType& psiMax,
590 const PsiMinType& psiMin,
591 const bool returnCorr
607 phiBDPf = phiPsiBf[patchi];
621 mesh.time().timeName(),
653 phiPsi = phiBD +
lambda*phiCorr;
674 const PsiMaxType& psiMax,
675 const PsiMinType& psiMin,
688 const scalar rDeltaT = 1.0/
mesh.time().deltaTValue();
694template<
class SurfaceScalarFieldList>
698 UPtrList<scalarField> phiPsiCorrsInternal(phiPsiCorrs.size());
707 const surfaceScalarField::Boundary& bfld =
708 phiPsiCorrs[0].boundaryField();
714 UPtrList<scalarField> phiPsiCorrsPatch(phiPsiCorrs.size());
730template<
class SurfaceScalarFieldList>
733 const SurfaceScalarFieldList& alphas,
734 SurfaceScalarFieldList& phiPsiCorrs,
735 const labelHashSet& fixed
739 UPtrList<const scalarField> alphasInternal(alphas.size());
744 UPtrList<scalarField> phiPsiCorrsInternal(phiPsiCorrs.size());
750 limitSum(alphasInternal, phiPsiCorrsInternal, fixed);
753 const surfaceScalarField::Boundary& bfld =
754 phiPsiCorrs[0].boundaryField();
760 UPtrList<const scalarField> alphasPatch(alphas.size());
766 &alphas[
phasei].boundaryField()[patchi]
769 UPtrList<scalarField> phiPsiCorrsPatch(phiPsiCorrs.size());
779 limitSum(alphasPatch, phiPsiCorrsPatch, fixed);
MULES: Multidimensional universal limiter for explicit solution.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
void size(const label n)
Older name for setAddressableSize.
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
Mesh data needed to do the Finite Volume discretisation.
virtual bool fixesValue() const
True if the patch field fixes a value.
virtual bool coupled() const
True if the patch field is coupled.
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
static const volScalarField & localRDeltaT(const fvMesh &mesh)
Return the reciprocal of the local time-step.
static bool enabled(const fvMesh &mesh)
Return true if LTS is enabled.
virtual bool coupled() const
True if the patch field is coupled.
virtual tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > flux(const GeometricField< Type, fvPatchField, volMesh > &) const
Return the interpolation weighting factors.
A class for managing temporary objects.
Upwind differencing scheme class.
A class representing the concept of a field of 0 used to avoid unnecessary manipulations for objects ...
const volScalarField & psi
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
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)
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)
void limitSum(UPtrList< scalarField > &phiPsiCorrs)
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
const dimensionSet dimless
Dimensionless.
IOstream & fixed(IOstream &io)
HashSet< label, Hash< label > > labelHashSet
A HashSet of labels, uses label hasher.
SlicedGeometricField< scalar, fvsPatchField, slicedFvsPatchField, surfaceMesh > slicedSurfaceScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
static constexpr const zero Zero
Global zero (0).
void Sp(fvMatrix< typename Expr::value_type > &m, const Expr2 &mult, const Expr &expression)
UList< label > labelUList
A UList of labels.
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
fvsPatchField< scalar > fvsPatchScalarField
fvPatchField< scalar > fvPatchScalarField
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
Us boundaryFieldRef().evaluateCoupled< coupledFaPatch >()
#define forAll(list, i)
Loop across all elements in list.