Go to the documentation of this file.00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026 #include "cloudAbsorptionEmission.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028
00029 #include <lagrangianIntermediate/thermoCloud.H>
00030
00031
00032
00033 namespace Foam
00034 {
00035 namespace radiation
00036 {
00037 defineTypeNameAndDebug(cloudAbsorptionEmission, 0);
00038
00039 addToRunTimeSelectionTable
00040 (
00041 absorptionEmissionModel,
00042 cloudAbsorptionEmission,
00043 dictionary
00044 );
00045 }
00046 }
00047
00048
00049
00050
00051 Foam::radiation::cloudAbsorptionEmission::cloudAbsorptionEmission
00052 (
00053 const dictionary& dict,
00054 const fvMesh& mesh
00055 )
00056 :
00057 absorptionEmissionModel(dict, mesh),
00058 coeffsDict_(dict.subDict(typeName + "Coeffs")),
00059 cloudNames_(coeffsDict_.lookup("cloudNames"))
00060 {}
00061
00062
00063
00064
00065 Foam::radiation::cloudAbsorptionEmission::~cloudAbsorptionEmission()
00066 {}
00067
00068
00069
00070
00071 Foam::tmp<Foam::volScalarField>
00072 Foam::radiation::cloudAbsorptionEmission::aDisp(const label) const
00073 {
00074 tmp<volScalarField> ta
00075 (
00076 new volScalarField
00077 (
00078 IOobject
00079 (
00080 "a",
00081 mesh_.time().timeName(),
00082 mesh_,
00083 IOobject::NO_READ,
00084 IOobject::NO_WRITE,
00085 false
00086 ),
00087 mesh_,
00088 dimensionedScalar("a", dimless/dimLength, 0.0)
00089 )
00090 );
00091
00092 forAll(cloudNames_, i)
00093 {
00094 const thermoCloud& tc
00095 (
00096 mesh_.objectRegistry::lookupObject<thermoCloud>(cloudNames_[i])
00097 );
00098
00099 ta() += tc.ap();
00100 }
00101
00102 return ta;
00103 }
00104
00105
00106 Foam::tmp<Foam::volScalarField>
00107 Foam::radiation::cloudAbsorptionEmission::eDisp(const label bandI) const
00108 {
00109 tmp<volScalarField> te
00110 (
00111 new volScalarField
00112 (
00113 IOobject
00114 (
00115 "e",
00116 mesh_.time().timeName(),
00117 mesh_,
00118 IOobject::NO_READ,
00119 IOobject::NO_WRITE,
00120 false
00121 ),
00122 mesh_,
00123 dimensionedScalar("e", dimless/dimLength, 0.0)
00124 )
00125 );
00126
00127 return te;
00128 }
00129
00130
00131 Foam::tmp<Foam::volScalarField>
00132 Foam::radiation::cloudAbsorptionEmission::EDisp(const label bandI) const
00133 {
00134 tmp<volScalarField> tE
00135 (
00136 new volScalarField
00137 (
00138 IOobject
00139 (
00140 "E",
00141 mesh_.time().timeName(),
00142 mesh_,
00143 IOobject::NO_READ,
00144 IOobject::NO_WRITE,
00145 false
00146 ),
00147 mesh_,
00148 dimensionedScalar("E", dimMass/dimLength/pow3(dimTime), 0.0)
00149 )
00150 );
00151
00152 forAll(cloudNames_, i)
00153 {
00154 const thermoCloud& tc
00155 (
00156 mesh_.objectRegistry::lookupObject<thermoCloud>(cloudNames_[i])
00157 );
00158
00159 tE() += tc.Ep();
00160 }
00161
00162 return tE;
00163 }
00164
00165
00166