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-2016 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"
32 
33 using namespace Foam::constant;
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  namespace radiation
40  {
41  defineTypeNameAndDebug(viewFactor, 0);
43  }
44 }
45 
46 
47 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
48 
49 void Foam::radiation::viewFactor::initialise()
50 {
51  const polyBoundaryMesh& coarsePatches = coarseMesh_.boundaryMesh();
52  const volScalarField::Boundary& Qrp = Qr_.boundaryField();
53 
54  label count = 0;
55  forAll(Qrp, patchi)
56  {
57  //const polyPatch& pp = mesh_.boundaryMesh()[patchi];
58  const fvPatchScalarField& QrPatchi = Qrp[patchi];
59 
60  if ((isA<fixedValueFvPatchScalarField>(QrPatchi)))
61  {
62  selectedPatches_[count] = QrPatchi.patch().index();
63  nLocalCoarseFaces_ += coarsePatches[patchi].size();
64  count++;
65  }
66  }
67 
68  selectedPatches_.resize(count--);
69 
70  if (debug)
71  {
72  Pout<< "radiation::viewFactor::initialise() Selected patches:"
73  << selectedPatches_ << endl;
74  Pout<< "radiation::viewFactor::initialise() Number of coarse faces:"
75  << nLocalCoarseFaces_ << endl;
76  }
77 
78  totalNCoarseFaces_ = nLocalCoarseFaces_;
79  reduce(totalNCoarseFaces_, sumOp<label>());
80 
81  if (debug && Pstream::master())
82  {
84  << "Total number of clusters : " << totalNCoarseFaces_ << endl;
85  }
86 
87  labelListIOList subMap
88  (
89  IOobject
90  (
91  "subMap",
92  mesh_.facesInstance(),
93  mesh_,
96  false
97  )
98  );
99 
100  labelListIOList constructMap
101  (
102  IOobject
103  (
104  "constructMap",
105  mesh_.facesInstance(),
106  mesh_,
109  false
110  )
111  );
112 
113  IOList<label> consMapDim
114  (
115  IOobject
116  (
117  "constructMapDim",
118  mesh_.facesInstance(),
119  mesh_,
122  false
123  )
124  );
125 
126  map_.reset
127  (
128  new mapDistribute
129  (
130  consMapDim[0],
131  Xfer<labelListList>(subMap),
132  Xfer<labelListList>(constructMap)
133  )
134  );
135 
136  scalarListIOList FmyProc
137  (
138  IOobject
139  (
140  "F",
141  mesh_.facesInstance(),
142  mesh_,
145  false
146  )
147  );
148 
149  labelListIOList globalFaceFaces
150  (
151  IOobject
152  (
153  "globalFaceFaces",
154  mesh_.facesInstance(),
155  mesh_,
158  false
159  )
160  );
161 
162  List<labelListList> globalFaceFacesProc(Pstream::nProcs());
163  globalFaceFacesProc[Pstream::myProcNo()] = globalFaceFaces;
164  Pstream::gatherList(globalFaceFacesProc);
165 
166  List<scalarListList> F(Pstream::nProcs());
167  F[Pstream::myProcNo()] = FmyProc;
169 
170  globalIndex globalNumbering(nLocalCoarseFaces_);
171 
172  if (Pstream::master())
173  {
174  Fmatrix_.reset
175  (
176  new scalarSquareMatrix(totalNCoarseFaces_, 0.0)
177  );
178 
179  if (debug)
180  {
182  << "Insert elements in the matrix..." << endl;
183  }
184 
185  for (label proci = 0; proci < Pstream::nProcs(); proci++)
186  {
187  insertMatrixElements
188  (
189  globalNumbering,
190  proci,
191  globalFaceFacesProc[proci],
192  F[proci],
193  Fmatrix_()
194  );
195  }
196 
197 
198  bool smoothing = readBool(coeffs_.lookup("smoothing"));
199  if (smoothing)
200  {
201  if (debug)
202  {
204  << "Smoothing the matrix..." << endl;
205  }
206 
207  for (label i=0; i<totalNCoarseFaces_; i++)
208  {
209  scalar sumF = 0.0;
210  for (label j=0; j<totalNCoarseFaces_; j++)
211  {
212  sumF += Fmatrix_()(i, j);
213  }
214  scalar delta = sumF - 1.0;
215  for (label j=0; j<totalNCoarseFaces_; j++)
216  {
217  Fmatrix_()(i, j) *= (1.0 - delta/(sumF + 0.001));
218  }
219  }
220  }
221 
222  constEmissivity_ = readBool(coeffs_.lookup("constantEmissivity"));
223  if (constEmissivity_)
224  {
225  CLU_.reset
226  (
227  new scalarSquareMatrix(totalNCoarseFaces_, 0.0)
228  );
229 
230  pivotIndices_.setSize(CLU_().m());
231  }
232  }
233 }
234 
235 
236 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
237 
238 Foam::radiation::viewFactor::viewFactor(const volScalarField& T)
239 :
240  radiationModel(typeName, T),
241  finalAgglom_
242  (
243  IOobject
244  (
245  "finalAgglom",
246  mesh_.facesInstance(),
247  mesh_,
248  IOobject::MUST_READ,
249  IOobject::NO_WRITE,
250  false
251  )
252  ),
253  map_(),
254  coarseMesh_
255  (
256  IOobject
257  (
258  mesh_.name(),
259  mesh_.polyMesh::instance(),
260  mesh_.time(),
261  IOobject::NO_READ,
262  IOobject::NO_WRITE
263  ),
264  mesh_,
265  finalAgglom_
266  ),
267  Qr_
268  (
269  IOobject
270  (
271  "Qr",
272  mesh_.time().timeName(),
273  mesh_,
274  IOobject::MUST_READ,
275  IOobject::AUTO_WRITE
276  ),
277  mesh_
278  ),
279  Fmatrix_(),
280  CLU_(),
281  selectedPatches_(mesh_.boundary().size(), -1),
282  totalNCoarseFaces_(0),
283  nLocalCoarseFaces_(0),
284  constEmissivity_(false),
285  iterCounter_(0),
286  pivotIndices_(0)
287 {
288  initialise();
289 }
290 
291 
292 Foam::radiation::viewFactor::viewFactor
293 (
294  const dictionary& dict,
295  const volScalarField& T
296 )
297 :
298  radiationModel(typeName, dict, T),
299  finalAgglom_
300  (
301  IOobject
302  (
303  "finalAgglom",
305  mesh_,
308  false
309  )
310  ),
311  map_(),
312  coarseMesh_
313  (
314  IOobject
315  (
316  mesh_.name(),
317  mesh_.polyMesh::instance(),
318  mesh_.time(),
321  ),
322  mesh_,
323  finalAgglom_
324  ),
325  Qr_
326  (
327  IOobject
328  (
329  "Qr",
330  mesh_.time().timeName(),
331  mesh_,
334  ),
335  mesh_
336  ),
337  Fmatrix_(),
338  CLU_(),
339  selectedPatches_(mesh_.boundary().size(), -1),
340  totalNCoarseFaces_(0),
341  nLocalCoarseFaces_(0),
342  constEmissivity_(false),
343  iterCounter_(0),
344  pivotIndices_(0)
345 {
346  initialise();
347 }
348 
349 
350 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
351 
353 {}
354 
355 
356 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
357 
359 {
360  if (radiationModel::read())
361  {
362  return true;
363  }
364  else
365  {
366  return false;
367  }
368 }
369 
370 
371 void Foam::radiation::viewFactor::insertMatrixElements
372 (
373  const globalIndex& globalNumbering,
374  const label proci,
375  const labelListList& globalFaceFaces,
376  const scalarListList& viewFactors,
377  scalarSquareMatrix& Fmatrix
378 )
379 {
380  forAll(viewFactors, facei)
381  {
382  const scalarList& vf = viewFactors[facei];
383  const labelList& globalFaces = globalFaceFaces[facei];
384 
385  label globalI = globalNumbering.toGlobal(proci, facei);
386  forAll(globalFaces, i)
387  {
388  Fmatrix[globalI][globalFaces[i]] = vf[i];
389  }
390  }
391 }
392 
393 
395 {
396  // Store previous iteration
397  Qr_.storePrevIter();
398 
399  scalarField compactCoarseT(map_->constructSize(), 0.0);
400  scalarField compactCoarseE(map_->constructSize(), 0.0);
401  scalarField compactCoarseHo(map_->constructSize(), 0.0);
402 
403  globalIndex globalNumbering(nLocalCoarseFaces_);
404 
405  // Fill local averaged(T), emissivity(E) and external heatFlux(Ho)
406  DynamicList<scalar> localCoarseTave(nLocalCoarseFaces_);
407  DynamicList<scalar> localCoarseEave(nLocalCoarseFaces_);
408  DynamicList<scalar> localCoarseHoave(nLocalCoarseFaces_);
409 
410  volScalarField::Boundary& QrBf = Qr_.boundaryFieldRef();
411 
412  forAll(selectedPatches_, i)
413  {
414  label patchID = selectedPatches_[i];
415 
416  const scalarField& Tp = T_.boundaryField()[patchID];
417  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
418 
419  fvPatchScalarField& QrPatch = QrBf[patchID];
420 
422  refCast
423  <
425  >(QrPatch);
426 
427  const scalarList eb = Qrp.emissivity();
428 
429  const scalarList& Hoi = Qrp.Qro();
430 
431  const polyPatch& pp = coarseMesh_.boundaryMesh()[patchID];
432  const labelList& coarsePatchFace = coarseMesh_.patchFaceMap()[patchID];
433 
434  scalarList Tave(pp.size(), 0.0);
435  scalarList Eave(Tave.size(), 0.0);
436  scalarList Hoiave(Tave.size(), 0.0);
437 
438  if (pp.size() > 0)
439  {
440  const labelList& agglom = finalAgglom_[patchID];
441  label nAgglom = max(agglom) + 1;
442 
443  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
444 
445  forAll(coarseToFine, coarseI)
446  {
447  const label coarseFaceID = coarsePatchFace[coarseI];
448  const labelList& fineFaces = coarseToFine[coarseFaceID];
449  UIndirectList<scalar> fineSf
450  (
451  sf,
452  fineFaces
453  );
454  scalar area = sum(fineSf());
455  // Temperature, emissivity and external flux area weighting
456  forAll(fineFaces, j)
457  {
458  label facei = fineFaces[j];
459  Tave[coarseI] += (Tp[facei]*sf[facei])/area;
460  Eave[coarseI] += (eb[facei]*sf[facei])/area;
461  Hoiave[coarseI] += (Hoi[facei]*sf[facei])/area;
462  }
463  }
464  }
465 
466  localCoarseTave.append(Tave);
467  localCoarseEave.append(Eave);
468  localCoarseHoave.append(Hoiave);
469  }
470 
471  // Fill the local values to distribute
472  SubList<scalar>(compactCoarseT,nLocalCoarseFaces_) = localCoarseTave;
473  SubList<scalar>(compactCoarseE,nLocalCoarseFaces_) = localCoarseEave;
474  SubList<scalar>(compactCoarseHo,nLocalCoarseFaces_) = localCoarseHoave;
475 
476  // Distribute data
477  map_->distribute(compactCoarseT);
478  map_->distribute(compactCoarseE);
479  map_->distribute(compactCoarseHo);
480 
481  // Distribute local global ID
482  labelList compactGlobalIds(map_->constructSize(), 0.0);
483 
484  labelList localGlobalIds(nLocalCoarseFaces_);
485 
486  for(label k = 0; k < nLocalCoarseFaces_; k++)
487  {
488  localGlobalIds[k] = globalNumbering.toGlobal(Pstream::myProcNo(), k);
489  }
490 
492  (
493  compactGlobalIds,
494  nLocalCoarseFaces_
495  ) = localGlobalIds;
496 
497  map_->distribute(compactGlobalIds);
498 
499  // Create global size vectors
500  scalarField T(totalNCoarseFaces_, 0.0);
501  scalarField E(totalNCoarseFaces_, 0.0);
502  scalarField QrExt(totalNCoarseFaces_, 0.0);
503 
504  // Fill lists from compact to global indexes.
505  forAll(compactCoarseT, i)
506  {
507  T[compactGlobalIds[i]] = compactCoarseT[i];
508  E[compactGlobalIds[i]] = compactCoarseE[i];
509  QrExt[compactGlobalIds[i]] = compactCoarseHo[i];
510  }
511 
515 
519 
520  // Net radiation
521  scalarField q(totalNCoarseFaces_, 0.0);
522 
523  if (Pstream::master())
524  {
525  // Variable emissivity
526  if (!constEmissivity_)
527  {
528  scalarSquareMatrix C(totalNCoarseFaces_, 0.0);
529 
530  for (label i=0; i<totalNCoarseFaces_; i++)
531  {
532  for (label j=0; j<totalNCoarseFaces_; j++)
533  {
534  scalar invEj = 1.0/E[j];
535  scalar sigmaT4 =
536  physicoChemical::sigma.value()*pow(T[j], 4.0);
537 
538  if (i==j)
539  {
540  C(i, j) = invEj - (invEj - 1.0)*Fmatrix_()(i, j);
541  q[i] += (Fmatrix_()(i, j) - 1.0)*sigmaT4 - QrExt[j];
542  }
543  else
544  {
545  C(i, j) = (1.0 - invEj)*Fmatrix_()(i, j);
546  q[i] += Fmatrix_()(i, j)*sigmaT4;
547  }
548 
549  }
550  }
551 
552  Info<< "\nSolving view factor equations..." << endl;
553 
554  // Negative coming into the fluid
555  LUsolve(C, q);
556  }
557  else //Constant emissivity
558  {
559  // Initial iter calculates CLU and chaches 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  scalar invEj = 1.0/E[j];
567  if (i==j)
568  {
569  CLU_()(i, j) = invEj-(invEj-1.0)*Fmatrix_()(i, j);
570  }
571  else
572  {
573  CLU_()(i, j) = (1.0 - invEj)*Fmatrix_()(i, j);
574  }
575  }
576  }
577 
578  if (debug)
579  {
581  << "\nDecomposing C matrix..." << endl;
582  }
583 
584  LUDecompose(CLU_(), pivotIndices_);
585  }
586 
587  for (label i=0; i<totalNCoarseFaces_; i++)
588  {
589  for (label j=0; j<totalNCoarseFaces_; j++)
590  {
591  scalar sigmaT4 =
593  *pow(T[j], 4.0);
594 
595  if (i==j)
596  {
597  q[i] += (Fmatrix_()(i, j) - 1.0)*sigmaT4 - QrExt[j];
598  }
599  else
600  {
601  q[i] += Fmatrix_()(i, j)*sigmaT4;
602  }
603  }
604  }
605 
606  if (debug)
607  {
609  << "\nLU Back substitute C matrix.." << endl;
610  }
611 
612  LUBacksubstitute(CLU_(), pivotIndices_, q);
613  iterCounter_ ++;
614  }
615  }
616 
617  // Scatter q and fill Qr
620 
621  label globCoarseId = 0;
622  forAll(selectedPatches_, i)
623  {
624  const label patchID = selectedPatches_[i];
625  const polyPatch& pp = mesh_.boundaryMesh()[patchID];
626  if (pp.size() > 0)
627  {
628  scalarField& Qrp = QrBf[patchID];
629  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
630  const labelList& agglom = finalAgglom_[patchID];
631  label nAgglom = max(agglom)+1;
632 
633  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
634 
635  const labelList& coarsePatchFace =
636  coarseMesh_.patchFaceMap()[patchID];
637 
638  scalar heatFlux = 0.0;
639  forAll(coarseToFine, coarseI)
640  {
641  label globalCoarse =
642  globalNumbering.toGlobal(Pstream::myProcNo(), globCoarseId);
643  const label coarseFaceID = coarsePatchFace[coarseI];
644  const labelList& fineFaces = coarseToFine[coarseFaceID];
645  forAll(fineFaces, k)
646  {
647  label facei = fineFaces[k];
648 
649  Qrp[facei] = q[globalCoarse];
650  heatFlux += Qrp[facei]*sf[facei];
651  }
652  globCoarseId ++;
653  }
654  }
655  }
656 
657  if (debug)
658  {
659  forAll(QrBf, patchID)
660  {
661  const scalarField& Qrp = QrBf[patchID];
662  const scalarField& magSf = mesh_.magSf().boundaryField()[patchID];
663  scalar heatFlux = gSum(Qrp*magSf);
664 
666  << "Total heat transfer rate at patch: "
667  << patchID << " "
668  << heatFlux << endl;
669  }
670  }
671 
672  // Relax Qr if necessary
673  Qr_.relax();
674 }
675 
676 
678 {
679  return tmp<volScalarField>
680  (
681  new volScalarField
682  (
683  IOobject
684  (
685  "Rp",
686  mesh_.time().timeName(),
687  mesh_,
690  false
691  ),
692  mesh_,
694  (
695  "zero",
697  0.0
698  )
699  )
700  );
701 }
702 
703 
706 {
708  (
710  (
711  IOobject
712  (
713  "Ru",
714  mesh_.time().timeName(),
715  mesh_,
718  false
719  ),
720  mesh_,
722  )
723  );
724 }
725 
726 // ************************************************************************* //
Collection of constants.
Foam::surfaceFields.
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
#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
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
IOList< labelList > labelListIOList
Label container classes.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:778
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
const volScalarField & T_
Reference to the temperature field.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:114
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:417
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:411
scalarField emissivity() const
Calculate corresponding emissivity field.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
const Type & value() const
Return const reference to value.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
label k
Boltzmann constant.
bool readBool(Istream &)
Definition: boolIO.C:60
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Macros for easy insertion into run-time selection tables.
const labelListList & patchFaceMap() const
From patchFace on this back to original mesh or agglomeration.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
A List obtained as a section of another List.
Definition: SubList.H:53
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: viewFactor.C:677
Type gSum(const FieldField< Field, Type > &f)
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
virtual bool read()=0
Read radiationProperties dictionary.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
fvPatchField< scalar > fvPatchScalarField
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
const fvMesh & mesh_
Reference to the mesh database.
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:292
Top level model for radiation modelling.
word timeName
Definition: getTimeIndex.H:3
faceListList boundary(nPatches)
void calculate()
Solve system of equation(s)
Definition: viewFactor.C:394
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
virtual ~viewFactor()
Destructor.
Definition: viewFactor.C:352
defineTypeNameAndDebug(combustionModel, 0)
bool read()
Read radiation properties dictionary.
Definition: viewFactor.C:358
volScalarField & C
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
volScalarField sf(fieldObject, mesh)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant)
Definition: viewFactor.C:705
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
This boundary condition provides a grey-diffuse condition for radiative heat flux, Qr, for use with the view factor model.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:399
label patchi
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
void relax(const scalar alpha)
Relax field (for steady-state solution).
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A List with indirect addressing.
Definition: fvMatrix.H:106
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensionedScalar pow4(const dimensionedScalar &ds)
IOList< scalarList > scalarListIOList
Scalar container classes.
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: PtrList.H:54
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
#define addToRadiationRunTimeSelectionTables(model)
SquareMatrix< scalar > scalarSquareMatrix
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:545
void storePrevIter() const
Store the field as the previous iteration value.
const word & name() const
Return reference to name.
Definition: fvMesh.H:257
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
Namespace for OpenFOAM.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
Solve the matrix using LU decomposition with pivoting.
void LUBacksubstitute(const scalarSquareMatrix &luMmatrix, const labelList &pivotIndices, List< Type > &source)
LU back-substitution with given source, returning the solution.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
#define InfoInFunction
Report an information message using Foam::Info.