Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
00051 #include <finiteVolume/fvCFD.H>
00052 #include <basicThermophysicalModels/basicPsiThermo.H>
00053 #include <finiteVolume/zeroGradientFvPatchFields.H>
00054 #include <rhoCentralBCs/fixedRhoFvPatchScalarField.H>
00055
00056
00057
00058 int main(int argc, char *argv[])
00059 {
00060 #include <OpenFOAM/setRootCase.H>
00061
00062 #include <OpenFOAM/createTime.H>
00063 #include <OpenFOAM/createMesh.H>
00064 #include "createFields.H"
00065 #include "readThermophysicalProperties.H"
00066 #include <finiteVolume/readTimeControls.H>
00067
00068
00069
00070 #include "readFluxScheme.H"
00071
00072 dimensionedScalar v_zero("v_zero",dimVolume/dimTime, 0.0);
00073
00074 Info<< "\nStarting time loop\n" << endl;
00075
00076 while (runTime.run())
00077 {
00078
00079
00080 surfaceScalarField rho_pos =
00081 fvc::interpolate(rho, pos, "reconstruct(rho)");
00082 surfaceScalarField rho_neg =
00083 fvc::interpolate(rho, neg, "reconstruct(rho)");
00084
00085 surfaceVectorField rhoU_pos =
00086 fvc::interpolate(rhoU, pos, "reconstruct(U)");
00087 surfaceVectorField rhoU_neg =
00088 fvc::interpolate(rhoU, neg, "reconstruct(U)");
00089
00090 volScalarField rPsi = 1.0/psi;
00091 surfaceScalarField rPsi_pos =
00092 fvc::interpolate(rPsi, pos, "reconstruct(T)");
00093 surfaceScalarField rPsi_neg =
00094 fvc::interpolate(rPsi, neg, "reconstruct(T)");
00095
00096 surfaceScalarField e_pos =
00097 fvc::interpolate(e, pos, "reconstruct(T)");
00098 surfaceScalarField e_neg =
00099 fvc::interpolate(e, neg, "reconstruct(T)");
00100
00101 surfaceVectorField U_pos = rhoU_pos/rho_pos;
00102 surfaceVectorField U_neg = rhoU_neg/rho_neg;
00103
00104 surfaceScalarField p_pos = rho_pos*rPsi_pos;
00105 surfaceScalarField p_neg = rho_neg*rPsi_neg;
00106
00107 surfaceScalarField phiv_pos = U_pos & mesh.Sf();
00108 surfaceScalarField phiv_neg = U_neg & mesh.Sf();
00109
00110 volScalarField c = sqrt(thermo.Cp()/thermo.Cv()*rPsi);
00111 surfaceScalarField cSf_pos = fvc::interpolate(c, pos, "reconstruct(T)")*mesh.magSf();
00112 surfaceScalarField cSf_neg = fvc::interpolate(c, neg, "reconstruct(T)")*mesh.magSf();
00113
00114 surfaceScalarField ap = max(max(phiv_pos + cSf_pos, phiv_neg + cSf_neg), v_zero);
00115 surfaceScalarField am = min(min(phiv_pos - cSf_pos, phiv_neg - cSf_neg), v_zero);
00116
00117 surfaceScalarField a_pos = ap/(ap - am);
00118
00119 surfaceScalarField amaxSf("amaxSf", max(mag(am), mag(ap)));
00120
00121 surfaceScalarField aSf = am*a_pos;
00122
00123 if (fluxScheme == "Tadmor")
00124 {
00125 aSf = -0.5*amaxSf;
00126 a_pos = 0.5;
00127 }
00128
00129 surfaceScalarField a_neg = (1.0 - a_pos);
00130
00131 phiv_pos *= a_pos;
00132 phiv_neg *= a_neg;
00133
00134 surfaceScalarField aphiv_pos = phiv_pos - aSf;
00135 surfaceScalarField aphiv_neg = phiv_neg + aSf;
00136
00137
00138
00139 amaxSf = max(mag(aphiv_pos), mag(aphiv_neg));
00140
00141 #include "compressibleCourantNo.H"
00142 #include <finiteVolume/readTimeControls.H>
00143 #include <finiteVolume/setDeltaT.H>
00144
00145 runTime++;
00146
00147 Info<< "Time = " << runTime.timeName() << nl << endl;
00148
00149 surfaceScalarField phi("phi", aphiv_pos*rho_pos + aphiv_neg*rho_neg);
00150
00151 surfaceVectorField phiUp =
00152 (aphiv_pos*rhoU_pos + aphiv_neg*rhoU_neg)
00153 + (a_pos*p_pos + a_neg*p_neg)*mesh.Sf();
00154
00155 surfaceScalarField phiEp =
00156 aphiv_pos*(rho_pos*(e_pos + 0.5*magSqr(U_pos)) + p_pos)
00157 + aphiv_neg*(rho_neg*(e_neg + 0.5*magSqr(U_neg)) + p_neg)
00158 + aSf*p_pos - aSf*p_neg;
00159
00160 volTensorField tauMC("tauMC", mu*dev2(fvc::grad(U)().T()));
00161
00162
00163 solve(fvm::ddt(rho) + fvc::div(phi));
00164
00165
00166 solve(fvm::ddt(rhoU) + fvc::div(phiUp));
00167
00168 U.dimensionedInternalField() =
00169 rhoU.dimensionedInternalField()
00170 /rho.dimensionedInternalField();
00171 U.correctBoundaryConditions();
00172 rhoU.boundaryField() = rho.boundaryField()*U.boundaryField();
00173
00174 volScalarField rhoBydt(rho/runTime.deltaT());
00175
00176 if (!inviscid)
00177 {
00178 solve
00179 (
00180 fvm::ddt(rho, U) - fvc::ddt(rho, U)
00181 - fvm::laplacian(mu, U)
00182 - fvc::div(tauMC)
00183 );
00184 rhoU = rho*U;
00185 }
00186
00187
00188 surfaceScalarField sigmaDotU =
00189 (
00190 (
00191 fvc::interpolate(mu)*mesh.magSf()*fvc::snGrad(U)
00192 + (mesh.Sf() & fvc::interpolate(tauMC))
00193 )
00194 & (a_pos*U_pos + a_neg*U_neg)
00195 );
00196
00197 solve
00198 (
00199 fvm::ddt(rhoE)
00200 + fvc::div(phiEp)
00201 - fvc::div(sigmaDotU)
00202 );
00203
00204 e = rhoE/rho - 0.5*magSqr(U);
00205 e.correctBoundaryConditions();
00206 thermo.correct();
00207 rhoE.boundaryField() =
00208 rho.boundaryField()*
00209 (
00210 e.boundaryField() + 0.5*magSqr(U.boundaryField())
00211 );
00212
00213 if (!inviscid)
00214 {
00215 volScalarField k("k", thermo.Cp()*mu/Pr);
00216 solve
00217 (
00218 fvm::ddt(rho, e) - fvc::ddt(rho, e)
00219 - fvm::laplacian(thermo.alpha(), e)
00220 + fvc::laplacian(thermo.alpha(), e)
00221 - fvc::laplacian(k, T)
00222 );
00223 thermo.correct();
00224 rhoE = rho*(e + 0.5*magSqr(U));
00225 }
00226
00227 p.dimensionedInternalField() =
00228 rho.dimensionedInternalField()
00229 /psi.dimensionedInternalField();
00230 p.correctBoundaryConditions();
00231 rho.boundaryField() = psi.boundaryField()*p.boundaryField();
00232
00233 runTime.write();
00234
00235 Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
00236 << " ClockTime = " << runTime.elapsedClockTime() << " s"
00237 << nl << endl;
00238 }
00239
00240 Info<< "End\n" << endl;
00241
00242 return 0;
00243 }
00244
00245