34 template<
class RdeltaTType,
class RhoType,
class SpType,
class SuType>
37 const RdeltaTType& rDeltaT,
59 )/(rho.field()*rDeltaT - Sp.field());
68 )/(rho.field()*rDeltaT - Sp.field());
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);
102 correct(rDeltaT, rho, psi, phiCorr, Sp, Su);
107 correct(rDeltaT, rho, psi, phiCorr, Sp, Su);
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);
166 correct(rDeltaT, rho, psi, phiCorr, Sp, Su);
185 correct(rDeltaT, rho, psi, phiCorr, Sp, Su);
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)
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 psiMaxn =
min(psiMaxn + extremaCoeff*(psiMax - psiMin), psiMax);
365 psiMinn =
max(psiMinn - extremaCoeff*(psiMax - psiMin), psiMin);
367 if (smoothLimiter > small)
370 min(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMaxn, psiMax);
372 max(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMinn, psiMin);
378 (rho.field()*rDeltaT - Sp.field())*psiMaxn
387 - (rho.field()*rDeltaT - Sp.field())*psiMinn
394 for (
int j=0; j<nLimiterIter; j++)
401 const label own = owner[facei];
402 const label nei = neighb[facei];
404 const scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
406 if (lambdaPhiCorrf > 0)
408 sumlPhip[own] += lambdaPhiCorrf;
409 mSumlPhim[nei] += lambdaPhiCorrf;
413 mSumlPhim[own] -= lambdaPhiCorrf;
414 sumlPhip[nei] -= lambdaPhiCorrf;
427 label pfCelli = pFaceCells[pFacei];
429 scalar lambdaPhiCorrf = lambdaPf[pFacei]*phiCorrfPf[pFacei];
431 if (lambdaPhiCorrf > 0)
433 sumlPhip[pfCelli] += lambdaPhiCorrf;
437 mSumlPhim[pfCelli] -= lambdaPhiCorrf;
447 (sumlPhip[celli] + psiMaxn[celli])
448 /(mSumPhim[celli] + rootVSmall),
455 (mSumlPhim[celli] + psiMinn[celli])
456 /(sumPhip[celli] + rootVSmall),
466 if (phiCorrIf[facei] > 0)
468 lambdaIf[facei] =
min 471 min(lambdap[owner[facei]], lambdam[neighb[facei]])
476 lambdaIf[facei] =
min 479 min(lambdam[owner[facei]], lambdap[neighb[facei]])
502 const label pfCelli = pFaceCells[pFacei];
504 if (phiCorrfPf[pFacei] > 0)
507 min(lambdaPf[pFacei], lambdap[pfCelli]);
512 min(lambdaPf[pFacei], lambdam[pfCelli]);
525 if ((phiPf[pFacei] + phiCorrfPf[pFacei]) > small*small)
527 const label pfCelli = pFaceCells[pFacei];
529 if (phiCorrfPf[pFacei] > 0)
532 min(lambdaPf[pFacei], lambdap[pfCelli]);
537 min(lambdaPf[pFacei], lambdam[pfCelli]);
547 surfaceScalarField::Internal::null(),
556 lambdaPf =
min(lambdaPf, lambdaNbrPf);
574 const RdeltaTType& rDeltaT,
581 const PsiMaxType& psiMax,
582 const PsiMinType& psiMin
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
#define forAll(list, i)
Loop across all elements in list.
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)
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...
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
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 Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
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.
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
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)
tmp< GeometricBoundaryField > boundaryNeighbourField() const
Return BoundaryField of the values on the other side of couples.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
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.
Mesh data needed to do the Finite Volume discretisation.
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 correctBoundaryConditions()
Correct boundary field.
A class for managing temporary objects.
virtual tmp< Field< Type > > patchNeighbourField(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking) const
Return patchField on the opposite patch of a coupled patch.
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...
CMULES: Multidimensional universal limiter for explicit corrected implicit solution.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.