47Foam::functionObjects::Curle::modeTypeNames_
49 { modeType::point,
"point" },
50 { modeType::surface,
"surface" },
64 writeFile(mesh_,
name),
70 surfaceWriterPtr_(nullptr)
87 dict.readIfPresent(
"p", pName_);
89 patchIDs_ =
pbm.patchSet(
dict.get<
wordRes>(
"patches")).sortedToc();
91 if (patchIDs_.
empty())
94 <<
"No patches defined" <<
endl;
98 const modeType inputMode = modeTypeNames_.
get(
"input",
dict);
102 case modeType::point:
104 observerPositions_ =
dict.get<List<point>>(
"observerPositions");
107 case modeType::surface:
109 const fileName fName = (
dict.get<fileName>(
"surface")).
expand();
110 inputSurface_ = MeshedSurface<face>(fName);
112 observerPositions_ = inputSurface_.Cf();
116 if (observerPositions_.empty())
119 <<
"No observer positions defined" <<
endl;
124 const auto outputMode = modeTypeNames_.get(
"output",
dict);
128 case modeType::point:
130 rawFilePtrs_.setSize(observerPositions_.size());
131 forAll(rawFilePtrs_, filei)
139 if (rawFilePtrs_.set(filei))
141 OFstream&
os = rawFilePtrs_[filei];
152 case modeType::surface:
154 if (inputMode != modeType::surface)
157 <<
"Surface output is only available when input mode is "
158 <<
"type '" << modeTypeNames_[modeType::surface] <<
"'"
162 const word writerType =
dict.get<word>(
"surfaceType");
171 surfaceWriterPtr_->useTimeDir(
true);
176 dict.readEntry(
"c0", c0_);
181 <<
"Reference speed of sound = " << c0_
182 <<
" cannot be negative or zero."
192 if (!foundObject<volScalarField>(pName_))
206 for (
const label patchi : patchIDs_)
213 forAll(observerPositions_, pointi)
215 const vectorField r(observerPositions_[pointi] - Cfp);
219 sum((
pp*
sqr(invMagR) + invMagR/c0_*dpdtp)*(Sfp & (r*invMagR)));
227 if (surfaceWriterPtr_)
232 surfaceWriterPtr_->beginTime(time_);
234 surfaceWriterPtr_->open
236 inputSurface_.points(),
237 inputSurface_.surfFaces(),
238 (baseFileDir()/
name()/
"surface"),
242 surfaceWriterPtr_->write(
"p(Curle)", pDash);
244 surfaceWriterPtr_->endTime();
245 surfaceWriterPtr_->clear();
250 forAll(observerPositions_, pointi)
252 if (rawFilePtrs_.set(pointi))
254 OFstream&
os = rawFilePtrs_[pointi];
255 writeCurrentTime(
os);
256 os << pDash[pointi] <<
endl;
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
EnumType get(const word &enumName) const
The enumeration corresponding to the given name.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
A surface geometry mesh with zone information, not to be confused with the similarly named surfaceMes...
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
bool empty() const noexcept
True if List is empty (ie, size() is zero).
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
A class for handling file names.
Abstract base-class for Time/database function objects.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
word scopedName(const word &name) const
Return a scoped (prefixed) name.
Computes the acoustic pressure based on Curle's analogy.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
Curle(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
virtual bool execute()
Execute the function-object operations.
virtual bool write()
Write the function-object results.
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
bool foundObject(const word &fieldName) const
Find object (eg, a field) in the (sub) objectRegistry.
const Time & time_
Reference to the time database.
fileName baseFileDir() const
Return the base directory for output.
virtual autoPtr< OFstream > newFileAtStartTime(const word &name) const
Return autoPtr to a new file using the simulation start time.
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
void writeHeaderValue(Ostream &os, const string &property, const Type &value) const
Write a (commented) header property and value pair.
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.
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.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
static autoPtr< surfaceWriter > New(const word &writeType)
Select construct a surfaceWriter.
static dictionary formatOptions(const dictionary &dict, const word &formatName, const word &entryName="formatOptions")
Same as fileFormats::getFormatOptions.
A List of wordRe with additional matching capabilities.
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.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
OBJstream os(runTime.globalPath()/outputName)
Calculate the first temporal derivative.
#define WarningInFunction
Report a warning using Foam::Warning.
constexpr scalar pi(M_PI)
Different types of constants.
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
string expand(const std::string &s, const HashTable< string > &mapping, const char sigil='$')
Expand occurrences of variables according to the mapping and return the expanded string.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Field< vector > vectorField
Specialisation of Field<T> for vector.
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
#define forAll(list, i)
Loop across all elements in list.