UcEqn.H
Go to the documentation of this file.
2 (
3  fvm::ddt(alphac, Uc) + fvm::div(alphaPhic, Uc)
4  - fvm::Sp(fvc::ddt(alphac) + fvc::div(alphaPhic), Uc)
5  + continuousPhaseTurbulence->divDevTau(Uc)
6  ==
7  (1.0/rhoc)*cloudSU
8 );
9 
10 UcEqn.relax();
11 
13 
14 volScalarField rAUc(1.0/UcEqn.A());
15 volScalarField rASpUc(1.0/(UcEqn.A() - cloudSUp/rhoc));
17 
19 (
20  fvc::flux(rASpUc*cloudSUu/rhoc)
21  + rASpUcf*(g & mesh.Sf())
22 );
24 (
25  fvc::interpolate(rASpUc*cloudSUp/rhoc)
26 );
27 
28 if (pimple.momentumPredictor())
29 {
30  solve
31  (
32  UcEqn
33  ==
35  (
37  - fvc::snGrad(p)*mesh.magSf()
38  )
39  + (1.0/rhoc)*(fvm::Sp(cloudSUp, Uc) - cloudSUp*Uc)
40  );
41 
43 }
pimpleNoLoopControl & pimple
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
volScalarField rASpUc(1.0/(UcEqn.A() - cloudSUp/rhoc))
surfaceScalarField rASpUcf("Dp", fvc::interpolate(rASpUc))
fvMesh & mesh
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
Foam::fvConstraints & fvConstraints
phic
Definition: correctPhic.H:2
volScalarField rhoc(IOobject(rhocValue.name(), runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, rhocValue)
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.
const tmp< volScalarField::Internal > & Sp
Definition: alphaSuSp.H:7
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Info<< "Creating field alphac\"<< endl;volScalarField alphac(IOobject(IOobject::groupName("alpha", continuousPhaseName), runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), mesh, dimensionedScalar(dimless, 0));Info<< "Constructing clouds"<< endl;parcelCloudList clouds(rhoc, Uc, muc, g);scalar alphacMin(1 - mesh.solution().solverDict(alphac.name()).lookup< scalar >"max"));alphac=max(1.0 - clouds.theta(), alphacMin);alphac.correctBoundaryConditions();surfaceScalarField alphacf("alphacf", fvc::interpolate(alphac));surfaceScalarField alphaPhic(IOobject::groupName("alphaPhi", continuousPhaseName), alphacf *phic);autoPtr< phaseIncompressible::momentumTransportModel > continuousPhaseTurbulence(phaseIncompressible::momentumTransportModel::New(alphac, Uc, alphaPhic, phic, continuousPhaseViscosity))
rhoEqn solve()
fvVectorMatrix UcEqn(fvm::ddt(alphac, Uc)+fvm::div(alphaPhic, Uc) - fvm::Sp(fvc::ddt(alphac)+fvc::div(alphaPhic), Uc)+continuousPhaseTurbulence->divDevTau(Uc)==(1.0/rhoc) *cloudSU)
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
surfaceScalarField phicSUSu(fvc::flux(rASpUc *cloudSUu/rhoc)+rASpUcf *(g &mesh.Sf()))
volScalarField rAUc(1.0/UcEqn.A())
volScalarField & p
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
const dimensionedVector & g
surfaceScalarField phicSUSp(fvc::interpolate(rASpUc *cloudSUp/rhoc))
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45