Loading...
Searching...
No Matches
createConstants.H
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ██╗████████╗██╗ ██╗ █████╗ ██████╗ █████╗ ███████╗██╗ ██╗
3 ██║╚══██╔══╝██║ ██║██╔══██╗██╔════╝██╔══██╗ ██╔════╝██║ ██║
4 ██║ ██║ ███████║███████║██║ ███████║█████╗█████╗ ██║ ██║
5 ██║ ██║ ██╔══██║██╔══██║██║ ██╔══██║╚════╝██╔══╝ ╚██╗ ██╔╝
6 ██║ ██║ ██║ ██║██║ ██║╚██████╗██║ ██║ ██║ ╚████╔╝
7 ╚═╝ ╚═╝ ╚═╝ ╚═╝╚═╝ ╚═╝ ╚═════╝╚═╝ ╚═╝ ╚═╝ ╚═══╝
8
9 * In real Time Highly Advanced Computational Applications for Finite Volumes
10 * Copyright (C) 2017 by the ITHACA-FV authors
11-------------------------------------------------------------------------------
12
13 License
14 This file is part of ITHACA-FV
15
16 ITHACA-FV is free software: you can redistribute it and/or modify
17 it under the terms of the GNU Lesser General Public License as published by
18 the Free Software Foundation, either version 3 of the License, or
19 (at your option) any later version.
20
21 ITHACA-FV is distributed in the hope that it will be useful,
22 but WITHOUT ANY WARRANTY; without even the implied warranty of
23 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
24 GNU Lesser General Public License for more details.
25
26 You should have received a copy of the GNU Lesser General Public License
27 along with ITHACA-FV. If not, see <http://www.gnu.org/licenses/>.
28
29\*---------------------------------------------------------------------------*/
30#include "readConstants.H"
31#include "readConstants_N.H"
32#include "readConstants_TH.H"
33
34
35_betaTot = autoPtr<dimensionedScalar>
36 (
37 new dimensionedScalar
38 (
40 )
41 );
42dimensionedScalar& betaTot = _betaTot();
43
44_decbetaTot = autoPtr<dimensionedScalar>
45 (
46 new dimensionedScalar
47 (
49 )
50 );
51dimensionedScalar& decbetaTot = _decbetaTot();
52
53_v = autoPtr<volScalarField>
54 (
55 new volScalarField
56 (
57 IOobject
58 (
59 "v",
60 runTime.timeName(),
61 mesh,
62 IOobject::NO_READ,
63 IOobject::AUTO_WRITE
64 ),
65 mesh,
66 dimensionedScalar("v", dimensionSet(-1, 3, 0, 0, 0, 0, 0), 1 / rhoRef.value())
67 )
68 );
69volScalarField& v = _v();
70
71volScalarField v0(v);
72
73_v0 = autoPtr<volScalarField>
74 (
75 new volScalarField(v0)
76 );
77
78_D = autoPtr<volScalarField>
79 (
80 new volScalarField
81 (
82 IOobject
83 (
84 "D",
85 runTime.timeName(),
86 mesh,
87 IOobject::NO_READ,
88 IOobject::AUTO_WRITE
89 ),
90 mesh,
91 dimensionedScalar("D", dimensionSet(0, 1, 0, 0, 0, 0, 0), D1_0.value())
92 )
93 );
94volScalarField& D = _D();
95
96volScalarField D0(D);
97
98_D0 = autoPtr<volScalarField>
99 (
100 new volScalarField(D0)
101 );
102
103_NSF = autoPtr<volScalarField>
104 (
105 new volScalarField
106 (
107 IOobject
108 (
109 "NSF",
110 runTime.timeName(),
111 mesh,
112 IOobject::NO_READ,
113 IOobject::AUTO_WRITE
114 ),
115 mesh,
116 dimensionedScalar("NSF", dimensionSet(0, -1, 0, 0, 0, 0, 0), NSF1_0.value())
117 )
118 );
119volScalarField& NSF = _NSF();
120
121volScalarField NSF0(NSF);
122
123_NSF0 = autoPtr<volScalarField>
124 (
125 new volScalarField(NSF0)
126 );
127
128_A = autoPtr<volScalarField>
129 (
130 new volScalarField
131 (
132 IOobject
133 (
134 "A",
135 runTime.timeName(),
136 mesh,
137 IOobject::NO_READ,
138 IOobject::AUTO_WRITE
139 ),
140 mesh,
141 dimensionedScalar("A", dimensionSet(0, -1, 0, 0, 0, 0, 0), A1_0.value())
142 )
143 );
144volScalarField& A = _A();
145
146volScalarField A0(A);
147
148_A0 = autoPtr<volScalarField>
149 (
150 new volScalarField(A0)
151 );
152
153
154_SP = autoPtr<volScalarField>
155 (
156 new volScalarField
157 (
158 IOobject
159 (
160 "SP",
161 runTime.timeName(),
162 mesh,
163 IOobject::NO_READ,
164 IOobject::AUTO_WRITE
165 ),
166 mesh,
167 dimensionedScalar("SP", dimensionSet(1, 1, -2, 0, 0, 0, 0), SP1_0.value())
168 )
169 );
170volScalarField& SP = _SP();
171
172volScalarField SP0(SP);
173
174_SP0 = autoPtr<volScalarField>
175 (
176 new volScalarField(SP0)
177 );
178
179_TXS = autoPtr<volScalarField>
180 (
181 new volScalarField
182 (
183 IOobject
184 (
185 "TXS",
186 runTime.timeName(),
187 mesh,
188 IOobject::NO_READ,
189 IOobject::AUTO_WRITE
190 ),
191 v* SP
192 )
193 );
194volScalarField& TXS = _TXS();
195
196volScalarField TXS0(TXS);
197
198_TXS0 = autoPtr<volScalarField>
199 (
200 new volScalarField(TXS0)
201 );
202
203_logT = autoPtr<volScalarField>
204 (
205 new volScalarField
206 (
207 IOobject
208 (
209 "logT",
210 runTime.timeName(),
211 mesh,
212 IOobject::NO_READ,
213 IOobject::AUTO_WRITE
214 ),
215 log(max(0.5, log(T / TrefXS)))
216 )
217 );
218volScalarField& logT = _logT();
219
220_alphat = autoPtr<volScalarField>
221 (
222 new volScalarField
223 (
224 IOobject
225 (
226 "alphat",
227 runTime.timeName(),
228 mesh,
229 IOobject::NO_READ,
230 IOobject::AUTO_WRITE
231 ),
232 turbulence->nut() / Prt
233 )
234 );
235volScalarField& alphat = _alphat();
236
237_difft = autoPtr<volScalarField>
238 (
239 new volScalarField
240 (
241 IOobject
242 (
243 "difft",
244 runTime.timeName(),
245 mesh,
246 IOobject::NO_READ,
247 IOobject::AUTO_WRITE
248 ),
249 turbulence->nut() / Sct
250 )
251 );
252volScalarField& difft = _difft();
Foam::fvMesh & mesh
Definition createMesh.H:47
Foam::Time & runTime
Definition createTime.H:33
volScalarField & T
Definition createT.H:46
_decbetaTot
volScalarField & logT
volScalarField & SP
volScalarField NSF0(NSF)
volScalarField D0(D)
dimensionedScalar & decbetaTot
volScalarField & NSF
volScalarField A0(A)
volScalarField & A
dimensionedScalar & betaTot
volScalarField & D
volScalarField v0(v)
volScalarField SP0(SP)
volScalarField & v
volScalarField & TXS
volScalarField & difft
volScalarField & alphat
_betaTot
volScalarField TXS0(TXS)
dimensionedScalar & rhoRef
dimensionedScalar & Prt
dimensionedScalar & Sct
dimensionedScalar & beta1
dimensionedScalar & beta3
dimensionedScalar & D1_0
dimensionedScalar & A1_0
dimensionedScalar & beta2
dimensionedScalar & beta5
dimensionedScalar & beta4
dimensionedScalar & beta8
dimensionedScalar & beta6
dimensionedScalar & NSF1_0
dimensionedScalar & TrefXS
dimensionedScalar & beta7
dimensionedScalar & SP1_0
dimensionedScalar & decBeta3
dimensionedScalar & decBeta1
dimensionedScalar & decBeta2
turbulence