initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / applications / solvers / compressible / rhoCentralFoam / rhoCentralFoam.C
bloba9146479317566e629ccf0658f3f3a6e926fbb76
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     rhoCentralFoam
28 Description
29     Density-based compressible flow solver based on central-upwind schemes of
30     Kurganov and Tadmor
32 \*---------------------------------------------------------------------------*/
34 #include "fvCFD.H"
35 #include "basicThermo.H"
36 #include "zeroGradientFvPatchFields.H"
37 #include "fixedRhoFvPatchScalarField.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 int main(int argc, char *argv[])
44 #   include "setRootCase.H"
46 #   include "createTime.H"
47 #   include "createMesh.H"
48 #   include "createFields.H"
49 #   include "readThermophysicalProperties.H"
50 #   include "readTimeControls.H"
52 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
54 #   include "readFluxScheme.H"
56     dimensionedScalar v_zero("v_zero",dimVolume/dimTime, 0.0);
58     Info<< "\nStarting time loop\n" << endl;
60     while (runTime.run())
61     {
62         // --- upwind interpolation of primitive fields on faces
64         surfaceScalarField rho_pos =
65             fvc::interpolate(rho, pos, "reconstruct(rho)");
66         surfaceScalarField rho_neg =
67             fvc::interpolate(rho, neg, "reconstruct(rho)");
69         surfaceVectorField rhoU_pos =
70             fvc::interpolate(rhoU, pos, "reconstruct(U)");
71         surfaceVectorField rhoU_neg =
72             fvc::interpolate(rhoU, neg, "reconstruct(U)");
74         volScalarField rPsi = 1.0/psi;
75         surfaceScalarField rPsi_pos =
76             fvc::interpolate(rPsi, pos, "reconstruct(T)");
77         surfaceScalarField rPsi_neg =
78             fvc::interpolate(rPsi, neg, "reconstruct(T)");
80         surfaceScalarField h_pos =
81             fvc::interpolate(h, pos, "reconstruct(T)");
82         surfaceScalarField h_neg =
83             fvc::interpolate(h, neg, "reconstruct(T)");
85         surfaceVectorField U_pos = rhoU_pos/rho_pos;
86         surfaceVectorField U_neg = rhoU_neg/rho_neg;
88         surfaceScalarField p_pos = rho_pos*rPsi_pos;
89         surfaceScalarField p_neg = rho_neg*rPsi_neg;
91         surfaceScalarField phiv_pos = U_pos & mesh.Sf();
92         surfaceScalarField phiv_neg = U_neg & mesh.Sf();
94         volScalarField c = sqrt(thermo->Cp()/thermo->Cv()*rPsi);
95         surfaceScalarField cSf_pos = fvc::interpolate(c, pos, "reconstruct(T)")*mesh.magSf();
96         surfaceScalarField cSf_neg = fvc::interpolate(c, neg, "reconstruct(T)")*mesh.magSf();
98         surfaceScalarField ap = max(max(phiv_pos + cSf_pos, phiv_neg + cSf_neg), v_zero);
99         surfaceScalarField am = min(min(phiv_pos - cSf_pos, phiv_neg - cSf_neg), v_zero);
101         surfaceScalarField a_pos = ap/(ap - am);
103         surfaceScalarField amaxSf("amaxSf", max(mag(am), mag(ap)));
105 #       include "compressibleCourantNo.H"
106 #       include "readTimeControls.H"
107 #       include "setDeltaT.H"
109         runTime++;
111         Info<< "Time = " << runTime.timeName() << nl << endl;
113         surfaceScalarField aSf = am*a_pos;
115         if (fluxScheme == "Tadmor")
116         {
117             aSf = -0.5*amaxSf;
118             a_pos = 0.5;
119         }
121         surfaceScalarField a_neg = (1.0 - a_pos);
123         phiv_pos *= a_pos;
124         phiv_neg *= a_neg;
126         surfaceScalarField aphiv_pos = phiv_pos - aSf;
127         surfaceScalarField aphiv_neg = phiv_neg + aSf;
129         surfaceScalarField phi("phi", aphiv_pos*rho_pos + aphiv_neg*rho_neg);
131         surfaceVectorField phiUp =
132             (aphiv_pos*rhoU_pos + aphiv_neg*rhoU_neg)
133           + (a_pos*p_pos + a_neg*p_neg)*mesh.Sf();
135         surfaceScalarField phiEp =
136             aphiv_pos*rho_pos*(h_pos + 0.5*magSqr(U_pos))
137           + aphiv_neg*rho_neg*(h_neg + 0.5*magSqr(U_neg))
138           + aSf*p_pos - aSf*p_neg;
140         volTensorField tauMC("tauMC", mu*dev2(fvc::grad(U)().T()));
142         // --- Solve density
143         solve(fvm::ddt(rho) + fvc::div(phi));
145         // --- Solve momentum
146         solve(fvm::ddt(rhoU) + fvc::div(phiUp));
148         U.dimensionedInternalField() =
149             rhoU.dimensionedInternalField()
150            /rho.dimensionedInternalField();
151         U.correctBoundaryConditions();
152         rhoU.boundaryField() = rho.boundaryField()*U.boundaryField();
154         volScalarField rhoBydt(rho/runTime.deltaT());
156         if (!inviscid)
157         {
158             solve
159             (
160                 fvm::ddt(rho, U) - fvc::ddt(rho,U)
161               - fvm::laplacian(mu, U)
162               - fvc::div(tauMC)
163             );
164             rhoU = rho*U;
165         }
167         // --- Solve energy
168         surfaceScalarField sigmaDotU =
169         (
170             (
171                 fvc::interpolate(mu)*mesh.magSf()*fvc::snGrad(U)
172               + (mesh.Sf() & fvc::interpolate(tauMC))
173             )
174             & (a_pos*U_pos + a_neg*U_neg)
175         );
177         solve
178         (
179             fvm::ddt(rhoE)
180           + fvc::div(phiEp)
181           - fvc::div(sigmaDotU)
182         );
184         h = (rhoE + p)/rho - 0.5*magSqr(U);
185         h.correctBoundaryConditions();
186         thermo->correct();
187         rhoE.boundaryField() =
188             rho.boundaryField()*
189             (
190                 h.boundaryField() + 0.5*magSqr(U.boundaryField())
191             )
192           - p.boundaryField();
194         if (!inviscid)
195         {
196             volScalarField k("k", thermo->Cp()*mu/Pr);
197             solve
198             (
199                 fvm::ddt(rho, h) - fvc::ddt(rho, h)
200               - fvm::laplacian(thermo->alpha(), h)
201               + fvc::laplacian(thermo->alpha(), h)
202               - fvc::laplacian(k, T)
203             );
204             thermo->correct();
205             rhoE = rho*(h + 0.5*magSqr(U)) - p;
206         }
208         p.dimensionedInternalField() =
209             rho.dimensionedInternalField()
210            /psi.dimensionedInternalField();
211         p.correctBoundaryConditions();
212         rho.boundaryField() = psi.boundaryField()*p.boundaryField();
214         runTime.write();
216         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
217             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
218             << nl << endl;
219     }
221     Info<< "End\n" << endl;
223     return(0);
226 // ************************************************************************* //