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"
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::GeometricBoundaryField& 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  {
83  InfoIn("radiation::viewFactor::initialise()")
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_, totalNCoarseFaces_, 0.0)
177  );
178 
179  if (debug)
180  {
181  InfoIn("radiation::viewFactor::initialise()")
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  {
203  InfoIn("radiation::viewFactor::initialise()")
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  (
228  (
229  totalNCoarseFaces_,
230  totalNCoarseFaces_,
231  0.0
232  )
233  );
234 
235  pivotIndices_.setSize(CLU_().n());
236  }
237  }
238 }
239 
240 
241 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
242 
243 Foam::radiation::viewFactor::viewFactor(const volScalarField& T)
244 :
245  radiationModel(typeName, T),
246  finalAgglom_
247  (
248  IOobject
249  (
250  "finalAgglom",
251  mesh_.facesInstance(),
252  mesh_,
253  IOobject::MUST_READ,
254  IOobject::NO_WRITE,
255  false
256  )
257  ),
258  map_(),
259  coarseMesh_
260  (
261  IOobject
262  (
263  mesh_.name(),
264  mesh_.polyMesh::instance(),
265  mesh_.time(),
266  IOobject::NO_READ,
267  IOobject::NO_WRITE
268  ),
269  mesh_,
270  finalAgglom_
271  ),
272  Qr_
273  (
274  IOobject
275  (
276  "Qr",
277  mesh_.time().timeName(),
278  mesh_,
279  IOobject::MUST_READ,
280  IOobject::AUTO_WRITE
281  ),
282  mesh_
283  ),
284  Fmatrix_(),
285  CLU_(),
286  selectedPatches_(mesh_.boundary().size(), -1),
287  totalNCoarseFaces_(0),
288  nLocalCoarseFaces_(0),
289  constEmissivity_(false),
290  iterCounter_(0),
291  pivotIndices_(0)
292 {
293  initialise();
294 }
295 
296 
297 Foam::radiation::viewFactor::viewFactor
298 (
299  const dictionary& dict,
300  const volScalarField& T
301 )
302 :
303  radiationModel(typeName, dict, T),
304  finalAgglom_
305  (
306  IOobject
307  (
308  "finalAgglom",
310  mesh_,
313  false
314  )
315  ),
316  map_(),
317  coarseMesh_
318  (
319  IOobject
320  (
321  mesh_.name(),
322  mesh_.polyMesh::instance(),
323  mesh_.time(),
326  ),
327  mesh_,
328  finalAgglom_
329  ),
330  Qr_
331  (
332  IOobject
333  (
334  "Qr",
335  mesh_.time().timeName(),
336  mesh_,
339  ),
340  mesh_
341  ),
342  Fmatrix_(),
343  CLU_(),
344  selectedPatches_(mesh_.boundary().size(), -1),
345  totalNCoarseFaces_(0),
346  nLocalCoarseFaces_(0),
347  constEmissivity_(false),
348  iterCounter_(0),
349  pivotIndices_(0)
350 {
351  initialise();
352 }
353 
354 
355 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
356 
358 {}
359 
360 
361 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
362 
364 {
365  if (radiationModel::read())
366  {
367  return true;
368  }
369  else
370  {
371  return false;
372  }
373 }
374 
375 
376 void Foam::radiation::viewFactor::insertMatrixElements
377 (
378  const globalIndex& globalNumbering,
379  const label procI,
380  const labelListList& globalFaceFaces,
381  const scalarListList& viewFactors,
382  scalarSquareMatrix& Fmatrix
383 )
384 {
385  forAll(viewFactors, faceI)
386  {
387  const scalarList& vf = viewFactors[faceI];
388  const labelList& globalFaces = globalFaceFaces[faceI];
389 
390  label globalI = globalNumbering.toGlobal(procI, faceI);
391  forAll(globalFaces, i)
392  {
393  Fmatrix[globalI][globalFaces[i]] = vf[i];
394  }
395  }
396 }
397 
398 
400 {
401  // Store previous iteration
402  Qr_.storePrevIter();
403 
404  scalarField compactCoarseT(map_->constructSize(), 0.0);
405  scalarField compactCoarseE(map_->constructSize(), 0.0);
406  scalarField compactCoarseHo(map_->constructSize(), 0.0);
407 
408  globalIndex globalNumbering(nLocalCoarseFaces_);
409 
410  // Fill local averaged(T), emissivity(E) and external heatFlux(Ho)
411  DynamicList<scalar> localCoarseTave(nLocalCoarseFaces_);
412  DynamicList<scalar> localCoarseEave(nLocalCoarseFaces_);
413  DynamicList<scalar> localCoarseHoave(nLocalCoarseFaces_);
414 
415  forAll(selectedPatches_, i)
416  {
417  label patchID = selectedPatches_[i];
418 
419  const scalarField& Tp = T_.boundaryField()[patchID];
420  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
421 
422  fvPatchScalarField& QrPatch = Qr_.boundaryField()[patchID];
423 
425  refCast
426  <
428  >(QrPatch);
429 
430  const scalarList eb = Qrp.emissivity();
431 
432  const scalarList& Hoi = Qrp.Qro();
433 
434  const polyPatch& pp = coarseMesh_.boundaryMesh()[patchID];
435  const labelList& coarsePatchFace = coarseMesh_.patchFaceMap()[patchID];
436 
437  scalarList Tave(pp.size(), 0.0);
438  scalarList Eave(Tave.size(), 0.0);
439  scalarList Hoiave(Tave.size(), 0.0);
440 
441  if (pp.size() > 0)
442  {
443  const labelList& agglom = finalAgglom_[patchID];
444  label nAgglom = max(agglom) + 1;
445 
446  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
447 
448  forAll(coarseToFine, coarseI)
449  {
450  const label coarseFaceID = coarsePatchFace[coarseI];
451  const labelList& fineFaces = coarseToFine[coarseFaceID];
452  UIndirectList<scalar> fineSf
453  (
454  sf,
455  fineFaces
456  );
457  scalar area = sum(fineSf());
458  // Temperature, emissivity and external flux area weighting
459  forAll(fineFaces, j)
460  {
461  label faceI = fineFaces[j];
462  Tave[coarseI] += (Tp[faceI]*sf[faceI])/area;
463  Eave[coarseI] += (eb[faceI]*sf[faceI])/area;
464  Hoiave[coarseI] += (Hoi[faceI]*sf[faceI])/area;
465  }
466  }
467  }
468 
469  localCoarseTave.append(Tave);
470  localCoarseEave.append(Eave);
471  localCoarseHoave.append(Hoiave);
472  }
473 
474  // Fill the local values to distribute
475  SubList<scalar>(compactCoarseT,nLocalCoarseFaces_).assign(localCoarseTave);
476  SubList<scalar>(compactCoarseE,nLocalCoarseFaces_).assign(localCoarseEave);
478  (compactCoarseHo,nLocalCoarseFaces_).assign(localCoarseHoave);
479 
480  // Distribute data
481  map_->distribute(compactCoarseT);
482  map_->distribute(compactCoarseE);
483  map_->distribute(compactCoarseHo);
484 
485  // Distribute local global ID
486  labelList compactGlobalIds(map_->constructSize(), 0.0);
487 
488  labelList localGlobalIds(nLocalCoarseFaces_);
489 
490  for(label k = 0; k < nLocalCoarseFaces_; k++)
491  {
492  localGlobalIds[k] = globalNumbering.toGlobal(Pstream::myProcNo(), k);
493  }
494 
496  (
497  compactGlobalIds,
498  nLocalCoarseFaces_
499  ).assign(localGlobalIds);
500 
501  map_->distribute(compactGlobalIds);
502 
503  // Create global size vectors
504  scalarField T(totalNCoarseFaces_, 0.0);
505  scalarField E(totalNCoarseFaces_, 0.0);
506  scalarField QrExt(totalNCoarseFaces_, 0.0);
507 
508  // Fill lists from compact to global indexes.
509  forAll(compactCoarseT, i)
510  {
511  T[compactGlobalIds[i]] = compactCoarseT[i];
512  E[compactGlobalIds[i]] = compactCoarseE[i];
513  QrExt[compactGlobalIds[i]] = compactCoarseHo[i];
514  }
515 
519 
523 
524  // Net radiation
525  scalarField q(totalNCoarseFaces_, 0.0);
526 
527  if (Pstream::master())
528  {
529  // Variable emissivity
530  if (!constEmissivity_)
531  {
532  scalarSquareMatrix C(totalNCoarseFaces_, totalNCoarseFaces_, 0.0);
533 
534  for (label i=0; i<totalNCoarseFaces_; i++)
535  {
536  for (label j=0; j<totalNCoarseFaces_; j++)
537  {
538  scalar invEj = 1.0/E[j];
539  scalar sigmaT4 =
540  physicoChemical::sigma.value()*pow(T[j], 4.0);
541 
542  if (i==j)
543  {
544  C[i][j] = invEj - (invEj - 1.0)*Fmatrix_()[i][j];
545  q[i] += (Fmatrix_()[i][j] - 1.0)*sigmaT4 - QrExt[j];
546  }
547  else
548  {
549  C[i][j] = (1.0 - invEj)*Fmatrix_()[i][j];
550  q[i] += Fmatrix_()[i][j]*sigmaT4;
551  }
552 
553  }
554  }
555 
556  Info<< "\nSolving view factor equations..." << endl;
557 
558  // Negative coming into the fluid
559  LUsolve(C, q);
560  }
561  else //Constant emissivity
562  {
563  // Initial iter calculates CLU and chaches it
564  if (iterCounter_ == 0)
565  {
566  for (label i=0; i<totalNCoarseFaces_; i++)
567  {
568  for (label j=0; j<totalNCoarseFaces_; j++)
569  {
570  scalar invEj = 1.0/E[j];
571  if (i==j)
572  {
573  CLU_()[i][j] = invEj-(invEj-1.0)*Fmatrix_()[i][j];
574  }
575  else
576  {
577  CLU_()[i][j] = (1.0 - invEj)*Fmatrix_()[i][j];
578  }
579  }
580  }
581 
582  if (debug)
583  {
584  InfoIn("radiation::viewFactor::initialise()")
585  << "\nDecomposing C matrix..." << endl;
586  }
587 
588  LUDecompose(CLU_(), pivotIndices_);
589  }
590 
591  for (label i=0; i<totalNCoarseFaces_; i++)
592  {
593  for (label j=0; j<totalNCoarseFaces_; j++)
594  {
595  scalar sigmaT4 =
597  *pow(T[j], 4.0);
598 
599  if (i==j)
600  {
601  q[i] += (Fmatrix_()[i][j] - 1.0)*sigmaT4 - QrExt[j];
602  }
603  else
604  {
605  q[i] += Fmatrix_()[i][j]*sigmaT4;
606  }
607  }
608  }
609 
610  if (debug)
611  {
612  InfoIn("radiation::viewFactor::initialise()")
613  << "\nLU Back substitute C matrix.." << endl;
614  }
615 
616  LUBacksubstitute(CLU_(), pivotIndices_, q);
617  iterCounter_ ++;
618  }
619  }
620 
621  // Scatter q and fill Qr
624 
625 
626  label globCoarseId = 0;
627  forAll(selectedPatches_, i)
628  {
629  const label patchID = selectedPatches_[i];
630  const polyPatch& pp = mesh_.boundaryMesh()[patchID];
631  if (pp.size() > 0)
632  {
633  scalarField& Qrp = Qr_.boundaryField()[patchID];
634  const scalarField& sf = mesh_.magSf().boundaryField()[patchID];
635  const labelList& agglom = finalAgglom_[patchID];
636  label nAgglom = max(agglom)+1;
637 
638  labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
639 
640  const labelList& coarsePatchFace =
641  coarseMesh_.patchFaceMap()[patchID];
642 
643  scalar heatFlux = 0.0;
644  forAll(coarseToFine, coarseI)
645  {
646  label globalCoarse =
647  globalNumbering.toGlobal(Pstream::myProcNo(), globCoarseId);
648  const label coarseFaceID = coarsePatchFace[coarseI];
649  const labelList& fineFaces = coarseToFine[coarseFaceID];
650  forAll(fineFaces, k)
651  {
652  label faceI = fineFaces[k];
653 
654  Qrp[faceI] = q[globalCoarse];
655  heatFlux += Qrp[faceI]*sf[faceI];
656  }
657  globCoarseId ++;
658  }
659  }
660  }
661 
662  if (debug)
663  {
664  forAll(Qr_.boundaryField(), patchID)
665  {
666  const scalarField& Qrp = Qr_.boundaryField()[patchID];
667  const scalarField& magSf = mesh_.magSf().boundaryField()[patchID];
668  scalar heatFlux = gSum(Qrp*magSf);
669 
670  InfoIn("radiation::viewFactor::initialise()")
671  << "Total heat transfer rate at patch: "
672  << patchID << " "
673  << heatFlux << endl;
674  }
675  }
676 
677  // Relax Qr if necessary
678  Qr_.relax();
679 }
680 
681 
683 {
684  return tmp<volScalarField>
685  (
686  new volScalarField
687  (
688  IOobject
689  (
690  "Rp",
691  mesh_.time().timeName(),
692  mesh_,
695  false
696  ),
697  mesh_,
699  (
700  "zero",
702  0.0
703  )
704  )
705  );
706 }
707 
708 
711 {
713  (
715  (
716  IOobject
717  (
718  "Ru",
719  mesh_.time().timeName(),
720  mesh_,
723  false
724  ),
725  mesh_,
727  )
728  );
729 }
730 
731 // ************************************************************************* //
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
dimensionedScalar pow3(const dimensionedScalar &ds)
bool readBool(Istream &)
Definition: boolIO.C:60
faceListList boundary(nPatches)
void LUBacksubstitute(const scalarSquareMatrix &luMmatrix, const labelList &pivotIndices, List< Type > &source)
LU back-substitution with given source, returning the solution.
Foam::surfaceFields.
const word & name() const
Return reference to name.
Definition: fvMesh.H:257
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
IOList< labelList > labelListIOList
Label container classes.
Top level model for radiation modelling.
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Collection of constants.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:67
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
This boundary condition provides a grey-diffuse condition for radiative heat flux, Qr, for use with the view factor model.
scalarField emissivity() const
Calculate corresponding emissivity field.
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
Definition: globalIndex.H:63
const volScalarField & T_
Reference to the temperature field.
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
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:310
messageStream Info
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
volScalarField sf(fieldObject, mesh)
virtual bool read()=0
Read radiationProperties dictionary.
Namespace for OpenFOAM.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Type gSum(const FieldField< Field, Type > &f)
label n
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:386
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
Definition: viewFactor.C:682
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
A List with indirect addressing.
Definition: fvMatrix.H:106
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
bool read()
Read radiation properties dictionary.
Definition: viewFactor.C:363
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
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)
Definition: UList.H:421
label k
Boltzmann constant.
#define addToRadiationRunTimeSelectionTables(model)
scalar delta
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Macros for easy insertion into run-time selection tables.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
#define InfoIn(functionName)
Report a information message using Foam::Info.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
fvPatchField< scalar > fvPatchScalarField
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:404
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
autoPtr< radiation::radiationModel > radiation(radiation::radiationModel::New(T))
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:824
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
label toGlobal(const label i) const
From local to global.
Definition: globalIndexI.H:82
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
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:131
SquareMatrix< scalar > scalarSquareMatrix
virtual tmp< DimensionedField< scalar, volMesh > > Ru() const
Source term component (constant)
Definition: viewFactor.C:710
IOList< scalarList > scalarListIOList
Scalar container classes.
dimensionedScalar pow4(const dimensionedScalar &ds)
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
virtual ~viewFactor()
Destructor.
Definition: viewFactor.C:357
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:398
A List obtained as a section of another List.
Definition: SubList.H:53
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
void calculate()
Solve system of equation(s)
Definition: viewFactor.C:399
To & refCast(From &r)
Reference type cast template function.
Definition: typeInfo.H:106
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
Solve the matrix using LU decomposition with pivoting.
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
void storePrevIter() const
Store the field as the previous iteration value.
const fvMesh & mesh_
Reference to the mesh database.
volScalarField & C
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:118
const labelListList & patchFaceMap() const
From patchFace on this back to original mesh or agglomeration.
void relax(const scalar alpha)
Relax field (for steady-state solution).
const Type & value() const
Return const reference to value.
defineTypeNameAndDebug(combustionModel, 0)
word timeName
Definition: getTimeIndex.H:3
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:552