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 | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Application
25  rhoCentralFoam
26 
27 Group
28  grpCompressibleSolvers
29 
30 Description
31  Density-based compressible flow solver based on central-upwind schemes of
32  Kurganov and Tadmor
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #include "fvCFD.H"
37 #include "psiThermo.H"
41 #include "directionInterpolate.H"
42 #include "localEulerDdtScheme.H"
43 #include "fvcSmooth.H"
44 
45 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
46 
47 int main(int argc, char *argv[])
48 {
49  #include "setRootCase.H"
50 
51  #include "createTime.H"
52  #include "createMesh.H"
53  #include "createFields.H"
54  #include "createTimeControls.H"
55  #include "createRDeltaT.H"
56 
57  turbulence->validate();
58 
59  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61  #include "readFluxScheme.H"
62 
63  dimensionedScalar v_zero("v_zero", dimVolume/dimTime, 0.0);
64 
65  // Courant numbers used to adjust the time-step
66  scalar CoNum = 0.0;
67  scalar meanCoNum = 0.0;
68 
69  Info<< "\nStarting time loop\n" << endl;
70 
71  while (runTime.run())
72  {
73  // --- Directed interpolation of primitive fields onto faces
74 
77 
78  surfaceVectorField rhoU_pos(interpolate(rhoU, pos, U.name()));
79  surfaceVectorField rhoU_neg(interpolate(rhoU, neg, U.name()));
80 
81  volScalarField rPsi("rPsi", 1.0/psi);
82  surfaceScalarField rPsi_pos(interpolate(rPsi, pos, T.name()));
83  surfaceScalarField rPsi_neg(interpolate(rPsi, neg, T.name()));
84 
85  surfaceScalarField e_pos(interpolate(e, pos, T.name()));
86  surfaceScalarField e_neg(interpolate(e, neg, T.name()));
87 
88  surfaceVectorField U_pos("U_pos", rhoU_pos/rho_pos);
89  surfaceVectorField U_neg("U_neg", rhoU_neg/rho_neg);
90 
91  surfaceScalarField p_pos("p_pos", rho_pos*rPsi_pos);
92  surfaceScalarField p_neg("p_neg", rho_neg*rPsi_neg);
93 
94  surfaceScalarField phiv_pos("phiv_pos", U_pos & mesh.Sf());
95  surfaceScalarField phiv_neg("phiv_neg", U_neg & mesh.Sf());
96 
97  volScalarField c("c", sqrt(thermo.Cp()/thermo.Cv()*rPsi));
98  surfaceScalarField cSf_pos
99  (
100  "cSf_pos",
101  interpolate(c, pos, T.name())*mesh.magSf()
102  );
103  surfaceScalarField cSf_neg
104  (
105  "cSf_neg",
106  interpolate(c, neg, T.name())*mesh.magSf()
107  );
108 
110  (
111  "ap",
112  max(max(phiv_pos + cSf_pos, phiv_neg + cSf_neg), v_zero)
113  );
115  (
116  "am",
117  min(min(phiv_pos - cSf_pos, phiv_neg - cSf_neg), v_zero)
118  );
119 
120  surfaceScalarField a_pos("a_pos", ap/(ap - am));
121 
122  surfaceScalarField amaxSf("amaxSf", max(mag(am), mag(ap)));
123 
124  surfaceScalarField aSf("aSf", am*a_pos);
125 
126  if (fluxScheme == "Tadmor")
127  {
128  aSf = -0.5*amaxSf;
129  a_pos = 0.5;
130  }
131 
132  surfaceScalarField a_neg("a_neg", 1.0 - a_pos);
133 
134  phiv_pos *= a_pos;
135  phiv_neg *= a_neg;
136 
137  surfaceScalarField aphiv_pos("aphiv_pos", phiv_pos - aSf);
138  surfaceScalarField aphiv_neg("aphiv_neg", phiv_neg + aSf);
139 
140  // Reuse amaxSf for the maximum positive and negative fluxes
141  // estimated by the central scheme
142  amaxSf = max(mag(aphiv_pos), mag(aphiv_neg));
143 
144  #include "centralCourantNo.H"
145  #include "readTimeControls.H"
146 
147  if (LTS)
148  {
149  #include "setRDeltaT.H"
150  }
151  else
152  {
153  #include "setDeltaT.H"
154  }
155 
156  runTime++;
157 
158  Info<< "Time = " << runTime.timeName() << nl << endl;
159 
160  phi = aphiv_pos*rho_pos + aphiv_neg*rho_neg;
161 
162  surfaceVectorField phiUp
163  (
164  (aphiv_pos*rhoU_pos + aphiv_neg*rhoU_neg)
165  + (a_pos*p_pos + a_neg*p_neg)*mesh.Sf()
166  );
167 
168  surfaceScalarField phiEp
169  (
170  "phiEp",
171  aphiv_pos*(rho_pos*(e_pos + 0.5*magSqr(U_pos)) + p_pos)
172  + aphiv_neg*(rho_neg*(e_neg + 0.5*magSqr(U_neg)) + p_neg)
173  + aSf*p_pos - aSf*p_neg
174  );
175 
176  volScalarField muEff("muEff", turbulence->muEff());
177  volTensorField tauMC("tauMC", muEff*dev2(Foam::T(fvc::grad(U))));
178 
179  // --- Solve density
181 
182  // --- Solve momentum
183  solve(fvm::ddt(rhoU) + fvc::div(phiUp));
184 
185  U.dimensionedInternalField() =
186  rhoU.dimensionedInternalField()
187  /rho.dimensionedInternalField();
188  U.correctBoundaryConditions();
189  rhoU.boundaryField() == rho.boundaryField()*U.boundaryField();
190 
191  if (!inviscid)
192  {
193  solve
194  (
195  fvm::ddt(rho, U) - fvc::ddt(rho, U)
196  - fvm::laplacian(muEff, U)
197  - fvc::div(tauMC)
198  );
199  rhoU = rho*U;
200  }
201 
202  // --- Solve energy
203  surfaceScalarField sigmaDotU
204  (
205  "sigmaDotU",
206  (
207  fvc::interpolate(muEff)*mesh.magSf()*fvc::snGrad(U)
208  + (mesh.Sf() & fvc::interpolate(tauMC))
209  )
210  & (a_pos*U_pos + a_neg*U_neg)
211  );
212 
213  solve
214  (
215  fvm::ddt(rhoE)
216  + fvc::div(phiEp)
217  - fvc::div(sigmaDotU)
218  );
219 
220  e = rhoE/rho - 0.5*magSqr(U);
221  e.correctBoundaryConditions();
222  thermo.correct();
223  rhoE.boundaryField() ==
224  rho.boundaryField()*
225  (
226  e.boundaryField() + 0.5*magSqr(U.boundaryField())
227  );
228 
229  if (!inviscid)
230  {
231  solve
232  (
233  fvm::ddt(rho, e) - fvc::ddt(rho, e)
234  - fvm::laplacian(turbulence->alphaEff(), e)
235  );
236  thermo.correct();
237  rhoE = rho*(e + 0.5*magSqr(U));
238  }
239 
240  p.dimensionedInternalField() =
241  rho.dimensionedInternalField()
242  /psi.dimensionedInternalField();
243  p.correctBoundaryConditions();
244  rho.boundaryField() == psi.boundaryField()*p.boundaryField();
245 
246  turbulence->correct();
247 
248  runTime.write();
249 
250  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
251  << " ClockTime = " << runTime.elapsedClockTime() << " s"
252  << nl << endl;
253  }
254 
255  Info<< "End\n" << endl;
256 
257  return 0;
258 }
259 
260 // ************************************************************************* //
Foam::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Foam::volTensorField
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:59
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
p
p
Definition: pEqn.H:62
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::MULES::interpolate
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Definition: IMULESTemplates.C:40
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
createRDeltaT.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::dev2
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:115
U
U
Definition: pEqn.H:46
inviscid
bool inviscid(true)
centralCourantNo.H
Calculates the mean and maximum wave speed based Courant Numbers.
directionInterpolate.H
Foam::fvc::laplacian
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
solve
rhoEqn solve()
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
localEulerDdtScheme.H
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
rhoE
volScalarField & rhoE
Definition: readMechanicalProperties.H:125
CoNum
scalar CoNum
Definition: compressibleCourantNo.H:29
fluxScheme
word fluxScheme("Kurganov")
fixedRhoFvPatchScalarField.H
LTS
bool LTS
Definition: createRDeltaT.H:1
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
rho
rho
Definition: pEqn.H:3
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
setRootCase.H
psiThermo.H
T
const volScalarField & T
Definition: createFields.H:25
psi
const volScalarField & psi
Definition: setRegionFluidFields.H:13
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
createTimeControls.H
Read the control parameters used by setDeltaT.
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
readTimeControls.H
Read the control parameters used by setDeltaT.
createMesh.H
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
meanCoNum
scalar meanCoNum
Definition: compressibleCourantNo.H:30
createTime.H
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
fvCFD.H
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
readFluxScheme.H
Foam::surfaceVectorField
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
Definition: surfaceFieldsFwd.H:55
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
zeroGradientFvPatchFields.H
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:201
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
turbulentFluidThermoModel.H
fvcSmooth.H
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190