rhoCentralDyMFoam.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  rhoCentralDyMFoam
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 "dynamicFvMesh.H"
35 #include "psiThermo.H"
38 #include "directionInterpolate.H"
39 #include "motionSolver.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 "createDynamicFvMesh.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  #include "readTimeControls.H"
72  #include "setDeltaT.H"
73 
74  runTime++;
75 
76  Info<< "Time = " << runTime.timeName() << nl << endl;
77 
78  // Do any mesh changes
79  mesh.update();
80 
81  // --- Directed interpolation of primitive fields onto faces
82 
85 
86  surfaceVectorField rhoU_pos(interpolate(rhoU, pos, U.name()));
87  surfaceVectorField rhoU_neg(interpolate(rhoU, neg, U.name()));
88 
89  volScalarField rPsi("rPsi", 1.0/psi);
90  surfaceScalarField rPsi_pos(interpolate(rPsi, pos, T.name()));
91  surfaceScalarField rPsi_neg(interpolate(rPsi, neg, T.name()));
92 
93  surfaceScalarField e_pos(interpolate(e, pos, T.name()));
94  surfaceScalarField e_neg(interpolate(e, neg, T.name()));
95 
96  surfaceVectorField U_pos("U_pos", rhoU_pos/rho_pos);
97  surfaceVectorField U_neg("U_neg", rhoU_neg/rho_neg);
98 
99  surfaceScalarField p_pos("p_pos", rho_pos*rPsi_pos);
100  surfaceScalarField p_neg("p_neg", rho_neg*rPsi_neg);
101 
102  surfaceScalarField phiv_pos("phiv_pos", U_pos & mesh.Sf());
103  surfaceScalarField phiv_neg("phiv_neg", U_neg & mesh.Sf());
104 
105  // Make fluxes relative to mesh-motion
106  if (mesh.moving())
107  {
108  phiv_pos -= mesh.phi();
109  phiv_neg -= mesh.phi();
110  }
111 
112  volScalarField c("c", sqrt(thermo.Cp()/thermo.Cv()*rPsi));
113  surfaceScalarField cSf_pos
114  (
115  "cSf_pos",
116  interpolate(c, pos, T.name())*mesh.magSf()
117  );
118  surfaceScalarField cSf_neg
119  (
120  "cSf_neg",
121  interpolate(c, neg, T.name())*mesh.magSf()
122  );
123 
125  (
126  "ap",
127  max(max(phiv_pos + cSf_pos, phiv_neg + cSf_neg), v_zero)
128  );
130  (
131  "am",
132  min(min(phiv_pos - cSf_pos, phiv_neg - cSf_neg), v_zero)
133  );
134 
135  surfaceScalarField a_pos("a_pos", ap/(ap - am));
136 
137  surfaceScalarField amaxSf("amaxSf", max(mag(am), mag(ap)));
138 
139  surfaceScalarField aSf("aSf", am*a_pos);
140 
141  if (fluxScheme == "Tadmor")
142  {
143  aSf = -0.5*amaxSf;
144  a_pos = 0.5;
145  }
146 
147  surfaceScalarField a_neg("a_neg", 1.0 - a_pos);
148 
149  phiv_pos *= a_pos;
150  phiv_neg *= a_neg;
151 
152  surfaceScalarField aphiv_pos("aphiv_pos", phiv_pos - aSf);
153  surfaceScalarField aphiv_neg("aphiv_neg", phiv_neg + aSf);
154 
155  // Reuse amaxSf for the maximum positive and negative fluxes
156  // estimated by the central scheme
157  amaxSf = max(mag(aphiv_pos), mag(aphiv_neg));
158 
159  #include "centralCourantNo.H"
160 
161  phi = aphiv_pos*rho_pos + aphiv_neg*rho_neg;
162 
163  surfaceVectorField phiUp
164  (
165  (aphiv_pos*rhoU_pos + aphiv_neg*rhoU_neg)
166  + (a_pos*p_pos + a_neg*p_neg)*mesh.Sf()
167  );
168 
169  surfaceScalarField phiEp
170  (
171  "phiEp",
172  aphiv_pos*(rho_pos*(e_pos + 0.5*magSqr(U_pos)) + p_pos)
173  + aphiv_neg*(rho_neg*(e_neg + 0.5*magSqr(U_neg)) + p_neg)
174  + aSf*p_pos - aSf*p_neg
175  );
176 
177  // Make flux for pressure-work absolute
178  if (mesh.moving())
179  {
180  phiEp += mesh.phi()*(a_pos*p_pos + a_neg*p_neg);
181  }
182 
183  volScalarField muEff("muEff", turbulence->muEff());
184  volTensorField tauMC("tauMC", muEff*dev2(Foam::T(fvc::grad(U))));
185 
186  // --- Solve density
188 
189  // --- Solve momentum
190  solve(fvm::ddt(rhoU) + fvc::div(phiUp));
191 
192  U.ref() =
193  rhoU()
194  /rho();
195  U.correctBoundaryConditions();
196  rhoU.boundaryFieldRef() == rho.boundaryField()*U.boundaryField();
197 
198  if (!inviscid)
199  {
200  solve
201  (
202  fvm::ddt(rho, U) - fvc::ddt(rho, U)
203  - fvm::laplacian(muEff, U)
204  - fvc::div(tauMC)
205  );
206  rhoU = rho*U;
207  }
208 
209  // --- Solve energy
210  surfaceScalarField sigmaDotU
211  (
212  "sigmaDotU",
213  (
214  fvc::interpolate(muEff)*mesh.magSf()*fvc::snGrad(U)
215  + fvc::dotInterpolate(mesh.Sf(), tauMC)
216  )
217  & (a_pos*U_pos + a_neg*U_neg)
218  );
219 
220  solve
221  (
222  fvm::ddt(rhoE)
223  + fvc::div(phiEp)
224  - fvc::div(sigmaDotU)
225  );
226 
227  e = rhoE/rho - 0.5*magSqr(U);
228  e.correctBoundaryConditions();
229  thermo.correct();
230  rhoE.boundaryFieldRef() ==
231  rho.boundaryField()*
232  (
233  e.boundaryField() + 0.5*magSqr(U.boundaryField())
234  );
235 
236  if (!inviscid)
237  {
238  solve
239  (
240  fvm::ddt(rho, e) - fvc::ddt(rho, e)
241  - fvm::laplacian(turbulence->alphaEff(), e)
242  );
243  thermo.correct();
244  rhoE = rho*(e + 0.5*magSqr(U));
245  }
246 
247  p.ref() =
248  rho()
249  /psi();
250  p.correctBoundaryConditions();
251  rho.boundaryFieldRef() == psi.boundaryField()*p.boundaryField();
252 
253  turbulence->correct();
254 
255  runTime.write();
256 
257  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
258  << " ClockTime = " << runTime.elapsedClockTime() << " s"
259  << nl << endl;
260  }
261 
262  Info<< "End\n" << endl;
263 
264  return 0;
265 }
266 
267 // ************************************************************************* //
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
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< surfaceScalarField > interpolate(const RhoType &rho)
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
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const volScalarField & T
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")
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