36#include "viscosityModel.H"
47 _args = autoPtr<argList>
49 new argList(argc, argv)
52 if (!
_args->checkRootCase())
54 Foam::FatalError.exit();
57 argList& args =
_args();
60 _simple = autoPtr<simpleControl>
70#include "initContinuityErrs.H"
87 middleExport =
para->ITHACAdict->lookupOrDefault<
bool>(
"middleExport",
true);
88 middleStep =
para->ITHACAdict->lookupOrDefault<label>(
"middleStep", 20);
99 volVectorField&
U =
_U();
101 volScalarField& E =
pThermo().he();
118 const volScalarField&
mu =
pThermo().mu();
119 volScalarField& mu_field =
const_cast<volScalarField&
>(
mu);
122 for (
int i = 0;
i < mu_field.boundaryFieldRef().size();
i++)
130 surfaceScalarField&
phi =
_phi();
131 fvVectorMatrix NLTerm = fvm::div(
phi,
U);
140 fvVectorMatrix viscTerm = - fvc::div((
rho * nuEff) * dev2(
T(fvc::grad(
141 U)))) - fvm::laplacian(
rho * nuEff,
U);
147 volVectorField gradP = fvc::grad(
p);
169 volScalarField& he =
pThermo().he();
170 surfaceScalarField&
phi =
_phi();
171 fvScalarMatrix fluxTerm = fvm::div(
phi, he);
178 surfaceScalarField&
phi =
_phi();
180 volScalarField kinEn = fvc::div(
phi, volScalarField(
"Ekp",
181 0.5 * magSqr(
U) +
p /
rho));
187 volScalarField& he =
pThermo().he();
188 fvScalarMatrix diffTerm = fvm::laplacian(
turbulence->alphaEff(), he);
195 volScalarField& he =
pThermo().he();
209 volVectorField&
U, volScalarField&
p)
212 volScalarField
rAU(1.0 /
214 HbyA.reset(
new volVectorField(constrainHbyA(
rAU * Ueqn.H(),
U,
216 phiHbyA.reset(
new surfaceScalarField(
"phiHbyA",
217 fvc::interpolate(
rho)*fvc::flux(
HbyA())));
222 volVectorField&
U, volScalarField&
p)
224 volScalarField divPhiHbyA = fvc::div(
getPhiHbyA(Ueqn,
U,
p));
231 volScalarField
rAU(1.0 /
233 rhorAUf.reset(
new surfaceScalarField(
"rhorAUf", fvc::interpolate(
rho *
rAU)));
240 fvScalarMatrix poissonTerm = fvm::laplacian(
getRhorAUf(Ueqn),
p);
245 volVectorField&
U, volScalarField&
p)
268 _mesh().objectRegistry::clear();
283 argList& args =
_args();
284 Info <<
"ReCreate time\n" << Foam::endl;
304 _simple = autoPtr<simpleControl>
314 fluidThermo::New(
mesh)
318 _p = autoPtr<volScalarField>
320 new volScalarField(
thermo.p())
322 volScalarField&
p =
thermo.p();
323 _rho = autoPtr<volScalarField>
332 IOobject::READ_IF_PRESENT,
339 _E = autoPtr<volScalarField>
341 new volScalarField(
thermo.he())
343 Info <<
"ReReading field U \n" << endl;
344 _U = autoPtr<volVectorField>
359 volVectorField&
U =
_U();
360 Info <<
"ReReading/calculating face flux field phi\n" << endl;
361 _phi = autoPtr<surfaceScalarField>
363 new surfaceScalarField
370 IOobject::READ_IF_PRESENT,
373 linearInterpolate(
rho *
U) &
mesh.Sf()
376 surfaceScalarField&
phi =
_phi();
381 mesh.setFluxRequired(
p.name());
382 Info <<
"ReCreating turbulence model\n" << endl;
383 turbulence = autoPtr<compressible::turbulenceModel>
385 compressible::turbulenceModel::New
395 new dimensionedScalar(fvc::domainIntegrate(
rho))
397 _psi = autoPtr<volScalarField>
399 new volScalarField(
thermo.psi())
402#include "initContinuityErrs.H"
Header file of the steadyNS class.
pressureControl & pressureControl
volScalarField getKinEnTerm(volVectorField &U, volScalarField &p)
volScalarField getDivPhiHbyA(fvVectorMatrix &Ueqn, volVectorField &U, volScalarField &p)
fvScalarMatrix getEmatrix(volVectorField &U, volScalarField &p)
autoPtr< fluidThermo > pThermo
CompressibleSteadyNS()
Null constructor.
volVectorField getGradP(volScalarField &p)
autoPtr< volScalarField > _p
autoPtr< dimensionedScalar > _initialMass
fvScalarMatrix getDiffTerm()
autoPtr< surfaceScalarField > phiHbyA
void changeViscosity(double mu_new)
Function to change the viscosity.
fvScalarMatrix getPmatrix(fvVectorMatrix &Ueqn, volVectorField &U, volScalarField &p)
bool middleExport
Export also intermediate fields.
surfaceScalarField getRhorAUf(fvVectorMatrix &Ueqn)
autoPtr< volScalarField > _E
autoPtr< fvScalarMatrix > Eeqn_global
fvScalarMatrix getFluxTerm()
autoPtr< compressible::turbulenceModel > turbulence
autoPtr< fvVectorMatrix > Ueqn_global
fvVectorMatrix getNLTerm(volVectorField &U)
label middleStep
Distancing between intermediate steps (for turbulent case only)
fvVectorMatrix getViscTerm(volVectorField &U)
autoPtr< pressureControl > _pressureControl
fvScalarMatrix getPoissonTerm(fvVectorMatrix &Ueqn, volScalarField &p)
autoPtr< volScalarField > _psi
autoPtr< volVectorField > HbyA
fvVectorMatrix getUmatrix(volVectorField &U)
autoPtr< volScalarField > _rho
surfaceScalarField getPhiHbyA(fvVectorMatrix &Ueqn, volVectorField &U, volScalarField &p)
void restart()
set all variables back to the values into the 0 folder
autoPtr< surfaceScalarField > rhorAUf
autoPtr< fvScalarMatrix > Peqn_global
Class for the definition of some general parameters, the parameters must be defined from the file ITH...
static ITHACAparameters * getInstance()
Gets an instance of ITHACAparameters, to be used if the instance is already existing.
autoPtr< volScalarField > _nut
Eddy viscosity field.
label counter
Counter used for the output of the full order solutions.
void assignBC(volVectorField &s, label BC_ind, Vector< double > &value)
Assign Boundary Condition to a volVectorField.
void assignIF(T &s, G &value)
Assign internal field condition.
bool offline
Boolean variable, it is 1 if the Offline phase has already been computed, else 0.
IOdictionary * ITHACAdict
dictionary to store input output infos
Eigen::MatrixXd mu
Row matrix of parameters.
autoPtr< argList > _args
argList
bool podex
Boolean variable, it is 1 if the POD has already been computed, else 0.
void truthSolve()
Perform a TruthSolve.
bool supex
Boolean variable to check the existence of the supremizer modes.
autoPtr< surfaceScalarField > _phi
Flux.
autoPtr< simpleControl > _simple
simpleControl
autoPtr< fv::options > _fvOptions
fvOptions
autoPtr< Time > _runTime
Time.
autoPtr< fvMesh > _mesh
Mesh.
autoPtr< volVectorField > _U
Velocity field.
bool check_pod()
Check if the POD data folder "./ITHACAoutput/POD" exists.
bool check_off()
Check if the offline data folder "./ITHACAoutput/Offline" exists.
bool check_sup()
Check if the supremizer folder exists.
volScalarField rAU(1.0/Ueqn.A())
simpleControl simple(mesh)