35 template<
class RdeltaTType,
class RhoType,
class SpType,
class SuType>
38 const RdeltaTType& rDeltaT,
61 )/(rho.field()*rDeltaT - Sp.field());
70 )/(rho.field()*rDeltaT - Sp.field());
77 template<
class RhoType,
class SpType,
class SuType>
92 if (fv::localEulerDdt::enabled(mesh))
94 const volScalarField& rDeltaT = fv::localEulerDdt::localRDeltaT(mesh);
108 correct(rDeltaT, rho, psi, phi, phiCorr, Sp, Su);
127 correct(rDeltaT, rho, psi, phi, phiCorr, Sp, Su);
132 template<
class RdeltaTType,
class RhoType,
class SpType,
class SuType>
136 const RdeltaTType& rDeltaT,
174 const surfaceScalarField::Boundary& phiBf =
178 const surfaceScalarField::Boundary& phiCorrBf =
199 surfaceScalarField::Boundary& lambdaBf =
210 label own = owner[facei];
211 label nei = neighb[facei];
213 psiMaxn[own] =
max(psiMaxn[own], psiIf[nei]);
214 psiMinn[own] =
min(psiMinn[own], psiIf[nei]);
216 psiMaxn[nei] =
max(psiMaxn[nei], psiIf[own]);
217 psiMinn[nei] =
min(psiMinn[nei], psiIf[own]);
219 scalar phiCorrf = phiCorrIf[facei];
223 sumPhip[own] += phiCorrf;
224 mSumPhim[nei] += phiCorrf;
228 mSumPhim[own] -= phiCorrf;
229 sumPhip[nei] -= phiCorrf;
246 label pfCelli = pFaceCells[pFacei];
248 psiMaxn[pfCelli] =
max(psiMaxn[pfCelli], psiPNf[pFacei]);
249 psiMinn[pfCelli] =
min(psiMinn[pfCelli], psiPNf[pFacei]);
256 label pfCelli = pFaceCells[pFacei];
258 psiMaxn[pfCelli] =
max(psiMaxn[pfCelli], psiPf[pFacei]);
259 psiMinn[pfCelli] =
min(psiMinn[pfCelli], psiPf[pFacei]);
265 label pfCelli = pFaceCells[pFacei];
267 scalar phiCorrf = phiCorrPf[pFacei];
271 sumPhip[pfCelli] += phiCorrf;
275 mSumPhim[pfCelli] -= phiCorrf;
280 psiMaxn =
min(psiMaxn + extremaCoeff*(psiMax - psiMin), psiMax);
281 psiMinn =
max(psiMinn - extremaCoeff*(psiMax - psiMin), psiMin);
283 if (smoothLimiter > SMALL)
286 min(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMaxn, psiMax);
288 max(smoothLimiter*psiIf + (1.0 - smoothLimiter)*psiMinn, psiMin);
294 (rho.field()*rDeltaT - Sp.field())*psiMaxn
303 - (rho.field()*rDeltaT - Sp.field())*psiMinn
310 for (
int j=0; j<nLimiterIter; j++)
317 label own = owner[facei];
318 label nei = neighb[facei];
320 scalar lambdaPhiCorrf = lambdaIf[facei]*phiCorrIf[facei];
322 if (lambdaPhiCorrf > 0.0)
324 sumlPhip[own] += lambdaPhiCorrf;
325 mSumlPhim[nei] += lambdaPhiCorrf;
329 mSumlPhim[own] -= lambdaPhiCorrf;
330 sumlPhip[nei] -= lambdaPhiCorrf;
343 label pfCelli = pFaceCells[pFacei];
345 scalar lambdaPhiCorrf = lambdaPf[pFacei]*phiCorrfPf[pFacei];
347 if (lambdaPhiCorrf > 0.0)
349 sumlPhip[pfCelli] += lambdaPhiCorrf;
353 mSumlPhim[pfCelli] -= lambdaPhiCorrf;
363 (sumlPhip[celli] + psiMaxn[celli])
364 /(mSumPhim[celli] - SMALL),
371 (mSumlPhim[celli] + psiMinn[celli])
372 /(sumPhip[celli] + SMALL),
382 if (phiCorrIf[facei] > 0.0)
384 lambdaIf[facei] =
min 387 min(lambdap[owner[facei]], lambdam[neighb[facei]])
392 lambdaIf[facei] =
min 395 min(lambdam[owner[facei]], lambdap[neighb[facei]])
418 label pfCelli = pFaceCells[pFacei];
420 if (phiCorrfPf[pFacei] > 0.0)
423 min(lambdaPf[pFacei], lambdap[pfCelli]);
428 min(lambdaPf[pFacei], lambdam[pfCelli]);
441 if (phiPf[pFacei] > SMALL*SMALL)
443 label pfCelli = pFaceCells[pFacei];
445 if (phiCorrfPf[pFacei] > 0.0)
448 min(lambdaPf[pFacei], lambdap[pfCelli]);
453 min(lambdaPf[pFacei], lambdam[pfCelli]);
465 template<
class RdeltaTType,
class RhoType,
class SpType,
class SuType>
468 const RdeltaTType& rDeltaT,
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField...
#define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
virtual bool coupled() const
Return true if this patch field is coupled.
void surfaceIntegrate(Field< Type > &ivf, const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
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 > &)
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...
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
bool moving() const
Is mesh moving.
virtual tmp< Field< Type > > patchNeighbourField() const
Return patchField on the opposite patch of a coupled patch.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const labelUList & neighbour() const
Internal face neighbour.
tmp< DimensionedField< scalar, volMesh > > Vsc() const
Return sub-cycle cell volumes.
label readLabel(Istream &is)
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.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil)-psi)*pSat, rhoMin);#1"/home/ubuntu/OpenFOAM-4.1/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()
scalar deltaTValue() const
Return time step value.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
void limitCorr(const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const scalar psiMax, const scalar psiMin)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const Mesh & mesh() const
Return mesh.
Mesh data needed to do the Finite Volume discretisation.
dimensionedScalar lambda(laminarTransport.lookup("lambda"))
void correctBoundaryConditions()
Correct boundary field.
const volScalarField & psi
const dictionary & solverDict(const word &name) const
Return the solver controls dictionary for the given field.
const labelUList & owner() const
Internal face owner.
void limiterCorr(scalarField &allLambda, const RdeltaTType &rDeltaT, const RhoType &rho, const volScalarField &psi, const surfaceScalarField &phi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su, const scalar psiMax, const scalar psiMin)
A class for managing temporary objects.
Specialization of GeometricField which holds slices of given complete fields in a form that they act ...
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
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.
const word & name() const
Return name.
const Time & time() const
Return the top-level database.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.