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-2016 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 Description
28  Density-based compressible flow solver based on central-upwind schemes of
29  Kurganov and Tadmor.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "fvCFD.H"
34 #include "psiThermo.H"
37 #include "directionInterpolate.H"
38 #include "localEulerDdtScheme.H"
39 #include "fvcSmooth.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 int main(int argc, char *argv[])
44 {
45  #define NO_CONTROL
46  #include "postProcess.H"
47 
48  #include "setRootCase.H"
49  #include "createTime.H"
50  #include "createMesh.H"
51  #include "createFields.H"
52  #include "createFieldRefs.H"
53  #include "createTimeControls.H"
54 
55  turbulence->validate();
56 
57  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
58 
59  #include "readFluxScheme.H"
60 
61  dimensionedScalar v_zero("v_zero", dimVolume/dimTime, 0.0);
62 
63  // Courant numbers used to adjust the time-step
64  scalar CoNum = 0.0;
65  scalar meanCoNum = 0.0;
66 
67  Info<< "\nStarting time loop\n" << endl;
68 
69  while (runTime.run())
70  {
71  // --- Directed interpolation of primitive fields onto faces
72 
75 
76  surfaceVectorField rhoU_pos(interpolate(rhoU, pos, U.name()));
77  surfaceVectorField rhoU_neg(interpolate(rhoU, neg, U.name()));
78 
79  volScalarField rPsi("rPsi", 1.0/psi);
80  surfaceScalarField rPsi_pos(interpolate(rPsi, pos, T.name()));
81  surfaceScalarField rPsi_neg(interpolate(rPsi, neg, T.name()));
82 
83  surfaceScalarField e_pos(interpolate(e, pos, T.name()));
84  surfaceScalarField e_neg(interpolate(e, neg, T.name()));
85 
86  surfaceVectorField U_pos("U_pos", rhoU_pos/rho_pos);
87  surfaceVectorField U_neg("U_neg", rhoU_neg/rho_neg);
88 
89  surfaceScalarField p_pos("p_pos", rho_pos*rPsi_pos);
90  surfaceScalarField p_neg("p_neg", rho_neg*rPsi_neg);
91 
92  surfaceScalarField phiv_pos("phiv_pos", U_pos & mesh.Sf());
93  surfaceScalarField phiv_neg("phiv_neg", U_neg & mesh.Sf());
94 
95  volScalarField c("c", sqrt(thermo.Cp()/thermo.Cv()*rPsi));
96  surfaceScalarField cSf_pos
97  (
98  "cSf_pos",
99  interpolate(c, pos, T.name())*mesh.magSf()
100  );
101  surfaceScalarField cSf_neg
102  (
103  "cSf_neg",
104  interpolate(c, neg, T.name())*mesh.magSf()
105  );
106 
108  (
109  "ap",
110  max(max(phiv_pos + cSf_pos, phiv_neg + cSf_neg), v_zero)
111  );
113  (
114  "am",
115  min(min(phiv_pos - cSf_pos, phiv_neg - cSf_neg), v_zero)
116  );
117 
118  surfaceScalarField a_pos("a_pos", ap/(ap - am));
119 
120  surfaceScalarField amaxSf("amaxSf", max(mag(am), mag(ap)));
121 
122  surfaceScalarField aSf("aSf", am*a_pos);
123 
124  if (fluxScheme == "Tadmor")
125  {
126  aSf = -0.5*amaxSf;
127  a_pos = 0.5;
128  }
129 
130  surfaceScalarField a_neg("a_neg", 1.0 - a_pos);
131 
132  phiv_pos *= a_pos;
133  phiv_neg *= a_neg;
134 
135  surfaceScalarField aphiv_pos("aphiv_pos", phiv_pos - aSf);
136  surfaceScalarField aphiv_neg("aphiv_neg", phiv_neg + aSf);
137 
138  // Reuse amaxSf for the maximum positive and negative fluxes
139  // estimated by the central scheme
140  amaxSf = max(mag(aphiv_pos), mag(aphiv_neg));
141 
142  #include "centralCourantNo.H"
143  #include "readTimeControls.H"
144 
145  if (LTS)
146  {
147  #include "setRDeltaT.H"
148  }
149  else
150  {
151  #include "setDeltaT.H"
152  }
153 
154  runTime++;
155 
156  Info<< "Time = " << runTime.timeName() << nl << endl;
157 
158  phi = aphiv_pos*rho_pos + aphiv_neg*rho_neg;
159 
160  surfaceVectorField phiUp
161  (
162  (aphiv_pos*rhoU_pos + aphiv_neg*rhoU_neg)
163  + (a_pos*p_pos + a_neg*p_neg)*mesh.Sf()
164  );
165 
166  surfaceScalarField phiEp
167  (
168  "phiEp",
169  aphiv_pos*(rho_pos*(e_pos + 0.5*magSqr(U_pos)) + p_pos)
170  + aphiv_neg*(rho_neg*(e_neg + 0.5*magSqr(U_neg)) + p_neg)
171  + aSf*p_pos - aSf*p_neg
172  );
173 
174  volScalarField muEff("muEff", turbulence->muEff());
175  volTensorField tauMC("tauMC", muEff*dev2(Foam::T(fvc::grad(U))));
176 
177  // --- Solve density
179 
180  // --- Solve momentum
181  solve(fvm::ddt(rhoU) + fvc::div(phiUp));
182 
183  U.ref() =
184  rhoU()
185  /rho();
186  U.correctBoundaryConditions();
187  rhoU.boundaryFieldRef() == rho.boundaryField()*U.boundaryField();
188 
189  if (!inviscid)
190  {
191  solve
192  (
193  fvm::ddt(rho, U) - fvc::ddt(rho, U)
194  - fvm::laplacian(muEff, U)
195  - fvc::div(tauMC)
196  );
197  rhoU = rho*U;
198  }
199 
200  // --- Solve energy
201  surfaceScalarField sigmaDotU
202  (
203  "sigmaDotU",
204  (
205  fvc::interpolate(muEff)*mesh.magSf()*fvc::snGrad(U)
206  + fvc::dotInterpolate(mesh.Sf(), tauMC)
207  )
208  & (a_pos*U_pos + a_neg*U_neg)
209  );
210 
211  solve
212  (
213  fvm::ddt(rhoE)
214  + fvc::div(phiEp)
215  - fvc::div(sigmaDotU)
216  );
217 
218  e = rhoE/rho - 0.5*magSqr(U);
219  e.correctBoundaryConditions();
220  thermo.correct();
221  rhoE.boundaryFieldRef() ==
222  rho.boundaryField()*
223  (
224  e.boundaryField() + 0.5*magSqr(U.boundaryField())
225  );
226 
227  if (!inviscid)
228  {
229  solve
230  (
231  fvm::ddt(rho, e) - fvc::ddt(rho, e)
232  - fvm::laplacian(turbulence->alphaEff(), e)
233  );
234  thermo.correct();
235  rhoE = rho*(e + 0.5*magSqr(U));
236  }
237 
238  p.ref() =
239  rho()
240  /psi();
241  p.correctBoundaryConditions();
242  rho.boundaryFieldRef() == psi.boundaryField()*p.boundaryField();
243 
244  turbulence->correct();
245 
246  runTime.write();
247 
248  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
249  << " ClockTime = " << runTime.elapsedClockTime() << " s"
250  << nl << endl;
251  }
252 
253  Info<< "End\n" << endl;
254 
255  return 0;
256 }
257 
258 // ************************************************************************* //
surfaceScalarField & phi
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
U
Definition: pEqn.H:83
volScalarField & rhoE
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:59
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
static tmp< GeometricField< typename innerProduct< vector, Type >::type, fvsPatchField, surfaceMesh > > dotInterpolate(const surfaceVectorField &Sf, const GeometricField< Type, fvPatchField, volMesh > &tvf)
Interpolate field onto faces.
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Read the control parameters used by setDeltaT.
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
dimensionedScalar neg(const dimensionedScalar &ds)
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
psiReactionThermo & thermo
Definition: createFields.H:31
dimensionedScalar pos(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:78
dynamicFvMesh & mesh
rhoEqn solve()
dimensioned< scalar > magSqr(const dimensioned< Type > &)
static const char nl
Definition: Ostream.H:262
bool LTS
Definition: createRDeltaT.H:1
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const volScalarField & T
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
bool inviscid(true)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionedScalar c
Speed of light in a vacuum.
Calculates the mean and maximum wave speed based Courant Numbers.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
const volScalarField & psi
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Execute application functionObjects to post-process existing results.
word fluxScheme("Kurganov")
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
Read the control parameters used by setDeltaT.
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45