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-2021 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;
34
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:
64 writeFile(mesh_, name),
65 pName_("p"),
66 patchSet_(),
67 observerPositions_(),
68 c0_(0),
69 rawFilePtrs_(),
70 inputSurface_(),
71 surfaceWriterPtr_(nullptr)
72{
73 read(dict);
74}
75
76
77// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
78
80{
82 {
83 return false;
84 }
85
86 dict.readIfPresent("p", pName_);
87
88 patchSet_ = mesh_.boundaryMesh().patchSet(dict.get<wordRes>("patches"));
89
90 if (patchSet_.empty())
91 {
93 << "No patches defined"
94 << endl;
95
96 return false;
97 }
98
99 const modeType inputMode = modeTypeNames_.get("input", dict);
100
101 switch (inputMode)
102 {
103 case modeType::point:
104 {
105 observerPositions_ = dict.get<List<point>>("observerPositions");
106 break;
107 }
108 case modeType::surface:
109 {
110 const fileName fName = (dict.get<fileName>("surface")).expand();
111 inputSurface_ = MeshedSurface<face>(fName);
112
113 observerPositions_ = inputSurface_.Cf();
114 }
115 }
116
117 if (observerPositions_.empty())
118 {
120 << "No observer positions defined"
121 << endl;
122
123 return false;
124 }
125
126 const auto outputMode = modeTypeNames_.get("output", dict);
127
128 switch (outputMode)
129 {
130 case modeType::point:
131 {
132 rawFilePtrs_.setSize(observerPositions_.size());
133 forAll(rawFilePtrs_, filei)
134 {
135 rawFilePtrs_.set
136 (
137 filei,
138 createFile("observer" + Foam::name(filei))
139 );
140
141 if (rawFilePtrs_.set(filei))
142 {
143 OFstream& os = rawFilePtrs_[filei];
144
145 writeHeaderValue(os, "Observer", filei);
146 writeHeaderValue(os, "Position", observerPositions_[filei]);
147 writeCommented(os, "Time");
148 writeTabbed(os, "p(Curle)");
149 os << endl;
150 }
151 }
152 break;
153 }
154 case modeType::surface:
155 {
156 if (inputMode != modeType::surface)
157 {
159 << "Surface output is only available when input mode is "
160 << "type '" << modeTypeNames_[modeType::surface] << "'"
162 }
163
164 const word surfaceType(dict.get<word>("surfaceType"));
165 surfaceWriterPtr_ = surfaceWriter::New
166 (
167 surfaceType,
168 dict.subOrEmptyDict("formatOptions").subOrEmptyDict(surfaceType)
169 );
170
171 // Use outputDir/TIME/surface-name
172 surfaceWriterPtr_->useTimeDir(true);
173 }
174 }
175
176 // Read the reference speed of sound
177 dict.readEntry("c0", c0_);
178
179 if (c0_ < VSMALL)
180 {
182 << "Reference speed of sound = " << c0_
183 << " cannot be negative or zero."
184 << abort(FatalError);
185 }
186
187 return true;
188}
189
190
192{
193 if (!foundObject<volScalarField>(pName_))
194 {
195 return false;
196 }
197
198 const volScalarField& p = lookupObject<volScalarField>(pName_);
199 const volScalarField::Boundary& pBf = p.boundaryField();
200 const volScalarField dpdt(scopedName("dpdt"), fvc::ddt(p));
201 const volScalarField::Boundary& dpdtBf = dpdt.boundaryField();
202 const surfaceVectorField::Boundary& SfBf = mesh_.Sf().boundaryField();
203 const surfaceVectorField::Boundary& CfBf = mesh_.Cf().boundaryField();
204
205 scalarField pDash(observerPositions_.size(), 0);
206
207 for (auto patchi : patchSet_)
208 {
209 const scalarField& pp = pBf[patchi];
210 const scalarField& dpdtp = dpdtBf[patchi];
211 const vectorField& Sfp = SfBf[patchi];
212 const vectorField& Cfp = CfBf[patchi];
213
214 forAll(observerPositions_, pointi)
215 {
216 const vectorField r(observerPositions_[pointi] - Cfp);
217 const scalarField invMagR(1/(mag(r) + ROOTVSMALL));
218
219 pDash[pointi] +=
220 sum((pp*sqr(invMagR) + invMagR/c0_*dpdtp)*(Sfp & (r*invMagR)));
221 }
222 }
223
224 pDash /= 4*mathematical::pi;
225
227
228 if (surfaceWriterPtr_)
229 {
230 if (Pstream::master())
231 {
232 // Time-aware, with time spliced into the output path
233 surfaceWriterPtr_->beginTime(time_);
234
235 surfaceWriterPtr_->open
236 (
237 inputSurface_.points(),
238 inputSurface_.surfFaces(),
239 (baseFileDir()/name()/"surface"),
240 false // serial - already merged
241 );
242
243 surfaceWriterPtr_->write("p(Curle)", pDash);
244
245 surfaceWriterPtr_->endTime();
246 surfaceWriterPtr_->clear();
247 }
248 }
249 else
250 {
251 forAll(observerPositions_, pointi)
252 {
253 if (rawFilePtrs_.set(pointi))
254 {
255 OFstream& os = rawFilePtrs_[pointi];
256 writeCurrentTime(os);
257 os << pDash[pointi] << endl;
258 }
259 }
260 }
261
262 return true;
263}
264
265
267{
268 return true;
269}
270
271
272// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
A surface geometry mesh with zone information, not to be confused with the similarly named surfaceMes...
Definition: MeshedSurface.H:99
Output to file stream, using an OSstream.
Definition: OFstream.H:57
static void listCombineAllGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
virtual bool read()
Re-read model coefficients if they have changed.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
A class for handling file names.
Definition: fileName.H:76
Abstract base-class for Time/database function objects.
Computes the acoustic pressure based on Curle's analogy.
Definition: Curle.H:252
virtual bool execute()
Called at each ++ or += of the time-loop.
Definition: Curle.C:191
virtual bool write()
Called at each ++ or += of the time-loop.
Definition: Curle.C:266
virtual bool read(const dictionary &)
Read the Curle data.
Definition: Curle.C:79
Watches for presence of the named trigger file in the case directory and signals a simulation stop (o...
Definition: abort.H:128
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Computes the magnitude of an input field.
Definition: mag.H:153
Base class for writing single files from the function objects.
Definition: writeFile.H:120
splitCell * master() const
Definition: splitCell.H:113
string & expand(const bool allowEmpty=false)
Definition: string.C:173
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
volScalarField & p
engineTime & runTime
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
volScalarField & dpdt
Calculate the first temporal derivative.
#define WarningInFunction
Report a warning using Foam::Warning.
constexpr scalar pi(M_PI)
Different types of constants.
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:47
Namespace for OpenFOAM.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
IOerror FatalIOError
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333