Backported rhoPorousMRFPimpleFoam, rhoPorousSimpleFoam, rhoSimpleFoam, rhoSonicFoam...
[foam-extend-4.0.git] / applications / solvers / compressible / rhoPorousSimpleFoam / createFields.H
blob0e1ba2436562984867ccd0ed09fc54d86d720452
1     Info<< "Reading thermophysical properties\n" << endl;
3     autoPtr<basicPsiThermo> pThermo
4     (
5         basicPsiThermo::New(mesh)
6     );
7     basicPsiThermo& thermo = pThermo();
9     volScalarField rho
10     (
11         IOobject
12         (
13             "rho",
14             runTime.timeName(),
15             mesh,
16             IOobject::READ_IF_PRESENT,
17             IOobject::AUTO_WRITE
18         ),
19         thermo.rho()
20     );
22     volScalarField& p = thermo.p();
23     volScalarField& h = thermo.h();
24     const volScalarField& psi = thermo.psi();
27     Info<< "Reading field U\n" << endl;
28     volVectorField U
29     (
30         IOobject
31         (
32             "U",
33             runTime.timeName(),
34             mesh,
35             IOobject::MUST_READ,
36             IOobject::AUTO_WRITE
37         ),
38         mesh
39     );
41     #include "compressibleCreatePhi.H"
44     label pRefCell = 0;
45     scalar pRefValue = 0.0;
46     setRefCell(p, mesh.solutionDict().subDict("SIMPLE"), pRefCell, pRefValue);
47     mesh.schemesDict().setFluxRequired(p.name());
49     dimensionedScalar pMin
50     (
51         mesh.solutionDict().subDict("SIMPLE").lookup("pMin")
52     );
54     Info<< "Creating turbulence model\n" << endl;
55     autoPtr<compressible::RASModel> turbulence
56     (
57         compressible::RASModel::New
58         (
59             rho,
60             U,
61             phi,
62             thermo
63         )
64     );
66     dimensionedScalar initialMass = fvc::domainIntegrate(rho);
68     porousZones pZones(mesh);
69     Switch pressureImplicitPorosity(false);
71     int nUCorr = 0;
72     if (pZones.size())
73     {
74         // nUCorrectors for pressureImplicitPorosity
75         if (mesh.solutionDict().subDict("SIMPLE").found("nUCorrectors"))
76         {
77             nUCorr = readInt
78             (
79                 mesh.solutionDict().subDict("SIMPLE").lookup("nUCorrectors")
80             );
81         }
83         if (nUCorr > 0)
84         {
85             pressureImplicitPorosity = true;
86         }
87     }