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