alphaEqn.H
Go to the documentation of this file.
1 {
2  word alphaScheme("div(phi,alpha)");
3  word alpharScheme("div(phirb,alpha)");
4 
5  // Set the off-centering coefficient according to ddt scheme
6  scalar ocCoeff = 0;
7  {
8  tmp<fv::ddtScheme<scalar>> tddtAlpha
9  (
11  (
12  mesh,
13  mesh.ddtScheme("ddt(alpha)")
14  )
15  );
16  const fv::ddtScheme<scalar>& ddtAlpha = tddtAlpha();
17 
18  if
19  (
20  isType<fv::EulerDdtScheme<scalar>>(ddtAlpha)
21  || isType<fv::localEulerDdtScheme<scalar>>(ddtAlpha)
22  )
23  {
24  ocCoeff = 0;
25  }
26  else if (isType<fv::CrankNicolsonDdtScheme<scalar>>(ddtAlpha))
27  {
28  if (nAlphaSubCycles > 1)
29  {
31  << "Sub-cycling is not supported "
32  "with the CrankNicolson ddt scheme"
33  << exit(FatalError);
34  }
35 
36  if
37  (
39  || mesh.time().timeIndex() > mesh.time().startTimeIndex() + 1
40  )
41  {
42  ocCoeff =
43  refCast<const fv::CrankNicolsonDdtScheme<scalar>>(ddtAlpha)
44  .ocCoeff();
45  }
46  }
47  else
48  {
50  << "Only Euler and CrankNicolson ddt schemes are supported"
51  << exit(FatalError);
52  }
53  }
54 
55  // Set the time blending factor, 1 for Euler
56  scalar cnCoeff = 1.0/(1.0 + ocCoeff);
57 
58  // Standard face-flux compression coefficient
59  surfaceScalarField phic(mixture.cAlpha()*mag(phi/mesh.magSf()));
60 
61  // Add the optional isotropic compression contribution
62  if (icAlpha > 0)
63  {
64  phic *= (1.0 - icAlpha);
65  phic += (mixture.cAlpha()*icAlpha)*fvc::interpolate(mag(U));
66  }
67 
68  // Add the optional shear compression contribution
69  if (scAlpha > 0)
70  {
71  phic +=
73  }
74 
75 
76  surfaceScalarField::Boundary& phicBf =
77  phic.boundaryFieldRef();
78 
79  // Do not compress interface at non-coupled boundary faces
80  // (inlets, outlets etc.)
81  forAll(phic.boundaryField(), patchi)
82  {
83  fvsPatchScalarField& phicp = phicBf[patchi];
84 
85  if (!phicp.coupled())
86  {
87  phicp == 0;
88  }
89  }
90 
91  tmp<surfaceScalarField> phiCN(phi);
92 
93  // Calculate the Crank-Nicolson off-centred volumetric flux
94  if (ocCoeff > 0)
95  {
96  phiCN = cnCoeff*phi + (1.0 - cnCoeff)*phi.oldTime();
97  }
98 
99  if (MULESCorr)
100  {
101  #include "alphaSuSp.H"
102 
103  fvScalarMatrix alpha1Eqn
104  (
105  (
106  LTS
107  ? fv::localEulerDdtScheme<scalar>(mesh).fvmDdt(alpha1)
108  : fv::EulerDdtScheme<scalar>(mesh).fvmDdt(alpha1)
109  )
110  + fv::gaussConvectionScheme<scalar>
111  (
112  mesh,
113  phiCN,
114  upwind<scalar>(mesh, phiCN)
115  ).fvmDiv(phiCN, alpha1)
116  // - fvm::Sp(fvc::ddt(dimensionedScalar("1", dimless, 1), mesh)
117  // + fvc::div(phiCN), alpha1)
118  ==
119  Su + fvm::Sp(Sp + divU, alpha1)
120  );
121 
122  alpha1Eqn.solve();
123 
124  Info<< "Phase-1 volume fraction = "
125  << alpha1.weightedAverage(mesh.Vsc()).value()
126  << " Min(" << alpha1.name() << ") = " << min(alpha1).value()
127  << " Max(" << alpha1.name() << ") = " << max(alpha1).value()
128  << endl;
129 
130  tmp<surfaceScalarField> talphaPhi1UD(alpha1Eqn.flux());
131  alphaPhi10 = talphaPhi1UD();
132 
133  if (alphaApplyPrevCorr && talphaPhi1Corr0.valid())
134  {
135  Info<< "Applying the previous iteration compression flux" << endl;
137  (
138  geometricOneField(),
139  alpha1,
140  alphaPhi10,
141  talphaPhi1Corr0.ref(),
142  oneField(),
143  zeroField()
144  );
145 
146  alphaPhi10 += talphaPhi1Corr0();
147  }
148 
149  // Cache the upwind-flux
150  talphaPhi1Corr0 = talphaPhi1UD;
151 
152  alpha2 = 1.0 - alpha1;
153 
154  mixture.correct();
155  }
156 
157 
158  for (int aCorr=0; aCorr<nAlphaCorr; aCorr++)
159  {
160  #include "alphaSuSp.H"
161 
163 
164  tmp<surfaceScalarField> talphaPhi1Un
165  (
166  fvc::flux
167  (
168  phiCN(),
169  cnCoeff*alpha1 + (1.0 - cnCoeff)*alpha1.oldTime(),
170  alphaScheme
171  )
172  + fvc::flux
173  (
175  alpha1,
177  )
178  );
179 
180  if (MULESCorr)
181  {
182  tmp<surfaceScalarField> talphaPhi1Corr(talphaPhi1Un() - alphaPhi10);
183  volScalarField alpha10("alpha10", alpha1);
184 
186  (
187  geometricOneField(),
188  alpha1,
189  talphaPhi1Un(),
190  talphaPhi1Corr.ref(),
191  Sp,
192  (-Sp*alpha1)(),
193  oneField(),
194  zeroField()
195  );
196 
197  // Under-relax the correction for all but the 1st corrector
198  if (aCorr == 0)
199  {
200  alphaPhi10 += talphaPhi1Corr();
201  }
202  else
203  {
204  alpha1 = 0.5*alpha1 + 0.5*alpha10;
205  alphaPhi10 += 0.5*talphaPhi1Corr();
206  }
207  }
208  else
209  {
210  alphaPhi10 = talphaPhi1Un;
211 
213  (
214  geometricOneField(),
215  alpha1,
216  phiCN,
217  alphaPhi10,
218  Sp,
219  (Su + divU*min(alpha1(), scalar(1)))(),
220  oneField(),
221  zeroField()
222  );
223  }
224 
225  alpha2 = 1.0 - alpha1;
226 
227  mixture.correct();
228  }
229 
231  {
233  talphaPhi1Corr0.ref().rename("alphaPhi1Corr0");
234  }
235  else
236  {
237  talphaPhi1Corr0.clear();
238  }
239 
240  #include "rhofs.H"
241 
242  if
243  (
244  word(mesh.ddtScheme("ddt(rho,U)"))
245  == fv::EulerDdtScheme<vector>::typeName
246  || word(mesh.ddtScheme("ddt(rho,U)"))
247  == fv::localEulerDdtScheme<vector>::typeName
248  )
249  {
250  rhoPhi = alphaPhi10*(rho1f - rho2f) + phiCN*rho2f;
251  }
252  else
253  {
254  if (ocCoeff > 0)
255  {
256  // Calculate the end-of-time-step alpha flux
257  alphaPhi10 =
258  (alphaPhi10 - (1.0 - cnCoeff)*alphaPhi10.oldTime())/cnCoeff;
259  }
260 
261  // Calculate the end-of-time-step mass flux
263  }
264 
265  Info<< "Phase-1 volume fraction = "
266  << alpha1.weightedAverage(mesh.Vsc()).value()
267  << " Min(" << alpha1.name() << ") = " << min(alpha1).value()
268  << " Max(" << alpha1.name() << ") = " << max(alpha1).value()
269  << endl;
270 }
forAll(phic.boundaryField(), patchi)
Definition: alphaEqn.H:81
fvsPatchField< scalar > fvsPatchScalarField
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
zeroField Su
Definition: alphaSuSp.H:1
surfaceScalarField & phi
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
surfaceScalarField rho1f(fvc::interpolate(rho1))
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
const bool alphaRestart
word alpharScheme("div(phirb,alpha)")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
scalar cnCoeff
Definition: alphaEqn.H:56
const fv::ddtScheme< scalar > & ddtAlpha
Definition: alphaEqn.H:16
surfaceScalarField alphaPhi10(alphaPhi10Header, phi *fvc::interpolate(alpha1))
label nAlphaCorr(readLabel(alphaControls.lookup("nAlphaCorr")))
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
scalar icAlpha(alphaControls.lookupOrDefault< scalar >("icAlpha", 0))
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
tmp< surfaceScalarField > phiCN(phi)
surfaceScalarField rho2f(fvc::interpolate(rho2))
dynamicFvMesh & mesh
rhoPhi
Definition: alphaEqn.H:116
volScalarField alpha10("alpha10", alpha1)
bool isType(const Type &t)
Check the typeid.
Definition: typeInfo.H:126
const volScalarField & alpha1
alpha2
Definition: alphaEqn.H:115
bool LTS
Definition: createRDeltaT.H:1
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-6/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()
Definition: pEqn.H:72
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Creating phaseChangeTwoPhaseMixture\"<< endl;autoPtr< phaseChangeTwoPhaseMixture > mixture
Definition: createFields.H:33
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
label patchi
U
Definition: pEqn.H:72
surfaceScalarField phic(mixture.cAlpha() *mag(phi/mesh.magSf()))
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
surfaceScalarField phir(IOobject("phir", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mixture.cAlpha() *mag(phi/mesh.magSf()) *mixture.nHatf())
bool alphaApplyPrevCorr(alphaControls.lookupOrDefault< Switch >("alphaApplyPrevCorr", false))
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
zeroField divU
Definition: alphaSuSp.H:3
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
void explicitSolve(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiPsi, const SpType &Sp, const SuType &Su)
bool MULESCorr(alphaControls.lookupOrDefault< Switch >("MULESCorr", false))
scalar scAlpha(alphaControls.lookupOrDefault< scalar >("scAlpha", 0))
label nAlphaSubCycles(readLabel(alphaControls.lookup("nAlphaSubCycles")))
zeroField Sp
Definition: alphaSuSp.H:2
scalar ocCoeff
Definition: alphaEqn.H:6
tmp< surfaceScalarField > talphaPhi1Corr0