-
Notifications
You must be signed in to change notification settings - Fork 0
/
rhoEnergyFoam.bak
263 lines (241 loc) · 8.69 KB
/
rhoEnergyFoam.bak
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
//// ========= |
// \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
// \\ / O peration |
// \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
// \\/ M anipulation |
//-------------------------------------------------------------------------------
//License
// This file is part of OpenFOAM.
//
// OpenFOAM is free software: you can redistribute it and/or modify it
// under the terms of the GNU General Public License as published by
// the Free Software Foundation, either version 3 of the License, or
// (at your option) any later version.
//
// OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
// FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
// for more details.
//
// You should have received a copy of the GNU General Public License
// along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
//
//Application
// rhoEnergyFoam
//
//Description
// Numerical solver for the solution of compressible shock-free flows.
// The convective terms are discretized using Pirozzoli's scheme(JCP 2010),
// density based splitting. The scheme allows conservation of the kinetic
// energy in the inviscid incompressible limit.
// midPoint interpolation must be selected in fvScheme dictionary.
// Viscous terms(Laplacians only) are evaluated directly, computing
// the face normal gradients.
// A third-order low-storage RK scheme is used for time integration.
// The OpenFOAM labrary for turbulence models is included.
//
// Author: Davide Modesti (davide.modesti@uniroma1.it)
// Last update 06/04/2017
// Reference
// D. Modesti and S. Pirozzoli A high-fidelity solver for turbulent compressible flows on unstructured meshes. Comput. & Fluids (2017)
//\*---------------------------------------------------------------------------*/
#include "fvCFD.H"
#include "psiThermo.H"
#include "turbulentFluidThermoModel.H"
#include "fixedRhoFvPatchScalarField.H"
//#include "directionInterpolate.H"
#include "localEulerDdtScheme.H"
#include "fvcSmooth.H"
//#include "Interpolation_IBM.C"
#include <fstream> // std::ofstream
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// AUSM useful functions
float m1 (float rm,float sgn ) // Return M_1
{
float r ;
r = 0.5*(rm + sgn*mag(rm)) ;
return r ;
}
///
float m2 (float rm,float sgn ) // Return M_2
{
float r ;
r = sgn*0.25*(rm + sgn)*(rm + sgn) ;
return r ;
}
///
float p5 (float rm,float sgn ,float alpha) // Return p_5
{
float r ;
if (abs(rm)<1.)
{
r = m2(rm,sgn)*( (sgn*2-rm) - sgn*16*alpha*rm*m2(rm,-sgn) ) ;
}
else
{
r = 1./rm*m1(rm,sgn) ;
}
return r ;
}
///
float m4 (float rm,float sgn ,float beta) // Return M_4
{
float r ;
if (abs(rm)<1)
{
r = m2(rm,sgn)*( 1 - sgn*16*beta*m2(rm,-sgn) ) ;
}
else
{
r = m1(rm,sgn) ;
}
return r ;
}
// Main
int main(int argc, char *argv[])
{
#define NO_CONTROL
#include "postProcess.H"
#include "addCheckCaseOptions.H"
#include "setRootCaseLists.H"
#include "createTime.H"
#include "createMesh.H"
#include "createFields.H"
#include "createFieldRefs.H"
#include "createTimeControls.H"
#include "readThermophysicalProperties.H"
#include "variables.H"
// #include "Set_Immersed_Boundary.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
turbulence->validate();
Info<< nl << "Starting time loop" << endl;
Info<< "Start Timing = " << runTime.clockTimeIncrement() << " s"
<< nl << endl;
// while (runTime.loop()) //Start time loop
while (runTime.run()) //Start time loop
{
Info<< "Time = " << runTime.timeName() << nl << endl;
runTime++;
// Saving quantities at preavious time step
rhoOld = rho;
rhoUOld = rhoU;
rhoEOld = rhoE;
// RK Time step
for (int cycle =0; cycle < rkCoeff.size(); cycle++)
{
// Speed of sound and Mach number
c = Foam::sqrt(thermo.Cp()/thermo.Cv()/psi);
Mach = U/c ;
// Interpolated quantities at cell faces
surfaceScalarField rhoave = fvc:: interpolate(rho) ;
surfaceVectorField Uave = fvc:: interpolate(U) ;
// Flux at the intercell
phi = fvc:: interpolate(U) & mesh.Sf() ;
phit = fvc:: interpolate(rhoU) & mesh.Sf() ; //flux in turbulent model
// Enthalpy
H = (rhoE + p)/rho ;
// Enthalpy at the intercell
surfaceScalarField Have = fvc:: interpolate(H) ;
// Pressure at the intercell
surfaceScalarField pave = fvc::interpolate(p) ;
#include "sensor.H" //Ducros sensor
if(pressArtDiff)
{
#include "AUSM.H" // AUSM+up dissipation on pressure term, add dissipation on pave
}
//
// Evaluate viscous terms
//
// Divergence of the velocity
vecDivU.component(0) = fvc::interpolate(fvc::div(U));
vecDivU.component(1) = fvc::interpolate(fvc::div(U));
vecDivU.component(2) = fvc::interpolate(fvc::div(U));
//
//
volScalarField muEff(turbulence->muEff());
// volScalarField mu = thermo.mu() ;
// volScalarField mut = muEff-mu ;
// muEff = mu;
volTensorField tauMC("tauMC", muEff*dev2(Foam::T(fvc::grad(U))));
// volTensorField tauMC("tauMC", mut*(fvc::grad(U) + Foam::T(fvc::grad(U))));
surfaceScalarField muave = fvc::interpolate(muEff);//mu at cell faces
//
volScalarField k("k", thermo.Cp()*muEff/Pr);//thermal diffusivity
// volScalarField k("k", thermo.Cp()*mut/Pr);//thermal diffusivity
//
surfaceScalarField kave=fvc::interpolate(k);//k at cell faces. alphaEff=muEff/Prt
//momentum viscous flux
surfaceVectorField momVisFlux = muave*(fvc::snGrad(U)*mesh.magSf());//- 2./3.*vecDivU*mesh.magSf()) ;
//energy viscous flux
surfaceScalarField heatFlux = kave*fvc::snGrad(T)*mesh.magSf();
surfaceScalarField visWork = (momVisFlux + fvc::dotInterpolate(mesh.Sf(), tauMC) ) & Uave;
// surfaceScalarField visWork = (momVisFlux ) & Uave;
enVisFlux = heatFlux + visWork ;
//
// Total fluxes, Eulerian + viscous
surfaceScalarField rhoFlux = rhoave*phi ;
momFlux = rhoave*Uave*phi + pave*mesh.Sf() - momVisFlux ;
enFlux = rhoave*Have*phi - enVisFlux ;
//
if(convArtDiff)
{
#include "AUSM_conv.H" // AUSM+up dissipation on convective terms, add dissipation on rhoFlux,momFlux,enFlux
}
//
//
volScalarField rhoFl = fvc::div(rhoFlux);
volVectorField momFl = fvc::div(momFlux) - fvc::div(tauMC);
volScalarField enFl = fvc::div(enFlux) ;
#include "pressureGrad.H" // forcing term to keep a constant mass flow rate in channel flow
// RK sub-step
rho = rhoOld + rkCoeff[cycle]*runTime.deltaT()*(
- rhoFl);
//
rhoU = rhoUOld + rkCoeff[cycle]*runTime.deltaT()*(
// - momFl );
- momFl + dpdx*ones_vec);
//
volScalarField work = U & ones_vec ;
rhoE = rhoEOld + rkCoeff[cycle]*runTime.deltaT()*(
// -enFl );
-enFl + dpdx*work);
//Update primitive variables and boundary conditions
U.ref() = rhoU() / rho();
U.correctBoundaryConditions();
/*Edited*/
rhoU.boundaryFieldRef() = rho.boundaryField()*U.boundaryField();
e = rhoE/rho - 0.5*magSqr(U);
e.correctBoundaryConditions();
//Thermodinamic library
thermo.correct();
/*Edited*/
rhoE.boundaryFieldRef() =
rho.boundaryField()*
(
e.boundaryField() + 0.5*magSqr(U.boundaryField())
);
//
p.ref() = rho() / psi();
p.correctBoundaryConditions();
/*Edited*/
rho.boundaryFieldRef() = psi.boundaryField()*p.boundaryField(); //psi=1/(R*T)
if (Tbulk_target>0.) {
#include "tbforce.H"
}
}//end of RK time integration
runTime.write();
turbulence->correct(); //turbulence model
//
#include "diagnostics.H" //print tke on diagnostics.dat
#include "step.H" //Evaluate Courant Number
#include "setDeltaT.H" //Andjust time step
//
//
}
runTime.write();
Info<< "Start Timing = " << runTime.clockTimeIncrement() << " s"
<< nl << endl;
Info<< "End\n" << endl;
return 0;
}
// ************************************************************************* //