44void Foam::functionObjects::radiometerProbes::writeFileHeader(
Ostream&
os)
51 for (label i = 0; i < szProbes_; ++i)
53 const vector& loc = locs[i];
57 <<
',' << loc.
x() <<
',' << loc.
y() <<
',' << loc.
z()
58 <<
',' <<
n.x() <<
',' <<
n.y() <<
',' <<
n.z()
63 for (
int i = 0; i < szProbes_; ++i)
110 if (!dom_.radiation())
113 <<
"The radiation model is inactive."
114 <<
"Skipping the function object " <<
type() <<
' ' <<
name()
123 szProbes_ = this->size();
129 if (!szProbes_ || !nRay_)
132 <<
"size(probe locations): " << szProbes_ <<
nl
133 <<
"size(rays): " << nRay_ <<
nl
134 <<
"The input size of probe locations and rays cannot be zero."
139 dict.readEntry(
"probeNormals", n_);
140 if (n_.size() != szProbes_)
143 <<
"size(probe locations): " << szProbes_ <<
nl
144 <<
"size(probe normals): " << n_.size() <<
nl
145 <<
"The input size of probe locations and normals must match."
153 n_dAve_.resize(nRay_);
156 for (label rayi = 0; rayi < nRay_; ++rayi)
158 const vector& dAvei = dom_.IRay(rayi).dAve();
164 Cray.resize(szProbes_,
false);
166 for (label
pi = 0;
pi < szProbes_; ++
pi)
168 n_dAveRay[
pi] = n_[
pi] & dAvei;
170 if (n_dAveRay[
pi] < 0)
177 qin_.resize(szProbes_);
197 if (!szProbes_ || !dom_.radiation() || !shouldCalcThisStep())
206 for (label rayi = 0; rayi < nRay_; ++rayi)
219 for (label
pi = 0;
pi < szProbes_; ++
pi)
223 qin_[
pi] += Ip[
pi]*n_dAveRay[
pi];
235 if (!szProbes_ || !dom_.radiation() || !shouldCalcThisStep())
244 Ostream&
os = file();
246 os << mesh_.time().timeOutputValue();
247 for (label
pi = 0;
pi < szProbes_; ++
pi)
249 os <<
',' << qin_[
pi];
constexpr scalar pi(M_PI)
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
void resize(const label len)
Adjust allocated size of list.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
const Cmpt & x() const noexcept
Access to the vector x component.
const Cmpt & z() const noexcept
Access to the vector z component.
const Cmpt & y() const noexcept
Access to the vector y component.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Abstract base-class for Time/database function objects.
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.
Probes the incident radiative heat flux, qin, at arbitrary points within a domain.
radiometerProbes(const radiometerProbes &)=delete
No copy construct.
virtual bool execute()
Execute the function object.
virtual bool write()
Write to data files/fields and to streams.
virtual bool read(const dictionary &)
Read the function object settings.
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
virtual bool read(const dictionary &dict)
Read optional controls.
Base class for writing single files from the function objects.
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 bool read(const dictionary &dict)
Read.
bool writtenHeader_
Flag to identify whether the header has been written.
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 resetFile(const word &name)
Reset internal file pointer to new file with new name.
virtual bool canWriteHeader() const
Flag to allow writing the header.
virtual bool canResetFile() const
Flag to allow resetting the file.
A utility class for probing field values at specified point locations within an fvMesh.
internalFieldProbe(const fvMesh &mesh, const dictionary &dict)
Construct from Time and dictionary.
tmp< Field< Type > > sample(const VolumeField< Type > &) const
Sample a volume field at all locations.
virtual bool read(const dictionary &)
Read the settings dictionary.
label size() const
Return number of probe locations.
const pointField & probeLocations() const
Return const reference to the probe locations.
Finite Volume Discrete Ordinates Method. Solves the RTE equation for n directions in a participating ...
A class for managing temporary objects.
const T & cref() const
Return const reference to the object or to the contents of a (non-null) managed pointer.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
A class for handling words, derived from Foam::string.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
OBJstream os(runTime.globalPath()/outputName)
#define WarningInFunction
Report a warning using Foam::Warning.
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for radiation modelling.
bool read(const char *buf, int32_t &val)
Same as readInt32.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
static const Identity< scalar > I
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
List< bool > boolList
A List of bools.
static constexpr const zero Zero
Global zero (0).
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
vectorField pointField
pointField is a vectorField.
errorManipArg< error, int > exit(error &err, const int errNo=1)
List< scalar > scalarList
List of scalar.
constexpr char nl
The newline '\n' character (0x0a).