00001
00002 PtrList<basicRhoThermo> thermoFluid(fluidRegions.size());
00003 PtrList<volScalarField> rhoFluid(fluidRegions.size());
00004 PtrList<volScalarField> KFluid(fluidRegions.size());
00005 PtrList<volVectorField> UFluid(fluidRegions.size());
00006 PtrList<surfaceScalarField> phiFluid(fluidRegions.size());
00007 PtrList<uniformDimensionedVectorField> gFluid(fluidRegions.size());
00008 PtrList<compressible::turbulenceModel> turbulence(fluidRegions.size());
00009 PtrList<volScalarField> p_rghFluid(fluidRegions.size());
00010 PtrList<volScalarField> ghFluid(fluidRegions.size());
00011 PtrList<surfaceScalarField> ghfFluid(fluidRegions.size());
00012 PtrList<volScalarField> DpDtFluid(fluidRegions.size());
00013
00014 List<scalar> initialMassFluid(fluidRegions.size());
00015
00016
00017 forAll(fluidRegions, i)
00018 {
00019 Info<< "*** Reading fluid mesh thermophysical properties for region "
00020 << fluidRegions[i].name() << nl << endl;
00021
00022 Info<< " Adding to thermoFluid\n" << endl;
00023 thermoFluid.set
00024 (
00025 i,
00026 basicRhoThermo::New(fluidRegions[i]).ptr()
00027 );
00028
00029 Info<< " Adding to rhoFluid\n" << endl;
00030 rhoFluid.set
00031 (
00032 i,
00033 new volScalarField
00034 (
00035 IOobject
00036 (
00037 "rho",
00038 runTime.timeName(),
00039 fluidRegions[i],
00040 IOobject::NO_READ,
00041 IOobject::AUTO_WRITE
00042 ),
00043 thermoFluid[i].rho()
00044 )
00045 );
00046
00047 Info<< " Adding to KFluid\n" << endl;
00048 KFluid.set
00049 (
00050 i,
00051 new volScalarField
00052 (
00053 IOobject
00054 (
00055 "Kcond",
00056 runTime.timeName(),
00057 fluidRegions[i],
00058 IOobject::NO_READ,
00059 IOobject::NO_WRITE
00060 ),
00061 thermoFluid[i].Cp()*thermoFluid[i].alpha()
00062 )
00063 );
00064
00065 Info<< " Adding to UFluid\n" << endl;
00066 UFluid.set
00067 (
00068 i,
00069 new volVectorField
00070 (
00071 IOobject
00072 (
00073 "U",
00074 runTime.timeName(),
00075 fluidRegions[i],
00076 IOobject::MUST_READ,
00077 IOobject::AUTO_WRITE
00078 ),
00079 fluidRegions[i]
00080 )
00081 );
00082
00083 Info<< " Adding to phiFluid\n" << endl;
00084 phiFluid.set
00085 (
00086 i,
00087 new surfaceScalarField
00088 (
00089 IOobject
00090 (
00091 "phi",
00092 runTime.timeName(),
00093 fluidRegions[i],
00094 IOobject::READ_IF_PRESENT,
00095 IOobject::AUTO_WRITE
00096 ),
00097 linearInterpolate(rhoFluid[i]*UFluid[i])
00098 & fluidRegions[i].Sf()
00099 )
00100 );
00101
00102 Info<< " Adding to gFluid\n" << endl;
00103 gFluid.set
00104 (
00105 i,
00106 new uniformDimensionedVectorField
00107 (
00108 IOobject
00109 (
00110 "g",
00111 runTime.constant(),
00112 fluidRegions[i],
00113 IOobject::MUST_READ,
00114 IOobject::NO_WRITE
00115 )
00116 )
00117 );
00118
00119 Info<< " Adding to turbulence\n" << endl;
00120 turbulence.set
00121 (
00122 i,
00123 autoPtr<compressible::turbulenceModel>
00124 (
00125 compressible::turbulenceModel::New
00126 (
00127 rhoFluid[i],
00128 UFluid[i],
00129 phiFluid[i],
00130 thermoFluid[i]
00131 )
00132 ).ptr()
00133 );
00134
00135 Info<< " Adding to ghFluid\n" << endl;
00136 ghFluid.set
00137 (
00138 i,
00139 new volScalarField("gh", gFluid[i] & fluidRegions[i].C())
00140 );
00141
00142 Info<< " Adding to ghfFluid\n" << endl;
00143 ghfFluid.set
00144 (
00145 i,
00146 new surfaceScalarField("ghf", gFluid[i] & fluidRegions[i].Cf())
00147 );
00148
00149 p_rghFluid.set
00150 (
00151 i,
00152 new volScalarField
00153 (
00154 IOobject
00155 (
00156 "p_rgh",
00157 runTime.timeName(),
00158 fluidRegions[i],
00159 IOobject::MUST_READ,
00160 IOobject::AUTO_WRITE
00161 ),
00162 fluidRegions[i]
00163 )
00164 );
00165
00166
00167 p_rghFluid[i] = thermoFluid[i].p() - rhoFluid[i]*ghFluid[i];
00168
00169 initialMassFluid[i] = fvc::domainIntegrate(rhoFluid[i]).value();
00170
00171 Info<< " Adding to DpDtFluid\n" << endl;
00172 DpDtFluid.set
00173 (
00174 i,
00175 new volScalarField
00176 (
00177 "DpDt",
00178 fvc::DDt
00179 (
00180 surfaceScalarField
00181 (
00182 "phiU",
00183 phiFluid[i]/fvc::interpolate(rhoFluid[i])
00184 ),
00185 thermoFluid[i].p()
00186 )
00187 )
00188 );
00189 }