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 "solidWallMixedTemperatureCoupledFvPatchScalarField.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <finiteVolume/fvPatchFieldMapper.H>
00029 #include <finiteVolume/volFields.H>
00030 #include <meshTools/directMappedPatchBase.H>
00031 #include <OpenFOAM/mapDistribute.H>
00032 #include <compressibleRASModels/regionProperties.H>
00033
00034
00035
00036 Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::
00037 solidWallMixedTemperatureCoupledFvPatchScalarField
00038 (
00039 const fvPatch& p,
00040 const DimensionedField<scalar, volMesh>& iF
00041 )
00042 :
00043 mixedFvPatchScalarField(p, iF),
00044 neighbourFieldName_("undefined-neighbourFieldName"),
00045 KName_("undefined-K")
00046 {
00047 this->refValue() = 0.0;
00048 this->refGrad() = 0.0;
00049 this->valueFraction() = 1.0;
00050 }
00051
00052
00053 Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::
00054 solidWallMixedTemperatureCoupledFvPatchScalarField
00055 (
00056 const solidWallMixedTemperatureCoupledFvPatchScalarField& ptf,
00057 const fvPatch& p,
00058 const DimensionedField<scalar, volMesh>& iF,
00059 const fvPatchFieldMapper& mapper
00060 )
00061 :
00062 mixedFvPatchScalarField(ptf, p, iF, mapper),
00063 neighbourFieldName_(ptf.neighbourFieldName_),
00064 KName_(ptf.KName_)
00065 {}
00066
00067
00068 Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::
00069 solidWallMixedTemperatureCoupledFvPatchScalarField
00070 (
00071 const fvPatch& p,
00072 const DimensionedField<scalar, volMesh>& iF,
00073 const dictionary& dict
00074 )
00075 :
00076 mixedFvPatchScalarField(p, iF),
00077 neighbourFieldName_(dict.lookup("neighbourFieldName")),
00078 KName_(dict.lookup("Kcond"))
00079 {
00080 if (!isA<directMappedPatchBase>(this->patch().patch()))
00081 {
00082 FatalErrorIn
00083 (
00084 "solidWallMixedTemperatureCoupledFvPatchScalarField::"
00085 "solidWallMixedTemperatureCoupledFvPatchScalarField\n"
00086 "(\n"
00087 " const fvPatch& p,\n"
00088 " const DimensionedField<scalar, volMesh>& iF,\n"
00089 " const dictionary& dict\n"
00090 ")\n"
00091 ) << "\n patch type '" << p.type()
00092 << "' not type '" << directMappedPatchBase::typeName << "'"
00093 << "\n for patch " << p.name()
00094 << " of field " << dimensionedInternalField().name()
00095 << " in file " << dimensionedInternalField().objectPath()
00096 << exit(FatalError);
00097 }
00098
00099 fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
00100
00101 if (dict.found("refValue"))
00102 {
00103
00104 refValue() = scalarField("refValue", dict, p.size());
00105 refGrad() = scalarField("refGradient", dict, p.size());
00106 valueFraction() = scalarField("valueFraction", dict, p.size());
00107 }
00108 else
00109 {
00110
00111 refValue() = *this;
00112 refGrad() = 0.0;
00113 valueFraction() = 1.0;
00114 }
00115 }
00116
00117
00118 Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::
00119 solidWallMixedTemperatureCoupledFvPatchScalarField
00120 (
00121 const solidWallMixedTemperatureCoupledFvPatchScalarField& wtcsf,
00122 const DimensionedField<scalar, volMesh>& iF
00123 )
00124 :
00125 mixedFvPatchScalarField(wtcsf, iF),
00126 neighbourFieldName_(wtcsf.neighbourFieldName_),
00127 KName_(wtcsf.KName_)
00128 {}
00129
00130
00131
00132
00133 const Foam::fvPatchScalarField&
00134 Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::K() const
00135 {
00136 return this->patch().lookupPatchField<volScalarField, scalar>(KName_);
00137 }
00138
00139
00140 void Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::updateCoeffs()
00141 {
00142 if (updated())
00143 {
00144 return;
00145 }
00146
00147
00148 const directMappedPatchBase& mpp = refCast<const directMappedPatchBase>
00149 (
00150 patch().patch()
00151 );
00152 const polyMesh& nbrMesh = mpp.sampleMesh();
00153 const fvPatch& nbrPatch = refCast<const fvMesh>
00154 (
00155 nbrMesh
00156 ).boundary()[mpp.samplePolyPatch().index()];
00157
00158
00159 const mapDistribute& distMap = mpp.map();
00160
00161 tmp<scalarField> intFld = patchInternalField();
00162
00163
00164 const solidWallMixedTemperatureCoupledFvPatchScalarField& nbrField =
00165 refCast<const solidWallMixedTemperatureCoupledFvPatchScalarField>
00166 (
00167 nbrPatch.lookupPatchField<volScalarField, scalar>
00168 (
00169 neighbourFieldName_
00170 )
00171 );
00172
00173
00174 scalarField nbrIntFld = nbrField.patchInternalField();
00175 mapDistribute::distribute
00176 (
00177 Pstream::defaultCommsType,
00178 distMap.schedule(),
00179 distMap.constructSize(),
00180 distMap.subMap(),
00181 distMap.constructMap(),
00182 nbrIntFld
00183 );
00184
00185
00186 scalarField nbrKDelta = nbrField.K()*nbrPatch.deltaCoeffs();
00187 mapDistribute::distribute
00188 (
00189 Pstream::defaultCommsType,
00190 distMap.schedule(),
00191 distMap.constructSize(),
00192 distMap.subMap(),
00193 distMap.constructMap(),
00194 nbrKDelta
00195 );
00196
00197 tmp<scalarField> myKDelta = K()*patch().deltaCoeffs();
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216 this->refValue() = nbrIntFld;
00217
00218 this->refGrad() = 0.0;
00219
00220 this->valueFraction() = nbrKDelta / (nbrKDelta + myKDelta());
00221
00222 mixedFvPatchScalarField::updateCoeffs();
00223
00224
00225 if (debug)
00226 {
00227 scalar Q = gSum(K()*patch().magSf()*snGrad());
00228
00229 Info<< patch().boundaryMesh().mesh().name() << ':'
00230 << patch().name() << ':'
00231 << this->dimensionedInternalField().name() << " <- "
00232 << nbrMesh.name() << ':'
00233 << nbrPatch.name() << ':'
00234 << this->dimensionedInternalField().name() << " :"
00235 << " heat[W]:" << Q
00236 << " walltemperature "
00237 << " min:" << gMin(*this)
00238 << " max:" << gMax(*this)
00239 << " avg:" << gAverage(*this)
00240 << endl;
00241 }
00242 }
00243
00244
00245 void Foam::solidWallMixedTemperatureCoupledFvPatchScalarField::write
00246 (
00247 Ostream& os
00248 ) const
00249 {
00250 mixedFvPatchScalarField::write(os);
00251 os.writeKeyword("neighbourFieldName")<< neighbourFieldName_
00252 << token::END_STATEMENT << nl;
00253 os.writeKeyword("Kcond") << KName_ << token::END_STATEMENT << nl;
00254 }
00255
00256
00257
00258
00259 namespace Foam
00260 {
00261
00262 makePatchTypeField
00263 (
00264 fvPatchScalarField,
00265 solidWallMixedTemperatureCoupledFvPatchScalarField
00266 );
00267
00268 }
00269
00270