39template<
class RdeltaTType,
class RhoType,
class SpType,
class SuType>
42 const RdeltaTType& rDeltaT,
58 fvc::surfaceIntegrate(psiIf, phiPsi);
68 )/(
rho.field()*rDeltaT -
Sp.field());
74 rho.oldTime().field()*psi0*rDeltaT
77 )/(
rho.field()*rDeltaT -
Sp.field());
80 psi.correctBoundaryConditions();
84template<
class RhoType>
96template<
class RhoType,
class SpType,
class SuType>
108 if (fv::localEulerDdt::enabled(
mesh))
115 const scalar rDeltaT = 1.0/
mesh.time().deltaTValue();
121template<
class RhoType,
class PsiMaxType,
class PsiMinType>
128 const PsiMaxType& psiMax,
129 const PsiMinType& psiMin
162 const PsiMaxType& psiMax,
163 const PsiMinType& psiMin
168 psi.correctBoundaryConditions();
170 if (fv::localEulerDdt::enabled(
mesh))
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 labelList& 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 labelList& 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]])
539 if (isA<wedgeFvPatch>(
mesh.boundary()[patchi]))
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,
681 if (fv::localEulerDdt::enabled(
mesh))
688 const scalar rDeltaT = 1.0/
mesh.time().deltaTValue();
694template<
class SurfaceScalarFieldList>
708 phiPsiCorrs[0].boundaryField();
720 &phiPsiCorrs[
phasei].boundaryFieldRef()[patchi]
730template<
class SurfaceScalarFieldList>
733 const SurfaceScalarFieldList& alphas,
734 SurfaceScalarFieldList& phiPsiCorrs,
754 phiPsiCorrs[0].boundaryField();
766 &alphas[
phasei].boundaryField()[patchi]
775 &phiPsiCorrs[
phasei].boundaryFieldRef()[patchi]
MULES: Multidimensional universal limiter for explicit solution.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Specialization of GeometricField which holds slices of given complete fields in a form that they act ...
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
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
Mesh data needed to do the Finite Volume discretisation.
virtual bool fixesValue() const
Return true if this patch field fixes a value.
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
virtual bool coupled() const
Return true if this patch field is coupled.
virtual bool coupled() const
Return true if this 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
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
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)
const dimensionSet dimless
Dimensionless.
IOstream & fixed(IOstream &io)
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
tmp< areaScalarField > limiter(const areaScalarField &phi)
dimensionedScalar lambda("lambda", dimTime/sqr(dimLength), laminarTransport)
#define forAll(list, i)
Loop across all elements in list.