1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright held by original author
7 -------------------------------------------------------------------------------
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
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
29 Density-based compressible flow solver based on central-upwind schemes of
32 \*---------------------------------------------------------------------------*/
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;
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")
112 surfaceScalarField a_neg = (1.0 - a_pos);
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"
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()));
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());
163 fvm::ddt(rho, U) - fvc::ddt(rho, U)
164 - fvm::laplacian(mu, U)
171 surfaceScalarField sigmaDotU =
174 fvc::interpolate(mu)*mesh.magSf()*fvc::snGrad(U)
175 + (mesh.Sf() & fvc::interpolate(tauMC))
177 & (a_pos*U_pos + a_neg*U_neg)
184 - fvc::div(sigmaDotU)
187 e = rhoE/rho - 0.5*magSqr(U);
188 e.correctBoundaryConditions();
190 rhoE.boundaryField() =
193 e.boundaryField() + 0.5*magSqr(U.boundaryField())
198 volScalarField k("k", thermo.Cp()*mu/Pr);
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)
207 rhoE = rho*(e + 0.5*magSqr(U));
210 p.dimensionedInternalField() =
211 rho.dimensionedInternalField()
212 /psi.dimensionedInternalField();
213 p.correctBoundaryConditions();
214 rho.boundaryField() = psi.boundaryField()*p.boundaryField();
218 Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
219 << " ClockTime = " << runTime.elapsedClockTime() << " s"
223 Info<< "End\n" << endl;
228 // ************************************************************************* //