FreeFOAM The Cross-Platform CFD Toolkit
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) 1991-2010 OpenCFD Ltd.
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 Usage
32  - sonicLiquidFoam [OPTION]
33 
34  @param -case <dir> \n
35  Specify the case directory
36 
37  @param -parallel \n
38  Run the case in parallel
39 
40  @param -help \n
41  Display short usage message
42 
43  @param -doc \n
44  Display Doxygen documentation page
45 
46  @param -srcDoc \n
47  Display source code
48 
49 \*---------------------------------------------------------------------------*/
50 
51 #include <finiteVolume/fvCFD.H>
52 
53 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 
55 int main(int argc, char *argv[])
56 {
57  #include <OpenFOAM/setRootCase.H>
58  #include <OpenFOAM/createTime.H>
59  #include <OpenFOAM/createMesh.H>
61  #include "readTransportProperties.H"
62  #include "createFields.H"
64 
65  // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66 
67  Info<< "\nStarting time loop\n" << endl;
68 
69  while (runTime.loop())
70  {
71  Info<< "Time = " << runTime.timeName() << nl << endl;
72 
75 
76  #include <finiteVolume/rhoEqn.H>
77 
79  (
80  fvm::ddt(rho, U)
81  + fvm::div(phi, U)
82  - fvm::laplacian(mu, U)
83  );
84 
85  solve(UEqn == -fvc::grad(p));
86 
87 
88  // --- PISO loop
89 
90  for (int corr=0; corr<nCorr; corr++)
91  {
92  volScalarField rUA = 1.0/UEqn.A();
93  U = rUA*UEqn.H();
94 
96  (
97  "phid",
98  psi
99  *(
100  (fvc::interpolate(U) & mesh.Sf())
101  + fvc::ddtPhiCorr(rUA, rho, U, phi)
102  )
103  );
104 
105  phi = (rhoO/psi)*phid;
106 
107  fvScalarMatrix pEqn
108  (
109  fvm::ddt(psi, p)
110  + fvc::div(phi)
111  + fvm::div(phid, p)
112  - fvm::laplacian(rho*rUA, p)
113  );
114 
115  pEqn.solve();
116 
117  phi += pEqn.flux();
118 
120 
121  U -= rUA*fvc::grad(p);
122  U.correctBoundaryConditions();
123  }
124 
125  rho = rhoO + psi*p;
126 
127  runTime.write();
128 
129  Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
130  << " ClockTime = " << runTime.elapsedClockTime() << " s"
131  << nl << endl;
132  }
133 
134  Info<< "End\n" << endl;
135 
136  return 0;
137 }
138 
139 
140 // ************************ vim: set sw=4 sts=4 et: ************************ //