turbulentHeatFluxTemperatureFvPatchScalarField.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"
29 #include "volFields.H"
30 
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36  // declare specialization within 'Foam' namespace
37  template<>
38  const char* NamedEnum
39  <
42  2
43  >::names[] =
44  {
45  "power",
46  "flux"
47  };
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 
53 
54 namespace Foam
55 {
56 
57 namespace compressible
58 {
59 
60 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
61 
62 const NamedEnum
63 <
65  2
66 > turbulentHeatFluxTemperatureFvPatchScalarField::heatSourceTypeNames_;
67 
68 
69 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
70 
73 (
74  const fvPatch& p,
76 )
77 :
78  fixedGradientFvPatchScalarField(p, iF),
79  temperatureCoupledBase(patch(), "undefined", "undefined", "undefined-K"),
80  heatSource_(hsPower),
81  q_(p.size(), 0.0),
82  QrName_("undefinedQr")
83 {}
84 
85 
88 (
90  const fvPatch& p,
92  const fvPatchFieldMapper& mapper
93 )
94 :
95  fixedGradientFvPatchScalarField(ptf, p, iF, mapper),
96  temperatureCoupledBase(patch(), ptf),
97  heatSource_(ptf.heatSource_),
98  q_(ptf.q_, mapper),
99  QrName_(ptf.QrName_)
100 {}
101 
102 
105 (
106  const fvPatch& p,
108  const dictionary& dict
109 )
110 :
111  fixedGradientFvPatchScalarField(p, iF),
112  temperatureCoupledBase(patch(), dict),
113  heatSource_(heatSourceTypeNames_.read(dict.lookup("heatSource"))),
114  q_("q", dict, p.size()),
115  QrName_(dict.lookupOrDefault<word>("Qr", "none"))
116 {
117  if (dict.found("value") && dict.found("gradient"))
118  {
120  gradient() = Field<scalar>("gradient", dict, p.size());
121  }
122  else
123  {
124  // Still reading so cannot yet evaluate. Make up a value.
125  fvPatchField<scalar>::operator=(patchInternalField());
126  gradient() = 0.0;
127  }
128 }
129 
130 
133 (
135 )
136 :
137  fixedGradientFvPatchScalarField(thftpsf),
138  temperatureCoupledBase(patch(), thftpsf),
139  heatSource_(thftpsf.heatSource_),
140  q_(thftpsf.q_),
141  QrName_(thftpsf.QrName_)
142 {}
143 
144 
147 (
150 )
151 :
152  fixedGradientFvPatchScalarField(thftpsf, iF),
153  temperatureCoupledBase(patch(), thftpsf),
154  heatSource_(thftpsf.heatSource_),
155  q_(thftpsf.q_),
156  QrName_(thftpsf.QrName_)
157 {}
158 
159 
160 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
161 
163 (
164  const fvPatchFieldMapper& m
165 )
166 {
167  fixedGradientFvPatchScalarField::autoMap(m);
168  q_.autoMap(m);
169 }
170 
171 
173 (
174  const fvPatchScalarField& ptf,
175  const labelList& addr
176 )
177 {
178  fixedGradientFvPatchScalarField::rmap(ptf, addr);
179 
181  refCast<const turbulentHeatFluxTemperatureFvPatchScalarField>
182  (
183  ptf
184  );
185 
186  q_.rmap(thftptf.q_, addr);
187 }
188 
189 
191 {
192  if (updated())
193  {
194  return;
195  }
196 
197  const scalarField& Tp = *this;
198 
199  scalarField qr(this->size(), 0.0);
200 
201  //- Qr is negative going into the domain
202  if (QrName_ != "none")
203  {
204  qr = patch().lookupPatchField<volScalarField, scalar>(QrName_);
205  }
206 
207  switch (heatSource_)
208  {
209  case hsPower:
210  {
211  const scalar Ap = gSum(patch().magSf());
212  gradient() = (q_/Ap + qr)/kappa(Tp);
213  break;
214  }
215  case hsFlux:
216  {
217  gradient() = (q_ + qr)/kappa(Tp);
218  break;
219  }
220  default:
221  {
223  << "Unknown heat source type. Valid types are: "
224  << heatSourceTypeNames_ << nl << exit(FatalError);
225  }
226  }
227 
228  fixedGradientFvPatchScalarField::updateCoeffs();
229 }
230 
231 
233 (
234  Ostream& os
235 ) const
236 {
238  os.writeKeyword("heatSource") << heatSourceTypeNames_[heatSource_]
239  << token::END_STATEMENT << nl;
241  q_.writeEntry("q", os);
242  os.writeKeyword("Qr")<< QrName_ << token::END_STATEMENT << nl;
243  writeEntry("value", os);
244 }
245 
246 
247 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
248 
250 (
253 );
254 
255 
256 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
257 
258 } // End namespace compressible
259 } // End namespace Foam
260 
261 
262 // ************************************************************************* //
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
dictionary dict
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
void write(Ostream &) const
Write.
Fixed heat boundary condition to specify temperature gradient. Input heat source either specified in ...
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
void writeEntry(const word &keyword, Ostream &os) const
Write the field as a dictionary entry.
Definition: Field.C:719
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
virtual void rmap(const fvPatchScalarField &, const labelList &)
Reverse map the given fvPatchField onto this fvPatchField.
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Type gSum(const FieldField< Field, Type > &f)
A class for handling words, derived from string.
Definition: word.H:59
Foam::fvPatchFieldMapper.
virtual label size() const
Return size.
Definition: fvPatch.H:161
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
virtual void autoMap(const fvPatchFieldMapper &)
Map (and resize as needed) from self given a mapping object.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:54
Common functions used in temperature coupled boundaries.
virtual void operator=(const UList< Type > &)
Definition: fvPatchField.C:396
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
tmp< scalarField > kappa(const scalarField &Tp) const
Given patch temperature calculate corresponding K field.
runTime write()
makePatchTypeField(fvPatchScalarField, thermalBaffleFvPatchScalarField)
turbulentHeatFluxTemperatureFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:451
bool compressible
Definition: pEqn.H:30