1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend 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 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
28 Density-based compressible flow solver based on central-upwind schemes of
31 \*---------------------------------------------------------------------------*/
34 #include "basicPsiThermo.H"
35 #include "zeroGradientFvPatchFields.H"
36 #include "fixedRhoFvPatchScalarField.H"
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 int main(int argc, char *argv[])
42 # include "setRootCase.H"
44 # include "createTime.H"
45 # include "createMesh.H"
46 # include "createFields.H"
47 # include "readThermophysicalProperties.H"
48 # include "readTimeControls.H"
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 # include "readFluxScheme.H"
54 dimensionedScalar v_zero("v_zero",dimVolume/dimTime, 0.0);
56 Info<< "\nStarting time loop\n" << endl;
60 // --- upwind interpolation of primitive fields on faces
62 surfaceScalarField rho_pos =
63 fvc::interpolate(rho, pos, "reconstruct(rho)");
64 surfaceScalarField rho_neg =
65 fvc::interpolate(rho, neg, "reconstruct(rho)");
67 surfaceVectorField rhoU_pos =
68 fvc::interpolate(rhoU, pos, "reconstruct(U)");
69 surfaceVectorField rhoU_neg =
70 fvc::interpolate(rhoU, neg, "reconstruct(U)");
72 volScalarField rPsi = 1.0/psi;
73 surfaceScalarField rPsi_pos =
74 fvc::interpolate(rPsi, pos, "reconstruct(T)");
75 surfaceScalarField rPsi_neg =
76 fvc::interpolate(rPsi, neg, "reconstruct(T)");
78 surfaceScalarField e_pos =
79 fvc::interpolate(e, pos, "reconstruct(T)");
80 surfaceScalarField e_neg =
81 fvc::interpolate(e, neg, "reconstruct(T)");
83 surfaceVectorField U_pos = rhoU_pos/rho_pos;
84 surfaceVectorField U_neg = rhoU_neg/rho_neg;
86 surfaceScalarField p_pos = rho_pos*rPsi_pos;
87 surfaceScalarField p_neg = rho_neg*rPsi_neg;
89 surfaceScalarField phiv_pos = U_pos & mesh.Sf();
90 surfaceScalarField phiv_neg = U_neg & mesh.Sf();
92 volScalarField c = sqrt(thermo.Cp()/thermo.Cv()*rPsi);
93 surfaceScalarField cSf_pos = fvc::interpolate(c, pos, "reconstruct(T)")*mesh.magSf();
94 surfaceScalarField cSf_neg = fvc::interpolate(c, neg, "reconstruct(T)")*mesh.magSf();
96 surfaceScalarField ap = max(max(phiv_pos + cSf_pos, phiv_neg + cSf_neg), v_zero);
97 surfaceScalarField am = min(min(phiv_pos - cSf_pos, phiv_neg - cSf_neg), v_zero);
99 surfaceScalarField a_pos = ap/(ap - am);
101 surfaceScalarField amaxSf("amaxSf", max(mag(am), mag(ap)));
103 surfaceScalarField aSf = am*a_pos;
105 if (fluxScheme == "Tadmor")
111 surfaceScalarField a_neg = (1.0 - a_pos);
116 surfaceScalarField aphiv_pos = phiv_pos - aSf;
117 surfaceScalarField aphiv_neg = phiv_neg + aSf;
119 // Reuse amaxSf for the maximum positive and negative fluxes
120 // estimated by the central scheme
121 amaxSf = max(mag(aphiv_pos), mag(aphiv_neg));
123 #include "compressibleCourantNo.H"
124 #include "readTimeControls.H"
125 #include "setDeltaT.H"
129 Info<< "Time = " << runTime.timeName() << nl << endl;
131 surfaceScalarField phi("phi", aphiv_pos*rho_pos + aphiv_neg*rho_neg);
133 surfaceVectorField phiUp =
134 (aphiv_pos*rhoU_pos + aphiv_neg*rhoU_neg)
135 + (a_pos*p_pos + a_neg*p_neg)*mesh.Sf();
137 surfaceScalarField phiEp =
138 aphiv_pos*(rho_pos*(e_pos + 0.5*magSqr(U_pos)) + p_pos)
139 + aphiv_neg*(rho_neg*(e_neg + 0.5*magSqr(U_neg)) + p_neg)
140 + aSf*p_pos - aSf*p_neg;
142 volTensorField tauMC("tauMC", mu*dev2(fvc::grad(U)().T()));
145 solve(fvm::ddt(rho) + fvc::div(phi));
147 // --- Solve momentum
148 solve(fvm::ddt(rhoU) + fvc::div(phiUp));
150 U.dimensionedInternalField() =
151 rhoU.dimensionedInternalField()
152 /rho.dimensionedInternalField();
153 U.correctBoundaryConditions();
154 rhoU.boundaryField() = rho.boundaryField()*U.boundaryField();
156 volScalarField rhoBydt(rho/runTime.deltaT());
162 fvm::ddt(rho, U) - fvc::ddt(rho, U)
163 - fvm::laplacian(mu, U)
170 surfaceScalarField sigmaDotU =
173 fvc::interpolate(mu)*mesh.magSf()*fvc::snGrad(U)
174 + (mesh.Sf() & fvc::interpolate(tauMC))
176 & (a_pos*U_pos + a_neg*U_neg)
183 - fvc::div(sigmaDotU)
186 e = rhoE/rho - 0.5*magSqr(U);
187 e.correctBoundaryConditions();
189 rhoE.boundaryField() =
192 e.boundaryField() + 0.5*magSqr(U.boundaryField())
197 volScalarField k("k", thermo.Cp()*mu/Pr);
200 fvm::ddt(rho, e) - fvc::ddt(rho, e)
201 - fvm::laplacian(thermo.alpha(), e)
202 + fvc::laplacian(thermo.alpha(), e)
203 - fvc::laplacian(k, T)
206 rhoE = rho*(e + 0.5*magSqr(U));
209 p.dimensionedInternalField() =
210 rho.dimensionedInternalField()
211 /psi.dimensionedInternalField();
212 p.correctBoundaryConditions();
213 rho.boundaryField() = psi.boundaryField()*p.boundaryField();
217 Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
218 << " ClockTime = " << runTime.elapsedClockTime() << " s"
222 Info<< "End\n" << endl;
227 // ************************************************************************* //