energyRegionCoupledFvPatchScalarField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2018 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
53 > Foam::energyRegionCoupledFvPatchScalarField::methodTypeNames_;
54 
55 
56 // * * * * * * * * * * * * * * * * Private members * * * * * * * * * * * * *//
57 
58 void Foam::energyRegionCoupledFvPatchScalarField::setMethod() const
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  {
80  nbrThermoPtr_ =
81  (
82  &regionCoupledPatch_.nbrMesh().lookupObject<basicThermo>
83  (
85  )
86  );
87  }
88 
89  if (!thermoPtr_)
90  {
91  thermoPtr_ =
92  (
93  &this->db().lookupObject<basicThermo>
94  (
96  )
97  );
98  }
99 }
100 
101 
102 Foam::tmp<Foam::scalarField> Foam::energyRegionCoupledFvPatchScalarField::
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 
148 Foam::tmp<Foam::scalarField> Foam::energyRegionCoupledFvPatchScalarField::
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.ref();
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_(nullptr),
220  thermoPtr_(nullptr)
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_(nullptr),
237  thermoPtr_(nullptr)
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_(nullptr),
253  thermoPtr_(nullptr)
254 {
255 
256  if (!isA<regionCoupledBase>(this->patch().patch()))
257  {
259  << "' not type '" << regionCoupledBase::typeName << "'"
260  << "\n for patch " << p.name()
261  << " of field " << internalField().name()
262  << " in file " << internalField().objectPath()
263  << exit(FatalError);
264  }
265 }
266 
267 
270 (
272 )
273 :
275  regionCoupledPatch_(ptf.regionCoupledPatch_),
276  method_(ptf.method_),
277  nbrThermoPtr_(nullptr),
278  thermoPtr_(nullptr)
279 {}
280 
281 
284 (
287 )
288 :
290  regionCoupledPatch_(ptf.regionCoupledPatch_),
291  method_(ptf.method_),
292  nbrThermoPtr_(nullptr),
293  thermoPtr_(nullptr)
294 {}
295 
296 
297 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
298 
300 snGrad() const
301 {
302  return
303  regionCoupledPatch_.patch().deltaCoeffs()
304  *(*this - patchInternalField());
305 }
306 
307 
309 snGrad(const scalarField&) const
310 {
311  return snGrad();
312 }
313 
314 
316 (
317  const Pstream::commsTypes
318 )
319 {
320  if (!updated())
321  {
322  updateCoeffs();
323  }
324 
325  label patchi = patch().index();
326  const scalarField& pp = thermoPtr_->p().boundaryField()[patchi];
327 
328  const scalarField lWeights(weights());
329 
331  (
332  thermoPtr_->he
333  (
334  pp,
335  lWeights*patchInternalTemperatureField()
336  + (1.0 - lWeights)*patchNeighbourTemperatureField(),
337  patchi
338  )
339  );
340 
342 }
343 
344 
348 {
349  const fvPatch& nbrPatch = regionCoupledPatch_.neighbFvPatch();
350 
351  const labelUList& nbrFaceCells = nbrPatch.faceCells();
352 
353  setMethod();
354 
355  const scalarField nbrIntT
356  (
357  nbrThermoPtr_->T().primitiveField(), nbrFaceCells
358  );
359 
360  scalarField intNbrT
361  (
362  regionCoupledPatch_.regionCoupledPatch().interpolate(nbrIntT)
363  );
364 
365  label patchi = patch().index();
366  const scalarField& pp = thermoPtr_->p().boundaryField()[patchi];
367  tmp<scalarField> tmyHE = thermoPtr_->he(pp, intNbrT, patchi);
368 
369  return tmyHE;
370 }
371 
372 
373 Foam::tmp<Foam::scalarField> Foam::energyRegionCoupledFvPatchScalarField::
374 patchNeighbourTemperatureField() const
375 {
376  const fvPatch& nbrPatch = regionCoupledPatch_.neighbFvPatch();
377 
378  const labelUList& nbrFaceCells = nbrPatch.faceCells();
379 
380  const scalarField nbrIntT
381  (
382  nbrThermoPtr_->T().primitiveField(), nbrFaceCells
383  );
384 
385  tmp<scalarField> tintNbrT =
386  regionCoupledPatch_.regionCoupledPatch().interpolate(nbrIntT);
387 
388  return tintNbrT;
389 }
390 
391 
392 Foam::tmp<Foam::scalarField> Foam::energyRegionCoupledFvPatchScalarField::
393 patchInternalTemperatureField() const
394 {
395  const labelUList& faceCells = regionCoupledPatch_.faceCells();
396 
397  tmp<scalarField> tintT
398  (
399  new scalarField(thermoPtr_->T().primitiveField(), faceCells)
400  );
401 
402  return tintT;
403 }
404 
405 
407 (
408  Field<scalar>& result,
409  const scalarField& psiInternal,
410  const scalarField& coeffs,
411  const direction cmpt,
412  const Pstream::commsTypes
413 ) const
414 {
415  setMethod();
416 
417  scalarField myHE(this->size());
418 
419  if (&psiInternal == &primitiveField())
420  {
421  label patchi = this->patch().index();
422  const scalarField& pp = thermoPtr_->p().boundaryField()[patchi];
423  const scalarField& Tp = thermoPtr_->T().boundaryField()[patchi];
424 
425  myHE = thermoPtr_->he(pp, Tp, patchi);
426  }
427  else
428  {
429  // NOTE: This is not correct for preconditioned solvers
430  // psiInternal is not the information needed of the slave
431  forAll(*this, facei)
432  {
433  myHE[facei] = psiInternal[regionCoupledPatch_.faceCells()[facei]];
434  }
435  }
436 
437  // Multiply the field by coefficients and add into the result
438  const labelUList& faceCells = regionCoupledPatch_.faceCells();
439 
440  forAll(faceCells, elemI)
441  {
442  result[faceCells[elemI]] -= coeffs[elemI]*myHE[elemI];
443  }
444 }
445 
446 
448 (
449  Field<scalar>& result,
450  const Field<scalar>& psiInternal,
451  const scalarField& coeffs,
452  const Pstream::commsTypes
453 ) const
454 {
456 }
457 
458 
460 {
462  this->writeEntry("value", os);
463 }
464 
465 
466 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
467 
468 namespace Foam
469 {
471  (
474  );
475 };
476 
477 
478 // ************************************************************************* //
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
virtual tmp< volScalarField > kappa() const =0
Thermal diffusivity for temperature of mixture [J/m/s/K].
Abstract base-class for fluid and solid thermodynamic properties.
Definition: basicThermo.H:52
virtual tmp< volScalarField > kappaEff() const
Effective thermal turbulent diffusivity for temperature.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
commsTypes
Types of communications.
Definition: UPstream.H:64
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
tmp< vectorField > nf() const
Return face normals.
Definition: fvPatch.C:124
uint8_t direction
Definition: direction.H:45
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
virtual void updateInterfaceMatrix(Field< scalar > &result, const scalarField &psiInternal, const scalarField &coeffs, const direction cmpt, const Pstream::commsTypes commsType) const
Update result field based on interface functionality.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:61
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:66
virtual void evaluate(const Pstream::commsTypes commsType)
Evaluate the patch field.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
rhoReactionThermo & thermo
Definition: createFields.H:28
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:362
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
const dimensionedScalar kappa
Coulomb constant: default SI units: [N.m2/C2].
Macros for easy insertion into run-time selection tables.
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:115
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
friend Ostream & operator(Ostream &, const Field< scalar > &)
virtual tmp< scalarField > patchNeighbourField() const
Return neighbour coupled internal cell data.
static const word propertiesName
Default name of the turbulence properties dictionary.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Energy region coupled implicit boundary condition. The fvPatch is treated as uncoupled from the delta...
virtual const labelUList & faceCells() const
Return faceCells.
Definition: fvPatch.C:93
const polyPatch & patch() const
Return the polyPatch.
Definition: fvPatch.H:143
Foam::fvPatchFieldMapper.
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:61
virtual tmp< scalarField > snGrad() const
Return patch-normal gradient.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
energyRegionCoupledFvPatchScalarField(const fvPatch &, const DimensionedField< scalar, volMesh > &)
Construct from patch and internal field.
virtual void evaluate(const Pstream::commsTypes commsType=Pstream::commsTypes::blocking)
Evaluate the patch field, sets Updated to false.
Definition: fvPatchField.C:331
label patchi
const GeometricField::Patch & lookupPatchField(const word &name, const GeometricField *=nullptr, const Type *=nullptr) const
Lookup and return the patchField of the named field from the.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
tmp< Field< Type > > T() const
Return the field transpose (only defined for second rank tensors)
Definition: Field.C:717
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:366
makePatchTypeField(fvPatchScalarField, atmBoundaryLayerInletEpsilonFvPatchScalarField)
virtual tmp< vectorField > delta() const
Return cell-centre to face-centre vector.
Definition: fvPatch.C:142
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > snGrad(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvcSnGrad.C:45
Namespace for OpenFOAM.