rhoCentralFoam.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2022 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 with support for mesh-motion and topology changes.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "fvCFD.H"
34 #include "psiThermo.H"
38 #include "directionInterpolate.H"
39 #include "localEulerDdtScheme.H"
40 #include "fvcSmooth.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 int main(int argc, char *argv[])
45 {
46  #define NO_CONTROL
47  #include "postProcess.H"
48 
49  #include "setRootCaseLists.H"
50  #include "createTime.H"
51  #include "createMesh.H"
52  #include "createFields.H"
53  #include "createFieldRefs.H"
54  #include "createTimeControls.H"
55 
56  turbulence->validate();
57 
58  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60  #include "readFluxScheme.H"
61 
62  dimensionedScalar v_zero("v_zero", dimVolume/dimTime, 0.0);
63 
64  // Courant numbers used to adjust the time-step
65  scalar CoNum = 0.0;
66  scalar meanCoNum = 0.0;
67 
68  Info<< "\nStarting time loop\n" << endl;
69 
70  while (runTime.run())
71  {
72  #include "readTimeControls.H"
73 
74  if (!LTS)
75  {
76  #include "setDeltaT.H"
77 
78  // Update the mesh for topology change, mesh to mesh mapping
79  mesh.update();
80 
81  runTime++;
82 
83  // Move the mesh
84  mesh.move();
85  }
86 
87  // --- Directed interpolation of primitive fields onto faces
88 
91 
92  surfaceVectorField rhoU_pos(interpolate(rhoU, pos, U.name()));
93  surfaceVectorField rhoU_neg(interpolate(rhoU, neg, U.name()));
94 
95  volScalarField rPsi("rPsi", 1.0/psi);
96  surfaceScalarField rPsi_pos(interpolate(rPsi, pos, T.name()));
97  surfaceScalarField rPsi_neg(interpolate(rPsi, neg, T.name()));
98 
99  surfaceScalarField e_pos(interpolate(e, pos, T.name()));
100  surfaceScalarField e_neg(interpolate(e, neg, T.name()));
101 
102  surfaceVectorField U_pos("U_pos", rhoU_pos/rho_pos);
103  surfaceVectorField U_neg("U_neg", rhoU_neg/rho_neg);
104 
105  surfaceScalarField p_pos("p_pos", rho_pos*rPsi_pos);
106  surfaceScalarField p_neg("p_neg", rho_neg*rPsi_neg);
107 
108  surfaceScalarField phiv_pos("phiv_pos", U_pos & mesh.Sf());
109  surfaceScalarField phiv_neg("phiv_neg", U_neg & mesh.Sf());
110 
111  // Make fluxes relative to mesh-motion
112  if (mesh.moving())
113  {
114  phiv_pos -= mesh.phi();
115  phiv_neg -= mesh.phi();
116  }
117 
118  volScalarField c("c", sqrt(thermo.Cp()/thermo.Cv()*rPsi));
119  surfaceScalarField cSf_pos
120  (
121  "cSf_pos",
122  interpolate(c, pos, T.name())*mesh.magSf()
123  );
124  surfaceScalarField cSf_neg
125  (
126  "cSf_neg",
127  interpolate(c, neg, T.name())*mesh.magSf()
128  );
129 
131  (
132  "ap",
133  max(max(phiv_pos + cSf_pos, phiv_neg + cSf_neg), v_zero)
134  );
136  (
137  "am",
138  min(min(phiv_pos - cSf_pos, phiv_neg - cSf_neg), v_zero)
139  );
140 
141  surfaceScalarField a_pos("a_pos", ap/(ap - am));
142 
143  surfaceScalarField amaxSf("amaxSf", max(mag(am), mag(ap)));
144 
145  surfaceScalarField aSf("aSf", am*a_pos);
146 
147  if (fluxScheme == "Tadmor")
148  {
149  aSf = -0.5*amaxSf;
150  a_pos = 0.5;
151  }
152 
153  surfaceScalarField a_neg("a_neg", 1.0 - a_pos);
154 
155  phiv_pos *= a_pos;
156  phiv_neg *= a_neg;
157 
158  surfaceScalarField aphiv_pos("aphiv_pos", phiv_pos - aSf);
159  surfaceScalarField aphiv_neg("aphiv_neg", phiv_neg + aSf);
160 
161  // Reuse amaxSf for the maximum positive and negative fluxes
162  // estimated by the central scheme
163  amaxSf = max(mag(aphiv_pos), mag(aphiv_neg));
164 
165  #include "centralCourantNo.H"
166 
167  if (LTS)
168  {
169  #include "setRDeltaT.H"
170  runTime++;
171  }
172 
173  Info<< "Time = " << runTime.userTimeName() << nl << endl;
174 
175  phi = aphiv_pos*rho_pos + aphiv_neg*rho_neg;
176 
177  surfaceVectorField phiUp
178  (
179  (aphiv_pos*rhoU_pos + aphiv_neg*rhoU_neg)
180  + (a_pos*p_pos + a_neg*p_neg)*mesh.Sf()
181  );
182 
183  surfaceScalarField phiEp
184  (
185  "phiEp",
186  aphiv_pos*(rho_pos*(e_pos + 0.5*magSqr(U_pos)) + p_pos)
187  + aphiv_neg*(rho_neg*(e_neg + 0.5*magSqr(U_neg)) + p_neg)
188  + aSf*p_pos - aSf*p_neg
189  );
190 
191  // Make flux for pressure-work absolute
192  if (mesh.moving())
193  {
194  phiEp += mesh.phi()*(a_pos*p_pos + a_neg*p_neg);
195  }
196 
197  volScalarField muEff("muEff", rho*turbulence->nuEff());
198  volTensorField tauMC("tauMC", muEff*dev2(Foam::T(fvc::grad(U))));
199 
200  // --- Solve density
202 
203  // --- Solve momentum
204  solve(fvm::ddt(rhoU) + fvc::div(phiUp));
205 
206  U.ref() =
207  rhoU()
208  /rho();
209  U.correctBoundaryConditions();
210  rhoU.boundaryFieldRef() == rho.boundaryField()*U.boundaryField();
211 
212  if (!inviscid)
213  {
214  solve
215  (
216  fvm::ddt(rho, U) - fvc::ddt(rho, U)
217  - fvm::laplacian(muEff, U)
218  - fvc::div(tauMC)
219  );
220  rhoU = rho*U;
221  }
222 
223  // --- Solve energy
224  surfaceScalarField sigmaDotU
225  (
226  "sigmaDotU",
227  (
228  fvc::interpolate(muEff)*mesh.magSf()*fvc::snGrad(U)
229  + fvc::dotInterpolate(mesh.Sf(), tauMC)
230  )
231  & (a_pos*U_pos + a_neg*U_neg)
232  );
233 
234  solve
235  (
236  fvm::ddt(rhoE)
237  + fvc::div(phiEp)
238  - fvc::div(sigmaDotU)
239  );
240 
241  e = rhoE/rho - 0.5*magSqr(U);
242  e.correctBoundaryConditions();
243  thermo.correct();
244  rhoE.boundaryFieldRef() ==
245  rho.boundaryField()*
246  (
247  e.boundaryField() + 0.5*magSqr(U.boundaryField())
248  );
249 
250  if (!inviscid)
251  {
252  solve
253  (
254  fvm::ddt(rho, e) - fvc::ddt(rho, e)
255  + thermophysicalTransport->divq(e)
256  );
257  thermo.correct();
258  rhoE = rho*(e + 0.5*magSqr(U));
259  }
260 
261  p.ref() =
262  rho()
263  /psi();
264  p.correctBoundaryConditions();
265  rho.boundaryFieldRef() == psi.boundaryField()*p.boundaryField();
266 
267  turbulence->correct();
268  thermophysicalTransport->correct();
269 
270  runTime.write();
271 
272  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
273  << " ClockTime = " << runTime.elapsedClockTime() << " s"
274  << nl << endl;
275  }
276 
277  Info<< "End\n" << endl;
278 
279  return 0;
280 }
281 
282 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
fluidReactionThermo & thermo
Definition: createFields.H:28
U
Definition: pEqn.H:72
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:63
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:251
Read the control parameters used by setDeltaT.
fvMesh & mesh
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
const dimensionedScalar c
Speed of light in a vacuum.
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:58
dimensionedScalar pos(const dimensionedScalar &ds)
const dimensionSet dimTime
const volScalarField & psi
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
phi
Definition: correctPhi.H:3
dimensioned< scalar > magSqr(const dimensioned< Type > &)
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(dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar(dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::momentumTransportModel > turbulence(compressible::momentumTransportModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
static const char nl
Definition: Ostream.H:260
bool LTS
Definition: createRDeltaT.H:1
fluidReactionThermophysicalTransportModel & thermophysicalTransport
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
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.
const volScalarField & T
rhoEqn solve()
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
bool inviscid(true)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Calculates the mean and maximum wave speed based Courant Numbers.
const dimensionSet dimVolume
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField & p
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Execute application functionObjects to post-process existing results.
word fluxScheme("Kurganov")
const dimensionedScalar e
Elementary charge.
Definition: doubleScalar.H:105
Provides functions smooth spread and sweep which use the FvFaceCellWave algorithm to smooth and redis...
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