46void Foam::functionObjects::viewFactorHeatFlux::initialise()
49 const labelList selectedPatches =
pbm.indices(
"viewFactorWall");
59 for (
const label patchi : selectedPatches)
72 const label patchi = selectedPatches[i];
73 auto slice = compactPatchID_.
slice(nFace,
pbm[patchi].size());
75 nFace += slice.
size();
88 mesh_.facesInstance(),
97 auto& mapDist = mapDistPtr_();
99 mapDist.distribute(compactPatchID_);
106 globalIndex gi(nFace);
108 SubList<label>(compactGlobalIds, nFace) =
identity(gi.range());
110 mapDist.distribute(compactGlobalIds);
112 const label nTotalFaces =
returnReduce(nFace, sumOp<label>());
114 const labelList globalToCompact(
invert(nTotalFaces, compactGlobalIds));
116 for (
labelList& visibleFaces : faceFaces_)
118 for (label& globalFacei : visibleFaces)
120 globalFacei = globalToCompact[globalFacei];
139 qrName_(
dict.getOrDefault<word>(
"qr",
"qr")),
140 mapDistPtr_(nullptr),
146 mesh_.facesInstance(),
150 IOobject::NO_REGISTER
158 mesh_.facesInstance(),
162 IOobject::NO_REGISTER
165 compactPatchID_(faceFaces_.size(),
Zero)
181 qrName_(
dict.getOrDefault<
word>(
"qr",
"qr")),
182 mapDistPtr_(nullptr),
188 mesh_.facesInstance(),
200 mesh_.facesInstance(),
207 compactPatchID_(faceFaces_.size(),
Zero)
219 dict.readIfPresent(
"qr", qrName_);
244 const auto& qr = *qrPtr;
252 mesh_.boundaryMesh().indices(
"viewFactorWall");
254 const auto&
pbm = mesh_.boundaryMesh();
255 const label nPatch = selectedPatches.
size();
263 forAll(selectedPatches, i)
265 const label patchi = selectedPatches[i];
266 const auto& qrp = qr.boundaryField()[patchi];
269 compactQr[compacti] = qrp[facei];
270 qrPatch[i] += qrp[facei];
275 reduce(qrPatch, sumOp<scalarList>());
280 mapDistPtr_->distribute(compactQr);
288 const scalar qr = compactQr[startFacei];
289 const labelList& visibleSlots = faceFaces_[startFacei];
290 const label i = compactPatchID_[startFacei];
292 forAll(visibleSlots, visSloti)
294 const label sloti = visibleSlots[visSloti];
295 const label j = compactPatchID_[sloti];
297 qrVf[i][j] += Fij_[startFacei][visSloti]*qr;
301 reduce(qrVf, sumOp<scalarSquareMatrix>());
305 Log <<
" Writing patch totals to " << file().name()
309 writeCurrentTime(file());
311 for (
const auto& qrp : qrPatch)
313 file() <<
tab << qrp;
320 auto osPtr = newFileAtTime(
"viewFactorQr", mesh_.time().value());
323 Log <<
" Writing view factor breakdown to " <<
os.name()
327 forAll(selectedPatches, i)
331 for (
const label patchj : selectedPatches)
338 const label patchi = selectedPatches[i];
340 os <<
pbm[patchi].name();;
342 forAll(selectedPatches, j)
344 os <<
tab << qrVf[i][j];
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const polyBoundaryMesh & pbm
@ NO_REGISTER
Do not request registration (bool: false).
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
SubList< T > slice(const label pos, label len=-1)
Return SubList slice (non-const access) - no range checking.
void size(const label n)
Older name for setAddressableSize.
static bool parRun(const bool on) noexcept
Set as parallel run on/off.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
static bool & parRun() noexcept
Test if this a parallel run.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Abstract base-class for Time/database function objects.
const word & name() const noexcept
Return the name of this functionObject.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
Reads fields from the time directories and adds them to the mesh database for further post-processing...
virtual const objectRegistry & obr() const
The region or sub-region registry being used.
Determines radiation heat flux between patches when using the viewFactor radiation model.
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
viewFactorHeatFlux(const word &name, const Time &runTime, const dictionary &dict, const bool readFields=true)
Construct from name, Time and dictionary.
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
virtual bool execute()
Execute the function-object operations (no-op).
virtual bool write()
Write the function-object results (per patch).
virtual bool read(const dictionary &)
Read the function-object dictionary.
Base class for writing single files from the function objects.
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
writeFile(const objectRegistry &obr, const fileName &prefix, const word &name="undefined", const bool writeToFile=true, const string &ext=".dat")
Construct from objectRegistry, prefix, fileName.
virtual void writeHeader(Ostream &os, const string &str) const
Write a commented header to stream.
virtual OFstream & file()
Return access to the file (if only 1).
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
virtual void writeCurrentTime(Ostream &os) const
Write the current time to stream.
virtual autoPtr< OFstream > newFileAtTime(const word &name, scalar timeValue) const
Return autoPtr to a new file for a given time.
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Registry of regIOobjects.
Mesh consisting of general polyhedral cells.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
OBJstream os(runTime.globalPath()/outputName)
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
List< label > labelList
A List of labels.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
labelList identity(const label len, label start=0)
Return an identity map of the given length with (map[i] == i), works like std::iota() but returning a...
SquareMatrix< scalar > scalarSquareMatrix
static constexpr const zero Zero
Global zero (0).
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
List< scalar > scalarList
List of scalar.
constexpr char nl
The newline '\n' character (0x0a).
constexpr char tab
The tab '\t' character(0x09).
#define forAll(list, i)
Loop across all elements in list.