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->divDevRhoReff(Uc)
6  ==
7  (1.0/rhoc)*cloudSU
8 );
9 
10 UcEqn.relax();
11 
12 volScalarField rAUc(1.0/UcEqn.A());
14 
16 (
17  fvc::flux(rAUc*cloudVolSUSu/rhoc) + rAUcf*(g & mesh.Sf())
18 );
19 
20 if (pimple.momentumPredictor())
21 {
22  solve
23  (
24  UcEqn
25  ==
27  (
28  phicForces/rAUcf - fvc::snGrad(p)*mesh.magSf()
29  )
30  );
31 }
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
const dictionary & pimple
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:67
tmp< surfaceScalarField > interpolate(const RhoType &rho)
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
surfaceScalarField rAUcf("Dp", fvc::interpolate(rAUc))
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
rhoEqn solve()
const dimensionedVector & g
Info<< "Reading field U\n"<< endl;volVectorField Uc(IOobject(IOobject::groupName("U", continuousPhaseName), runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field p\n"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading/calculating continuous-phase face flux field phic\n"<< endl;surfaceScalarField phic(IOobject(IOobject::groupName("phi", continuousPhaseName), runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), linearInterpolate(Uc)&mesh.Sf());label pRefCell=0;scalar pRefValue=0.0;setRefCell(p, pimple.dict(), pRefCell, pRefValue);mesh.setFluxRequired(p.name());Info<< "Creating turbulence model\n"<< endl;singlePhaseTransportModel continuousPhaseTransport(Uc, phic);dimensionedScalar rhocValue(IOobject::groupName("rho", continuousPhaseName), dimDensity, continuousPhaseTransport.lookup(IOobject::groupName("rho", continuousPhaseName)));volScalarField rhoc(IOobject(rhocValue.name(), runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, rhocValue);volScalarField muc(IOobject(IOobject::groupName("mu", continuousPhaseName), runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), rhoc *continuousPhaseTransport.nu());Info<< "Creating field alphac\n"<< endl;volScalarField alphac(IOobject(IOobject::groupName("alpha", continuousPhaseName), runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), mesh, dimensionedScalar("0", dimless, 0));word kinematicCloudName("kinematicCloud");args.optionReadIfPresent("cloudName", kinematicCloudName);Info<< "Constructing kinematicCloud "<< kinematicCloudName<< endl;basicKinematicTypeCloud kinematicCloud(kinematicCloudName, rhoc, Uc, muc, g);scalar alphacMin(1.0-readScalar(kinematicCloud.particleProperties().subDict("constantProperties").lookup("alphaMax")));alphac=max(1.0-kinematicCloud.theta(), alphacMin);alphac.correctBoundaryConditions();surfaceScalarField alphacf("alphacf", fvc::interpolate(alphac));surfaceScalarField alphaPhic("alphaPhic", alphacf *phic);autoPtr< PhaseIncompressibleTurbulenceModel< singlePhaseTransportModel > > continuousPhaseTurbulence(PhaseIncompressibleTurbulenceModel< singlePhaseTransportModel >::New(alphac, Uc, alphaPhic, phic, continuousPhaseTransport))
Definition: createFields.H:158
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> reconstruct(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
volScalarField rAUc(1.0/UcEqn.A())
fvVectorMatrix UcEqn(fvm::ddt(alphac, Uc)+fvm::div(alphaPhic, Uc)-fvm::Sp(fvc::ddt(alphac)+fvc::div(alphaPhic), Uc)+continuousPhaseTurbulence->divDevRhoReff(Uc)==(1.0/rhoc)*cloudSU)
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
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
surfaceScalarField phicForces(fvc::flux(rAUc *cloudVolSUSu/rhoc)+rAUcf *(g &mesh.Sf()))