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 | Copyright (C) 2011-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 
26 #include "viewFactor.H"
27 #include "surfaceFields.H"
28 #include "constants.H"
30 #include "typeInfo.H"
33 
34 using namespace Foam::constant;
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  namespace radiation
41  {
42  defineTypeNameAndDebug(viewFactor, 0);
44  }
45 }
46 
48  = "viewFactorWall";
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
53 {
54  const polyBoundaryMesh& coarsePatches = coarseMesh_.boundaryMesh();
55 
56  selectedPatches_ = mesh_.boundaryMesh().findIndices(viewFactorWalls);
57  forAll(selectedPatches_, i)
58  {
59  const label patchI = selectedPatches_[i];
60  nLocalCoarseFaces_ += coarsePatches[patchI].size();
61  }
62 
63  if (debug)
64  {
65  Pout<< "radiation::viewFactor::initialise() Selected patches:"
66  << selectedPatches_ << endl;
67  Pout<< "radiation::viewFactor::initialise() Number of coarse faces:"
68  << nLocalCoarseFaces_ << endl;
69  }
70 
71  totalNCoarseFaces_ = nLocalCoarseFaces_;
72  reduce(totalNCoarseFaces_, sumOp<label>());
73 
74  if (debug && Pstream::master())
75  {
77  << "Total number of clusters : " << totalNCoarseFaces_ << endl;
78  }
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_, totalNCoarseFaces_, 0.0)
137  );
138 
139  if (debug)
140  {
142  << "Insert elements in the matrix..." << endl;
143  }
144 
145  for (label procI = 0; procI < Pstream::nProcs(); procI++)
146  {
147  insertMatrixElements
148  (
149  globalNumbering,
150  procI,
151  globalFaceFacesProc[procI],
152  F[procI],
153  Fmatrix_()
154  );
155  }
156 
157 
158  bool smoothing = readBool(coeffs_.lookup("smoothing"));
159  if (smoothing)
160  {
161  if (debug)
162  {
164  << "Smoothing the matrix..." << endl;
165  }
166 
167  for (label i=0; i<totalNCoarseFaces_; i++)
168  {
169  scalar sumF = 0.0;
170  for (label j=0; j<totalNCoarseFaces_; j++)
171  {
172  sumF += Fmatrix_()[i][j];
173  }
174  scalar delta = sumF - 1.0;
175  for (label j=0; j<totalNCoarseFaces_; j++)
176  {
177  Fmatrix_()[i][j] *= (1.0 - delta/(sumF + 0.001));
178  }
179  }
180  }
181 
182  constEmissivity_ = readBool(coeffs_.lookup("constantEmissivity"));
183  if (constEmissivity_)
184  {
185  CLU_.reset
186  (
188  (
189  totalNCoarseFaces_,
190  totalNCoarseFaces_,
191  0.0
192  )
193  );
194 
195  pivotIndices_.setSize(CLU_().n());
196  }
197  }
198 
199  if (this->found("useSolarLoad"))
200  {
201  this->lookup("useSolarLoad") >> useSolarLoad_;
202  }
203 
204  if (useSolarLoad_)
205  {
206  const dictionary& solarDict = this->subDict("solarLoarCoeffs");
207  solarLoad_.reset
208  (
209  new solarLoad(solarDict, T_, externalRadHeatFieldName_)
210  );
211 
212  if (solarLoad_->nBands() > 1)
213  {
215  << "Requested solar radiation with fvDOM. Using "
216  << "more thant one band for the solar load is not allowed"
217  << abort(FatalError);
218  }
219 
220  Info<< "Creating Solar Load Model " << nl;
221  }
222 }
223 
224 
225 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
226 
228 :
229  radiationModel(typeName, T),
230  finalAgglom_
231  (
232  IOobject
233  (
234  "finalAgglom",
235  mesh_.facesInstance(),
236  mesh_,
237  IOobject::MUST_READ,
238  IOobject::NO_WRITE,
239  false
240  )
241  ),
242  map_(),
243  coarseMesh_
244  (
245  IOobject
246  (
247  "coarse:" + mesh_.name(),
248  mesh_.polyMesh::instance(),
249  mesh_.time(),
250  IOobject::NO_READ,
251  IOobject::NO_WRITE
252  ),
253  mesh_,
254  finalAgglom_
255  ),
256  Qr_
257  (
258  IOobject
259  (
260  "Qr",
261  mesh_.time().timeName(),
262  mesh_,
263  IOobject::MUST_READ,
264  IOobject::AUTO_WRITE
265  ),
266  mesh_
267  ),
268  Fmatrix_(),
269  CLU_(),
270  selectedPatches_(mesh_.boundary().size(), -1),
271  totalNCoarseFaces_(0),
272  nLocalCoarseFaces_(0),
273  constEmissivity_(false),
274  iterCounter_(0),
275  pivotIndices_(0),
276  useSolarLoad_(false),
277  solarLoad_()
278 {
279  initialise();
280 }
281 
282 
284 (
285  const dictionary& dict,
286  const volScalarField& T
287 )
288 :
289  radiationModel(typeName, dict, T),
290  finalAgglom_
291  (
292  IOobject
293  (
294  "finalAgglom",
295  mesh_.facesInstance(),
296  mesh_,
299  false
300  )
301  ),
302  map_(),
303  coarseMesh_
304  (
305  IOobject
306  (
307  "coarse:" + mesh_.name(),
308  mesh_.polyMesh::instance(),
309  mesh_.time(),
312  ),
313  mesh_,
314  finalAgglom_
315  ),
316  Qr_
317  (
318  IOobject
319  (
320  "Qr",
321  mesh_.time().timeName(),
322  mesh_,
325  ),
326  mesh_
327  ),
328  Fmatrix_(),
329  CLU_(),
330  selectedPatches_(mesh_.boundary().size(), -1),
331  totalNCoarseFaces_(0),
332  nLocalCoarseFaces_(0),
333  constEmissivity_(false),
334  iterCounter_(0),
335  pivotIndices_(0),
336  useSolarLoad_(false),
337  solarLoad_()
338 {
339  initialise();
340 }
341 
342 
343 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
344 
346 {}
347 
348 
349 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
350 
352 {
353  if (radiationModel::read())
354  {
355  return true;
356  }
357  else
358  {
359  return false;
360  }
361 }
362 
363 
365 (
366  const globalIndex& globalNumbering,
367  const label procI,
368  const labelListList& globalFaceFaces,
369  const scalarListList& viewFactors,
370  scalarSquareMatrix& Fmatrix
371 )
372 {
373  forAll(viewFactors, faceI)
374  {
375  const scalarList& vf = viewFactors[faceI];
376  const labelList& globalFaces = globalFaceFaces[faceI];
377 
378  label globalI = globalNumbering.toGlobal(procI, faceI);
379  forAll(globalFaces, i)
380  {
381  Fmatrix[globalI][globalFaces[i]] = vf[i];
382  }
383  }
384 }
385 
386 
388 {
389  // Store previous iteration
390  Qr_.storePrevIter();
391 
392  if (useSolarLoad_)
393  {
394  solarLoad_->calculate();
395  }
396 
397  scalarField compactCoarseT(map_->constructSize(), 0.0);
398  scalarField compactCoarseE(map_->constructSize(), 0.0);
399  scalarField compactCoarseHo(map_->constructSize(), 0.0);
400 
401  globalIndex globalNumbering(nLocalCoarseFaces_);
402 
403  // Fill local averaged(T), emissivity(E) and external heatFlux(Ho)
404  DynamicList<scalar> localCoarseTave(nLocalCoarseFaces_);
405  DynamicList<scalar> localCoarseEave(nLocalCoarseFaces_);
406  DynamicList<scalar> localCoarseHoave(nLocalCoarseFaces_);
407 
408  const boundaryRadiationProperties& boundaryRadiation =
410 
411  forAll(selectedPatches_, i)
412  {
413  label patchID = selectedPatches_[i];
414 
415  const scalarField& Tp = T_.boundaryField()[patchID];
416  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
417 
418  fvPatchScalarField& QrPatch = Qr_.boundaryField()[patchID];
419 
421  refCast
422  <
424  >(QrPatch);
425 
426  const tmp<scalarField> teb = boundaryRadiation.emissivity(patchID);
427  const scalarField& eb = teb();
428 
429  const tmp<scalarField> tHoi = Qrp.Qro();
430  const scalarField& Hoi = tHoi();
431 
432  const polyPatch& pp = coarseMesh_.boundaryMesh()[patchID];
433  const labelList& coarsePatchFace = coarseMesh_.patchFaceMap()[patchID];
434 
435  scalarList Tave(pp.size(), 0.0);
436  scalarList Eave(Tave.size(), 0.0);
437  scalarList Hoiave(Tave.size(), 0.0);
438 
439  if (pp.size() > 0)
440  {
441  const labelList& agglom = finalAgglom_[patchID];
442  label nAgglom = max(agglom) + 1;
443 
444  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
445 
446  forAll(coarseToFine, coarseI)
447  {
448  const label coarseFaceID = coarsePatchFace[coarseI];
449  const labelList& fineFaces = coarseToFine[coarseFaceID];
450  UIndirectList<scalar> fineSf
451  (
452  sf,
453  fineFaces
454  );
455  scalar area = sum(fineSf());
456  // Temperature, emissivity and external flux area weighting
457  forAll(fineFaces, j)
458  {
459  label faceI = fineFaces[j];
460  Tave[coarseI] += (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  localCoarseTave.append(Tave);
468  localCoarseEave.append(Eave);
469  localCoarseHoave.append(Hoiave);
470  }
471 
472  // Fill the local values to distribute
473  SubList<scalar>(compactCoarseT,nLocalCoarseFaces_).assign(localCoarseTave);
474  SubList<scalar>(compactCoarseE,nLocalCoarseFaces_).assign(localCoarseEave);
476  (compactCoarseHo,nLocalCoarseFaces_).assign(localCoarseHoave);
477 
478  // Distribute data
479  map_->distribute(compactCoarseT);
480  map_->distribute(compactCoarseE);
481  map_->distribute(compactCoarseHo);
482 
483  // Distribute local global ID
484  labelList compactGlobalIds(map_->constructSize(), 0.0);
485 
486  labelList localGlobalIds(nLocalCoarseFaces_);
487 
488  for(label k = 0; k < nLocalCoarseFaces_; k++)
489  {
490  localGlobalIds[k] = globalNumbering.toGlobal(Pstream::myProcNo(), k);
491  }
492 
494  (
495  compactGlobalIds,
496  nLocalCoarseFaces_
497  ).assign(localGlobalIds);
498 
499  map_->distribute(compactGlobalIds);
500 
501  // Create global size vectors
502  scalarField T(totalNCoarseFaces_, 0.0);
503  scalarField E(totalNCoarseFaces_, 0.0);
504  scalarField QrExt(totalNCoarseFaces_, 0.0);
505 
506  // Fill lists from compact to global indexes.
507  forAll(compactCoarseT, i)
508  {
509  T[compactGlobalIds[i]] = compactCoarseT[i];
510  E[compactGlobalIds[i]] = compactCoarseE[i];
511  QrExt[compactGlobalIds[i]] = compactCoarseHo[i];
512  }
513 
517 
521 
522  // Net radiation
523  scalarField q(totalNCoarseFaces_, 0.0);
524 
525  if (Pstream::master())
526  {
527  // Variable emissivity
528  if (!constEmissivity_)
529  {
530  scalarSquareMatrix C(totalNCoarseFaces_, totalNCoarseFaces_, 0.0);
531 
532  for (label i=0; i<totalNCoarseFaces_; i++)
533  {
534  for (label j=0; j<totalNCoarseFaces_; j++)
535  {
536  scalar invEj = 1.0/E[j];
537  scalar sigmaT4 =
538  physicoChemical::sigma.value()*pow(T[j], 4.0);
539 
540  if (i==j)
541  {
542  C[i][j] = invEj - (invEj - 1.0)*Fmatrix_()[i][j];
543  q[i] += (Fmatrix_()[i][j] - 1.0)*sigmaT4 - QrExt[j];
544  }
545  else
546  {
547  C[i][j] = (1.0 - invEj)*Fmatrix_()[i][j];
548  q[i] += Fmatrix_()[i][j]*sigmaT4;
549  }
550 
551  }
552  }
553 
554  Info<< "\nSolving view factor equations..." << endl;
555 
556  // Negative coming into the fluid
557  LUsolve(C, q);
558  }
559  else //Constant emissivity
560  {
561  // Initial iter calculates CLU and chaches it
562  if (iterCounter_ == 0)
563  {
564  for (label i=0; i<totalNCoarseFaces_; i++)
565  {
566  for (label j=0; j<totalNCoarseFaces_; j++)
567  {
568  scalar invEj = 1.0/E[j];
569  if (i==j)
570  {
571  CLU_()[i][j] = invEj-(invEj-1.0)*Fmatrix_()[i][j];
572  }
573  else
574  {
575  CLU_()[i][j] = (1.0 - invEj)*Fmatrix_()[i][j];
576  }
577  }
578  }
579 
580  if (debug)
581  {
583  << "\nDecomposing C matrix..." << endl;
584  }
585 
586  LUDecompose(CLU_(), pivotIndices_);
587  }
588 
589  for (label i=0; i<totalNCoarseFaces_; i++)
590  {
591  for (label j=0; j<totalNCoarseFaces_; j++)
592  {
593  scalar sigmaT4 =
595  *pow(T[j], 4.0);
596 
597  if (i==j)
598  {
599  q[i] += (Fmatrix_()[i][j] - 1.0)*sigmaT4 - QrExt[j];
600  }
601  else
602  {
603  q[i] += Fmatrix_()[i][j]*sigmaT4;
604  }
605  }
606  }
607 
608  if (debug)
609  {
611  << "\nLU Back substitute C matrix.." << endl;
612  }
613 
614  LUBacksubstitute(CLU_(), pivotIndices_, q);
615  iterCounter_ ++;
616  }
617  }
618 
619  // Scatter q and fill Qr
622 
623 
624  label globCoarseId = 0;
625  forAll(selectedPatches_, i)
626  {
627  const label patchID = selectedPatches_[i];
628  const polyPatch& pp = mesh_.boundaryMesh()[patchID];
629  if (pp.size() > 0)
630  {
631  scalarField& Qrp = Qr_.boundaryField()[patchID];
632  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
633  const labelList& agglom = finalAgglom_[patchID];
634  label nAgglom = max(agglom)+1;
635 
636  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
637 
638  const labelList& coarsePatchFace =
639  coarseMesh_.patchFaceMap()[patchID];
640 
641  scalar heatFlux = 0.0;
642  forAll(coarseToFine, coarseI)
643  {
644  label globalCoarse =
645  globalNumbering.toGlobal(Pstream::myProcNo(), globCoarseId);
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(Qr_.boundaryField(), patchID)
663  {
664  const scalarField& Qrp = Qr_.boundaryField()[patchID];
665  const scalarField& magSf = mesh_.magSf().boundaryField()[patchID];
666  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 {
682  return tmp<volScalarField>
683  (
684  new volScalarField
685  (
686  IOobject
687  (
688  "Rp",
689  mesh_.time().timeName(),
690  mesh_,
693  false
694  ),
695  mesh_,
697  (
698  "zero",
700  0.0
701  )
702  )
703  );
704 }
705 
706 
709 {
711  (
713  (
714  IOobject
715  (
716  "Ru",
717  mesh_.time().timeName(),
718  mesh_,
721  false
722  ),
723  mesh_,
725  )
726  );
727 }
728 
729 // ************************************************************************* //
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::radiation::boundaryRadiationProperties::emissivity
tmp< scalarField > emissivity(const label patchId, const label bandI=0) const
Access boundary emissivity on patch.
Definition: boundaryRadiationProperties.C:110
Foam::LUsolve
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
Solve the matrix using LU decomposition with pivoting.
Definition: scalarMatricesTemplates.C:215
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:365
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
InfoInFunction
#define InfoInFunction
Report a information message using Foam::Info.
Definition: messageStream.H:281
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
typeInfo.H
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
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::DynamicList
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:56
Foam::invertOneToMany
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
Foam::constant
Collection of constants.
Definition: atomicConstants.C:37
Foam::IOmapDistribute
IOmapDistribute is derived from mapDistribute and IOobject to give the mapDistribute automatic IO fun...
Definition: IOmapDistribute.H:51
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:387
boundary
faceListList boundary(nPatches)
Foam::maxEqOp
Definition: ops.H:77
Foam::radiation::viewFactor::viewFactor
viewFactor(const viewFactor &)
Disallow default bitwise copy construct.
Foam::radiation::viewFactor::Rp
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: viewFactor.C:680
greyDiffusiveViewFactorFixedValueFvPatchScalarField.H
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::dimensioned::value
const Type & value() const
Return const reference to value.
Definition: dimensionedType.C:261
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::radiation::viewFactor::~viewFactor
virtual ~viewFactor()
Destructor.
Definition: viewFactor.C:345
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:118
Foam::MeshObject< fvMesh, Foam::GeometricMeshObject, boundaryRadiationProperties >::New
static const boundaryRadiationProperties & New(const fvMesh &mesh)
Definition: MeshObject.C:44
Foam::constant::physicoChemical::F
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
C
volScalarField & C
Definition: readThermalProperties.H:104
Foam::radiation::viewFactor::initialise
void initialise()
Initialise.
Definition: viewFactor.C:52
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:98
Foam::radiation::radiationModel::read
virtual bool read()=0
Read radiationProperties dictionary.
Definition: radiationModel.C:197
Foam::polyBoundaryMesh::findIndices
labelList findIndices(const keyType &, const bool usePatchGroups=true) const
Return patch indices for all matches. Optionally matches patchGroups.
Definition: polyBoundaryMesh.C:573
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::polyPatch::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:302
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:87
Foam::IOobject::NO_READ
@ NO_READ
Definition: IOobject.H:111
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField::Qro
tmp< scalarField > Qro() const
Return external + solar load radiative heat flux.
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.C:145
Foam::UList::assign
void assign(const UList< T > &)
Assign elements to those from UList.
Definition: UList.C:37
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:55
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:73
radiation
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Foam::radiation::boundaryRadiationProperties
Definition: boundaryRadiationProperties.H:56
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::scalarSquareMatrix
SquareMatrix< scalar > scalarSquareMatrix
Definition: scalarMatrices.H:55
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
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::radiation::viewFactor::calculate
void calculate()
Solve system of equation(s)
Definition: viewFactor.C:387
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
sf
volScalarField sf(fieldObject, mesh)
Foam::globalIndex
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
Foam::refCast
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
Foam::LUBacksubstitute
void LUBacksubstitute(const scalarSquareMatrix &luMmatrix, const labelList &pivotIndices, List< Type > &source)
LU back-substitution with given source, returning the solution.
Definition: scalarMatricesTemplates.C:120
constants.H
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::LUDecompose
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
Definition: scalarMatrices.C:32
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:399
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:282
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::SquareMatrix< scalar >
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
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:49
Foam::radiation::solarLoad
The solar load radiation model includes Sun primary hits, their reflective fluxes and diffusive sky r...
Definition: solarLoad.H:80
Foam::sumOp
Definition: ops.H:162
Foam::radiation::greyDiffusiveViewFactorFixedValueFvPatchScalarField
This boundary condition provides a grey-diffuse condition for radiative heat flux,...
Definition: greyDiffusiveViewFactorFixedValueFvPatchScalarField.H:86
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::List< labelListList >
Foam::radiation::radiationModel
Top level model for radiation modelling.
Definition: radiationModel.H:73
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::IOList
A List of objects of type <T> with automated input and output.
Definition: IOList.H:50
Foam::PtrList::size
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:333
timeName
word timeName
Definition: getTimeIndex.H:3
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:52
Foam::radiation::viewFactor::viewFactorWalls
static const word viewFactorWalls
Static name for view factor walls.
Definition: viewFactor.H:76
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::radiation::viewFactor::Ru
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant)
Definition: viewFactor.C:708
Foam::radiation::viewFactor::read
bool read()
Read radiation properties dictionary.
Definition: viewFactor.C:351
Foam::readBool
bool readBool(Istream &)
Definition: boolIO.C:60
addToRadiationRunTimeSelectionTables
#define addToRadiationRunTimeSelectionTables(model)
Definition: radiationModel.H:263
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
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:417
Foam::C
Graphite solid properties.
Definition: C.H:57
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
boundaryRadiationProperties.H
viewFactor.H
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::globalIndex::toGlobal
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress