intersection with triangle plane for miss
[OpenFOAM-1.5.x.git] / applications / solvers / compressible / sonicFoamAutoMotion / sonicFoamAutoMotion.C
blob445f88d578254a0cd9f9d65f7921dc6a111030a2
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     sonicFoamAutoMotion
28 Description
29     Transient solver for trans-sonic/supersonic, laminar flow of a 
30     compressible gas with mesh motion..
32 \*---------------------------------------------------------------------------*/
34 #include "fvCFD.H"
35 #include "motionSolver.H"
37 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 int main(int argc, char *argv[])
42 #   include "setRootCase.H"
43 #   include "createTime.H"
44 #   include "createMesh.H"
45 #   include "readThermodynamicProperties.H"
46 #   include "readTransportProperties.H"
47 #   include "createFields.H"
48 #   include "initContinuityErrs.H"
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52     Info<< "\nStarting time loop\n" << endl;
54     autoPtr<Foam::motionSolver> motionPtr = motionSolver::New(mesh);
56     for (runTime++; !runTime.end(); runTime++)
57     {
58         Info<< "Time = " << runTime.timeName() << nl << endl;
60 #       include "readPISOControls.H"
61 #       include "compressibleCourantNo.H"
63         mesh.movePoints(motionPtr->newPoints());
65 #       include "rhoEqn.H"
67         fvVectorMatrix UEqn
68         (
69             fvm::ddt(rho, U)
70           + fvm::div(phi, U)
71           - fvm::laplacian(mu, U)
72         );
74         solve(UEqn == -fvc::grad(p));
76         solve
77         (
78             fvm::ddt(rho, e)
79           + fvm::div(phi, e)
80           - fvm::laplacian(mu, e)
81           ==
82           - p*fvc::div(phi/fvc::interpolate(rho) + fvc::meshPhi(rho, U))
83           + mu*magSqr(symm(fvc::grad(U)))
84         );
86         T = e/Cv;
87         psi = 1.0/(R*T);
89         // --- PISO loop
91         for (int corr=0; corr<nCorr; corr++)
92         {
93             U = UEqn.H()/UEqn.A();
95             surfaceScalarField phid
96             (
97                 "phid",
98                 fvc::interpolate(psi)*
99                 (
100                     (fvc::interpolate(U) & mesh.Sf()) - fvc::meshPhi(rho, U)
101                 )
102             );
104             for (int nonOrth=0; nonOrth<=nNonOrthCorr; nonOrth++)
105             {
106                 fvScalarMatrix pEqn
107                 (
108                     fvm::ddt(psi, p)
109                   + fvm::div(phid, p)
110                   - fvm::laplacian(rho/UEqn.A(), p)
111                 );
113                 pEqn.solve();
115                 phi = pEqn.flux();
116             }
118 #           include "compressibleContinuityErrs.H"
120             U -= fvc::grad(p)/UEqn.A();
121             U.correctBoundaryConditions();
122         }
124         rho = psi*p;
126         runTime.write();
128         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
129             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
130             << nl << endl;
131     }
133     Info<< "End\n" << endl;
135     return(0);
139 // ************************************************************************* //