39 namespace radiationModels
49 void Foam::radiationModels::viewFactor::initialise()
51 const polyBoundaryMesh& coarsePatches = coarseMesh_.boundaryMesh();
60 if ((isA<fixedValueFvPatchScalarField>(qrPatchi)))
62 selectedPatches_[
count] = qrPatchi.patch().index();
63 nLocalCoarseFaces_ += coarsePatches[
patchi].
size();
68 selectedPatches_.resize(count--);
72 Pout<<
"radiationModels::viewFactor::initialise() " 73 <<
"Selected patches:" << selectedPatches_ <<
endl;
74 Pout<<
"radiationModels::viewFactor::initialise() " 75 <<
"Number of coarse faces:" << nLocalCoarseFaces_ <<
endl;
78 totalNCoarseFaces_ = nLocalCoarseFaces_;
79 reduce(totalNCoarseFaces_, sumOp<label>());
84 <<
"Total number of clusters : " << totalNCoarseFaces_ <<
endl;
92 mesh_.facesInstance(),
105 mesh_.facesInstance(),
113 IOList<label> consMapDim
118 mesh_.facesInstance(),
141 mesh_.facesInstance(),
154 mesh_.facesInstance(),
170 globalIndex globalNumbering(nLocalCoarseFaces_);
182 <<
"Insert elements in the matrix..." <<
endl;
191 globalFaceFacesProc[proci],
198 bool smoothing =
readBool(coeffs_.lookup(
"smoothing"));
204 <<
"Smoothing the matrix..." <<
endl;
207 for (
label i=0; i<totalNCoarseFaces_; i++)
210 for (
label j=0; j<totalNCoarseFaces_; j++)
212 sumF += Fmatrix_()(i, j);
215 const scalar delta = sumF - 1.0;
216 for (
label j=0; j<totalNCoarseFaces_; j++)
218 Fmatrix_()(i, j) *= (1.0 - delta/(sumF + 0.001));
223 constEmissivity_ =
readBool(coeffs_.lookup(
"constantEmissivity"));
224 if (constEmissivity_)
231 pivotIndices_.setSize(CLU_().m());
247 mesh_.facesInstance(),
282 selectedPatches_(mesh_.
boundary().size(), -1),
283 totalNCoarseFaces_(0),
284 nLocalCoarseFaces_(0),
285 constEmissivity_(false),
318 mesh_.polyMesh::instance(),
341 totalNCoarseFaces_(0),
342 nLocalCoarseFaces_(0),
343 constEmissivity_(
false),
372 void Foam::radiationModels::viewFactor::insertMatrixElements
381 forAll(viewFactors, facei)
384 const labelList& globalFaces = globalFaceFaces[facei];
389 Fmatrix[globalI][globalFaces[i]] = vf[i];
400 scalarField compactCoarseT4(map_->constructSize(), 0.0);
401 scalarField compactCoarseE(map_->constructSize(), 0.0);
402 scalarField compactCoarseHo(map_->constructSize(), 0.0);
413 forAll(selectedPatches_, i)
415 label patchID = selectedPatches_[i];
441 const labelList& agglom = finalAgglom_[patchID];
446 forAll(coarseToFine, coarseI)
448 const label coarseFaceID = coarsePatchFace[coarseI];
449 const labelList& fineFaces = coarseToFine[coarseFaceID];
456 const scalar area =
sum(fineSf());
461 label facei = fineFaces[j];
462 T4ave[coarseI] += (
pow4(Tp[facei])*sf[facei])/area;
463 Eave[coarseI] += (eb[facei]*sf[facei])/area;
464 Hoiave[coarseI] += (Hoi[facei]*sf[facei])/area;
469 localCoarseT4ave.
append(T4ave);
470 localCoarseEave.
append(Eave);
471 localCoarseHoave.
append(Hoiave);
475 SubList<scalar>(compactCoarseT4, nLocalCoarseFaces_) = localCoarseT4ave;
477 SubList<scalar>(compactCoarseHo, nLocalCoarseFaces_) = localCoarseHoave;
480 map_->distribute(compactCoarseT4);
481 map_->distribute(compactCoarseE);
482 map_->distribute(compactCoarseHo);
485 labelList compactGlobalIds(map_->constructSize(), 0.0);
487 labelList localGlobalIds(nLocalCoarseFaces_);
489 for(
label k = 0;
k < nLocalCoarseFaces_;
k++)
500 map_->distribute(compactGlobalIds);
508 forAll(compactCoarseT4, i)
510 T4[compactGlobalIds[i]] = compactCoarseT4[i];
511 E[compactGlobalIds[i]] = compactCoarseE[i];
512 qrExt[compactGlobalIds[i]] = compactCoarseHo[i];
529 if (!constEmissivity_)
533 for (
label i=0; i<totalNCoarseFaces_; i++)
535 for (
label j=0; j<totalNCoarseFaces_; j++)
537 const scalar invEj = 1.0/E[j];
542 C(i, j) = invEj - (invEj - 1.0)*Fmatrix_()(i, j);
543 q[i] += (Fmatrix_()(i, j) - 1.0)*sigmaT4 - qrExt[j];
547 C(i, j) = (1.0 - invEj)*Fmatrix_()(i, j);
548 q[i] += Fmatrix_()(i, j)*sigmaT4;
554 Info<<
"\nSolving view factor equations..." <<
endl;
562 if (iterCounter_ == 0)
564 for (
label i=0; i<totalNCoarseFaces_; i++)
566 for (
label j=0; j<totalNCoarseFaces_; j++)
568 const scalar invEj = 1.0/E[j];
571 CLU_()(i, j) = invEj-(invEj-1.0)*Fmatrix_()(i, j);
575 CLU_()(i, j) = (1.0 - invEj)*Fmatrix_()(i, j);
583 <<
"\nDecomposing C matrix..." <<
endl;
589 for (
label i=0; i<totalNCoarseFaces_; i++)
591 for (
label j=0; j<totalNCoarseFaces_; j++)
593 const scalar sigmaT4 =
598 q[i] += (Fmatrix_()(i, j) - 1.0)*sigmaT4 - qrExt[j];
602 q[i] += Fmatrix_()(i, j)*sigmaT4;
610 <<
"\nLU Back substitute C matrix.." <<
endl;
622 label globCoarseId = 0;
623 forAll(selectedPatches_, i)
625 const label patchID = selectedPatches_[i];
631 const labelList& agglom = finalAgglom_[patchID];
639 scalar heatFlux = 0.0;
640 forAll(coarseToFine, coarseI)
644 const label coarseFaceID = coarsePatchFace[coarseI];
645 const labelList& fineFaces = coarseToFine[coarseFaceID];
648 label facei = fineFaces[
k];
650 qrp[facei] = q[globalCoarse];
651 heatFlux += qrp[facei]*sf[facei];
664 const scalar heatFlux =
gSum(qrp*magSf);
667 <<
"Total heat transfer rate at patch: "
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
This boundary condition provides a grey-diffuse condition for radiative heat flux, qr, for use with the view factor model.
Graphite solid properties.
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.
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &)
Return a temporary field constructed from name, mesh.
const labelListList & patchFaceMap() const
From patchFace on this back to original mesh or agglomeration.
const dimensionedScalar F
Faraday constant: default SI units: [C/mol].
const fileName & facesInstance() const
Return the current instance directory for faces.
IOList< labelList > labelListIOList
Label container classes.
A list of keyword definitions, which are a keyword followed by any number of values (e...
const Boundary & boundaryField() const
Return const-reference to the boundary field.
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.
To & refCast(From &r)
Reference type cast template function.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< scalar >> &)
Return a temporary field constructed from name,.
Ostream & endl(Ostream &os)
Add newline and flush stream.
static bool master(const label communicator=0)
Am I the master process.
label count(const ListType &l, typename ListType::const_reference x)
Count the number of occurrences of a value in a list.
virtual bool read()=0
Read radiationProperties dictionary.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
scalarField emissivity() const
Calculate corresponding emissivity field.
label k
Boltzmann constant.
const Time & time() const
Return the top-level database.
void calculate()
Solve system of equation(s)
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const volScalarField & T_
Reference to the temperature field.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionSet dimTime
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
A List obtained as a section of another List.
const scalarList & qro()
Return external radiative heat flux.
Type gSum(const FieldField< Field, Type > &f)
fvPatchField< scalar > fvPatchScalarField
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Top level model for radiation modelling.
virtual tmp< volScalarField::Internal > Ru() const
Source term component (constant)
virtual ~viewFactor()
Destructor.
bool read()
Read radiation properties dictionary.
const Type & value() const
Return const reference to value.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
faceListList boundary(nPatches)
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
defineTypeNameAndDebug(combustionModel, 0)
const dimensionSet dimMass
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.
label size() const
Return the number of elements in the UPtrList.
viewFactor(const volScalarField &T)
Construct from components.
dimensionedScalar pow3(const dimensionedScalar &ds)
label toGlobal(const label i) const
From local to global.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
void relax(const scalar alpha)
Relax field (for steady-state solution).
const word & name() const
Return reference to name.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
prefixOSstream Pout(cout, "Pout")
radiationModel(const volScalarField &T)
Null constructor.
dimensionedScalar pow4(const dimensionedScalar &ds)
IOList< scalarList > scalarListIOList
Scalar container classes.
const fvMesh & mesh_
Reference to the mesh database.
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
void storePrevIter() const
Store the field as the previous iteration value.
Mesh consisting of general polyhedral cells.
A class for managing temporary objects.
A patch is a list of labels that address the faces in the global face list.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
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 dimensionSet dimTemperature
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.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
#define InfoInFunction
Report an information message using Foam::Info.