Loading...
Searching...
No Matches
Curle.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2017-2022 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "Curle.H"
29#include "fvcDdt.H"
32
33using namespace Foam::constant;
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace functionObjects
40{
43}
44}
45
47Foam::functionObjects::Curle::modeTypeNames_
48({
49 { modeType::point, "point" },
50 { modeType::surface, "surface" },
51});
52
53
54// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55
57(
58 const word& name,
59 const Time& runTime,
60 const dictionary& dict
61)
62:
63 fvMeshFunctionObject(name, runTime, dict),
64 writeFile(mesh_, name),
65 pName_("p"),
66 observerPositions_(),
67 c0_(0),
68 rawFilePtrs_(),
69 inputSurface_(),
70 surfaceWriterPtr_(nullptr)
72 read(dict);
73}
74
75
76// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77
79{
81
83 {
84 return false;
85 }
86
87 dict.readIfPresent("p", pName_);
88
89 patchIDs_ = pbm.patchSet(dict.get<wordRes>("patches")).sortedToc();
90
91 if (patchIDs_.empty())
92 {
94 << "No patches defined" << endl;
95 return false;
96 }
97
98 const modeType inputMode = modeTypeNames_.get("input", dict);
99
100 switch (inputMode)
101 {
102 case modeType::point:
103 {
104 observerPositions_ = dict.get<List<point>>("observerPositions");
105 break;
106 }
107 case modeType::surface:
108 {
109 const fileName fName = (dict.get<fileName>("surface")).expand();
110 inputSurface_ = MeshedSurface<face>(fName);
111
112 observerPositions_ = inputSurface_.Cf();
113 }
114 }
115
116 if (observerPositions_.empty())
117 {
119 << "No observer positions defined" << endl;
120
121 return false;
122 }
123
124 const auto outputMode = modeTypeNames_.get("output", dict);
125
126 switch (outputMode)
127 {
128 case modeType::point:
129 {
130 rawFilePtrs_.setSize(observerPositions_.size());
131 forAll(rawFilePtrs_, filei)
132 {
133 rawFilePtrs_.set
134 (
135 filei,
136 newFileAtStartTime("observer" + Foam::name(filei))
137 );
138
139 if (rawFilePtrs_.set(filei))
140 {
141 OFstream& os = rawFilePtrs_[filei];
142
143 writeHeaderValue(os, "Observer", filei);
144 writeHeaderValue(os, "Position", observerPositions_[filei]);
145 writeCommented(os, "Time");
146 writeTabbed(os, "p(Curle)");
147 os << endl;
148 }
149 }
150 break;
151 }
152 case modeType::surface:
153 {
154 if (inputMode != modeType::surface)
155 {
157 << "Surface output is only available when input mode is "
158 << "type '" << modeTypeNames_[modeType::surface] << "'"
160 }
161
162 const word writerType = dict.get<word>("surfaceType");
163
164 surfaceWriterPtr_ = surfaceWriter::New
165 (
166 writerType,
168 );
169
170 // Use outputDir/TIME/surface-name
171 surfaceWriterPtr_->useTimeDir(true);
172 }
173 }
174
175 // Read the reference speed of sound
176 dict.readEntry("c0", c0_);
177
178 if (c0_ < VSMALL)
179 {
181 << "Reference speed of sound = " << c0_
182 << " cannot be negative or zero."
184 }
185
186 return true;
187}
188
189
191{
192 if (!foundObject<volScalarField>(pName_))
193 {
194 return false;
195 }
196
197 const volScalarField& p = lookupObject<volScalarField>(pName_);
198 const volScalarField::Boundary& pBf = p.boundaryField();
199 const volScalarField dpdt(scopedName("dpdt"), fvc::ddt(p));
200 const volScalarField::Boundary& dpdtBf = dpdt.boundaryField();
201 const surfaceVectorField::Boundary& SfBf = mesh_.Sf().boundaryField();
202 const surfaceVectorField::Boundary& CfBf = mesh_.Cf().boundaryField();
203
204 scalarField pDash(observerPositions_.size(), 0);
205
206 for (const label patchi : patchIDs_)
207 {
208 const scalarField& pp = pBf[patchi];
209 const scalarField& dpdtp = dpdtBf[patchi];
210 const vectorField& Sfp = SfBf[patchi];
211 const vectorField& Cfp = CfBf[patchi];
212
213 forAll(observerPositions_, pointi)
214 {
215 const vectorField r(observerPositions_[pointi] - Cfp);
216 const scalarField invMagR(1/(mag(r) + ROOTVSMALL));
217
218 pDash[pointi] +=
219 sum((pp*sqr(invMagR) + invMagR/c0_*dpdtp)*(Sfp & (r*invMagR)));
220 }
221 }
222
223 pDash /= 4*mathematical::pi;
224
225 Pstream::listReduce(pDash, sumOp<scalar>());
226
227 if (surfaceWriterPtr_)
228 {
229 if (UPstream::master())
230 {
231 // Time-aware, with time spliced into the output path
232 surfaceWriterPtr_->beginTime(time_);
233
234 surfaceWriterPtr_->open
235 (
236 inputSurface_.points(),
237 inputSurface_.surfFaces(),
238 (baseFileDir()/name()/"surface"),
239 false // serial - already merged
240 );
241
242 surfaceWriterPtr_->write("p(Curle)", pDash);
243
244 surfaceWriterPtr_->endTime();
245 surfaceWriterPtr_->clear();
246 }
247 }
248 else
249 {
250 forAll(observerPositions_, pointi)
251 {
252 if (rawFilePtrs_.set(pointi))
253 {
254 OFstream& os = rawFilePtrs_[pointi];
255 writeCurrentTime(os);
256 os << pDash[pointi] << endl;
257 }
259 }
260
261 return true;
262}
263
264
266{
267 return true;
268}
269
270
271// ************************************************************************* //
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...
Definition Enum.H:57
EnumType get(const word &enumName) const
The enumeration corresponding to the given name.
Definition Enum.C:68
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...
Definition List.H:72
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.
Definition OFstream.H:75
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.
Definition Time.H:75
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A class for handling file names.
Definition fileName.H:75
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.
Definition Curle.H:234
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
Definition Curle.C:71
Curle(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
Definition Curle.C:50
virtual bool execute()
Execute the function-object operations.
Definition Curle.C:183
virtual bool write()
Write the function-object results.
Definition Curle.C:258
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.
Definition writeFile.C:43
virtual autoPtr< OFstream > newFileAtStartTime(const word &name) const
Return autoPtr to a new file using the simulation start time.
Definition writeFile.C:156
virtual void writeTabbed(Ostream &os, const string &str) const
Write a tabbed string to stream.
Definition writeFile.C:334
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.
Definition writeFile.C:200
virtual bool read(const dictionary &dict)
Read.
Definition writeFile.C:240
virtual void writeCommented(Ostream &os, const string &str) const
Write a commented string to stream.
Definition writeFile.C:318
virtual void writeCurrentTime(Ostream &os) const
Write the current time to stream.
Definition writeFile.C:354
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
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.
Definition wordRes.H:56
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
volScalarField & p
engineTime & runTime
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
volScalarField & dpdt
Calculate the first temporal derivative.
auto & name
#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)
Definition fvcDdt.C:40
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.
Namespace for OpenFOAM.
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.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
errorManip< error > abort(error &err)
Definition errorManip.H:139
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.
Definition exprTraits.C:127
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299