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-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  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 "postProcess.H"
41 
42  #include "setRootCase.H"
43  #include "createTime.H"
44  #include "createMesh.H"
45  #include "createControl.H"
46  #include "createFields.H"
47  #include "initContinuityErrs.H"
48 
49  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
50 
51  Info<< "\nStarting time loop\n" << endl;
52 
53  while (runTime.loop())
54  {
55  Info<< "Time = " << runTime.timeName() << nl << endl;
56 
57  #include "compressibleCourantNo.H"
58 
60 
61  // --- Pressure-velocity PIMPLE corrector loop
62  while (pimple.loop())
63  {
65  (
66  fvm::ddt(rho, U)
67  + fvm::div(phi, U)
68  - fvm::laplacian(mu, U)
69  );
70 
71  solve(UEqn == -fvc::grad(p));
72 
73  // --- Pressure corrector loop
74  while (pimple.correct())
75  {
76  volScalarField rAU("rAU", 1.0/UEqn.A());
78  (
79  "rhorAUf",
81  );
82 
83  U = rAU*UEqn.H();
84 
86  (
87  "phid",
88  psi
89  *(
90  fvc::flux(U)
92  )
93  );
94 
95  phi = (rhoO/psi)*phid;
96 
97  fvScalarMatrix pEqn
98  (
99  fvm::ddt(psi, p)
100  + fvc::div(phi)
101  + fvm::div(phid, p)
102  - fvm::laplacian(rhorAUf, p)
103  );
104 
105  pEqn.solve();
106 
107  phi += pEqn.flux();
108 
110  #include "compressibleContinuityErrs.H"
111 
112  U -= rAU*fvc::grad(p);
113  U.correctBoundaryConditions();
114  }
115  }
116 
117  rho = rhoO + psi*p;
118 
119  runTime.write();
120 
121  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
122  << " ClockTime = " << runTime.elapsedClockTime() << " s"
123  << nl << endl;
124  }
125 
126  Info<< "End\n" << endl;
127 
128  return 0;
129 }
130 
131 
132 // ************************************************************************* //
surfaceScalarField & phi
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
U
Definition: pEqn.H:83
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
tmp< GeometricField< typename flux< Type >::type, fvsPatchField, surfaceMesh > > ddtCorr(const GeometricField< Type, fvPatchField, volMesh > &U, const GeometricField< Type, fvsPatchField, surfaceMesh > &Uf)
Definition: fvcDdt.C:170
const dictionary & pimple
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcLaplacian.C:45
surfaceScalarField rhorAUf("rhorAUf", fvc::interpolate(rho *rAU))
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
rhoEqn solve()
static const char nl
Definition: Ostream.H:262
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.
surfaceScalarField phid("phid", fvc::interpolate(psi) *(fvc::flux(HbyA)+rhorAUf *fvc::ddtCorr(rho, U, phi)/fvc::interpolate(rho)))
const dimensionedScalar mu
Atomic mass unit.
fvVectorMatrix & UEqn
Definition: UEqn.H:13
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:45
messageStream Info
const volScalarField & psi
volScalarField & p
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
Definition: fvcFlux.C:32
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Execute application functionObjects to post-process existing results.
volScalarField rAU(1.0/UEqn.A())