BUGFIX: Uninitialised member variables
[foam-extend-3.2.git] / applications / solvers / compressible / rhoCentralFoam / rhoCentralFoam.C
blob5b7d57f4334aab3d69152d409aac4b6e61b1fa55
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright held by original author
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 "basicPsiThermo.H"
36 #include "zeroGradientFvPatchFields.H"
37 #include "fixedRhoFvPatchScalarField.H"
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 int main(int argc, char *argv[])
43 #   include "setRootCase.H"
45 #   include "createTime.H"
46 #   include "createMesh.H"
47 #   include "createFields.H"
48 #   include "readThermophysicalProperties.H"
49 #   include "readTimeControls.H"
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 #   include "readFluxScheme.H"
55     dimensionedScalar v_zero("v_zero",dimVolume/dimTime, 0.0);
57     Info<< "\nStarting time loop\n" << endl;
59     while (runTime.run())
60     {
61         // --- upwind interpolation of primitive fields on faces
63         surfaceScalarField rho_pos =
64             fvc::interpolate(rho, pos, "reconstruct(rho)");
65         surfaceScalarField rho_neg =
66             fvc::interpolate(rho, neg, "reconstruct(rho)");
68         surfaceVectorField rhoU_pos =
69             fvc::interpolate(rhoU, pos, "reconstruct(U)");
70         surfaceVectorField rhoU_neg =
71             fvc::interpolate(rhoU, neg, "reconstruct(U)");
73         volScalarField rPsi = 1.0/psi;
74         surfaceScalarField rPsi_pos =
75             fvc::interpolate(rPsi, pos, "reconstruct(T)");
76         surfaceScalarField rPsi_neg =
77             fvc::interpolate(rPsi, neg, "reconstruct(T)");
79         surfaceScalarField e_pos =
80             fvc::interpolate(e, pos, "reconstruct(T)");
81         surfaceScalarField e_neg =
82             fvc::interpolate(e, neg, "reconstruct(T)");
84         surfaceVectorField U_pos = rhoU_pos/rho_pos;
85         surfaceVectorField U_neg = rhoU_neg/rho_neg;
87         surfaceScalarField p_pos = rho_pos*rPsi_pos;
88         surfaceScalarField p_neg = rho_neg*rPsi_neg;
90         surfaceScalarField phiv_pos = U_pos & mesh.Sf();
91         surfaceScalarField phiv_neg = U_neg & mesh.Sf();
93         volScalarField c = sqrt(thermo.Cp()/thermo.Cv()*rPsi);
94         surfaceScalarField cSf_pos = fvc::interpolate(c, pos, "reconstruct(T)")*mesh.magSf();
95         surfaceScalarField cSf_neg = fvc::interpolate(c, neg, "reconstruct(T)")*mesh.magSf();
97         surfaceScalarField ap = max(max(phiv_pos + cSf_pos, phiv_neg + cSf_neg), v_zero);
98         surfaceScalarField am = min(min(phiv_pos - cSf_pos, phiv_neg - cSf_neg), v_zero);
100         surfaceScalarField a_pos = ap/(ap - am);
102         surfaceScalarField amaxSf("amaxSf", max(mag(am), mag(ap)));
104         surfaceScalarField aSf = am*a_pos;
106         if (fluxScheme == "Tadmor")
107         {
108             aSf = -0.5*amaxSf;
109             a_pos = 0.5;
110         }
112         surfaceScalarField a_neg = (1.0 - a_pos);
114         phiv_pos *= a_pos;
115         phiv_neg *= a_neg;
117         surfaceScalarField aphiv_pos = phiv_pos - aSf;
118         surfaceScalarField aphiv_neg = phiv_neg + aSf;
120         // Reuse amaxSf for the maximum positive and negative fluxes
121         // estimated by the central scheme
122         amaxSf = max(mag(aphiv_pos), mag(aphiv_neg));
124         #include "compressibleCourantNo.H"
125         #include "readTimeControls.H"
126         #include "setDeltaT.H"
128         runTime++;
130         Info<< "Time = " << runTime.timeName() << nl << endl;
132         surfaceScalarField phi("phi", aphiv_pos*rho_pos + aphiv_neg*rho_neg);
134         surfaceVectorField phiUp =
135             (aphiv_pos*rhoU_pos + aphiv_neg*rhoU_neg)
136           + (a_pos*p_pos + a_neg*p_neg)*mesh.Sf();
138         surfaceScalarField phiEp =
139             aphiv_pos*(rho_pos*(e_pos + 0.5*magSqr(U_pos)) + p_pos)
140           + aphiv_neg*(rho_neg*(e_neg + 0.5*magSqr(U_neg)) + p_neg)
141           + aSf*p_pos - aSf*p_neg;
143         volTensorField tauMC("tauMC", mu*dev2(fvc::grad(U)().T()));
145         // --- Solve density
146         solve(fvm::ddt(rho) + fvc::div(phi));
148         // --- Solve momentum
149         solve(fvm::ddt(rhoU) + fvc::div(phiUp));
151         U.dimensionedInternalField() =
152             rhoU.dimensionedInternalField()
153            /rho.dimensionedInternalField();
154         U.correctBoundaryConditions();
155         rhoU.boundaryField() = rho.boundaryField()*U.boundaryField();
157         volScalarField rhoBydt(rho/runTime.deltaT());
159         if (!inviscid)
160         {
161             solve
162             (
163                 fvm::ddt(rho, U) - fvc::ddt(rho, U)
164               - fvm::laplacian(mu, U)
165               - fvc::div(tauMC)
166             );
167             rhoU = rho*U;
168         }
170         // --- Solve energy
171         surfaceScalarField sigmaDotU =
172         (
173             (
174                 fvc::interpolate(mu)*mesh.magSf()*fvc::snGrad(U)
175               + (mesh.Sf() & fvc::interpolate(tauMC))
176             )
177             & (a_pos*U_pos + a_neg*U_neg)
178         );
180         solve
181         (
182             fvm::ddt(rhoE)
183           + fvc::div(phiEp)
184           - fvc::div(sigmaDotU)
185         );
187         e = rhoE/rho - 0.5*magSqr(U);
188         e.correctBoundaryConditions();
189         thermo.correct();
190         rhoE.boundaryField() =
191             rho.boundaryField()*
192             (
193                 e.boundaryField() + 0.5*magSqr(U.boundaryField())
194             );
196         if (!inviscid)
197         {
198             volScalarField k("k", thermo.Cp()*mu/Pr);
199             solve
200             (
201                 fvm::ddt(rho, e) - fvc::ddt(rho, e)
202               - fvm::laplacian(thermo.alpha(), e)
203               + fvc::laplacian(thermo.alpha(), e)
204               - fvc::laplacian(k, T)
205             );
206             thermo.correct();
207             rhoE = rho*(e + 0.5*magSqr(U));
208         }
210         p.dimensionedInternalField() =
211             rho.dimensionedInternalField()
212            /psi.dimensionedInternalField();
213         p.correctBoundaryConditions();
214         rho.boundaryField() = psi.boundaryField()*p.boundaryField();
216         runTime.write();
218         Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
219             << "  ClockTime = " << runTime.elapsedClockTime() << " s"
220             << nl << endl;
221     }
223     Info<< "End\n" << endl;
225     return 0;
228 // ************************************************************************* //