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-2019 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 
77  if (debug && Pstream::master())
78  {
80  << "Total number of clusters : " << totalNCoarseFaces_ << endl;
81  }
82 
83  map_.reset
84  (
85  new IOmapDistribute
86  (
87  IOobject
88  (
89  "mapDist",
90  mesh_.facesInstance(),
91  mesh_,
94  false
95  )
96  )
97  );
98 
99  scalarListIOList FmyProc
100  (
101  IOobject
102  (
103  "F",
104  mesh_.facesInstance(),
105  mesh_,
108  false
109  )
110  );
111 
112  labelListIOList globalFaceFaces
113  (
114  IOobject
115  (
116  "globalFaceFaces",
117  mesh_.facesInstance(),
118  mesh_,
121  false
122  )
123  );
124 
125  List<labelListList> globalFaceFacesProc(Pstream::nProcs());
126  globalFaceFacesProc[Pstream::myProcNo()] = globalFaceFaces;
127  Pstream::gatherList(globalFaceFacesProc);
128 
130  F[Pstream::myProcNo()] = FmyProc;
132 
133  globalIndex globalNumbering(nLocalCoarseFaces_);
134 
135  if (Pstream::master())
136  {
137  Fmatrix_.reset
138  (
139  new scalarSquareMatrix(totalNCoarseFaces_, Zero)
140  );
141 
142  if (debug)
143  {
145  << "Insert elements in the matrix..." << endl;
146  }
147 
148  for (label procI = 0; procI < Pstream::nProcs(); procI++)
149  {
150  insertMatrixElements
151  (
152  globalNumbering,
153  procI,
154  globalFaceFacesProc[procI],
155  F[procI],
156  Fmatrix_()
157  );
158  }
159 
160 
161  if (coeffs_.get<bool>("smoothing"))
162  {
163  if (debug)
164  {
166  << "Smoothing the matrix..." << endl;
167  }
168 
169  for (label i=0; i<totalNCoarseFaces_; i++)
170  {
171  scalar sumF = 0.0;
172  for (label j=0; j<totalNCoarseFaces_; j++)
173  {
174  sumF += Fmatrix_()(i, j);
175  }
176 
177  const scalar delta = sumF - 1.0;
178  for (label j=0; j<totalNCoarseFaces_; j++)
179  {
180  Fmatrix_()(i, j) *= (1.0 - delta/(sumF + 0.001));
181  }
182  }
183  }
184 
185  coeffs_.readEntry("constantEmissivity", constEmissivity_);
186  if (constEmissivity_)
187  {
188  CLU_.reset
189  (
190  new scalarSquareMatrix(totalNCoarseFaces_, Zero)
191  );
192 
193  pivotIndices_.setSize(CLU_().m());
194  }
195  }
196 
197  coeffs_.readIfPresent("useSolarLoad", useSolarLoad_);
198 
199  if (useSolarLoad_)
200  {
201  const dictionary& solarDict = this->subDict("solarLoadCoeffs");
202  solarLoad_.reset(new solarLoad(solarDict, T_));
203 
204  if (solarLoad_->nBands() != nBands_)
205  {
207  << "Solar radiation and view factor band numbers "
208  << "are different"
209  << abort(FatalError);
210  }
211 
212  Info<< "Creating Solar Load Model " << nl;
213  }
214 }
215 
216 
217 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
218 
219 Foam::radiation::viewFactor::viewFactor(const volScalarField& T)
220 :
221  radiationModel(typeName, T),
222  finalAgglom_
223  (
224  IOobject
225  (
226  "finalAgglom",
227  mesh_.facesInstance(),
228  mesh_,
229  IOobject::MUST_READ,
230  IOobject::NO_WRITE,
231  false
232  )
233  ),
234  map_(),
235  coarseMesh_
236  (
237  IOobject
238  (
239  "coarse:" + mesh_.name(),
240  mesh_.polyMesh::instance(),
241  mesh_.time(),
242  IOobject::NO_READ,
243  IOobject::NO_WRITE,
244  false
245  ),
246  mesh_,
247  finalAgglom_
248  ),
249  qr_
250  (
251  IOobject
252  (
253  "qr",
254  mesh_.time().timeName(),
255  mesh_,
256  IOobject::MUST_READ,
257  IOobject::AUTO_WRITE
258  ),
259  mesh_
260  ),
261  Fmatrix_(),
262  CLU_(),
263  selectedPatches_(mesh_.boundary().size(), -1),
264  totalNCoarseFaces_(0),
265  nLocalCoarseFaces_(0),
266  constEmissivity_(false),
267  iterCounter_(0),
268  pivotIndices_(0),
269  useSolarLoad_(false),
270  solarLoad_(),
271  nBands_(coeffs_.lookupOrDefault<label>("nBands", 1))
272 {
273  initialise();
274 }
275 
276 
277 Foam::radiation::viewFactor::viewFactor
278 (
279  const dictionary& dict,
280  const volScalarField& T
281 )
282 :
283  radiationModel(typeName, dict, T),
284  finalAgglom_
285  (
286  IOobject
287  (
288  "finalAgglom",
289  mesh_.facesInstance(),
290  mesh_,
293  false
294  )
295  ),
296  map_(),
297  coarseMesh_
298  (
299  IOobject
300  (
301  "coarse:" + mesh_.name(),
302  mesh_.polyMesh::instance(),
303  mesh_.time(),
306  false
307  ),
308  mesh_,
309  finalAgglom_
310  ),
311  qr_
312  (
313  IOobject
314  (
315  "qr",
316  mesh_.time().timeName(),
317  mesh_,
320  ),
321  mesh_
322  ),
323  Fmatrix_(),
324  CLU_(),
325  selectedPatches_(mesh_.boundary().size(), -1),
326  totalNCoarseFaces_(0),
327  nLocalCoarseFaces_(0),
328  constEmissivity_(false),
329  iterCounter_(0),
330  pivotIndices_(0),
331  useSolarLoad_(false),
332  solarLoad_(),
333  nBands_(coeffs_.lookupOrDefault<label>("nBands", 1))
334 {
335  initialise();
336 }
337 
338 
339 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
340 
342 {
343  if (radiationModel::read())
344  {
345  return true;
346  }
347 
348  return false;
349 }
350 
351 
353 (
354  const globalIndex& globalNumbering,
355  const label procI,
356  const labelListList& globalFaceFaces,
357  const scalarListList& viewFactors,
358  scalarSquareMatrix& Fmatrix
359 )
360 {
361  forAll(viewFactors, faceI)
362  {
363  const scalarList& vf = viewFactors[faceI];
364  const labelList& globalFaces = globalFaceFaces[faceI];
365 
366  label globalI = globalNumbering.toGlobal(procI, faceI);
367  forAll(globalFaces, i)
368  {
369  Fmatrix[globalI][globalFaces[i]] = vf[i];
370  }
371  }
372 }
373 
374 
376 {
377  // Store previous iteration
378  qr_.storePrevIter();
379 
380  if (useSolarLoad_)
381  {
382  solarLoad_->calculate();
383  }
384 
385  // Net radiation
386  scalarField q(totalNCoarseFaces_, 0.0);
387  volScalarField::Boundary& qrBf = qr_.boundaryFieldRef();
388 
389  globalIndex globalNumbering(nLocalCoarseFaces_);
390 
391  const boundaryRadiationProperties& boundaryRadiation =
393 
394  for (label bandI = 0; bandI < nBands_; bandI++)
395  {
396  scalarField compactCoarseT4(map_->constructSize(), 0.0);
397  scalarField compactCoarseE(map_->constructSize(), 0.0);
398  scalarField compactCoarseHo(map_->constructSize(), 0.0);
399 
400  // Fill local averaged(T), emissivity(E) and external heatFlux(Ho)
401  DynamicList<scalar> localCoarseT4ave(nLocalCoarseFaces_);
402  DynamicList<scalar> localCoarseEave(nLocalCoarseFaces_);
403  DynamicList<scalar> localCoarseHoave(nLocalCoarseFaces_);
404 
405  forAll(selectedPatches_, i)
406  {
407  label patchID = selectedPatches_[i];
408 
409  const scalarField& Tp = T_.boundaryField()[patchID];
410  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
411 
412  fvPatchScalarField& qrPatch = qrBf[patchID];
413 
415  refCast
416  <
418  >(qrPatch);
419 
420  const tmp<scalarField> teb =
421  boundaryRadiation.emissivity(patchID, bandI);
422 
423  const scalarField& eb = teb();
424 
425  const tmp<scalarField> tHoi = qrp.qro(bandI);
426  const scalarField& Hoi = tHoi();
427 
428  const polyPatch& pp = coarseMesh_.boundaryMesh()[patchID];
429 
430  const labelList& coarsePatchFace =
431  coarseMesh_.patchFaceMap()[patchID];
432 
433  scalarList T4ave(pp.size(), 0.0);
434  scalarList Eave(pp.size(), 0.0);
435  scalarList Hoiave(pp.size(), 0.0);
436 
437  if (pp.size() > 0)
438  {
439  const labelList& agglom = finalAgglom_[patchID];
440  label nAgglom = max(agglom) + 1;
441 
442  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
443 
444  forAll(coarseToFine, coarseI)
445  {
446  const label coarseFaceID = coarsePatchFace[coarseI];
447  const labelList& fineFaces = coarseToFine[coarseFaceID];
448  UIndirectList<scalar> fineSf
449  (
450  sf,
451  fineFaces
452  );
453 
454  const scalar area = sum(fineSf());
455 
456  // Temperature, emissivity and external flux area weighting
457  forAll(fineFaces, j)
458  {
459  label facei = fineFaces[j];
460  T4ave[coarseI] += (pow4(Tp[facei])*sf[facei])/area;
461  Eave[coarseI] += (eb[facei]*sf[facei])/area;
462  Hoiave[coarseI] += (Hoi[facei]*sf[facei])/area;
463  }
464  }
465  }
466 
467  localCoarseT4ave.append(T4ave);
468  localCoarseEave.append(Eave);
469  localCoarseHoave.append(Hoiave);
470 
471  }
472 
473  // Fill the local values to distribute
474  SubList<scalar>(compactCoarseT4, nLocalCoarseFaces_) =
475  localCoarseT4ave;
476  SubList<scalar>(compactCoarseE, nLocalCoarseFaces_) = localCoarseEave;
477  SubList<scalar>(compactCoarseHo, nLocalCoarseFaces_) =
478  localCoarseHoave;
479 
480  // Distribute data
481  map_->distribute(compactCoarseT4);
482  map_->distribute(compactCoarseE);
483  map_->distribute(compactCoarseHo);
484 
485  // Distribute local global ID
486  labelList compactGlobalIds(map_->constructSize(), Zero);
487 
489  (
490  compactGlobalIds,
491  nLocalCoarseFaces_
492  ) = identity
493  (
494  globalNumbering.localSize(),
495  globalNumbering.localStart()
496  );
497 
498  map_->distribute(compactGlobalIds);
499 
500  // Create global size vectors
501  scalarField T4(totalNCoarseFaces_, 0.0);
502  scalarField E(totalNCoarseFaces_, 0.0);
503  scalarField qrExt(totalNCoarseFaces_, 0.0);
504 
505  // Fill lists from compact to global indexes.
506  forAll(compactCoarseT4, i)
507  {
508  T4[compactGlobalIds[i]] = compactCoarseT4[i];
509  E[compactGlobalIds[i]] = compactCoarseE[i];
510  qrExt[compactGlobalIds[i]] = compactCoarseHo[i];
511  }
512 
516 
520 
521  if (Pstream::master())
522  {
523  // Variable emissivity
524  if (!constEmissivity_)
525  {
526  scalarSquareMatrix C(totalNCoarseFaces_, 0.0);
527 
528  for (label i=0; i<totalNCoarseFaces_; i++)
529  {
530  for (label j=0; j<totalNCoarseFaces_; j++)
531  {
532  const scalar invEj = 1.0/E[j];
533  const scalar sigmaT4 =
534  physicoChemical::sigma.value()*T4[j];
535 
536  if (i==j)
537  {
538  C(i, j) = invEj - (invEj - 1.0)*Fmatrix_()(i, j);
539  q[i] +=
540  (Fmatrix_()(i, j) - 1.0)*sigmaT4 + qrExt[j];
541  }
542  else
543  {
544  C(i, j) = (1.0 - invEj)*Fmatrix_()(i, j);
545  q[i] += Fmatrix_()(i, j)*sigmaT4;
546  }
547 
548  }
549  }
550 
551  Info<< "Solving view factor equations for band :"
552  << bandI << endl;
553 
554  // Negative coming into the fluid
555  LUsolve(C, q);
556  }
557  else //Constant emissivity
558  {
559  // Initial iter calculates CLU and caches it
560  if (iterCounter_ == 0)
561  {
562  for (label i=0; i<totalNCoarseFaces_; i++)
563  {
564  for (label j=0; j<totalNCoarseFaces_; j++)
565  {
566  const scalar invEj = 1.0/E[j];
567  if (i==j)
568  {
569  CLU_()(i, j) =
570  invEj-(invEj-1.0)*Fmatrix_()(i, j);
571  }
572  else
573  {
574  CLU_()(i, j) = (1.0 - invEj)*Fmatrix_()(i, j);
575  }
576  }
577  }
578 
579  if (debug)
580  {
582  << "\nDecomposing C matrix..." << endl;
583  }
584 
585  LUDecompose(CLU_(), pivotIndices_);
586  }
587 
588  for (label i=0; i<totalNCoarseFaces_; i++)
589  {
590  for (label j=0; j<totalNCoarseFaces_; j++)
591  {
592  const scalar sigmaT4 =
594 
595  if (i==j)
596  {
597  q[i] +=
598  (Fmatrix_()(i, j) - 1.0)*sigmaT4 + qrExt[j];
599  }
600  else
601  {
602  q[i] += Fmatrix_()(i, j)*sigmaT4;
603  }
604  }
605  }
606 
607 
608  Info<< "Solving view factor equations for band : "
609  << bandI << endl;
610 
611 
612  LUBacksubstitute(CLU_(), pivotIndices_, q);
613  iterCounter_ ++;
614  }
615  }
616 
617  }
618  // Scatter q and fill qr
621 
622  label globCoarseId = 0;
623  for (const label patchID : selectedPatches_)
624  {
625  const polyPatch& pp = mesh_.boundaryMesh()[patchID];
626 
627  if (pp.size() > 0)
628  {
629  scalarField& qrp = qrBf[patchID];
630  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
631  const labelList& agglom = finalAgglom_[patchID];
632  label nAgglom = max(agglom)+1;
633 
634  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
635 
636  const labelList& coarsePatchFace =
637  coarseMesh_.patchFaceMap()[patchID];
638 
639  scalar heatFlux = 0.0;
640  forAll(coarseToFine, coarseI)
641  {
642  label globalCoarse =
643  globalNumbering.toGlobal
644  (Pstream::myProcNo(), globCoarseId);
645 
646  const label coarseFaceID = coarsePatchFace[coarseI];
647  const labelList& fineFaces = coarseToFine[coarseFaceID];
648  forAll(fineFaces, k)
649  {
650  label faceI = fineFaces[k];
651 
652  qrp[faceI] = q[globalCoarse];
653  heatFlux += qrp[faceI]*sf[faceI];
654  }
655  globCoarseId ++;
656  }
657  }
658  }
659 
660  if (debug)
661  {
662  forAll(qrBf, patchID)
663  {
664  const scalarField& qrp = qrBf[patchID];
665  const scalarField& magSf = mesh_.magSf().boundaryField()[patchID];
666  const scalar heatFlux = gSum(qrp*magSf);
667 
669  << "Total heat transfer rate at patch: "
670  << patchID << " "
671  << heatFlux << endl;
672  }
673  }
674 
675  // Relax qr if necessary
676  qr_.relax();
677 }
678 
679 
681 {
683  (
684  IOobject
685  (
686  "Rp",
687  mesh_.time().timeName(),
688  mesh_,
691  false
692  ),
693  mesh_,
695  (
697  )
698  );
699 }
700 
701 
704 {
706  (
707  IOobject
708  (
709  "Ru",
710  mesh_.time().timeName(),
711  mesh_,
714  false
715  ),
716  mesh_,
718  );
719 }
720 
721 
722 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::fvPatchField< scalar >
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
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:104
InfoInFunction
#define InfoInFunction
Report an information message using Foam::Info.
Definition: messageStream.H:316
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:129
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
typeInfo.H
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::globalIndex::localStart
label localStart() const
My local start.
Definition: globalIndexI.H:112
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
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:53
Foam::UPstream::nProcs
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:426
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::radiation::viewFactor::Rp
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: viewFactor.C:680
F
volVectorField F(fluid.F())
Foam::globalIndex::localSize
label localSize() const
My local size.
Definition: globalIndexI.H:124
greyDiffusiveViewFactorFixedValueFvPatchScalarField.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
surfaceFields.H
Foam::surfaceFields.
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:404
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
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:353
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:290
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:54
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::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:62
Foam::Field< scalar >
Foam::polyPatch::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:305
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::DynamicList::append
DynamicList< T, SizeMin > & append(const T &val)
Append an element to the end of this list.
Definition: DynamicListI.H:472
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:43
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:121
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:375
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
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::polyBoundaryMesh::indices
labelList indices(const keyType &key, const bool useGroups=true) const
Return patch indices for all matches.
Definition: polyBoundaryMesh.C:669
Foam::LUDecompose
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
Definition: scalarMatrices.C:34
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:444
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:438
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:355
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::radiation::solarLoad
The solar load radiation model includes Sun primary hits, their reflective fluxes and diffusive sky r...
Definition: solarLoad.H:79
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField
This boundary condition provides a grey-diffuse condition for radiative heat flux,...
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.H:88
Foam::List< labelListList >
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:75
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
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:55
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::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:109
Foam::radiation::viewFactor::Ru
virtual tmp< volScalarField::Internal > Ru() const
Source term component (constant)
Definition: viewFactor.C:703
Foam::radiation::viewFactor::read
bool read()
Read radiation properties dictionary.
Definition: viewFactor.C:341
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:123
Foam::invertOneToMany
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:113
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
boundary
faceListList boundary
Definition: createBlockMesh.H:4
boundaryRadiationProperties.H
viewFactor.H
Foam::globalIndex::toGlobal
label toGlobal(const label i) const
From local to global index.
Definition: globalIndexI.H:164
Foam::IOobject::MUST_READ
Definition: IOobject.H:120