This repository has been archived on 2023-11-20. You can view files and clone it, but cannot push or open issues or pull requests.
foam-extend4.1-coherent-io/applications/solvers/heatTransfer/chtMultiRegionFoam/fluid/createFluidFields.H

153 lines
4.2 KiB
C
Raw Normal View History

// Initialise fluid field pointer lists
2010-08-26 14:22:03 +00:00
PtrList<basicPsiThermo> thermoFluid(fluidRegions.size());
PtrList<volScalarField> rhoFluid(fluidRegions.size());
2015-08-05 11:34:57 +00:00
PtrList<volScalarField> KappaFluid(fluidRegions.size());
2010-08-26 14:22:03 +00:00
PtrList<volVectorField> UFluid(fluidRegions.size());
PtrList<surfaceScalarField> phiFluid(fluidRegions.size());
PtrList<uniformDimensionedVectorField> gFluid(fluidRegions.size());
PtrList<compressible::turbulenceModel> turbulence(fluidRegions.size());
PtrList<volScalarField> DpDtFluid(fluidRegions.size());
2010-08-26 14:22:03 +00:00
List<scalar> initialMassFluid(fluidRegions.size());
// Populate fluid field pointer lists
forAll(fluidRegions, i)
{
Info<< "*** Reading fluid mesh thermophysical properties for region "
<< fluidRegions[i].name() << nl << endl;
2010-08-26 14:22:03 +00:00
Info<< " Adding to thermoFluid\n" << endl;
thermoFluid.set
(
i,
2010-08-26 14:22:03 +00:00
basicPsiThermo::New(fluidRegions[i]).ptr()
);
2010-08-26 14:22:03 +00:00
Info<< " Adding to rhoFluid\n" << endl;
rhoFluid.set
(
i,
new volScalarField
(
IOobject
(
"rho",
runTime.timeName(),
fluidRegions[i],
IOobject::NO_READ,
IOobject::AUTO_WRITE
),
2010-08-26 14:22:03 +00:00
thermoFluid[i].rho()
)
);
2015-08-05 11:34:57 +00:00
Info<< " Adding to KappaFluid\n" << endl;
KappaFluid.set
(
i,
new volScalarField
(
IOobject
(
2015-08-05 11:34:57 +00:00
"Kappa",
runTime.timeName(),
fluidRegions[i],
IOobject::NO_READ,
IOobject::NO_WRITE
),
2010-08-26 14:22:03 +00:00
thermoFluid[i].Cp()*thermoFluid[i].alpha()
)
);
2010-08-26 14:22:03 +00:00
Info<< " Adding to UFluid\n" << endl;
UFluid.set
(
i,
new volVectorField
(
IOobject
(
"U",
runTime.timeName(),
fluidRegions[i],
IOobject::MUST_READ,
IOobject::AUTO_WRITE
),
fluidRegions[i]
)
);
2010-08-26 14:22:03 +00:00
Info<< " Adding to phiFluid\n" << endl;
phiFluid.set
(
i,
new surfaceScalarField
(
IOobject
(
"phi",
runTime.timeName(),
fluidRegions[i],
IOobject::READ_IF_PRESENT,
IOobject::AUTO_WRITE
),
2010-08-26 14:22:03 +00:00
linearInterpolate(rhoFluid[i]*UFluid[i])
& fluidRegions[i].Sf()
)
);
2010-08-26 14:22:03 +00:00
Info<< " Adding to gFluid\n" << endl;
gFluid.set
(
i,
new uniformDimensionedVectorField
(
IOobject
(
"g",
runTime.constant(),
fluidRegions[i],
IOobject::MUST_READ,
IOobject::NO_WRITE
)
)
);
Info<< " Adding to turbulence\n" << endl;
turbulence.set
(
i,
2010-08-26 14:22:03 +00:00
autoPtr<compressible::turbulenceModel>
(
2010-08-26 14:22:03 +00:00
compressible::turbulenceModel::New
(
2010-08-26 14:22:03 +00:00
rhoFluid[i],
UFluid[i],
phiFluid[i],
thermoFluid[i]
)
).ptr()
);
2010-08-26 14:22:03 +00:00
Info<< " Adding to DpDtFluid\n" << endl;
DpDtFluid.set
(
i,
new volScalarField
(
2010-08-26 14:22:03 +00:00
"DpDt",
fvc::DDt
(
surfaceScalarField
(
"phiU",
2010-08-26 14:22:03 +00:00
phiFluid[i]/fvc::interpolate(rhoFluid[i])
),
2010-08-26 14:22:03 +00:00
thermoFluid[i].p()
)
)
);
2010-08-26 14:22:03 +00:00
initialMassFluid[i] = fvc::domainIntegrate(rhoFluid[i]).value();
}