openfoam/applications/solvers/compressible/rhoSimpleFoam/createFields.H
Henry Weller 798ac98aef BUG: simpleFoam: moved createFvOptions.H into createFields.H for -postProcess option
To unsure fvOptions are instantiated for post-processing createFvOptions.H must
be included in createFields.H rather than in the solver directly.

Resolves bug-report https://bugs.openfoam.org/view.php?id=2733

BUG: porousSimpleFoam: moved createFvOptions.H into createFields.H for -postProcess option

Resolves bug-report https://bugs.openfoam.org/view.php?id=2733

BUG: solvers: Moved fvOption construction into createFields.H for post-processing

This ensures that the fvOptions are constructed for the -postProcessing option
so that functionObjects which process fvOption data operate correctly in this
mode.
2017-10-23 22:20:52 +01:00

61 lines
1.0 KiB
C

Info<< "Reading thermophysical properties\n" << endl;
autoPtr<fluidThermo> pThermo
(
fluidThermo::New(mesh)
);
fluidThermo& thermo = pThermo();
thermo.validate(args.executable(), "h", "e");
volScalarField& p = thermo.p();
volScalarField rho
(
IOobject
(
"rho",
runTime.timeName(),
mesh,
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
thermo.rho()
);
Info<< "Reading field U\n" << endl;
volVectorField U
(
IOobject
(
"U",
runTime.timeName(),
mesh,
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
mesh
);
#include "compressibleCreatePhi.H"
pressureControl pressureControl(p, rho, simple.dict());
mesh.setFluxRequired(p.name());
Info<< "Creating turbulence model\n" << endl;
autoPtr<compressible::turbulenceModel> turbulence
(
compressible::turbulenceModel::New
(
rho,
U,
phi,
thermo
)
);
dimensionedScalar initialMass = fvc::domainIntegrate(rho);
#include "createMRF.H"
#include "createFvOptions.H"