energyRegionCoupledFvPatchScalarField.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) 2012-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 "Time.H"
30 
31 // * * * * * * * * * * * * * Static Member Data * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35  template<>
36  const char* Foam::NamedEnum
37  <
39  3
40  >::names[] =
41  {
42  "solid",
43  "fluid",
44  "undefined"
45  };
46 }
47 
48 
49 const Foam::NamedEnum
50 <
52  3
54 
55 
56 // * * * * * * * * * * * * * * * * Private members * * * * * * * * * * * * *//
57 
59 {
60  if (method_ == UNDEFINED)
61  {
62  if
63  (
64  this->db().foundObject<compressible::turbulenceModel>
65  (
67  )
68  )
69  {
70  method_ = FLUID;
71  }
72  else
73  {
74  method_ = SOLID;
75  }
76  }
77 
78  if (!nbrThermoPtr_)
79  {
81  (
83  (
85  )
86  );
87  }
88 
89  if (!thermoPtr_)
90  {
91  thermoPtr_ =
92  (
93  &this->db().lookupObject<basicThermo>
94  (
96  )
97  );
98  }
99 }
100 
101 
103 kappa() const
104 {
105  switch (method_)
106  {
107  case FLUID:
108  {
109  const compressible::turbulenceModel& turbModel =
110  this->db().lookupObject<compressible::turbulenceModel>
111  (
113  );
114 
115  return turbModel.kappaEff(patch().index());
116  }
117  break;
118 
119  case SOLID:
120  {
121  const basicThermo& thermo =
122  this->db().lookupObject<basicThermo>
123  (
125  );
126 
127  return thermo.kappa(patch().index());
128  }
129  break;
130 
131  case UNDEFINED:
132  {
134  << " on mesh " << this->db().name() << " patch "
135  << patch().name()
136  << " could not find a method in. Methods are: "
137  << methodTypeNames_.toc()
138  << " Not turbulenceModel or thermophysicalProperties"
139  << " were found"
140  << exit(FatalError);
141  }
142  break;
143  }
144  return scalarField(0);
145 }
146 
147 
149 weights() const
150 {
151  const fvPatch& patch = regionCoupledPatch_.patch();
152 
153  const scalarField deltas
154  (
155  patch.nf() & patch.delta()
156  );
157 
158  const scalarField alphaDelta(kappa()/deltas);
159 
160  const fvPatch& nbrPatch = regionCoupledPatch_.neighbFvPatch();
161 
162  const energyRegionCoupledFvPatchScalarField& nbrField =
163  refCast
164  <
166  >
167  (
168  nbrPatch.lookupPatchField<volScalarField, scalar>("T")
169  );
170 
171  // Needed in the first calculation of weights
172  nbrField.setMethod();
173 
174  const scalarField nbrAlpha
175  (
176  regionCoupledPatch_.regionCoupledPatch().interpolate
177  (
178  nbrField.kappa()
179  )
180  );
181 
182  const scalarField nbrDeltas
183  (
184  regionCoupledPatch_.regionCoupledPatch().interpolate
185  (
186  nbrPatch.nf() & nbrPatch.delta()
187  )
188  );
189 
190  const scalarField nbrAlphaDelta(nbrAlpha/nbrDeltas);
191 
192  tmp<scalarField> tw(new scalarField(deltas.size()));
193  scalarField& w = tw();
194 
195  forAll(alphaDelta, faceI)
196  {
197  scalar di = alphaDelta[faceI];
198  scalar dni = nbrAlphaDelta[faceI];
199 
200  w[faceI] = di/(di + dni);
201  }
202 
203  return tw;
204 }
205 
206 
207 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
208 
211 (
212  const fvPatch& p,
214 )
215 :
217  regionCoupledPatch_(refCast<const regionCoupledBaseFvPatch>(p)),
218  method_(UNDEFINED),
219  nbrThermoPtr_(NULL),
220  thermoPtr_(NULL)
221 {}
222 
223 
226 (
228  const fvPatch& p,
230  const fvPatchFieldMapper& mapper
231 )
232 :
233  coupledFvPatchField<scalar>(ptf, p, iF, mapper),
234  regionCoupledPatch_(refCast<const regionCoupledBaseFvPatch>(p)),
235  method_(ptf.method_),
236  nbrThermoPtr_(NULL),
237  thermoPtr_(NULL)
238 {}
239 
240 
243 (
244  const fvPatch& p,
246  const dictionary& dict
247 )
248 :
250  regionCoupledPatch_(refCast<const regionCoupledBaseFvPatch>(p)),
251  method_(UNDEFINED),
252  nbrThermoPtr_(NULL),
253  thermoPtr_(NULL)
254 {
255 
256  if (!isA<regionCoupledBase>(this->patch().patch()))
257  {
259  << "' not type '" << regionCoupledBase::typeName << "'"
260  << "\n for patch " << p.name()
261  << " of field " << dimensionedInternalField().name()
262  << " in file " << dimensionedInternalField().objectPath()
263  << exit(FatalError);
264  }
265 }
266 
267 
270 (
272 )
273 :
275  regionCoupledPatch_(ptf.regionCoupledPatch_),
276  method_(ptf.method_),
277  nbrThermoPtr_(NULL),
278  thermoPtr_(NULL)
279 {}
280 
281 
284 (
287 )
288 :
290  regionCoupledPatch_(ptf.regionCoupledPatch_),
291  method_(ptf.method_),
292  nbrThermoPtr_(NULL),
293  thermoPtr_(NULL)
294 {}
295 
296 
297 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
298 
300 snGrad() const
301 {
302  Debug("snGrad");
303  return
304  regionCoupledPatch_.patch().deltaCoeffs()
305  *(*this - patchInternalField());
306 }
307 
308 
310 snGrad(const scalarField&) const
311 {
312  Debug("snGrad");
313  return snGrad();
314 }
315 
316 
318 (
319  const Pstream::commsTypes
320 )
321 {
322  if (!updated())
323  {
324  updateCoeffs();
325  }
326 
327  label patchi = patch().index();
328  const scalarField& pp = thermoPtr_->p().boundaryField()[patchi];
329 
330  const scalarField lWeights(weights());
331 
333  (
334  thermoPtr_->he
335  (
336  pp,
337  lWeights*patchInternalTemperatureField()
338  + (1.0 - lWeights)*patchNeighbourTemperatureField(),
339  patchi
340  )
341  );
342 
344 }
345 
346 
350 {
351  const fvPatch& nbrPatch = regionCoupledPatch_.neighbFvPatch();
352 
353  const labelUList& nbrFaceCells = nbrPatch.faceCells();
354 
355  setMethod();
356 
357  const scalarField nbrIntT
358  (
359  nbrThermoPtr_->T().internalField(), nbrFaceCells
360  );
361 
362  scalarField intNbrT
363  (
364  regionCoupledPatch_.regionCoupledPatch().interpolate(nbrIntT)
365  );
366 
367  label patchi = patch().index();
368  const scalarField& pp = thermoPtr_->p().boundaryField()[patchi];
369  tmp<scalarField> tmyHE = thermoPtr_->he(pp, intNbrT, patchi);
370 
371  return tmyHE;
372 }
373 
374 
377 {
378  const fvPatch& nbrPatch = regionCoupledPatch_.neighbFvPatch();
379 
380  const labelUList& nbrFaceCells = nbrPatch.faceCells();
381 
382  const scalarField nbrIntT
383  (
384  nbrThermoPtr_->T().internalField(), nbrFaceCells
385  );
386 
387  tmp<scalarField> tintNbrT =
388  regionCoupledPatch_.regionCoupledPatch().interpolate(nbrIntT);
389 
390  return tintNbrT;
391 }
392 
393 
396 {
397  const labelUList& faceCells = regionCoupledPatch_.faceCells();
398 
399  tmp<scalarField> tintT
400  (
401  new scalarField(thermoPtr_->T().internalField(), faceCells)
402  );
403 
404  return tintT;
405 }
406 
407 
409 (
410  Field<scalar>& result,
411  const scalarField& psiInternal,
412  const scalarField& coeffs,
413  const direction cmpt,
414  const Pstream::commsTypes
415 ) const
416 {
417  setMethod();
418 
419  scalarField myHE(this->size());
420 
421  if (&psiInternal == &internalField())
422  {
423  label patchi = this->patch().index();
424  const scalarField& pp = thermoPtr_->p().boundaryField()[patchi];
425  const scalarField& Tp = thermoPtr_->T().boundaryField()[patchi];
426 
427  myHE = thermoPtr_->he(pp, Tp, patchi);
428  }
429  else
430  {
431  //NOTE: This is not correct for preconditioned solvers
432  // psiInternal is not the information needed of the slave
433  forAll(*this, facei)
434  {
435  myHE[facei] = psiInternal[regionCoupledPatch_.faceCells()[facei]];
436  }
437  }
438 
439  // Multiply the field by coefficients and add into the result
440  const labelUList& faceCells = regionCoupledPatch_.faceCells();
441 
442  forAll(faceCells, elemI)
443  {
444  result[faceCells[elemI]] -= coeffs[elemI]*myHE[elemI];
445  }
446 }
447 
448 
450 (
451  Field<scalar>& result,
452  const Field<scalar>& psiInternal,
453  const scalarField& coeffs,
454  const Pstream::commsTypes
455 ) const
456 {
458 }
459 
460 
462 {
464  this->writeEntry("value", os);
465 }
466 
467 
468 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
469 
470 namespace Foam
471 {
473  (
476  );
477 };
478 
479 
480 // ************************************************************************* //
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::fvc::snGrad
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:48
Foam::fvPatchField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:318
w
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
Foam::regionCoupledBaseFvPatch::nbrMesh
virtual const polyMesh & nbrMesh() const
Returns neighbour polyMesh.
Definition: regionCoupledBaseFvPatch.H:147
p
p
Definition: pEqn.H:62
Foam::energyRegionCoupledFvPatchScalarField::methodTypeNames_
static const NamedEnum< kappaMethodType, 3 > methodTypeNames_
Methof to extract kappa.
Definition: energyRegionCoupledFvPatchScalarField.H:80
Foam::energyRegionCoupledFvPatchScalarField::energyRegionCoupledFvPatchScalarField
energyRegionCoupledFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
Definition: energyRegionCoupledFvPatchScalarField.C:211
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::energyRegionCoupledFvPatchScalarField::method_
kappaMethodType method_
How to get K.
Definition: energyRegionCoupledFvPatchScalarField.H:83
Foam::basicThermo::dictName
static const word dictName
Definition: basicThermo.H:176
Foam::fvPatch::delta
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: fvPatch.C:142
dimensionedInternalField
rDeltaT dimensionedInternalField()
Foam::energyRegionCoupledFvPatchScalarField::patchNeighbourTemperatureField
tmp< scalarField > patchNeighbourTemperatureField() const
Return nbr temperature internal field.
Definition: energyRegionCoupledFvPatchScalarField.C:376
Foam::Field::operator
friend Ostream & operator(Ostream &, const Field< Type > &)
Foam::turbulenceModel::propertiesName
static const word propertiesName
Default name of the turbulence properties dictionary.
Definition: turbulenceModel.H:97
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::energyRegionCoupledFvPatchScalarField::evaluate
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
Definition: energyRegionCoupledFvPatchScalarField.C:318
Foam::basicThermo
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:52
Foam::energyRegionCoupledFvPatchScalarField::thermoPtr_
const basicThermo * thermoPtr_
AutoPtr to my thermo.
Definition: energyRegionCoupledFvPatchScalarField.H:89
Foam::energyRegionCoupledFvPatchScalarField::nbrThermoPtr_
const basicThermo * nbrThermoPtr_
AutoPtr to nbr thermo.
Definition: energyRegionCoupledFvPatchScalarField.H:86
Foam::energyRegionCoupledFvPatchScalarField::patchInternalTemperatureField
tmp< scalarField > patchInternalTemperatureField() const
Return local temperature internal field.
Definition: energyRegionCoupledFvPatchScalarField.C:395
Foam::constant::electromagnetic::kappa
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Debug
#define Debug(var)
Report a variable name and value.
Definition: messageStream.H:298
energyRegionCoupledFvPatchScalarField.H
Foam::energyRegionCoupledFvPatchScalarField::UNDEFINED
@ UNDEFINED
Definition: energyRegionCoupledFvPatchScalarField.H:68
Foam::fvPatch::nf
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:124
NotImplemented
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:365
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::Field::T
tmp< Field< Type > > T() const
Return the field transpose (only defined for second rank tensors)
Definition: Field.C:691
Foam::fvPatch
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Foam::energyRegionCoupledFvPatchScalarField::weights
tmp< scalarField > weights() const
Local weight for this coupled field.
Definition: energyRegionCoupledFvPatchScalarField.C:149
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Foam::energyRegionCoupledFvPatchScalarField::SOLID
@ SOLID
Definition: energyRegionCoupledFvPatchScalarField.H:66
Foam::LduInterfaceField< scalar >::updateInterfaceMatrix
virtual void updateInterfaceMatrix(scalarField &, const scalarField &, const scalarField &, const direction, const Pstream::commsTypes commsType) const =0
Inherit updateInterfaceMatrix from lduInterfaceField.
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
Foam::energyRegionCoupledFvPatchScalarField::regionCoupledPatch_
const regionCoupledBaseFvPatch & regionCoupledPatch_
Local reference to region couple patch.
Definition: energyRegionCoupledFvPatchScalarField.H:77
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::makePatchTypeField
makePatchTypeField(fvPatchVectorField, SRFFreestreamVelocityFvPatchVectorField)
Foam::refCast
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
internalField
conserve internalField()+
Foam::ThermalDiffusivity::kappaEff
virtual tmp< volScalarField > kappaEff() const
Return the effective turbulent thermal diffusivity for temperature.
Definition: ThermalDiffusivity.H:144
Foam::UPstream::commsTypes
commsTypes
Types of communications.
Definition: UPstream.H:64
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::energyRegionCoupledFvPatchScalarField::FLUID
@ FLUID
Definition: energyRegionCoupledFvPatchScalarField.H:67
Foam::energyRegionCoupledFvPatchScalarField::kappa
tmp< scalarField > kappa() const
Return kappa.
Definition: energyRegionCoupledFvPatchScalarField.C:103
Foam::fvPatchField< scalar >::db
const objectRegistry & db() const
Return local objectRegistry.
Definition: fvPatchField.C:181
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: ThermalDiffusivity.H:48
Foam::coupledFvPatchField< scalar >
Foam::energyRegionCoupledFvPatchScalarField::setMethod
void setMethod() const
Set method.
Definition: energyRegionCoupledFvPatchScalarField.C:58
Foam::energyRegionCoupledFvPatchScalarField
Energy region coupled implicit boundary condition. The fvPatch is treated as uncoupled from the delta...
Definition: energyRegionCoupledFvPatchScalarField.H:57
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::direction
unsigned char direction
Definition: direction.H:43
Foam::fvPatch::lookupPatchField
const GeometricField::PatchFieldType & lookupPatchField(const word &name, const GeometricField *=NULL, const Type *=NULL) const
Lookup and return the patchField of the named field from the.
Definition: fvPatchFvMeshTemplates.C:32
Foam::fvPatch::faceCells
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
Foam::fvPatch::patch
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:143
Foam::fvPatchFieldMapper
Foam::fvPatchFieldMapper.
Definition: fvPatchFieldMapper.H:45
Foam::energyRegionCoupledFvPatchScalarField::write
virtual void write(Ostream &) const
Write.
Definition: energyRegionCoupledFvPatchScalarField.C:461
Foam::energyRegionCoupledFvPatchScalarField::patchNeighbourField
virtual tmp< scalarField > patchNeighbourField() const
Return neighbour coupled internal cell data.
Definition: energyRegionCoupledFvPatchScalarField.C:349
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Definition: objectRegistryTemplates.C:165
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::energyRegionCoupledFvPatchScalarField::kappaMethodType
kappaMethodType
Definition: energyRegionCoupledFvPatchScalarField.H:64
turbulentFluidThermoModel.H
Foam::NamedEnum
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
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::energyRegionCoupledFvPatchScalarField::snGrad
virtual tmp< scalarField > snGrad() const
Return patch-normal gradient.
Definition: energyRegionCoupledFvPatchScalarField.C:300