Loading...
Searching...
No Matches
createFields.H
Go to the documentation of this file.
2
3Info<< "Reading field p_rgh\n" << endl;
4volScalarField p_rgh
5(
6 IOobject
7 (
8 "p_rgh",
9 runTime.timeName(),
10 mesh,
11 IOobject::MUST_READ,
12 IOobject::AUTO_WRITE
13 ),
14 mesh
15);
16
17Info<< "Reading field U\n" << endl;
18volVectorField U
19(
20 IOobject
21 (
22 "U",
23 runTime.timeName(),
24 mesh,
25 IOobject::MUST_READ,
26 IOobject::AUTO_WRITE
27 ),
28 mesh
29);
30
31#include "createPhi.H"
32
33//- Overset specific
34
35// Add solver-specific interpolations
36{
37 wordHashSet& nonInt =
38 const_cast<wordHashSet&>(Stencil::New(mesh).nonInterpolatedFields());
39
40 nonInt.insert("HbyA");
41 nonInt.insert("grad(p_rgh)");
42 nonInt.insert("nHat");
43 nonInt.insert("surfaceIntegrate(nHatf)");
44 nonInt.insert("surfaceIntegrate(phi+meshPhi)");
45 nonInt.insert("surfaceIntegrate(phi)");
46 nonInt.insert("surfaceIntegrate(phiHbyA)");
47 nonInt.insert("surfaceSum(((S|magSf)*S)");
48 nonInt.insert("surfaceIntegrate(((rAUf*magSf)*snGradCorr(p_rgh)))");
49 nonInt.insert("cellMask");
50 nonInt.insert("cellDisplacement");
51 nonInt.insert("interpolatedCells");
52 nonInt.insert("cellInterpolationWeight");
53 nonInt.insert("pcorr");
56
57// Mask field for zeroing out contributions on hole cells
58#include "createCellMask.H"
59
60// Create bool field with interpolated cells
62
63
64Info<< "Reading transportProperties\n" << endl;
65immiscibleIncompressibleTwoPhaseMixture mixture(U, phi);
66
67volScalarField& alpha1(mixture.alpha1());
68volScalarField& alpha2(mixture.alpha2());
69
70const dimensionedScalar& rho1 = mixture.rho1();
71const dimensionedScalar& rho2 = mixture.rho2();
72
73
74// Need to store rho for ddt(rho, U)
75volScalarField rho
76(
77 IOobject
78 (
79 "rho",
80 runTime.timeName(),
81 mesh,
82 IOobject::READ_IF_PRESENT
83 ),
85);
86rho.oldTime();
87
88
89// Mass flux
90surfaceScalarField rhoPhi
91(
92 IOobject
93 (
94 "rhoPhi",
95 runTime.timeName(),
96 mesh,
97 IOobject::NO_READ,
98 IOobject::NO_WRITE
99 ),
100 fvc::interpolate(rho)*phi
101);
102
103
104// Construct incompressible turbulence model
105autoPtr<incompressible::turbulenceModel> turbulence
106(
107 incompressible::turbulenceModel::New(U, phi, mixture)
108);
109
110
111#include "readGravitationalAcceleration.H"
112#include "readhRef.H"
113#include "gh.H"
114
115
116volScalarField p
117(
118 IOobject
119 (
120 "p",
121 runTime.timeName(),
122 mesh,
123 IOobject::NO_READ,
124 IOobject::AUTO_WRITE
125 ),
126 p_rgh + rho*gh
127);
128
129label pRefCell = 0;
130scalar pRefValue = 0.0;
134 p_rgh,
135 pimple.dict(),
136 pRefCell,
138);
139
140if (p_rgh.needReference())
141{
142 p += dimensionedScalar
143 (
144 "p",
145 p.dimensions(),
146 pRefValue - getRefCellValue(p, pRefCell)
147 );
148 p_rgh = p - rho*gh;
149}
150
151mesh.setFluxRequired(p_rgh.name());
152mesh.setFluxRequired(alpha1.name());
154#include "createMRF.H"
rhoPhi
Definition rhoEqn.H:10
volScalarField & p_rgh
const scalar pRefValue
const label pRefCell
const volScalarField & gh
pimpleControl & pimple
const volScalarField & alpha1
volScalarField & rho2
const volScalarField & alpha2
volScalarField & rho1
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition HashSet.H:229
U
Definition pEqn.H:72
volScalarField & p
Creates mask for blocked out cells.
dynamicFvMesh & mesh
engineTime & runTime
Creates mask for interpolated cells.
compressible::turbulenceModel & turbulence
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
setRefCell(p, pimple.dict(), pRefCell, pRefValue)
Info<< "Creating temperaturePhaseChangeTwoPhaseMixture\n"<< endl;autoPtr< temperaturePhaseChangeTwoPhaseMixture > mixture