Loading...
Searching...
No Matches
rhoCentralFoam.C
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) 2021-2025 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
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
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27Application
28 rhoCentralFoam
29
30Group
31 grpCompressibleSolvers
32
33Description
34 Density-based compressible flow solver based on
35 central-upwind schemes of Kurganov and Tadmor with
36 support for mesh-motion and topology changes.
37
38\*---------------------------------------------------------------------------*/
39
40#include "fvCFD.H"
41#include "dynamicFvMesh.H"
42#include "psiThermo.H"
46#include "localEulerDdtScheme.H"
47#include "fvcSmooth.H"
48
49// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50
51int main(int argc, char *argv[])
52{
53 argList::addNote
54 (
55 "Density-based compressible flow solver based on"
56 " central-upwind schemes of Kurganov and Tadmor with"
57 " support for mesh-motion and topology changes."
58 );
59
60 #define NO_CONTROL
61 #include "postProcess.H"
62
63 #include "addCheckCaseOptions.H"
64 #include "setRootCaseLists.H"
65 #include "createTime.H"
66 #include "createDynamicFvMesh.H"
67 #include "createFields.H"
68 #include "createFieldRefs.H"
69 #include "createTimeControls.H"
70
71 turbulence->validate();
72
73 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
74
75 #include "readFluxScheme.H"
76
77 const dimensionedScalar v_zero(dimVolume/dimTime, Zero);
78
79 // Courant numbers used to adjust the time-step
80 scalar CoNum = 0.0;
81 scalar meanCoNum = 0.0;
82 scalar minCoNum = 0.0;
83
84 Info<< "\nStarting time loop\n" << endl;
85
86 while (runTime.run())
87 {
88 #include "readTimeControls.H"
89
90 if (!LTS)
91 {
92 #include "setDeltaT.H"
93
94 ++runTime;
95
96 // Do any mesh changes
97 mesh.update();
98 }
99
100 // --- Directed interpolation of primitive fields onto faces
101
102 surfaceScalarField rho_pos(interpolate(rho, pos));
103 surfaceScalarField rho_neg(interpolate(rho, neg));
104
105 surfaceVectorField rhoU_pos(interpolate(rhoU, pos, U.name()));
106 surfaceVectorField rhoU_neg(interpolate(rhoU, neg, U.name()));
107
108 volScalarField rPsi("rPsi", 1.0/psi);
109 surfaceScalarField rPsi_pos(interpolate(rPsi, pos, T.name()));
110 surfaceScalarField rPsi_neg(interpolate(rPsi, neg, T.name()));
111
112 surfaceScalarField e_pos(interpolate(e, pos, T.name()));
113 surfaceScalarField e_neg(interpolate(e, neg, T.name()));
114
115 surfaceVectorField U_pos("U_pos", rhoU_pos/rho_pos);
116 surfaceVectorField U_neg("U_neg", rhoU_neg/rho_neg);
117
118 surfaceScalarField p_pos("p_pos", rho_pos*rPsi_pos);
119 surfaceScalarField p_neg("p_neg", rho_neg*rPsi_neg);
120
121 surfaceScalarField phiv_pos("phiv_pos", U_pos & mesh.Sf());
122 // Note: extracted out the orientation so becomes unoriented
123 phiv_pos.setOriented(false);
124 surfaceScalarField phiv_neg("phiv_neg", U_neg & mesh.Sf());
125 phiv_neg.setOriented(false);
126
127 // Make fluxes relative to mesh-motion
128 if (mesh.moving())
129 {
131 meshPhi.setOriented(false);
132 phiv_pos -= meshPhi;
133 phiv_neg -= meshPhi;
134 }
135
136 volScalarField c("c", sqrt(thermo.Cp()/thermo.Cv()*rPsi));
137 surfaceScalarField cSf_pos
138 (
139 "cSf_pos",
140 interpolate(c, pos, T.name())*mesh.magSf()
141 );
142
143 surfaceScalarField cSf_neg
144 (
145 "cSf_neg",
146 interpolate(c, neg, T.name())*mesh.magSf()
147 );
148
150 (
151 "ap",
152 max(max(phiv_pos + cSf_pos, phiv_neg + cSf_neg), v_zero)
153 );
154
156 (
157 "am",
158 min(min(phiv_pos - cSf_pos, phiv_neg - cSf_neg), v_zero)
159 );
160
161 surfaceScalarField a_pos("a_pos", ap/(ap - am));
162
163 surfaceScalarField amaxSf("amaxSf", max(mag(am), mag(ap)));
164
165 surfaceScalarField aSf("aSf", am*a_pos);
166
167 if (fluxScheme == "Tadmor")
168 {
169 aSf = -0.5*amaxSf;
170 a_pos = 0.5;
171 }
172
173 surfaceScalarField a_neg("a_neg", 1.0 - a_pos);
174
175 phiv_pos *= a_pos;
176 phiv_neg *= a_neg;
177
178 surfaceScalarField aphiv_pos("aphiv_pos", phiv_pos - aSf);
179 surfaceScalarField aphiv_neg("aphiv_neg", phiv_neg + aSf);
180
181 // Reuse amaxSf for the maximum positive and negative fluxes
182 // estimated by the central scheme
183 amaxSf = max(mag(aphiv_pos), mag(aphiv_neg));
184
185 if (LTS)
186 {
187 #include "setRDeltaT.H"
188
189 ++runTime;
190 }
191
192 #include "centralCourantNo.H"
193
194 Info<< "Time = " << runTime.timeName() << nl << endl;
195
196 phi = aphiv_pos*rho_pos + aphiv_neg*rho_neg;
197
198 surfaceVectorField phiU(aphiv_pos*rhoU_pos + aphiv_neg*rhoU_neg);
199 // Note: reassembled orientation from the pos and neg parts so becomes
200 // oriented
201 phiU.setOriented(true);
202
203 surfaceVectorField phiUp(phiU + (a_pos*p_pos + a_neg*p_neg)*mesh.Sf());
204
206 (
207 "phiEp",
208 aphiv_pos*(rho_pos*(e_pos + 0.5*magSqr(U_pos)) + p_pos)
209 + aphiv_neg*(rho_neg*(e_neg + 0.5*magSqr(U_neg)) + p_neg)
210 + aSf*p_pos - aSf*p_neg
211 );
212
213 // Make flux for pressure-work absolute
214 if (mesh.moving())
215 {
217 meshPhi.setOriented(false);
218 phiEp += meshPhi*(a_pos*p_pos + a_neg*p_neg);
219 }
220
221 volScalarField muEff("muEff", turbulence->muEff());
222 volTensorField tauMC("tauMC", muEff*dev2(Foam::T(fvc::grad(U))));
223
224 // --- Solve density
225 solve(fvm::ddt(rho) + fvc::div(phi));
226
227 // --- Solve momentum
228 solve(fvm::ddt(rhoU) + fvc::div(phiUp));
229
230 U.ref() =
231 rhoU()
232 /rho();
233 U.correctBoundaryConditions();
234 rhoU.boundaryFieldRef() == rho.boundaryField()*U.boundaryField();
235
236 if (!inviscid)
237 {
238 solve
239 (
240 fvm::ddt(rho, U) - fvc::ddt(rho, U)
241 - fvm::laplacian(muEff, U)
242 - fvc::div(tauMC)
243 );
244 rhoU = rho*U;
245 }
246
247 // --- Solve energy
248 surfaceScalarField sigmaDotU
249 (
250 "sigmaDotU",
251 (
252 fvc::interpolate(muEff)*mesh.magSf()*fvc::snGrad(U)
253 + fvc::dotInterpolate(mesh.Sf(), tauMC)
254 )
255 & (a_pos*U_pos + a_neg*U_neg)
256 );
257
258 solve
259 (
260 fvm::ddt(rhoE)
261 + fvc::div(phiEp)
262 - fvc::div(sigmaDotU)
263 );
264
265 e = rhoE/rho - 0.5*magSqr(U);
266 e.correctBoundaryConditions();
267 thermo.correct();
268 rhoE.boundaryFieldRef() ==
269 rho.boundaryField()*
270 (
271 e.boundaryField() + 0.5*magSqr(U.boundaryField())
272 );
273
274 if (!inviscid)
275 {
276 solve
277 (
278 fvm::ddt(rho, e) - fvc::ddt(rho, e)
279 - fvm::laplacian(turbulence->alphaEff(), e)
280 );
281 thermo.correct();
282 rhoE = rho*(e + 0.5*magSqr(U));
283 }
284
285 p.ref() =
286 rho()
287 /psi();
288 p.correctBoundaryConditions();
289 rho.boundaryFieldRef() == psi.boundaryField()*p.boundaryField();
290
291 turbulence->correct();
292
293 runTime.write();
294
295 runTime.printExecutionTime(Info);
296 }
297
298 Info<< "End\n" << endl;
299
300 return 0;
301}
302
303// ************************************************************************* //
Required Classes.
Calculates the mean and maximum wave speed based Courant Numbers.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
U
Definition pEqn.H:72
volScalarField & p
const volScalarField & psi
const volScalarField & T
bool inviscid(true)
mesh interpolate(rAU)
dynamicFvMesh & mesh
engineTime & runTime
bool LTS
Read the control parameters used by setDeltaT.
compressible::turbulenceModel & turbulence
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
const dimensionedScalar c
Speed of light in a vacuum.
tmp< surfaceScalarField > meshPhi(const volVectorField &U)
Definition fvcMeshPhi.C:29
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
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
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Execute application functionObjects to post-process existing results.
word fluxScheme("Kurganov")
volScalarField & rhoE
Read the control parameters used by setDeltaT.
CEqn solve()
volScalarField & e
scalar CoNum
scalar meanCoNum