FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "turbulentTemperatureCoupledBaffleMixedFvPatchScalarField.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <finiteVolume/fvPatchFieldMapper.H>
00029 #include <finiteVolume/volFields.H>
00030 #include <meshTools/directMappedPatchBase.H>
00031 #include <compressibleRASModels/regionProperties.H>
00032 #include <basicThermophysicalModels/basicThermo.H>
00033 #include <compressibleRASModels/RASModel.H>
00034 #include <OpenFOAM/mapDistribute.H>
00035 
00036 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00037 
00038 namespace Foam
00039 {
00040 namespace compressible
00041 {
00042 
00043 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00044 
00045 
00046 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00047 
00048 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::
00049 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
00050 (
00051     const fvPatch& p,
00052     const DimensionedField<scalar, volMesh>& iF
00053 )
00054 :
00055     mixedFvPatchScalarField(p, iF),
00056     neighbourFieldName_("undefined-neighbourFieldName"),
00057     KName_("undefined-K")
00058 {
00059     this->refValue() = 0.0;
00060     this->refGrad() = 0.0;
00061     this->valueFraction() = 1.0;
00062 }
00063 
00064 
00065 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::
00066 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
00067 (
00068     const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField& ptf,
00069     const fvPatch& p,
00070     const DimensionedField<scalar, volMesh>& iF,
00071     const fvPatchFieldMapper& mapper
00072 )
00073 :
00074     mixedFvPatchScalarField(ptf, p, iF, mapper),
00075     neighbourFieldName_(ptf.neighbourFieldName_),
00076     KName_(ptf.KName_)
00077 {}
00078 
00079 
00080 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::
00081 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
00082 (
00083     const fvPatch& p,
00084     const DimensionedField<scalar, volMesh>& iF,
00085     const dictionary& dict
00086 )
00087 :
00088     mixedFvPatchScalarField(p, iF),
00089     neighbourFieldName_(dict.lookup("neighbourFieldName")),
00090     KName_(dict.lookup("Kcond"))
00091 {
00092     if (!isA<directMappedPatchBase>(this->patch().patch()))
00093     {
00094         FatalErrorIn
00095         (
00096             "turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::"
00097             "turbulentTemperatureCoupledBaffleMixedFvPatchScalarField\n"
00098             "(\n"
00099             "    const fvPatch& p,\n"
00100             "    const DimensionedField<scalar, volMesh>& iF,\n"
00101             "    const dictionary& dict\n"
00102             ")\n"
00103         )   << "\n    patch type '" << p.type()
00104             << "' not type '" << directMappedPatchBase::typeName << "'"
00105             << "\n    for patch " << p.name()
00106             << " of field " << dimensionedInternalField().name()
00107             << " in file " << dimensionedInternalField().objectPath()
00108             << exit(FatalError);
00109     }
00110 
00111     fvPatchScalarField::operator=(scalarField("value", dict, p.size()));
00112 
00113     if (dict.found("refValue"))
00114     {
00115         // Full restart
00116         refValue() = scalarField("refValue", dict, p.size());
00117         refGrad() = scalarField("refGradient", dict, p.size());
00118         valueFraction() = scalarField("valueFraction", dict, p.size());
00119     }
00120     else
00121     {
00122         // Start from user entered data. Assume fixedValue.
00123         refValue() = *this;
00124         refGrad() = 0.0;
00125         valueFraction() = 1.0;
00126     }
00127 }
00128 
00129 
00130 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::
00131 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
00132 (
00133     const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField& wtcsf,
00134     const DimensionedField<scalar, volMesh>& iF
00135 )
00136 :
00137     mixedFvPatchScalarField(wtcsf, iF),
00138     neighbourFieldName_(wtcsf.neighbourFieldName_),
00139     KName_(wtcsf.KName_)
00140 {}
00141 
00142 
00143 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00144 
00145 tmp<scalarField>
00146 turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::K() const
00147 {
00148     const fvMesh& mesh = patch().boundaryMesh().mesh();
00149 
00150     if (KName_ == "none")
00151     {
00152         const compressible::RASModel& model =
00153             db().lookupObject<compressible::RASModel>("RASProperties");
00154 
00155         const basicThermo& thermo =
00156             db().lookupObject<basicThermo>("thermophysicalProperties");
00157 
00158         return
00159             model.alphaEff()().boundaryField()[patch().index()]
00160            *thermo.Cp()().boundaryField()[patch().index()];
00161     }
00162     else if (mesh.objectRegistry::foundObject<volScalarField>(KName_))
00163     {
00164         return patch().lookupPatchField<volScalarField, scalar>(KName_);
00165     }
00166     else if (mesh.objectRegistry::foundObject<volSymmTensorField>(KName_))
00167     {
00168         const symmTensorField& KWall =
00169             patch().lookupPatchField<volSymmTensorField, scalar>(KName_);
00170 
00171         vectorField n = patch().nf();
00172 
00173         return n & KWall & n;
00174     }
00175     else
00176     {
00177         FatalErrorIn
00178         (
00179             "turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::K()"
00180             " const"
00181         )   << "Did not find field " << KName_
00182             << " on mesh " << mesh.name() << " patch " << patch().name()
00183             << endl
00184             << "Please set 'K' to 'none', a valid volScalarField"
00185             << " or a valid volSymmTensorField." << exit(FatalError);
00186 
00187         return scalarField(0);
00188     }
00189 }
00190 
00191 
00192 void turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::updateCoeffs()
00193 {
00194     if (updated())
00195     {
00196         return;
00197     }
00198 
00199     // Get the coupling information from the directMappedPatchBase
00200     const directMappedPatchBase& mpp = refCast<const directMappedPatchBase>
00201     (
00202         patch().patch()
00203     );
00204     const polyMesh& nbrMesh = mpp.sampleMesh();
00205     const fvPatch& nbrPatch = refCast<const fvMesh>
00206     (
00207         nbrMesh
00208     ).boundary()[mpp.samplePolyPatch().index()];
00209 
00210     // Force recalculation of mapping and schedule
00211     const mapDistribute& distMap = mpp.map();
00212 
00213     tmp<scalarField> intFld = patchInternalField();
00214 
00215 
00216     // Calculate the temperature by harmonic averaging
00217     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
00218 
00219     const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField& nbrField =
00220     refCast
00221     <
00222         const turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
00223     >
00224     (
00225         nbrPatch.lookupPatchField<volScalarField, scalar>
00226         (
00227             neighbourFieldName_
00228         )
00229     );
00230 
00231     // Swap to obtain full local values of neighbour internal field
00232     scalarField nbrIntFld = nbrField.patchInternalField();
00233     mapDistribute::distribute
00234     (
00235         Pstream::defaultCommsType,
00236         distMap.schedule(),
00237         distMap.constructSize(),
00238         distMap.subMap(),           // what to send
00239         distMap.constructMap(),     // what to receive
00240         nbrIntFld
00241     );
00242 
00243     // Swap to obtain full local values of neighbour K*delta
00244     scalarField nbrKDelta = nbrField.K()*nbrPatch.deltaCoeffs();
00245     mapDistribute::distribute
00246     (
00247         Pstream::defaultCommsType,
00248         distMap.schedule(),
00249         distMap.constructSize(),
00250         distMap.subMap(),           // what to send
00251         distMap.constructMap(),     // what to receive
00252         nbrKDelta
00253     );
00254 
00255     tmp<scalarField> myKDelta = K()*patch().deltaCoeffs();
00256 
00257 
00258     // Both sides agree on
00259     // - temperature : (myKDelta*fld + nbrKDelta*nbrFld)/(myKDelta+nbrKDelta)
00260     // - gradient    : (temperature-fld)*delta
00261     // We've got a degree of freedom in how to implement this in a mixed bc.
00262     // (what gradient, what fixedValue and mixing coefficient)
00263     // Two reasonable choices:
00264     // 1. specify above temperature on one side (preferentially the high side)
00265     //    and above gradient on the other. So this will switch between pure
00266     //    fixedvalue and pure fixedgradient
00267     // 2. specify gradient and temperature such that the equations are the
00268     //    same on both sides. This leads to the choice of
00269     //    - refGradient = zero gradient
00270     //    - refValue = neighbour value
00271     //    - mixFraction = nbrKDelta / (nbrKDelta + myKDelta())
00272 
00273 
00274     this->refValue() = nbrIntFld;
00275 
00276     this->refGrad() = 0.0;
00277 
00278     this->valueFraction() = nbrKDelta / (nbrKDelta + myKDelta());
00279 
00280     mixedFvPatchScalarField::updateCoeffs();
00281 
00282 
00283     if (debug)
00284     {
00285         scalar Q = gSum(K()*patch().magSf()*snGrad());
00286 
00287         Info<< patch().boundaryMesh().mesh().name() << ':'
00288             << patch().name() << ':'
00289             << this->dimensionedInternalField().name() << " <- "
00290             << nbrMesh.name() << ':'
00291             << nbrPatch.name() << ':'
00292             << this->dimensionedInternalField().name() << " :"
00293             << " heat[W]:" << Q
00294             << " walltemperature "
00295             << " min:" << gMin(*this)
00296             << " max:" << gMax(*this)
00297             << " avg:" << gAverage(*this)
00298             << endl;
00299     }
00300 }
00301 
00302 
00303 void turbulentTemperatureCoupledBaffleMixedFvPatchScalarField::write
00304 (
00305     Ostream& os
00306 ) const
00307 {
00308     mixedFvPatchScalarField::write(os);
00309     os.writeKeyword("neighbourFieldName")<< neighbourFieldName_
00310         << token::END_STATEMENT << nl;
00311     os.writeKeyword("Kcond") << KName_ << token::END_STATEMENT << nl;
00312 }
00313 
00314 
00315 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00316 
00317 makePatchTypeField
00318 (
00319     fvPatchScalarField,
00320     turbulentTemperatureCoupledBaffleMixedFvPatchScalarField
00321 );
00322 
00323 
00324 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00325 
00326 } // End namespace compressible
00327 } // End namespace Foam
00328 
00329 
00330 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines