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 <compressibleTurbulenceModel/turbulenceModel.H>
00054 #include <finiteVolume/zeroGradientFvPatchFields.H>
00055 #include <rhoCentralBCs/fixedRhoFvPatchScalarField.H>
00056 #include <dynamicMesh/motionSolver.H>
00057
00058
00059
00060 int main(int argc, char *argv[])
00061 {
00062 #include <OpenFOAM/setRootCase.H>
00063
00064 #include <OpenFOAM/createTime.H>
00065 #include <OpenFOAM/createMesh.H>
00066 #include "createFields.H"
00067 #include "../readThermophysicalProperties.H"
00068 #include <finiteVolume/readTimeControls.H>
00069
00070
00071
00072 #include "../readFluxScheme.H"
00073
00074 dimensionedScalar v_zero("v_zero", dimVolume/dimTime, 0.0);
00075
00076 Info<< "\nStarting time loop\n" << endl;
00077
00078 autoPtr<Foam::motionSolver> motionPtr = motionSolver::New(mesh);
00079
00080 while (runTime.run())
00081 {
00082
00083
00084 surfaceScalarField rho_pos =
00085 fvc::interpolate(rho, pos, "reconstruct(rho)");
00086 surfaceScalarField rho_neg =
00087 fvc::interpolate(rho, neg, "reconstruct(rho)");
00088
00089 surfaceVectorField rhoU_pos =
00090 fvc::interpolate(rhoU, pos, "reconstruct(U)");
00091 surfaceVectorField rhoU_neg =
00092 fvc::interpolate(rhoU, neg, "reconstruct(U)");
00093
00094 volScalarField rPsi = 1.0/psi;
00095 surfaceScalarField rPsi_pos =
00096 fvc::interpolate(rPsi, pos, "reconstruct(T)");
00097 surfaceScalarField rPsi_neg =
00098 fvc::interpolate(rPsi, neg, "reconstruct(T)");
00099
00100 surfaceScalarField e_pos =
00101 fvc::interpolate(e, pos, "reconstruct(T)");
00102 surfaceScalarField e_neg =
00103 fvc::interpolate(e, neg, "reconstruct(T)");
00104
00105 surfaceVectorField U_pos = rhoU_pos/rho_pos;
00106 surfaceVectorField U_neg = rhoU_neg/rho_neg;
00107
00108 surfaceScalarField p_pos = rho_pos*rPsi_pos;
00109 surfaceScalarField p_neg = rho_neg*rPsi_neg;
00110
00111 surfaceScalarField phiv_pos = U_pos & mesh.Sf();
00112 surfaceScalarField phiv_neg = U_neg & mesh.Sf();
00113
00114 volScalarField c = sqrt(thermo.Cp()/thermo.Cv()*rPsi);
00115 surfaceScalarField cSf_pos =
00116 fvc::interpolate(c, pos, "reconstruct(T)")*mesh.magSf();
00117 surfaceScalarField cSf_neg =
00118 fvc::interpolate(c, neg, "reconstruct(T)")*mesh.magSf();
00119
00120 surfaceScalarField ap =
00121 max(max(phiv_pos + cSf_pos, phiv_neg + cSf_neg), v_zero);
00122 surfaceScalarField am =
00123 min(min(phiv_pos - cSf_pos, phiv_neg - cSf_neg), v_zero);
00124
00125 surfaceScalarField a_pos = ap/(ap - am);
00126
00127 surfaceScalarField amaxSf("amaxSf", max(mag(am), mag(ap)));
00128
00129 surfaceScalarField aSf = am*a_pos;
00130
00131 if (fluxScheme == "Tadmor")
00132 {
00133 aSf = -0.5*amaxSf;
00134 a_pos = 0.5;
00135 }
00136
00137 surfaceScalarField a_neg = (1.0 - a_pos);
00138
00139 phiv_pos *= a_pos;
00140 phiv_neg *= a_neg;
00141
00142 surfaceScalarField aphiv_pos = phiv_pos - aSf;
00143 surfaceScalarField aphiv_neg = phiv_neg + aSf;
00144
00145
00146
00147 amaxSf = max(mag(aphiv_pos), mag(aphiv_neg));
00148
00149 #include "../compressibleCourantNo.H"
00150 #include <finiteVolume/readTimeControls.H>
00151 #include <finiteVolume/setDeltaT.H>
00152
00153 runTime++;
00154
00155 Info<< "Time = " << runTime.timeName() << nl << endl;
00156
00157 mesh.movePoints(motionPtr->newPoints());
00158 phiv_pos = U_pos & mesh.Sf();
00159 phiv_neg = U_neg & mesh.Sf();
00160 fvc::makeRelative(phiv_pos, U);
00161 fvc::makeRelative(phiv_neg, U);
00162 phiv_neg -= mesh.phi();
00163 phiv_pos *= a_pos;
00164 phiv_neg *= a_neg;
00165 aphiv_pos = phiv_pos - aSf;
00166 aphiv_neg = phiv_neg + aSf;
00167
00168 phi = aphiv_pos*rho_pos + aphiv_neg*rho_neg;
00169
00170 surfaceVectorField phiUp =
00171 (aphiv_pos*rhoU_pos + aphiv_neg*rhoU_neg)
00172 + (a_pos*p_pos + a_neg*p_neg)*mesh.Sf();
00173
00174 surfaceScalarField phiEp =
00175 aphiv_pos*(rho_pos*(e_pos + 0.5*magSqr(U_pos)) + p_pos)
00176 + aphiv_neg*(rho_neg*(e_neg + 0.5*magSqr(U_neg)) + p_neg)
00177 + aSf*p_pos - aSf*p_neg;
00178
00179 volScalarField muEff = turbulence->muEff();
00180 volTensorField tauMC("tauMC", muEff*dev2(fvc::grad(U)().T()));
00181
00182
00183 solve(fvm::ddt(rho) + fvc::div(phi));
00184
00185
00186 solve(fvm::ddt(rhoU) + fvc::div(phiUp));
00187
00188 U.dimensionedInternalField() =
00189 rhoU.dimensionedInternalField()
00190 /rho.dimensionedInternalField();
00191 U.correctBoundaryConditions();
00192 rhoU.boundaryField() = rho.boundaryField()*U.boundaryField();
00193
00194 if (!inviscid)
00195 {
00196 solve
00197 (
00198 fvm::ddt(rho, U) - fvc::ddt(rho, U)
00199 - fvm::laplacian(muEff, U)
00200 - fvc::div(tauMC)
00201 );
00202 rhoU = rho*U;
00203 }
00204
00205
00206 surfaceScalarField sigmaDotU =
00207 (
00208 (
00209 fvc::interpolate(muEff)*mesh.magSf()*fvc::snGrad(U)
00210 + (mesh.Sf() & fvc::interpolate(tauMC))
00211 )
00212 & (a_pos*U_pos + a_neg*U_neg)
00213 );
00214
00215 solve
00216 (
00217 fvm::ddt(rhoE)
00218 + fvc::div(phiEp)
00219 - fvc::div(sigmaDotU)
00220 );
00221
00222 e = rhoE/rho - 0.5*magSqr(U);
00223 e.correctBoundaryConditions();
00224 thermo.correct();
00225 rhoE.boundaryField() =
00226 rho.boundaryField()*
00227 (
00228 e.boundaryField() + 0.5*magSqr(U.boundaryField())
00229 );
00230
00231 if (!inviscid)
00232 {
00233 volScalarField k("k", thermo.Cp()*muEff/Pr);
00234 solve
00235 (
00236 fvm::ddt(rho, e) - fvc::ddt(rho, e)
00237 - fvm::laplacian(turbulence->alphaEff(), e)
00238 + fvc::laplacian(turbulence->alpha(), e)
00239 - fvc::laplacian(k, T)
00240 );
00241 thermo.correct();
00242 rhoE = rho*(e + 0.5*magSqr(U));
00243 }
00244
00245 p.dimensionedInternalField() =
00246 rho.dimensionedInternalField()
00247 /psi.dimensionedInternalField();
00248 p.correctBoundaryConditions();
00249 rho.boundaryField() = psi.boundaryField()*p.boundaryField();
00250
00251 turbulence->correct();
00252
00253 runTime.write();
00254
00255 Info<< "ExecutionTime = " << runTime.elapsedCpuTime() << " s"
00256 << " ClockTime = " << runTime.elapsedClockTime() << " s"
00257 << nl << endl;
00258 }
00259
00260 Info<< "End\n" << endl;
00261
00262 return 0;
00263 }
00264
00265