All Classes Namespaces Files Functions Variables Typedefs Friends Macros Groups
createFields.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10  Copyright (C) 2016-2019 ISP RAS (www.ispras.ru) UniCFD Group (www.unicfd.ru)
11 -------------------------------------------------------------------------------
12 License
13  This file is part of QGDsolver library, based on OpenFOAM+.
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22  You should have received a copy of the GNU General Public License
23  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
24 Global
25  createFields
26 Description
27  Creating fields for calculations
28 SourceFile
29  mulesQHDFoam.C
30 \*---------------------------------------------------------------------------*/
31 
32 Info<< "\nReading gravitationalProperties" << endl;
33 
34 IOdictionary gravitationalProperties
35 (
36  IOobject
37  (
38  "gravitationalProperties",
39  runTime.constant(),
40  mesh,
41  IOobject::MUST_READ_IF_MODIFIED,
42  IOobject::NO_WRITE
43  )
44 );
45 
46 const dimensionedVector g(gravitationalProperties.lookup("g"));
47 
48 #include "readhRef.H"
49 #include "gh.H"
50 
51 Info << "Reading thermophysical properties\n" << endl;
52 
53 autoPtr<rhoQGDThermo> pThermo
54 (
55  rhoQGDThermo::New(mesh)
56 );
57 rhoQGDThermo& thermo = pThermo();
58 thermo.correct();
59 
60 volScalarField& e = thermo.he();
61 
62 volScalarField& p = thermo.p();
63 volScalarField& T = const_cast<volScalarField&>(thermo.T());
64 const surfaceScalarField& hQGDf = thermo.hQGDf();
65 const surfaceScalarField& tauQGDf = thermo.tauQGDf();
66 
67 Info << "Thermo corrected" << endl;
68 
69 autoPtr<compressible::turbulenceModel> turbulence;
70 
71 volVectorField U
72 (
73  IOobject
74  (
75  "U",
76  runTime.timeName(),
77  mesh,
78  IOobject::MUST_READ,
79  IOobject::AUTO_WRITE
80  ),
81  mesh
82 );
83 
84 volScalarField T0
85 (
86  "T0",
87  T
88 );
89 
90 volScalarField rho
91 (
92  IOobject
93  (
94  "rho",
95  runTime.timeName(),
96  mesh,
97  IOobject::NO_READ,
98  IOobject::AUTO_WRITE
99  ),
100  thermo.rho()
101 );
102 
103 volVectorField W
104 (
105  IOobject
106  (
107  "W",
108  runTime.timeName(),
109  mesh,
110  IOobject::NO_READ,
111  IOobject::NO_WRITE
112  ),
113  U
114 );
115 
116 dimensionedScalar beta
117 (
118  "beta",
119  dimless / dimTemperature,
120  thermo.subDict("mixture").subDict("transport")
121 );
122 
123 
124 surfaceScalarField phiu
125 (
126  "phiu",
127  mesh.Sf() & linearInterpolate(U)
128 );
129 
130 surfaceScalarField phiwo
131 (
132  "phiwStar",
133  mesh.Sf() & linearInterpolate(W)
134 );
135 
136 surfaceScalarField phi
137 (
138  "phi",
139  mesh.Sf() & (linearInterpolate(U) - linearInterpolate(W))
140 );
141 
142 surfaceScalarField phiRhof
143 (
144  "phiRhof",
145  linearInterpolate(rho)*phi
146 );
147 
148 
149 volVectorField BdFrc
150 (
151  "BdFrc",
152  T* g* beta
153 );
154 
155 Switch implicitDiffusion (thermo.implicitDiffusion());
156 
157 Info << "Creating turbulence model\n" << endl;
158 turbulence.reset
159 (
160  compressible::turbulenceModel::New
161  (
162  rho,
163  U,
164  phiRhof,
165  thermo
166  ).ptr()
167 );
168 rho.oldTime();
169 
170 label pRefCell = 0;
171 scalar pRefValue = 0.0;
172 setRefCell(p, thermo.subDict("QGD"), pRefCell, pRefValue);
173 
174 #include "createIcoZeroSources.H"
175 
176 //
177 //END-OF-FILE
178 //
179 
volScalarField T0("T0", T)
volVectorField W(IOobject("W", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), U)
const surfaceScalarField & tauQGDf
Definition: createFields.H:55
const surfaceScalarField & hQGDf
Definition: createFields.H:54
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh)
Switch implicitDiffusion(thermo.implicitDiffusion())
dimensionedScalar beta("beta", dimless/dimTemperature, thermo.subDict("mixture").subDict("transport"))
surfaceScalarField phiRhof("phiRhof", linearInterpolate(rho)*phi)
BdFrc
Definition: updateFields.H:35
volScalarField & e
Definition: createFields.H:50
rhoQGDThermo & thermo
Definition: createFields.H:47
Info<< "\nReading gravitationalProperties"<< endl;IOdictionary gravitationalProperties(IOobject("gravitationalProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE));const dimensionedVector g(gravitationalProperties.lookup("g"));Info<< "Reading thermophysical properties\n"<< endl;autoPtr< rhoQGDThermo > pThermo(rhoQGDThermo::New(mesh))
volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho())
phiu
——–Start———
Definition: updateFluxes.H:33
phiwo
Definition: updateFluxes.H:8
Info<< "Thermo corrected"<< endl;autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:59
volScalarField & T
Definition: createFields.H:53
volScalarField & p
Definition: createFields.H:52