createFields.H
Go to the documentation of this file.
2 (
3  IOdictionary
4  (
5  IOobject
6  (
7  "transportProperties",
8  runTime.constant(),
9  mesh,
10  IOobject::MUST_READ
11  )
12  ).lookup("continuousPhaseName")
13 );
14 
15 Info<< "Reading field U\n" << endl;
17 (
18  IOobject
19  (
20  IOobject::groupName("U", continuousPhaseName),
21  runTime.timeName(),
22  mesh,
23  IOobject::MUST_READ,
24  IOobject::AUTO_WRITE
25  ),
26  mesh
27 );
28 
29 Info<< "Reading field p\n" << endl;
31 (
32  IOobject
33  (
34  "p",
35  runTime.timeName(),
36  mesh,
37  IOobject::MUST_READ,
38  IOobject::AUTO_WRITE
39  ),
40  mesh
41 );
42 
43 
44 Info<< "Reading/calculating continuous-phase face flux field phic\n"
45  << endl;
46 
48 (
49  IOobject
50  (
51  IOobject::groupName("phi", continuousPhaseName),
52  runTime.timeName(),
53  mesh,
54  IOobject::READ_IF_PRESENT,
55  IOobject::AUTO_WRITE
56  ),
57  linearInterpolate(Uc) & mesh.Sf()
58 );
59 
60 label pRefCell = 0;
61 scalar pRefValue = 0.0;
63 mesh.setFluxRequired(p.name());
64 
65 Info<< "Creating turbulence model\n" << endl;
66 
67 singlePhaseTransportModel continuousPhaseTransport(Uc, phic);
68 
69 dimensionedScalar rhocValue
70 (
71  IOobject::groupName("rho", continuousPhaseName),
72  dimDensity,
73  continuousPhaseTransport.lookup
74  (
75  IOobject::groupName("rho", continuousPhaseName)
76  )
77 );
78 
79 volScalarField rhoc
80 (
81  IOobject
82  (
83  rhocValue.name(),
84  runTime.timeName(),
85  mesh,
86  IOobject::NO_READ,
87  IOobject::AUTO_WRITE
88  ),
89  mesh,
90  rhocValue
91 );
92 
94 (
95  IOobject
96  (
97  IOobject::groupName("mu", continuousPhaseName),
98  runTime.timeName(),
99  mesh,
100  IOobject::NO_READ,
101  IOobject::AUTO_WRITE
102  ),
103  rhoc*continuousPhaseTransport.nu()
104 );
105 
106 Info << "Creating field alphac\n" << endl;
107 // alphac must be constructed before the cloud
108 // so that the drag-models can find it
109 volScalarField alphac
110 (
111  IOobject
112  (
113  IOobject::groupName("alpha", continuousPhaseName),
114  runTime.timeName(),
115  mesh,
116  IOobject::READ_IF_PRESENT,
117  IOobject::AUTO_WRITE
118  ),
119  mesh,
120  dimensionedScalar("0", dimless, 0)
121 );
122 
123 word kinematicCloudName("kinematicCloud");
125 
126 Info<< "Constructing kinematicCloud " << kinematicCloudName << endl;
127 basicKinematicTypeCloud kinematicCloud
128 (
130  rhoc,
131  Uc,
132  muc,
133  g
134 );
135 
136 // Particle fraction upper limit
137 scalar alphacMin
138 (
139  1.0
140  - readScalar
141  (
142  kinematicCloud.particleProperties().subDict("constantProperties")
143  .lookup("alphaMax")
144  )
145 );
146 
147 // Update alphac from the particle locations
148 alphac = max(1.0 - kinematicCloud.theta(), alphacMin);
149 alphac.correctBoundaryConditions();
150 
151 surfaceScalarField alphacf("alphacf", fvc::interpolate(alphac));
152 surfaceScalarField alphaPhic("alphaPhic", alphacf*phic);
153 
154 autoPtr<PhaseIncompressibleTurbulenceModel<singlePhaseTransportModel> >
156 (
158  (
159  alphac,
160  Uc,
161  alphaPhic,
162  phic,
163  continuousPhaseTransport
164  )
165 );
pRefCell
label pRefCell
Definition: createFields.H:112
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
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
Foam::dimDensity
const dimensionSet dimDensity
setRefCell
setRefCell(p, p_rgh, pimple.dict(), pRefCell, pRefValue)
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
g
const dimensionedVector & g
Definition: setRegionFluidFields.H:33
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
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
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
continuousPhaseTurbulence
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:156
kinematicCloudName
word kinematicCloudName("kinematicCloud")
p
volScalarField & p
Definition: createFields.H:51
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
phic
phic
Definition: alphaEqnsSubCycle.H:5
readScalar
#define readScalar
Definition: doubleScalar.C:38
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::linearInterpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: linear.H:107
continuousPhaseName
word continuousPhaseName(IOdictionary(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ)).lookup("continuousPhaseName"))
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
pimple
const dictionary & pimple
Definition: readFluidMultiRegionPIMPLEControls.H:1
args
Foam::argList args(argc, argv)
Foam::argList::optionReadIfPresent
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198