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),
305 mesh_.facesInstance(),
340 selectedPatches_(mesh_.
boundary().size(), -1),
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];
417 const scalarField& Tp = T_.boundaryField()[patchID];
418 const scalarField&
sf = mesh_.magSf().boundaryField()[patchID];
433 const labelList& coarsePatchFace = coarseMesh_.patchFaceMap()[patchID];
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<<
nl <<
"Solving 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];
630 const labelList& agglom = finalAgglom_[patchID];
636 coarseMesh_.patchFaceMap()[patchID];
638 forAll(coarseToFine, coarseI)
642 const label coarseFaceID = coarsePatchFace[coarseI];
643 const labelList& fineFaces = coarseToFine[coarseFaceID];
646 label facei = fineFaces[
k];
648 qrp[facei] = q[globalCoarse];
660 const scalarField& magSf = mesh_.magSf().boundaryField()[patchID];
661 const scalar heatFlux =
gSum(qrp*magSf);
664 <<
"Total heat transfer rate at patch: "
static const Foam::dimensionedScalar C("C", Foam::dimTemperature, 234.5)
#define forAll(list, i)
Loop across all elements in list.
Macros for easy insertion into run-time selection tables.
Graphite solid properties.
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &, const Field< Type > &)
Return a temporary field constructed from name, mesh,.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Generic GeometricBoundaryField class.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
A List obtained as a section of another List.
label size() const
Return the number of elements in the UList.
static bool master(const label communicator=0)
Am I the master process.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
label size() const
Return the number of elements in the UPtrList.
A list of keyword definitions, which are a keyword followed by any number of values (e....
const Type & value() const
Return const reference to value.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Calculates a unique integer (label so might not have enough room - 2G max) for processor + local inde...
label toGlobal(const label i) const
From local to global.
This boundary condition provides a grey-diffuse condition for radiative heat flux,...
const scalarList & qro()
Return external radiative heat flux.
Mesh consisting of general polyhedral cells.
A patch is a list of labels that address the faces in the global face list.
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
tmp< scalarField > emissivity() const
Calculate corresponding emissivity field.
Top level model for radiation modelling.
virtual bool read()=0
Read radiationProperties dictionary.
View factor radiation model. The system solved is: C q = b where: Cij = deltaij/Ej - (1/Ej - 1)Fij q ...
virtual ~viewFactor()
Destructor.
virtual tmp< volScalarField::Internal > Ru() const
Source term component (constant)
virtual tmp< volScalarField > Rp() const
Source term component (for power of T^4)
viewFactor(const volScalarField &T)
Construct from components.
bool read()
Read radiation properties dictionary.
void calculate()
Solve system of equation(s)
A class for managing temporary objects.
volScalarField sf(fieldObject, mesh)
#define InfoInFunction
Report an information message using Foam::Info.
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
const dimensionedScalar F
Faraday constant: default SI units: [C/kmol].
Type gSum(const FieldField< Field, Type > &f)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
To & refCast(From &r)
Reference type cast template function.
Ostream & endl(Ostream &os)
Add newline and flush stream.
IOList< labelList > labelListIOList
Label container classes.
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
void LUBacksubstitute(const scalarSquareMatrix &luMmatrix, const labelList &pivotIndices, List< Type > &source)
LU back-substitution with given source, returning the solution.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
void LUDecompose(scalarSquareMatrix &matrix, labelList &pivotIndices)
LU decompose the matrix with pivoting.
const dimensionSet dimLength
const dimensionSet dimTemperature
labelListList invertOneToMany(const label len, const labelUList &)
Invert one-to-many map. Unmapped elements will be size 0.
const dimensionSet dimTime
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow4(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
SquareMatrix< scalar > scalarSquareMatrix
const dimensionSet dimMass
prefixOSstream Pout(cout, "Pout")
IOList< scalarList > scalarListIOList
Scalar container classes.
label count(const ListType &l, typename ListType::const_reference x)
Count the number of occurrences of a value in a list.
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
Solve the matrix using LU decomposition with pivoting.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fvPatchField< scalar > fvPatchScalarField
#define addToRadiationRunTimeSelectionTables(model)
faceListList boundary(nPatches)
Basic run-time type information using word as the type's name. Used to enhance the standard RTTI to c...