wallHeatTransferFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
28 #include "fvPatchFieldMapper.H"
31 
32 namespace Foam
33 {
34 
35 namespace incompressible
36 {
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
41 (
42  const fvPatch& p,
44 )
45 :
46  mixedFvPatchScalarField(p, iF),
47  Tinf_(p.size(), 0.0),
48  alphaWall_(p.size(), 0.0)
49 {
50  refValue() = 0.0;
51  refGrad() = 0.0;
52  valueFraction() = 0.0;
53 }
54 
55 
57 (
59  const fvPatch& p,
61  const fvPatchFieldMapper& mapper
62 )
63 :
64  mixedFvPatchScalarField(ptf, p, iF, mapper),
65  Tinf_(ptf.Tinf_, mapper),
66  alphaWall_(ptf.alphaWall_, mapper)
67 {}
68 
69 
71 (
72  const fvPatch& p,
74  const dictionary& dict
75 )
76 :
77  mixedFvPatchScalarField(p, iF),
78  Tinf_("Tinf", dict, p.size()),
79  alphaWall_("alphaWall", dict, p.size())
80 {
81  refValue() = Tinf_;
82  refGrad() = 0.0;
83  valueFraction() = 0.0;
84 
85  if (dict.found("value"))
86  {
88  (
89  scalarField("value", dict, p.size())
90  );
91  }
92  else
93  {
94  evaluate();
95  }
96 }
97 
98 
100 (
102 )
103 :
104  mixedFvPatchScalarField(tppsf),
105  Tinf_(tppsf.Tinf_),
106  alphaWall_(tppsf.alphaWall_)
107 {}
108 
109 
111 (
114 )
115 :
116  mixedFvPatchScalarField(tppsf, iF),
117  Tinf_(tppsf.Tinf_),
118  alphaWall_(tppsf.alphaWall_)
119 {}
120 
121 
122 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
123 
125 (
126  const fvPatchFieldMapper& m
127 )
128 {
130  Tinf_.autoMap(m);
131  alphaWall_.autoMap(m);
132 }
133 
134 
136 (
137  const fvPatchScalarField& ptf,
138  const labelList& addr
139 )
140 {
141  mixedFvPatchScalarField::rmap(ptf, addr);
142 
144  refCast<const wallHeatTransferFvPatchScalarField>(ptf);
145 
146  Tinf_.rmap(tiptf.Tinf_, addr);
147  alphaWall_.rmap(tiptf.alphaWall_, addr);
148 }
149 
150 
152 {
153  if (updated())
154  {
155  return;
156  }
157 
158  const scalarField& alphat =
159  patch().lookupPatchField<volScalarField, scalar>("alphat");
160 
161  typedef incompressible::turbulenceModel iTModel;
162  const iTModel& tModel =
163  db().lookupObject<iTModel>
164  (
166  (
169  )
170  );
171 
172  // need kappaEff = kappa(laminar) + kappat, where kappal = nu/Pr
173  // See applications/solvers/heatTransfer/buoyantBoussinesqSimpleFoam/TEqn.H
174  const label patchI = patch().index();
175  const scalarField& nu = tModel.nu()->boundaryField()[patchI];
176  const singlePhaseTransportModel & transport =
177  static_cast<const singlePhaseTransportModel &>(tModel.transport());
178  const scalar
179  Pr(dimensionedScalar(transport.lookup("Pr")).value());
180  const scalar Cp0(dimensionedScalar(transport.lookup("Cp")).value());
181  const scalar rho(dimensionedScalar(transport.lookup("rho")).value());
182 
183  const scalarField kappaEff(rho * Cp0 * (nu/Pr + alphat));
184 
185  valueFraction() =
186  1.0/
187  (
188  1.0
189  + kappaEff*patch().deltaCoeffs()/alphaWall_
190  );
191 
192  mixedFvPatchScalarField::updateCoeffs();
193 }
194 
195 
197 {
199  Tinf_.writeEntry("Tinf", os);
200  alphaWall_.writeEntry("alphaWall", os);
201  writeEntry("value", os);
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
206 
208 (
211  );
212 
213 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
214 
215 } // End namespace incompressible
216 } // End namespace Foam
217 
218 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Foam::fvPatchField::write
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:349
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::incompressible::wallHeatTransferFvPatchScalarField
This boundary condition provides an enthalpy condition for wall heat transfer.
Definition: wallHeatTransferFvPatchScalarField.H:87
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
p
p
Definition: pEqn.H:62
Foam::incompressible::wallHeatTransferFvPatchScalarField::Tinf_
scalarField Tinf_
Temperature at the wall.
Definition: wallHeatTransferFvPatchScalarField.H:94
Foam::incompressible::wallHeatTransferFvPatchScalarField::wallHeatTransferFvPatchScalarField
wallHeatTransferFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: wallHeatTransferFvPatchScalarField.C:41
wallHeatTransferFvPatchScalarField.H
Foam::Field::autoMap
void autoMap(const FieldMapper &map, const bool applyFlip=true)
Map from self.
Definition: Field.C:508
singlePhaseTransportModel.H
dimensionedInternalField
rDeltaT dimensionedInternalField()
alphat
alphat
Definition: TEqn.H:2
turbulentTransportModel.H
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
Foam::constant::atomic::group
const char *const group
Group name for atomic constants.
Definition: atomicConstants.C:40
fvPatchFieldMapper.H
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::singlePhaseTransportModel
A simple single-phase transport model based on viscosityModel.
Definition: singlePhaseTransportModel.H:55
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::incompressible::wallHeatTransferFvPatchScalarField::rmap
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Definition: wallHeatTransferFvPatchScalarField.C:136
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::incompressible::wallHeatTransferFvPatchScalarField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: wallHeatTransferFvPatchScalarField.C:151
Pr
dimensionedScalar Pr("Pr", dimless, laminarTransport)
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::incompressible::wallHeatTransferFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: wallHeatTransferFvPatchScalarField.C:196
Foam::Field::writeEntry
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:700
rho
rho
Definition: pEqn.H:3
Foam::incompressible::wallHeatTransferFvPatchScalarField::alphaWall_
scalarField alphaWall_
Thermal diffusivity at the wall.
Definition: wallHeatTransferFvPatchScalarField.H:97
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::incompressible::wallHeatTransferFvPatchScalarField::autoMap
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
Definition: wallHeatTransferFvPatchScalarField.C:125
Foam::IncompressibleTurbulenceModel
Templated abstract base class for single-phase incompressible turbulence models.
Definition: IncompressibleTurbulenceModel.H:52
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:51
Foam::incompressible::makePatchTypeField
makePatchTypeField(fvPatchScalarField, turbulentBoundaryCoupledFvPatchScalarField)