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 #include "OFstream.H"
32 
33 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
34 
36 (
37  OFstream& os
38 ) const
39 {
40  os << "# Values: magSf value qDot htc" << endl;
41 }
42 
43 
44 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
45 
48 (
49  const fvPatch& p,
51 )
52 :
54 {}
55 
56 
59 (
61  const fvPatch& p,
63  const fvPatchFieldMapper& mapper
64 )
65 :
67 {}
68 
69 
72 (
73  const fvPatch& p,
75  const dictionary& dict
76 )
77 :
79 {}
80 
81 
84 (
86 )
87 :
89 {}
90 
91 
94 (
97 )
98 :
100 {}
101 
102 
103 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
104 
107 {}
108 
109 
110 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111 
113 (
114  OFstream& os
115 ) const
116 {
117  if (log())
118  {
119  Info<< type() << ": " << this->patch().name()
120  << ": writing data to " << os.name()
121  << endl;
122  }
123 
124  const label patchI = patch().index();
125 
126  // heat flux [W/m2]
127  scalarField qDot(this->patch().size(), 0.0);
128 
129  typedef compressible::turbulenceModel cmpTurbModelType;
130 
131  static word turbName
132  (
134  (
137  )
138  );
139 
140  static word thermoName(basicThermo::dictName);
141 
142  if (db().foundObject<cmpTurbModelType>(turbName))
143  {
144  const cmpTurbModelType& turbModel =
145  db().lookupObject<cmpTurbModelType>(turbName);
146 
147  const basicThermo& thermo = turbModel.transport();
148 
149  const fvPatchScalarField& hep = thermo.he().boundaryField()[patchI];
150 
151  qDot = turbModel.alphaEff(patchI)*hep.snGrad();
152  }
153  else if (db().foundObject<basicThermo>(thermoName))
154  {
155  const basicThermo& thermo = db().lookupObject<basicThermo>(thermoName);
156 
157  const fvPatchScalarField& hep = thermo.he().boundaryField()[patchI];
158 
159  qDot = thermo.alpha().boundaryField()[patchI]*hep.snGrad();
160  }
161  else
162  {
164  (
165  "void Foam::externalCoupledTemperatureMixedFvPatchScalarField::"
166  "transferData"
167  "("
168  "OFstream&"
169  ") const"
170  ) << "Condition requires either compressible turbulence and/or "
171  << "thermo model to be available" << exit(FatalError);
172  }
173 
174  // patch temperature [K]
175  const scalarField Tp(*this);
176 
177  // near wall cell temperature [K]
178  const scalarField Tc(patchInternalField());
179 
180  // heat transfer coefficient [W/m2/K]
181  const scalarField htc(qDot/(Tp - Tc + ROOTVSMALL));
182 
183  if (Pstream::parRun())
184  {
185  int tag = Pstream::msgType() + 1;
186 
188  magSfs[Pstream::myProcNo()].setSize(this->patch().size());
189  magSfs[Pstream::myProcNo()] = this->patch().magSf();
190  Pstream::gatherList(magSfs, tag);
191 
193  values[Pstream::myProcNo()].setSize(this->patch().size());
194  values[Pstream::myProcNo()] = Tp;
195  Pstream::gatherList(values, tag);
196 
198  qDots[Pstream::myProcNo()].setSize(this->patch().size());
199  qDots[Pstream::myProcNo()] = qDot;
200  Pstream::gatherList(qDots, tag);
201 
203  htcs[Pstream::myProcNo()].setSize(this->patch().size());
204  htcs[Pstream::myProcNo()] = htc;
205  Pstream::gatherList(htcs, tag);
206 
207  if (Pstream::master())
208  {
209  forAll(values, procI)
210  {
211  const Field<scalar>& magSf = magSfs[procI];
212  const Field<scalar>& value = values[procI];
213  const Field<scalar>& qDot = qDots[procI];
214  const Field<scalar>& htc = htcs[procI];
215 
216  forAll(magSf, faceI)
217  {
218  os << magSf[faceI] << token::SPACE
219  << value[faceI] << token::SPACE
220  << qDot[faceI] << token::SPACE
221  << htc[faceI] << token::SPACE
222  << nl;
223  }
224  }
225 
226  os.flush();
227  }
228  }
229  else
230  {
231  const Field<scalar>& magSf(this->patch().magSf());
232 
233  forAll(patch(), faceI)
234  {
235  os << magSf[faceI] << token::SPACE
236  << Tp[faceI] << token::SPACE
237  << qDot[faceI] << token::SPACE
238  << htc[faceI] << token::SPACE
239  << nl;
240  }
241 
242  os.flush();
243  }
244 }
245 
246 
248 (
249  const Pstream::commsTypes comms
250 )
251 {
253 }
254 
255 
257 (
258  Ostream& os
259 ) const
260 {
262 }
263 
264 
265 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
266 
267 namespace Foam
268 {
270  (
273  );
274 }
275 
276 
277 // ************************************************************************* //
virtual void transferData(OFstream &os) const
Transfer data for external source.
Output to file stream.
Definition: OFstream.H:81
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:380
virtual tmp< Field< Type > > snGrad() const
Return patch-normal gradient.
Definition: fvPatchField.C:207
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
This boundary condition provides a temperatue interface to an external application. Values are transferred as plain text files, where OpenFOAM data is written as:
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
const char *const group
Group name for atomic constants.
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:52
A class for handling words, derived from string.
Definition: word.H:59
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::fvPatchFieldMapper.
messageStream Info
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Namespace for OpenFOAM.
const word & name() const
Return name.
Definition: fvPatch.H:149
static word groupName(Name name, const word &group)
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
dictionary dict
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:386
static const char nl
Definition: Ostream.H:260
void setSize(const label)
Reset size of List.
Definition: List.C:318
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual void flush()
Flush stream.
Definition: OSstream.C:246
virtual tmp< Field< Type > > patchInternalField() const
Return internal field next to patch as patch field.
Definition: fvPatchField.C:215
volScalarField & p
Definition: createFields.H:51
const fileName & name() const
Return the name of the stream.
Definition: OFstream.H:118
commsTypes
Types of communications.
Definition: UPstream.H:64
#define forAll(list, i)
Definition: UList.H:421
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
static int & msgType()
Message tag of standard messages.
Definition: UPstream.H:451
Macros for easy insertion into run-time selection tables.
externalCoupledTemperatureMixedFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:404
label index() const
Return the index of this patch in the fvBoundaryMesh.
Definition: fvPatch.H:179
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field.
error FatalError
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:188
virtual const volScalarField & alpha() const
Thermal diffusivity for enthalpy of mixture [kg/m/s].
Definition: basicThermo.C:490
label size() const
Return the number of elements in the UList.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const word propertiesName
Default name of the turbulence properties dictionary.
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:398
const DimensionedField< Type, volMesh > & dimensionedInternalField() const
Return dimensioned internal field reference.
Definition: fvPatchField.H:307
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:115
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:136
const fvPatch & patch() const
Return patch.
Definition: fvPatchField.H:300
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
virtual volScalarField & he()=0
Enthalpy/Internal energy [J/kg].
virtual void writeHeader(OFstream &os) const
Write header to transfer file.