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 -------------------------------------------------------------------------------
10 License
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"
30 #include "mathematicalConstants.H"
32 
33 using namespace Foam::constant;
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace functionObjects
40 {
41  defineTypeNameAndDebug(Curle, 0);
42  addToRunTimeSelectionTable(functionObject, Curle, dictionary);
43 }
44 }
45 
47 Foam::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] << "'"
161  << abort(FatalIOError);
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 
228 
229  if (surfaceWriterPtr_)
230  {
231  if (Pstream::master())
232  {
233  // Time-aware, with time spliced into the output path
234  surfaceWriterPtr_->beginTime(time_);
235 
236  surfaceWriterPtr_->open
237  (
238  inputSurface_.points(),
239  inputSurface_.surfFaces(),
240  (baseFileDir()/name()/"surface"),
241  false // serial - already merged
242  );
243 
244  surfaceWriterPtr_->write("p(Curle)", pDash);
245 
246  surfaceWriterPtr_->endTime();
247  surfaceWriterPtr_->clear();
248  }
249  }
250  else
251  {
252  forAll(observerPositions_, pointi)
253  {
254  if (rawFilePtrs_.set(pointi))
255  {
256  OFstream& os = rawFilePtrs_[pointi];
257  writeCurrentTime(os);
258  os << pDash[pointi] << endl;
259  }
260  }
261  }
262 
263  return true;
264 }
265 
266 
268 {
269  return true;
270 }
271 
272 
273 // ************************************************************************* //
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(decompositionMethod, kahipDecomp, dictionary)
Foam::functionObjects::Curle::write
virtual bool write()
Called at each ++ or += of the time-loop.
Definition: Curle.C:267
runTime
engineTime & runTime
Definition: createEngineTime.H:13
mathematicalConstants.H
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
Foam::constant
Different types of constants.
Definition: atomicConstants.C:38
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::blockMeshTools::read
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:57
Foam::Field< scalar >
Foam::functionObjects::Curle::execute
virtual bool execute()
Called at each ++ or += of the time-loop.
Definition: Curle.C:191
Foam::functionObjects::Curle::Curle
Curle(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: Curle.C:57
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::functionObjects::Curle::read
virtual bool read(const dictionary &)
Read the Curle data.
Definition: Curle.C:79
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::surfaceWriter::New
static autoPtr< surfaceWriter > New(const word &writeType)
Return a reference to the selected surfaceWriter.
Definition: surfaceWriter.C:64
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:63
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::string::expand
string & expand(const bool allowEmpty=false)
Definition: string.C:173
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::functionObjects::writeFile
Base class for writing single files from the function objects.
Definition: writeFile.H:119
fvcDdt.H
Calculate the first temporal derivative.
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::fvc::ddt
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:47
dpdt
volScalarField & dpdt
Definition: setRegionFluidFields.H:32
Foam::plusEqOp
Definition: ops.H:72
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:432
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::MeshedSurface< face >
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Curle.H