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 "solidWallHeatFluxTemperatureFvPatchScalarField.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <finiteVolume/fvPatchFieldMapper.H>
00029 #include <finiteVolume/volFields.H>
00030
00031
00032
00033 Foam::solidWallHeatFluxTemperatureFvPatchScalarField::
00034 solidWallHeatFluxTemperatureFvPatchScalarField
00035 (
00036 const fvPatch& p,
00037 const DimensionedField<scalar, volMesh>& iF
00038 )
00039 :
00040 fixedGradientFvPatchScalarField(p, iF),
00041 q_(p.size(), 0.0),
00042 KName_("undefined-K")
00043 {}
00044
00045
00046 Foam::solidWallHeatFluxTemperatureFvPatchScalarField::
00047 solidWallHeatFluxTemperatureFvPatchScalarField
00048 (
00049 const solidWallHeatFluxTemperatureFvPatchScalarField& ptf,
00050 const fvPatch& p,
00051 const DimensionedField<scalar, volMesh>& iF,
00052 const fvPatchFieldMapper& mapper
00053 )
00054 :
00055 fixedGradientFvPatchScalarField(ptf, p, iF, mapper),
00056 q_(ptf.q_, mapper),
00057 KName_(ptf.KName_)
00058 {}
00059
00060
00061 Foam::solidWallHeatFluxTemperatureFvPatchScalarField::
00062 solidWallHeatFluxTemperatureFvPatchScalarField
00063 (
00064 const fvPatch& p,
00065 const DimensionedField<scalar, volMesh>& iF,
00066 const dictionary& dict
00067 )
00068 :
00069 fixedGradientFvPatchScalarField(p, iF, dict),
00070 q_("q", dict, p.size()),
00071 KName_(dict.lookup("Kcond"))
00072 {}
00073
00074
00075 Foam::solidWallHeatFluxTemperatureFvPatchScalarField::
00076 solidWallHeatFluxTemperatureFvPatchScalarField
00077 (
00078 const solidWallHeatFluxTemperatureFvPatchScalarField& tppsf
00079 )
00080 :
00081 fixedGradientFvPatchScalarField(tppsf),
00082 q_(tppsf.q_),
00083 KName_(tppsf.KName_)
00084 {}
00085
00086
00087 Foam::solidWallHeatFluxTemperatureFvPatchScalarField::
00088 solidWallHeatFluxTemperatureFvPatchScalarField
00089 (
00090 const solidWallHeatFluxTemperatureFvPatchScalarField& tppsf,
00091 const DimensionedField<scalar, volMesh>& iF
00092 )
00093 :
00094 fixedGradientFvPatchScalarField(tppsf, iF),
00095 q_(tppsf.q_),
00096 KName_(tppsf.KName_)
00097 {}
00098
00099
00100
00101
00102 void Foam::solidWallHeatFluxTemperatureFvPatchScalarField::autoMap
00103 (
00104 const fvPatchFieldMapper& m
00105 )
00106 {
00107 fixedGradientFvPatchScalarField::autoMap(m);
00108 q_.autoMap(m);
00109 }
00110
00111
00112 void Foam::solidWallHeatFluxTemperatureFvPatchScalarField::rmap
00113 (
00114 const fvPatchScalarField& ptf,
00115 const labelList& addr
00116 )
00117 {
00118 fixedGradientFvPatchScalarField::rmap(ptf, addr);
00119
00120 const solidWallHeatFluxTemperatureFvPatchScalarField& hfptf =
00121 refCast<const solidWallHeatFluxTemperatureFvPatchScalarField>(ptf);
00122
00123 q_.rmap(hfptf.q_, addr);
00124 }
00125
00126
00127 Foam::tmp<Foam::scalarField>
00128 Foam::solidWallHeatFluxTemperatureFvPatchScalarField::K() const
00129 {
00130 const fvMesh& mesh = patch().boundaryMesh().mesh();
00131
00132 if (mesh.objectRegistry::foundObject<volScalarField>(KName_))
00133 {
00134 return patch().lookupPatchField<volScalarField, scalar>(KName_);
00135 }
00136 else if (mesh.objectRegistry::foundObject<volSymmTensorField>(KName_))
00137 {
00138 const symmTensorField& KWall =
00139 patch().lookupPatchField<volSymmTensorField, scalar>(KName_);
00140
00141 vectorField n = patch().nf();
00142
00143 return n & KWall & n;
00144 }
00145 else
00146 {
00147 FatalErrorIn
00148 (
00149 "solidWallHeatFluxTemperatureFvPatchScalarField::K()"
00150 " const"
00151 ) << "Did not find field " << KName_
00152 << " on mesh " << mesh.name() << " patch " << patch().name()
00153 << endl
00154 << "Please set 'K' to a valid volScalarField"
00155 << " or a valid volSymmTensorField." << exit(FatalError);
00156
00157 return scalarField(0);
00158 }
00159 }
00160
00161
00162 void Foam::solidWallHeatFluxTemperatureFvPatchScalarField::updateCoeffs()
00163 {
00164 if (updated())
00165 {
00166 return;
00167 }
00168
00169 gradient() = q_/K();
00170
00171 fixedGradientFvPatchScalarField::updateCoeffs();
00172
00173 if (debug)
00174 {
00175 scalar Q = gSum(K()*patch().magSf()*snGrad());
00176
00177 Info<< patch().boundaryMesh().mesh().name() << ':'
00178 << patch().name() << ':'
00179 << this->dimensionedInternalField().name() << " :"
00180 << " heatFlux:" << Q
00181 << " walltemperature "
00182 << " min:" << gMin(*this)
00183 << " max:" << gMax(*this)
00184 << " avg:" << gAverage(*this)
00185 << endl;
00186 }
00187 }
00188
00189
00190 void Foam::solidWallHeatFluxTemperatureFvPatchScalarField::write
00191 (
00192 Ostream& os
00193 ) const
00194 {
00195 fixedGradientFvPatchScalarField::write(os);
00196 q_.writeEntry("q", os);
00197 os.writeKeyword("Kcond") << KName_ << token::END_STATEMENT << nl;
00198 this->writeEntry("value", os);
00199 }
00200
00201
00202
00203
00204 namespace Foam
00205 {
00206 makePatchTypeField
00207 (
00208 fvPatchScalarField,
00209 solidWallHeatFluxTemperatureFvPatchScalarField
00210 );
00211 }
00212
00213