Go to the documentation of this file.
47 int main(
int argc,
char *argv[])
53 #include "createFields.H"
69 Info<<
"\nStarting time loop\n" <<
endl;
112 max(
max(phiv_pos + cSf_pos, phiv_neg + cSf_neg), v_zero)
117 min(
min(phiv_pos - cSf_pos, phiv_neg - cSf_neg), v_zero)
142 amaxSf =
max(
mag(aphiv_pos),
mag(aphiv_neg));
149 #include "setRDeltaT.H"
153 #include "setDeltaT.H"
158 Info<<
"Time = " << runTime.timeName() <<
nl <<
endl;
160 phi = aphiv_pos*rho_pos + aphiv_neg*rho_neg;
164 (aphiv_pos*rhoU_pos + aphiv_neg*rhoU_neg)
165 + (a_pos*p_pos + a_neg*p_neg)*
mesh.Sf()
171 aphiv_pos*(rho_pos*(e_pos + 0.5*
magSqr(U_pos)) + p_pos)
172 + aphiv_neg*(rho_neg*(e_neg + 0.5*
magSqr(U_neg)) + p_neg)
173 + aSf*p_pos - aSf*p_neg
185 U.dimensionedInternalField() =
186 rhoU.dimensionedInternalField()
187 /
rho.dimensionedInternalField();
188 U.correctBoundaryConditions();
189 rhoU.boundaryField() ==
rho.boundaryField()*
U.boundaryField();
210 & (a_pos*U_pos + a_neg*U_neg)
221 e.correctBoundaryConditions();
223 rhoE.boundaryField() ==
226 e.boundaryField() + 0.5*
magSqr(
U.boundaryField())
240 p.dimensionedInternalField() =
241 rho.dimensionedInternalField()
242 /
psi.dimensionedInternalField();
243 p.correctBoundaryConditions();
244 rho.boundaryField() ==
psi.boundaryField()*
p.boundaryField();
250 Info<<
"ExecutionTime = " << runTime.elapsedCpuTime() <<
" s"
251 <<
" ClockTime = " << runTime.elapsedClockTime() <<
" s"
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
GeometricField< tensor, fvPatchField, volMesh > volTensorField
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
tmp< surfaceScalarField > interpolate(const RhoType &rho)
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
Calculates the mean and maximum wave speed based Courant Numbers.
tmp< GeometricField< Type, fvPatchField, volMesh > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
autoPtr< compressible::turbulenceModel > turbulence
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
int main(int argc, char *argv[])
GeometricField< scalar, fvPatchField, volMesh > volScalarField
word fluxScheme("Kurganov")
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
const volScalarField & psi
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Read the control parameters used by setDeltaT.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Read the control parameters used by setDeltaT.
const dimensionedScalar e
Elementary charge.
const dimensionedScalar c
Speed of light in a vacuum.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar neg(const dimensionedScalar &ds)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Provides functions smooth spread and sweep which use the FaceCellWave algorithm to smooth and redistr...
dimensionedScalar pos(const dimensionedScalar &ds)