sonicLiquidFoam.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  sonicLiquidFoam
26 
27 Description
28  Transient solver for trans-sonic/supersonic, laminar flow of a
29  compressible liquid.
30 
31 \*---------------------------------------------------------------------------*/
32 
33 #include "fvCFD.H"
34 #include "pimpleControl.H"
35 
36 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37 
38 int main(int argc, char *argv[])
39 {
40  #include "setRootCase.H"
41  #include "createTime.H"
42  #include "createMesh.H"
43 
44  pimpleControl pimple(mesh);
45 
46  #include "readThermodynamicProperties.H"
47  #include "readTransportProperties.H"
48  #include "createFields.H"
49  #include "initContinuityErrs.H"
50 
51  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53  Info<< "\nStarting time loop\n" << endl;
54 
55  while (runTime.loop())
56  {
57  Info<< "Time = " << runTime.timeName() << nl << endl;
58 
59  #include "compressibleCourantNo.H"
60 
62 
63  // --- Pressure-velocity PIMPLE corrector loop
64  while (pimple.loop())
65  {
67  (
68  fvm::ddt(rho, U)
69  + fvm::div(phi, U)
70  - fvm::laplacian(mu, U)
71  );
72 
73  solve(UEqn == -fvc::grad(p));
74 
75  // --- Pressure corrector loop
76  while (pimple.correct())
77  {
78  volScalarField rAU("rAU", 1.0/UEqn.A());
80  (
81  "rhorAUf",
83  );
84 
85  U = rAU*UEqn.H();
86 
88  (
89  "phid",
90  psi
91  *(
92  (fvc::interpolate(U) & mesh.Sf())
94  )
95  );
96 
97  phi = (rhoO/psi)*phid;
98 
99  fvScalarMatrix pEqn
100  (
101  fvm::ddt(psi, p)
102  + fvc::div(phi)
103  + fvm::div(phid, p)
104  - fvm::laplacian(rhorAUf, p)
105  );
106 
107  pEqn.solve();
108 
109  phi += pEqn.flux();
110 
112  #include "compressibleContinuityErrs.H"
113 
114  U -= rAU*fvc::grad(p);
115  U.correctBoundaryConditions();
116  }
117  }
118 
119  rho = rhoO + psi*p;
120 
121  runTime.write();
122 
123  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
124  << " ClockTime = " << runTime.elapsedClockTime() << " s"
125  << nl << endl;
126  }
127 
128  Info<< "End\n" << endl;
129 
130  return 0;
131 }
132 
133 
134 // ************************************************************************* //
const dimensionedScalar mu
Atomic mass unit.
Definition: createFields.H:13
rhoEqn solve()
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:155
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
messageStream Info
tmp< surfaceScalarField > interpolate(const RhoType &rho)
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:45
dynamicFvMesh & mesh
const dictionary & pimple
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
fvVectorMatrix UEqn(fvm::ddt(rho, U)+fvm::div(phi, U)+turbulence->divDevRhoReff(U)==fvOptions(rho, U))
surfaceScalarField phid("phid", fvc::interpolate(psi)*( (mesh.Sf()&fvc::interpolate(HbyA)) +rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho) ))
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
volScalarField & p
Definition: createFields.H:51
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
const volScalarField & psi
Definition: createFields.H:24
surfaceScalarField & phi
volScalarField rAU(1.0/UEqn.A())
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
int main(int argc, char *argv[])
Definition: postCalc.C:54
U
Definition: pEqn.H:82