viewFactor.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "viewFactor.H"
30 #include "surfaceFields.H"
31 #include "constants.H"
33 #include "typeInfo.H"
36 
37 using namespace Foam::constant;
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  namespace radiation
44  {
45  defineTypeNameAndDebug(viewFactor, 0);
47  }
48 }
49 
51  = "viewFactorWall";
52 
53 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
54 
56 {
57  const polyBoundaryMesh& coarsePatches = coarseMesh_.boundaryMesh();
58 
59  selectedPatches_ = mesh_.boundaryMesh().indices(viewFactorWalls);
60 
61  for (const label patchi : selectedPatches_)
62  {
63  nLocalCoarseFaces_ += coarsePatches[patchi].size();
64  }
65 
66  if (debug)
67  {
68  Pout<< "radiation::viewFactor::initialise() Selected patches:"
69  << selectedPatches_ << endl;
70  Pout<< "radiation::viewFactor::initialise() Number of coarse faces:"
71  << nLocalCoarseFaces_ << endl;
72  }
73 
74  totalNCoarseFaces_ = nLocalCoarseFaces_;
75  reduce(totalNCoarseFaces_, sumOp<label>());
76 
78  << "Total number of clusters : " << totalNCoarseFaces_ << endl;
79 
80  map_.reset
81  (
82  new IOmapDistribute
83  (
84  IOobject
85  (
86  "mapDist",
87  mesh_.facesInstance(),
88  mesh_,
91  false
92  )
93  )
94  );
95 
96  scalarListIOList FmyProc
97  (
98  IOobject
99  (
100  "F",
101  mesh_.facesInstance(),
102  mesh_,
105  false
106  )
107  );
108 
109  labelListIOList globalFaceFaces
110  (
111  IOobject
112  (
113  "globalFaceFaces",
114  mesh_.facesInstance(),
115  mesh_,
118  false
119  )
120  );
121 
122  List<labelListList> globalFaceFacesProc(Pstream::nProcs());
123  globalFaceFacesProc[Pstream::myProcNo()] = globalFaceFaces;
124  Pstream::gatherList(globalFaceFacesProc);
125 
127  F[Pstream::myProcNo()] = FmyProc;
129 
130  globalIndex globalNumbering(nLocalCoarseFaces_);
131 
132  if (Pstream::master())
133  {
134  Fmatrix_.reset
135  (
136  new scalarSquareMatrix(totalNCoarseFaces_, Zero)
137  );
138 
140  << "Insert elements in the matrix..." << endl;
141 
142  for (const int procI : Pstream::allProcs())
143  {
144  insertMatrixElements
145  (
146  globalNumbering,
147  procI,
148  globalFaceFacesProc[procI],
149  F[procI],
150  Fmatrix_()
151  );
152  }
153 
154 
155  if (coeffs_.get<bool>("smoothing"))
156  {
157  DebugInFunction << "Smoothing the matrix..." << endl;
158 
159  for (label i=0; i<totalNCoarseFaces_; i++)
160  {
161  scalar sumF = 0.0;
162  for (label j=0; j<totalNCoarseFaces_; j++)
163  {
164  sumF += Fmatrix_()(i, j);
165  }
166 
167  const scalar delta = sumF - 1.0;
168  for (label j=0; j<totalNCoarseFaces_; j++)
169  {
170  Fmatrix_()(i, j) *= (1.0 - delta/(sumF + 0.001));
171  }
172  }
173  }
174 
175  coeffs_.readEntry("constantEmissivity", constEmissivity_);
176  if (constEmissivity_)
177  {
178  CLU_.reset
179  (
180  new scalarSquareMatrix(totalNCoarseFaces_, Zero)
181  );
182 
183  pivotIndices_.setSize(CLU_().m());
184  }
185  }
186 
187  coeffs_.readIfPresent("useSolarLoad", useSolarLoad_);
188 
189  if (useSolarLoad_)
190  {
191  const dictionary& solarDict = this->subDict("solarLoadCoeffs");
192  solarLoad_.reset(new solarLoad(solarDict, T_));
193 
194  if (solarLoad_->nBands() != nBands_)
195  {
197  << "Solar radiation and view factor band numbers "
198  << "are different"
199  << abort(FatalError);
200  }
201 
202  Info<< "Creating Solar Load Model " << nl;
203  }
204 }
205 
206 
207 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
208 
209 Foam::radiation::viewFactor::viewFactor(const volScalarField& T)
210 :
211  radiationModel(typeName, T),
212  finalAgglom_
213  (
214  IOobject
215  (
216  "finalAgglom",
217  mesh_.facesInstance(),
218  mesh_,
219  IOobject::MUST_READ,
220  IOobject::NO_WRITE,
221  false
222  )
223  ),
224  map_(),
225  coarseMesh_
226  (
227  IOobject
228  (
229  "coarse:" + mesh_.name(),
230  mesh_.polyMesh::instance(),
231  mesh_.time(),
232  IOobject::NO_READ,
233  IOobject::NO_WRITE,
234  false
235  ),
236  mesh_,
237  finalAgglom_
238  ),
239  qr_
240  (
241  IOobject
242  (
243  "qr",
244  mesh_.time().timeName(),
245  mesh_,
246  IOobject::MUST_READ,
247  IOobject::AUTO_WRITE
248  ),
249  mesh_
250  ),
251  Fmatrix_(),
252  CLU_(),
253  selectedPatches_(mesh_.boundary().size(), -1),
254  totalNCoarseFaces_(0),
255  nLocalCoarseFaces_(0),
256  constEmissivity_(false),
257  iterCounter_(0),
258  pivotIndices_(0),
259  useSolarLoad_(false),
260  solarLoad_(),
261  nBands_(coeffs_.getOrDefault<label>("nBands", 1))
262 {
263  initialise();
264 }
265 
266 
267 Foam::radiation::viewFactor::viewFactor
268 (
269  const dictionary& dict,
270  const volScalarField& T
271 )
272 :
273  radiationModel(typeName, dict, T),
274  finalAgglom_
275  (
276  IOobject
277  (
278  "finalAgglom",
279  mesh_.facesInstance(),
280  mesh_,
283  false
284  )
285  ),
286  map_(),
287  coarseMesh_
288  (
289  IOobject
290  (
291  "coarse:" + mesh_.name(),
292  mesh_.polyMesh::instance(),
293  mesh_.time(),
296  false
297  ),
298  mesh_,
299  finalAgglom_
300  ),
301  qr_
302  (
303  IOobject
304  (
305  "qr",
306  mesh_.time().timeName(),
307  mesh_,
310  ),
311  mesh_
312  ),
313  Fmatrix_(),
314  CLU_(),
315  selectedPatches_(mesh_.boundary().size(), -1),
316  totalNCoarseFaces_(0),
317  nLocalCoarseFaces_(0),
318  constEmissivity_(false),
319  iterCounter_(0),
320  pivotIndices_(0),
321  useSolarLoad_(false),
322  solarLoad_(),
323  nBands_(coeffs_.getOrDefault<label>("nBands", 1))
324 {
325  initialise();
326 }
327 
328 
329 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
330 
332 {
333  if (radiationModel::read())
334  {
335  return true;
336  }
337 
338  return false;
339 }
340 
341 
343 (
344  const globalIndex& globalNumbering,
345  const label procI,
346  const labelListList& globalFaceFaces,
347  const scalarListList& viewFactors,
348  scalarSquareMatrix& Fmatrix
349 )
350 {
351  forAll(viewFactors, faceI)
352  {
353  const scalarList& vf = viewFactors[faceI];
354  const labelList& globalFaces = globalFaceFaces[faceI];
355 
356  label globalI = globalNumbering.toGlobal(procI, faceI);
357  forAll(globalFaces, i)
358  {
359  Fmatrix[globalI][globalFaces[i]] = vf[i];
360  }
361  }
362 }
363 
364 
366 {
367  // Store previous iteration
368  qr_.storePrevIter();
369 
370  if (useSolarLoad_)
371  {
372  solarLoad_->calculate();
373  }
374 
375  // Net radiation
376  scalarField q(totalNCoarseFaces_, 0.0);
377  volScalarField::Boundary& qrBf = qr_.boundaryFieldRef();
378 
379  globalIndex globalNumbering(nLocalCoarseFaces_);
380 
381  const boundaryRadiationProperties& boundaryRadiation =
383 
384  for (label bandI = 0; bandI < nBands_; bandI++)
385  {
386  scalarField compactCoarseT4(map_->constructSize(), 0.0);
387  scalarField compactCoarseE(map_->constructSize(), 0.0);
388  scalarField compactCoarseHo(map_->constructSize(), 0.0);
389 
390  // Fill local averaged(T), emissivity(E) and external heatFlux(Ho)
391  DynamicList<scalar> localCoarseT4ave(nLocalCoarseFaces_);
392  DynamicList<scalar> localCoarseEave(nLocalCoarseFaces_);
393  DynamicList<scalar> localCoarseHoave(nLocalCoarseFaces_);
394 
395  forAll(selectedPatches_, i)
396  {
397  label patchID = selectedPatches_[i];
398 
399  const scalarField& Tp = T_.boundaryField()[patchID];
400  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
401 
402  fvPatchScalarField& qrPatch = qrBf[patchID];
403 
405  refCast
406  <
408  >(qrPatch);
409 
410  const tmp<scalarField> teb =
411  boundaryRadiation.emissivity(patchID, bandI);
412 
413  const scalarField& eb = teb();
414 
415  const tmp<scalarField> tHoi = qrp.qro(bandI);
416  const scalarField& Hoi = tHoi();
417 
418  const polyPatch& pp = coarseMesh_.boundaryMesh()[patchID];
419 
420  const labelList& coarsePatchFace =
421  coarseMesh_.patchFaceMap()[patchID];
422 
423  scalarList T4ave(pp.size(), 0.0);
424  scalarList Eave(pp.size(), 0.0);
425  scalarList Hoiave(pp.size(), 0.0);
426 
427  if (pp.size() > 0)
428  {
429  const labelList& agglom = finalAgglom_[patchID];
430  label nAgglom = max(agglom) + 1;
431 
432  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
433 
434  forAll(coarseToFine, coarseI)
435  {
436  const label coarseFaceID = coarsePatchFace[coarseI];
437  const labelList& fineFaces = coarseToFine[coarseFaceID];
438  UIndirectList<scalar> fineSf
439  (
440  sf,
441  fineFaces
442  );
443 
444  const scalar area = sum(fineSf());
445 
446  // Temperature, emissivity and external flux area weighting
447  forAll(fineFaces, j)
448  {
449  label facei = fineFaces[j];
450  T4ave[coarseI] += (pow4(Tp[facei])*sf[facei])/area;
451  Eave[coarseI] += (eb[facei]*sf[facei])/area;
452  Hoiave[coarseI] += (Hoi[facei]*sf[facei])/area;
453  }
454  }
455  }
456 
457  localCoarseT4ave.append(T4ave);
458  localCoarseEave.append(Eave);
459  localCoarseHoave.append(Hoiave);
460 
461  }
462 
463  // Fill the local values to distribute
464  SubList<scalar>(compactCoarseT4, nLocalCoarseFaces_) =
465  localCoarseT4ave;
466  SubList<scalar>(compactCoarseE, nLocalCoarseFaces_) = localCoarseEave;
467  SubList<scalar>(compactCoarseHo, nLocalCoarseFaces_) =
468  localCoarseHoave;
469 
470  // Distribute data
471  map_->distribute(compactCoarseT4);
472  map_->distribute(compactCoarseE);
473  map_->distribute(compactCoarseHo);
474 
475  // Distribute local global ID
476  labelList compactGlobalIds(map_->constructSize(), Zero);
477 
479  (
480  compactGlobalIds,
481  nLocalCoarseFaces_
482  ) = identity
483  (
484  globalNumbering.localSize(),
485  globalNumbering.localStart()
486  );
487 
488  map_->distribute(compactGlobalIds);
489 
490  // Create global size vectors
491  scalarField T4(totalNCoarseFaces_, 0.0);
492  scalarField E(totalNCoarseFaces_, 0.0);
493  scalarField qrExt(totalNCoarseFaces_, 0.0);
494 
495  // Fill lists from compact to global indexes.
496  forAll(compactCoarseT4, i)
497  {
498  T4[compactGlobalIds[i]] = compactCoarseT4[i];
499  E[compactGlobalIds[i]] = compactCoarseE[i];
500  qrExt[compactGlobalIds[i]] = compactCoarseHo[i];
501  }
502 
506 
510 
511  if (Pstream::master())
512  {
513  // Variable emissivity
514  if (!constEmissivity_)
515  {
516  scalarSquareMatrix C(totalNCoarseFaces_, 0.0);
517 
518  for (label i=0; i<totalNCoarseFaces_; i++)
519  {
520  for (label j=0; j<totalNCoarseFaces_; j++)
521  {
522  const scalar invEj = 1.0/E[j];
523  const scalar sigmaT4 =
524  physicoChemical::sigma.value()*T4[j];
525 
526  if (i==j)
527  {
528  C(i, j) = invEj - (invEj - 1.0)*Fmatrix_()(i, j);
529  q[i] +=
530  (Fmatrix_()(i, j) - 1.0)*sigmaT4 + qrExt[j];
531  }
532  else
533  {
534  C(i, j) = (1.0 - invEj)*Fmatrix_()(i, j);
535  q[i] += Fmatrix_()(i, j)*sigmaT4;
536  }
537 
538  }
539  }
540 
541  Info<< "Solving view factor equations for band :"
542  << bandI << endl;
543 
544  // Negative coming into the fluid
545  LUsolve(C, q);
546  }
547  else //Constant emissivity
548  {
549  // Initial iter calculates CLU and caches it
550  if (iterCounter_ == 0)
551  {
552  for (label i=0; i<totalNCoarseFaces_; i++)
553  {
554  for (label j=0; j<totalNCoarseFaces_; j++)
555  {
556  const scalar invEj = 1.0/E[j];
557  if (i==j)
558  {
559  CLU_()(i, j) =
560  invEj-(invEj-1.0)*Fmatrix_()(i, j);
561  }
562  else
563  {
564  CLU_()(i, j) = (1.0 - invEj)*Fmatrix_()(i, j);
565  }
566  }
567  }
568 
569  DebugInFunction << "\nDecomposing C matrix..." << endl;
570 
571  LUDecompose(CLU_(), pivotIndices_);
572  }
573 
574  for (label i=0; i<totalNCoarseFaces_; i++)
575  {
576  for (label j=0; j<totalNCoarseFaces_; j++)
577  {
578  const scalar sigmaT4 =
580 
581  if (i==j)
582  {
583  q[i] +=
584  (Fmatrix_()(i, j) - 1.0)*sigmaT4 + qrExt[j];
585  }
586  else
587  {
588  q[i] += Fmatrix_()(i, j)*sigmaT4;
589  }
590  }
591  }
592 
593 
594  Info<< "Solving view factor equations for band : "
595  << bandI << endl;
596 
597 
598  LUBacksubstitute(CLU_(), pivotIndices_, q);
599  iterCounter_ ++;
600  }
601  }
602 
603  }
604  // Scatter q and fill qr
607 
608  label globCoarseId = 0;
609  for (const label patchID : selectedPatches_)
610  {
611  const polyPatch& pp = mesh_.boundaryMesh()[patchID];
612 
613  if (pp.size() > 0)
614  {
615  scalarField& qrp = qrBf[patchID];
616  //const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
617  const labelList& agglom = finalAgglom_[patchID];
618  label nAgglom = max(agglom)+1;
619 
620  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
621 
622  const labelList& coarsePatchFace =
623  coarseMesh_.patchFaceMap()[patchID];
624 
626  forAll(coarseToFine, coarseI)
627  {
628  label globalCoarse =
629  globalNumbering.toGlobal
630  (Pstream::myProcNo(), globCoarseId);
631 
632  const label coarseFaceID = coarsePatchFace[coarseI];
633  const labelList& fineFaces = coarseToFine[coarseFaceID];
634 
635  for (const label facei : fineFaces)
636  {
637  qrp[facei] = q[globalCoarse];
639  }
640  globCoarseId ++;
641  }
642  }
643  }
644 
645  if (debug)
646  {
647  forAll(qrBf, patchID)
648  {
649  const scalarField& qrp = qrBf[patchID];
650  const scalarField& magSf = mesh_.magSf().boundaryField()[patchID];
651  const scalar heatFlux = gSum(qrp*magSf);
652 
654  << "Total heat transfer rate at patch: "
655  << patchID << " "
656  << heatFlux << endl;
657  }
658  }
659 
660  // Relax qr if necessary
661  qr_.relax();
662 }
663 
664 
666 {
668  (
669  IOobject
670  (
671  "Rp",
672  mesh_.time().timeName(),
673  mesh_,
676  false
677  ),
678  mesh_,
680  (
682  )
683  );
684 }
685 
686 
689 {
691  (
692  IOobject
693  (
694  "Ru",
695  mesh_.time().timeName(),
696  mesh_,
699  false
700  ),
701  mesh_,
703  );
704 }
705 
706 
707 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField< scalar >
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
Foam::LUsolve
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
Definition: scalarMatricesTemplates.C:215
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:350
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
typeInfo.H
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::globalIndex::localStart
label localStart() const
My local start.
Definition: globalIndexI.H:175
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList< scalar >
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::IOmapDistribute
IOmapDistribute is derived from mapDistribute and IOobject to give the mapDistribute automatic IO fun...
Definition: IOmapDistribute.H:54
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:54
Foam::MeshObject< fvMesh, Foam::GeometricMeshObject, boundaryRadiationProperties >::New
static const boundaryRadiationProperties & New(const fvMesh &mesh, Args &&... args)
Get existing or create a new MeshObject.
Definition: MeshObject.C:48
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::radiation::viewFactor::Rp
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: viewFactor.C:665
F
volVectorField F(fluid.F())
Foam::globalIndex::localSize
label localSize() const
My local size.
Definition: globalIndexI.H:187
greyDiffusiveViewFactorFixedValueFvPatchScalarField.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
surfaceFields.H
Foam::surfaceFields.
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:434
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::radiation::viewFactor::insertMatrixElements
void insertMatrixElements(const globalIndex &index, const label fromProci, const labelListList &globalFaceFaces, const scalarListList &viewFactors, scalarSquareMatrix &matrix)
Insert view factors into main matrix.
Definition: viewFactor.C:343
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField::qro
tmp< scalarField > qro(label bandI=0) const
Return external + solar load radiative heat flux.
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.C:167
C
volScalarField & C
Definition: readThermalProperties.H:102
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::radiation::boundaryRadiationProperties::emissivity
tmp< scalarField > emissivity(const label patchI, const label bandI=0, vectorField *incomingDirection=nullptr, scalarField *T=nullptr) const
Access boundary emissivity on patch.
Definition: boundaryRadiationProperties.C:111
Foam::radiation::viewFactor::initialise
void initialise()
Initialise.
Definition: viewFactor.C:55
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::radiation::radiationModel::read
virtual bool read()=0
Read radiationProperties dictionary.
Definition: radiationModel.C:210
Foam::Field< scalar >
Foam::polyPatch::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:315
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
DebugInFunction
#define DebugInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:388
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:511
patchID
label patchID
Definition: boundaryProcessorFaPatchPoints.H:5
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:58
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::radiation::boundaryRadiationProperties
Boundary radiation properties holder.
Definition: boundaryRadiationProperties.H:56
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::scalarSquareMatrix
SquareMatrix< scalar > scalarSquareMatrix
Definition: scalarMatrices.H:57
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
Foam::radiation::viewFactor::calculate
void calculate()
Solve system of equation(s)
Definition: viewFactor.C:365
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:68
Foam::refCast
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:131
Foam::maxEqOp
Definition: ops.H:80
Foam::LUBacksubstitute
void LUBacksubstitute(const scalarSquareMatrix &luMmatrix, const labelList &pivotIndices, List< Type > &source)
Definition: scalarMatricesTemplates.C:120
constants.H
Foam::LUDecompose
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
Definition: scalarMatrices.C:34
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
Foam::SquareMatrix< scalar >
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:52
Foam::UPstream::allProcs
static rangeType allProcs(const label communicator=worldComm)
Range of process indices for all processes.
Definition: UPstream.H:508
Foam::radiation::solarLoad
The solarLoad radiation model includes Sun primary hits, their reflective fluxes and diffusive sky ra...
Definition: solarLoad.H:162
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField
This boundary condition provides a grey-diffuse condition for radiative heat flux,...
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.H:88
Foam::polyBoundaryMesh::indices
labelList indices(const wordRe &matcher, const bool useGroups=true) const
Return (sorted) patch indices for all matches.
Definition: polyBoundaryMesh.C:642
Foam::List< labelListList >
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:75
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::IOList
A List of objects of type <T> with automated input and output.
Definition: IOList.H:53
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::dimTemperature
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::radiation::viewFactor::viewFactorWalls
static const word viewFactorWalls
Static name for view factor walls.
Definition: viewFactor.H:93
Foam::fieldTypes::area
const wordList area
Standard area field types (scalar, vector, tensor, etc)
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::UIndirectList
A List with indirect addressing.
Definition: faMatrix.H:60
Foam::radiation::viewFactor::Ru
virtual tmp< volScalarField::Internal > Ru() const
Source term component (constant)
Definition: viewFactor.C:688
Foam::radiation::viewFactor::read
bool read()
Read radiation properties dictionary.
Definition: viewFactor.C:331
addToRadiationRunTimeSelectionTables
#define addToRadiationRunTimeSelectionTables(model)
Definition: radiationModel.H:288
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:432
Foam::C
Graphite solid properties.
Definition: C.H:50
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::IOobject::NO_READ
Definition: IOobject.H:188
Foam::invertOneToMany
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:114
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
boundary
faceListList boundary
Definition: createBlockMesh.H:4
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:445
boundaryRadiationProperties.H
viewFactor.H
Foam::globalIndex::toGlobal
label toGlobal(const label i) const
From local to global index.
Definition: globalIndexI.H:240
Foam::IOobject::MUST_READ
Definition: IOobject.H:185