externalCoupledTemperatureMixedFvPatchScalarField.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) 2013-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 
29 #include "fvPatchFieldMapper.H"
30 #include "volFields.H"
31 
32 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
33 
35 (
36  Ostream& os
37 ) const
38 {
39  os << "# Values: magSf T qDot htc" << endl;
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
47 (
48  const fvPatch& p,
50 )
51 :
53 {}
54 
55 
58 (
60  const fvPatch& p,
62  const fvPatchFieldMapper& mapper
63 )
64 :
66 {}
67 
68 
71 (
72  const fvPatch& p,
74  const dictionary& dict
75 )
76 :
77  //externalCoupledMixedFvPatchField<scalar>(p, iF, dict)
79 {
80  if (dict.found("refValue"))
81  {
82  // Initialise same way as mixed
83  this->refValue() = scalarField("refValue", dict, p.size());
84  this->refGrad() = scalarField("refGradient", dict, p.size());
85  this->valueFraction() = scalarField("valueFraction", dict, p.size());
86 
87  evaluate();
88  }
89  else
90  {
91  // For convenience: initialise as fixedValue with either read value
92  // or extrapolated value
93  if (dict.found("value"))
94  {
96  (
97  scalarField("value", dict, p.size())
98  );
99  }
100  else
101  {
102  fvPatchField<scalar>::operator=(this->patchInternalField());
103  }
104 
105  // Initialise as a fixed value
106  this->refValue() = *this;
107  this->refGrad() = pTraits<scalar>::zero;
108  this->valueFraction() = 1.0;
109  }
110 }
111 
112 
115 (
117 )
118 :
120 {}
121 
122 
125 (
128 )
129 :
131 {}
132 
133 
134 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
135 
138 {}
139 
140 
141 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
142 
144 (
145  Ostream& os
146 ) const
147 {
148  const label patchI = patch().index();
149 
150  // Heat flux [W/m2]
151  scalarField qDot(this->patch().size(), 0.0);
152 
153  typedef compressible::turbulenceModel cmpTurbModelType;
154 
155  static word turbName
156  (
158  (
161  )
162  );
163 
164  static word thermoName("thermophysicalProperties");
165 
166  if (db().foundObject<cmpTurbModelType>(turbName))
167  {
168  const cmpTurbModelType& turbModel =
169  db().lookupObject<cmpTurbModelType>(turbName);
170 
171  const basicThermo& thermo = turbModel.transport();
172 
173  const fvPatchScalarField& hep = thermo.he().boundaryField()[patchI];
174 
175  qDot = turbModel.alphaEff(patchI)*hep.snGrad();
176  }
177  else if (db().foundObject<basicThermo>(thermoName))
178  {
179  const basicThermo& thermo = db().lookupObject<basicThermo>(thermoName);
180 
181  const fvPatchScalarField& hep = thermo.he().boundaryField()[patchI];
182 
183  qDot = thermo.alpha().boundaryField()[patchI]*hep.snGrad();
184  }
185  else
186  {
188  << "Condition requires either compressible turbulence and/or "
189  << "thermo model to be available" << exit(FatalError);
190  }
191 
192  // Patch temperature [K]
193  const scalarField& Tp(*this);
194 
195  // Near wall cell temperature [K]
196  const scalarField Tc(patchInternalField());
197 
198  // Heat transfer coefficient [W/m2/K]
199  const scalarField htc(qDot/(Tp - Tc + ROOTVSMALL));
200 
201  const Field<scalar>& magSf(this->patch().magSf());
202 
203  forAll(patch(), faceI)
204  {
205  os << magSf[faceI] << token::SPACE
206  << Tp[faceI] << token::SPACE
207  << qDot[faceI] << token::SPACE
208  << htc[faceI] << token::SPACE
209  << nl;
210  }
211 }
212 
213 
215 (
216  Istream& is
217 )
218 {
219  // Assume generic input stream so we can do line-based format and skip
220  // unused columns
221  ISstream& iss = dynamic_cast<ISstream&>(is);
222 
223  string line;
224 
225  forAll(*this, faceI)
226  {
227  iss.getLine(line);
228  IStringStream lineStr(line);
229 
230  lineStr
231  >> this->refValue()[faceI]
232  >> this->refGrad()[faceI]
233  >> this->valueFraction()[faceI];
234  }
235 }
236 
237 
238 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
239 
240 namespace Foam
241 {
243  (
246  );
247 }
248 
249 
250 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
volFields.H
Foam::fvPatchField::snGrad
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:200
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
Foam::ISstream
Generic input stream.
Definition: ISstream.H:51
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::basicThermo
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:52
Foam::constant::atomic::group
const char *const group
Group name for atomic constants.
Definition: atomicConstants.C:40
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
fvPatchFieldMapper.H
Foam::externalCoupledMixedFvPatchField< scalar >
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::Istream
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:57
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::ISstream::getLine
ISstream & getLine(string &)
Raw, low-level getline into a string function.
Definition: ISstreamI.H:77
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::externalCoupledTemperatureMixedFvPatchScalarField::~externalCoupledTemperatureMixedFvPatchScalarField
virtual ~externalCoupledTemperatureMixedFvPatchScalarField()
Destructor.
Definition: externalCoupledTemperatureMixedFvPatchScalarField.C:137
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
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::externalCoupledTemperatureMixedFvPatchScalarField::writeData
virtual void writeData(Ostream &) const
Write data.
Definition: externalCoupledTemperatureMixedFvPatchScalarField.C:144
Foam::IStringStream
Input from memory buffer stream.
Definition: IStringStream.H:49
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::makePatchTypeField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
externalCoupledTemperatureMixedFvPatchScalarField.H
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: ThermalDiffusivity.H:48
Foam::externalCoupledTemperatureMixedFvPatchScalarField::externalCoupledTemperatureMixedFvPatchScalarField
externalCoupledTemperatureMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: externalCoupledTemperatureMixedFvPatchScalarField.C:47
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
scalarField
volScalarField scalarField(fieldObject, mesh)
Foam::line
A line primitive.
Definition: line.H:56
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::externalCoupledTemperatureMixedFvPatchScalarField::writeHeader
virtual void writeHeader(Ostream &) const
Write header.
Definition: externalCoupledTemperatureMixedFvPatchScalarField.C:35
Foam::externalCoupledTemperatureMixedFvPatchScalarField
This boundary condition provides a temperatue interface to an external application....
Definition: externalCoupledTemperatureMixedFvPatchScalarField.H:98
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::externalCoupledTemperatureMixedFvPatchScalarField::readData
virtual void readData(Istream &)
Read data.
Definition: externalCoupledTemperatureMixedFvPatchScalarField.C:215
turbulentFluidThermoModel.H
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::token::SPACE
@ SPACE
Definition: token.H:95