createFields.H
Go to the documentation of this file.
1 Info<< "Reading field p_rgh\n" << endl;
3 (
4  IOobject
5  (
6  "p_rgh",
7  runTime.timeName(),
8  mesh,
9  IOobject::MUST_READ,
10  IOobject::AUTO_WRITE
11  ),
12  mesh
13 );
14 
16 (
17  IOobject
18  (
19  "U",
20  runTime.timeName(),
21  mesh,
22  IOobject::NO_READ,
23  IOobject::AUTO_WRITE
24  ),
25  mesh,
26  dimensionedVector("U", dimVelocity, vector::zero)
27 );
28 
30 (
31  IOobject
32  (
33  "phi",
34  runTime.timeName(),
35  mesh,
36  IOobject::NO_READ,
37  IOobject::AUTO_WRITE
38  ),
39  mesh,
41 );
42 
43 multiphaseSystem fluid(U, phi);
44 
45 forAllIter(PtrDictionary<phaseModel>, fluid.phases(), iter)
46 {
47  phaseModel& phase = iter();
48  const volScalarField& alpha = phase;
49 
50  U += alpha*phase.U();
51  phi += fvc::interpolate(alpha)*phase.phi();
52 }
53 
54 scalar slamDampCoeff
55 (
56  fluid.lookupOrDefault<scalar>("slamDampCoeff", 1)
57 );
58 
59 dimensionedScalar maxSlamVelocity
60 (
61  "maxSlamVelocity",
63  fluid.lookupOrDefault<scalar>("maxSlamVelocity", GREAT)
64 );
65 
66 
68 (
69  IOobject
70  (
71  "rho",
72  runTime.timeName(),
73  mesh,
74  IOobject::NO_READ,
75  IOobject::AUTO_WRITE
76  ),
77  fluid.rho()
78 );
79 
80 
81 // Construct incompressible turbulence model
82 autoPtr<incompressible::turbulenceModel> turbulence
83 (
85 );
86 
87 
88 #include "readGravitationalAcceleration.H"
89 #include "readhRef.H"
90 #include "gh.H"
91 
92 
94 (
95  IOobject
96  (
97  "p",
98  runTime.timeName(),
99  mesh,
100  IOobject::NO_READ,
101  IOobject::AUTO_WRITE
102  ),
103  p_rgh + rho*gh
104 );
105 
107 scalar pRefValue = 0.0;
109 (
110  p,
111  p_rgh,
112  pimple.dict(),
113  pRefCell,
114  pRefValue
115 );
116 mesh.setFluxRequired(p_rgh.name());
pRefCell
label pRefCell
Definition: createFields.H:112
gh.H
pRefValue
scalar pRefValue
Definition: createFields.H:113
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
forAllIter
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:431
Foam::dimVelocity
const dimensionSet dimVelocity
fluid
multiphaseSystem & fluid
Definition: createFields.H:10
setRefCell
setRefCell(p, p_rgh, pimple.dict(), pRefCell, pRefValue)
phi
surfaceScalarField & phi
Definition: createFields.H:13
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
rho
rho
Definition: createFields.H:79
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:48
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::Info
messageStream Info
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
p
volScalarField & p
Definition: createFields.H:51
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
gh
const volScalarField & gh
Definition: setRegionFluidFields.H:34
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
readhRef.H
alpha
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
pimple
const dictionary & pimple
Definition: readFluidMultiRegionPIMPLEControls.H:1
U
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), Urel+SRF->U())
p_rgh
volScalarField & p_rgh
Definition: setRegionFluidFields.H:31