Loading...
Searching...
No Matches
momentum.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) 2018-2024 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 "momentum.H"
29#include "fvMesh.H"
30#include "volFields.H"
31#include "cellSet.H"
32#include "cylindricalRotation.H"
33#include "emptyPolyPatch.H"
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37
38namespace Foam
39{
40namespace functionObjects
41{
44}
45}
46
47
48// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
49
50void Foam::functionObjects::momentum::purgeFields()
51{
52 obr_.checkOut(scopedName("momentum"));
53 obr_.checkOut(scopedName("angularMomentum"));
54 obr_.checkOut(scopedName("angularVelocity"));
55}
56
57
58template<class GeoField>
60Foam::functionObjects::momentum::newField
61(
62 const word& baseName,
63 const dimensionSet& dims,
64 bool registerObject
65) const
66{
67 return
69 (
71 (
72 scopedName(baseName),
73 time_.timeName(),
74 mesh_.thisDb(),
77 registerObject
78 ),
79 mesh_,
80 Foam::zero{}, // value
81 dims
82 );
83}
84
85
86void Foam::functionObjects::momentum::calc()
87{
88 initialise();
89
90 // Ensure volRegion is properly up-to-date.
91 // Purge old fields if anything has changed (eg, mesh size etc)
93 {
94 purgeFields();
95 }
96
97 // When field writing is not enabled we need our local storage
98 // for the momentum and angular velocity fields
99 autoPtr<volVectorField> tmomentum, tAngularMom, tAngularVel;
100
101
102 // The base fields required
103 const auto& U = lookupObject<volVectorField>(UName_);
104 const auto* rhoPtr = findObject<volScalarField>(rhoName_);
105
106 const dimensionedScalar rhoRef("rho", dimDensity, rhoRef_);
107
108 // For quantities such as the mass-averaged angular velocity,
109 // we would calculate the mass per-cell here.
110
111 // tmp<volScalarField::Internal> tmass =
112 // (
113 // rhoPtr
114 // ? (mesh_.V() * (*rhoPtr))
115 // : (mesh_.V() * rhoRef)
116 // );
117
118
119 // Linear momentum
120 // ~~~~~~~~~~~~~~~
121
122 auto* pmomentum = getObjectPtr<volVectorField>(scopedName("momentum"));
123
124 if (!pmomentum)
125 {
126 tmomentum = newField<volVectorField>("momentum", dimVelocity*dimMass);
127 pmomentum = tmomentum.get(); // get(), not release()
128 }
129 auto& momentum = *pmomentum;
130
131 if (rhoPtr)
132 {
133 momentum.ref() = (U * mesh_.V() * (*rhoPtr));
134 }
135 else
136 {
137 momentum.ref() = (U * mesh_.V() * rhoRef);
138 }
139 momentum.correctBoundaryConditions();
140
141
142 // Angular momentum
143 // ~~~~~~~~~~~~~~~~
144
145 auto* pAngularMom =
146 getObjectPtr<volVectorField>(scopedName("angularMomentum"));
147
148 if (hasCsys_ && !pAngularMom)
149 {
150 tAngularMom =
151 newField<volVectorField>("angularMomentum", dimVelocity*dimMass);
152 pAngularMom = tAngularMom.get(); // get(), not release()
153 }
154 else if (!pAngularMom)
155 {
156 // Angular momentum not requested, but alias to normal momentum
157 // to simplify logic when calculating the summations
158 pAngularMom = pmomentum;
159 }
160 auto& angularMom = *pAngularMom;
161
162
163 // Angular velocity
164 // ~~~~~~~~~~~~~~~~
165
166 auto* pAngularVel =
167 getObjectPtr<volVectorField>(scopedName("angularVelocity"));
168
169 if (hasCsys_)
170 {
171 if (!pAngularVel)
172 {
173 tAngularVel =
174 newField<volVectorField>("angularVelocity", dimVelocity);
175 pAngularVel = tAngularVel.get(); // get(), not release()
176 }
177 auto& angularVel = *pAngularVel;
178
179
180 // Global to local
181
182 angularVel.primitiveFieldRef() =
183 csys_.invTransform(mesh_.cellCentres(), U.internalField());
184
185 angularVel.correctBoundaryConditions();
186
187 if (rhoPtr)
188 {
189 angularMom.ref() = (angularVel * mesh_.V() * (*rhoPtr));
190 }
191 else
192 {
193 angularMom.ref() = (angularVel * mesh_.V() * rhoRef);
194 }
195
196 angularMom.correctBoundaryConditions();
197 }
198
199
200 // Integrate the selection
201
202 sumMomentum_ = Zero;
203 sumAngularMom_ = Zero;
204
206 {
207 for (label celli=0; celli < mesh_.nCells(); ++celli)
208 {
209 sumMomentum_ += momentum[celli];
210 sumAngularMom_ += angularMom[celli];
211 }
212 }
213 else
214 {
215 for (const label celli : cellIDs())
216 {
217 sumMomentum_ += momentum[celli];
218 sumAngularMom_ += angularMom[celli];
219 }
220 }
221
224}
225
226
227// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
228
230{
231 if (!writeToFile() || writtenHeader_)
232 {
233 return;
234 }
235
236 if (hasCsys_)
237 {
238 writeHeader(os, "Momentum, Angular Momentum");
239 writeHeaderValue(os, "origin", csys_.origin());
240 writeHeaderValue(os, "axis", csys_.e3());
241 }
242 else
243 {
244 writeHeader(os, "Momentum");
245 }
246
248 {
250 (
251 os,
252 "Selection " + regionTypeNames_[regionType_]
253 + " = " + regionName_
254 );
255 }
256
257 writeHeader(os, "");
258 writeCommented(os, "Time");
259 writeTabbed(os, "(momentum_x momentum_y momentum_z)");
260
261 if (hasCsys_)
262 {
263 writeTabbed(os, "(momentum_r momentum_rtheta momentum_axis)");
264 }
265
266 writeTabbed(os, "volume");
267 os << endl;
268
269 writtenHeader_ = true;
270}
271
272
274{
275 if (initialised_)
276 {
277 return;
278 }
279
280 if (!foundObject<volVectorField>(UName_))
281 {
283 << "Could not find U: " << UName_ << " in database"
284 << exit(FatalError);
285 }
286
287
288 const auto* pPtr = cfindObject<volScalarField>(pName_);
289
290 if (pPtr && pPtr->dimensions() == dimPressure)
291 {
292 // Compressible - rho is mandatory
293
294 if (!foundObject<volScalarField>(rhoName_))
295 {
297 << "Could not find rho:" << rhoName_
298 << exit(FatalError);
300 }
301
302 initialised_ = true;
303}
304
305
307{
308 if (log)
309 {
310 Info<< type() << " " << name() << " write:" << nl;
311
312 Info<< " Sum of Momentum";
313
315 {
316 Info<< ' ' << regionTypeNames_[regionType_]
317 << ' ' << regionName_;
318 }
319
320 Info<< " (volume " << volRegion::V() << ')' << nl
321 << " linear : " << sumMomentum_ << nl;
322
323 if (hasCsys_)
324 {
325 Info<< " angular : " << sumAngularMom_ << nl;
326 }
327
328 Info<< endl;
329 }
330
331 if (writeToFile())
332 {
333 writeCurrentTime(os);
334
335 os << tab << sumMomentum_;
336
337 if (hasCsys_)
338 {
339 os << tab << sumAngularMom_;
340 }
341
343 }
344}
345
346
347// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
348
350(
351 const word& name,
352 const Time& runTime,
353 const dictionary& dict,
354 bool readFields
355)
356:
357 fvMeshFunctionObject(name, runTime, dict),
358 volRegion(fvMeshFunctionObject::mesh_, dict),
359 writeFile(mesh_, name, typeName, dict),
360 sumMomentum_(Zero),
361 sumAngularMom_(Zero),
362 UName_(),
363 pName_(),
364 rhoName_(),
365 rhoRef_(1.0),
366 csys_(),
367 hasCsys_(false),
368 writeMomentum_(false),
369 writeVelocity_(false),
370 writePosition_(false),
371 initialised_(false)
372{
373 if (readFields)
375 read(dict);
376 Log << endl;
377 }
378}
379
380
382(
383 const word& name,
384 const objectRegistry& obr,
385 const dictionary& dict,
386 bool readFields
387)
388:
389 fvMeshFunctionObject(name, obr, dict),
390 volRegion(fvMeshFunctionObject::mesh_, dict),
391 writeFile(mesh_, name, typeName, dict),
392 sumMomentum_(Zero),
393 sumAngularMom_(Zero),
394 UName_(),
395 pName_(),
396 rhoName_(),
397 rhoRef_(1.0),
398 csys_(),
399 hasCsys_(false),
400 writeMomentum_(false),
401 writeVelocity_(false),
402 writePosition_(false),
403 initialised_(false)
404{
405 if (readFields)
406 {
407 read(dict);
409 }
410}
411
412
413// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
414
416{
420
421 initialised_ = false;
422
423 Info<< type() << " " << name() << ":" << nl;
424
425 // Optional entries U and p
426 UName_ = dict.getOrDefault<word>("U", "U");
427 pName_ = dict.getOrDefault<word>("p", "p");
428 rhoName_ = dict.getOrDefault<word>("rho", "rho");
429 rhoRef_ = dict.getOrDefault<scalar>("rhoRef", 1.0);
430 hasCsys_ = dict.getOrDefault("cylindrical", false);
431
432 if (hasCsys_)
433 {
435 }
436
437 writeMomentum_ = dict.getOrDefault("writeMomentum", false);
438 writeVelocity_ = dict.getOrDefault("writeVelocity", false);
439 writePosition_ = dict.getOrDefault("writePosition", false);
440
441 Info<<"Integrating for selection: "
442 << regionTypeNames_[regionType_]
443 << " (" << regionName_ << ")" << nl;
444
445 if (writeMomentum_)
446 {
447 Info<< " Momentum fields will be written" << endl;
448
450 (
451 newField<volVectorField>("momentum", dimVelocity*dimMass)
452 );
453
454 if (hasCsys_)
455 {
457 (
458 newField<volVectorField>("angularMomentum", dimVelocity*dimMass)
459 );
460 }
461 }
462
463 if (hasCsys_)
464 {
465 if (writeVelocity_)
466 {
467 Info<< " Angular velocity will be written" << endl;
468
470 (
471 newField<volVectorField>("angularVelocity", dimVelocity)
472 );
473 }
474
475 if (writePosition_)
476 {
477 Info<< " Angular position will be written" << endl;
479 }
480
481 return true;
482}
483
484
486{
487 calc();
488
489 if (Pstream::master())
490 {
491 writeFileHeader(file());
492
493 writeValues(file());
494
495 Log << endl;
496 }
497
498 // Write state/results information
499 setResult("momentum_x", sumMomentum_[0]);
500 setResult("momentum_y", sumMomentum_[1]);
501 setResult("momentum_z", sumMomentum_[2]);
502
503 setResult("momentum_r", sumAngularMom_[0]);
504 setResult("momentum_rtheta", sumAngularMom_[1]);
505 setResult("momentum_axis", sumAngularMom_[2]);
506
507 return true;
508}
509
510
512{
513 if (writeMomentum_ || (hasCsys_ && (writeVelocity_ || writePosition_)))
514 {
515 Log << "Writing fields" << nl;
516
517 const volVectorField* fieldPtr;
518
519 fieldPtr = findObject<volVectorField>(scopedName("momentum"));
520 if (fieldPtr) fieldPtr->write();
521
522 fieldPtr = findObject<volVectorField>(scopedName("angularMomentum"));
523 if (fieldPtr) fieldPtr->write();
524
525 fieldPtr = findObject<volVectorField>(scopedName("angularVelocity"));
526 if (fieldPtr) fieldPtr->write();
527
528 if (hasCsys_ && writePosition_)
529 {
530 // Clunky, but currently no simple means of handling
531 // component-wise conversion and output
532
533 auto cyl_r = newField<volScalarField>("cyl_r", dimLength);
534 auto cyl_t = newField<volScalarField>("cyl_theta", dimless);
535 auto cyl_z = newField<volScalarField>("cyl_z", dimLength);
536
537 // Internal
538 {
539 const auto& pts = mesh_.cellCentres();
540 const label len = pts.size();
541
542 UList<scalar>& r = cyl_r->primitiveFieldRef(false);
543 UList<scalar>& t = cyl_t->primitiveFieldRef(false);
544 UList<scalar>& z = cyl_z->primitiveFieldRef(false);
545
546 for (label i=0; i < len; ++i)
547 {
548 point p(csys_.localPosition(pts[i]));
549
550 r[i] = p.x();
551 t[i] = p.y();
552 z[i] = p.z();
553 }
554 }
555
556 // Boundary
557 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
558
559 forAll(pbm, patchi)
560 {
561 if (isA<emptyPolyPatch>(pbm[patchi]))
562 {
563 continue;
564 }
565
566 const auto& pts = pbm[patchi].faceCentres();
567 const label len = pts.size();
568
569 UList<scalar>& r = cyl_r->boundaryFieldRef(false)[patchi];
570 UList<scalar>& t = cyl_t->boundaryFieldRef(false)[patchi];
571 UList<scalar>& z = cyl_z->boundaryFieldRef(false)[patchi];
572
573 for (label i=0; i < len; ++i)
574 {
575 point p(csys_.localPosition(pts[i]));
576
577 r[i] = p.x();
578 t[i] = p.y();
579 z[i] = p.z();
580 }
581 }
582
583 cyl_r->write();
584 cyl_t->write();
585 cyl_z->write();
587 }
588
589 return true;
590}
591
592
594{
596 purgeFields(); // Mesh motion makes calculated fields dubious
597}
598
599
601{
603 purgeFields(); // Mesh motion makes calculated fields dubious
604}
605
606
607// ************************************************************************* //
#define Log
Definition PDRblock.C:28
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
const polyBoundaryMesh & pbm
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition UList.H:89
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition UPstream.H:1714
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
T * get() noexcept
Return pointer to managed object without nullptr checking.
Definition autoPtr.H:216
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition autoPtr.H:178
A cylindrical coordinate system (r-theta-z). The coordinate system angle theta is always in radians.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Dimension set for the base types, which can be used to implement rigorous dimension checking for alge...
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.
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
Computes the natural logarithm of an input volScalarField.
Definition log.H:212
Computes linear/angular momentum, reporting integral values and optionally writing the fields.
Definition momentum.H:242
word pName_
The pressure field name (optional).
Definition momentum.H:294
scalar rhoRef_
Reference density (for incompressible).
Definition momentum.H:304
void initialise()
Initialise the fields.
Definition momentum.C:266
word UName_
The velocity field name (optional).
Definition momentum.H:287
bool hasCsys_
Are we using the cylindrical coordinate system?
Definition momentum.H:314
bool writePosition_
Write fields flag.
Definition momentum.H:329
bool writeMomentum_
Write fields flag.
Definition momentum.H:319
bool initialised_
Initialised flag.
Definition momentum.H:334
void writeValues(Ostream &os)
Write momentum data.
Definition momentum.C:299
momentum(const word &name, const Time &runTime, const dictionary &dict, const bool readFields=true)
Construct from name, Time and dictionary.
Definition momentum.C:343
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
Definition momentum.C:408
vector sumMomentum_
Integral (linear) momentum.
Definition momentum.H:274
word rhoName_
The density field name (optional).
Definition momentum.H:299
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition momentum.C:593
vector sumAngularMom_
Integral angular momentum.
Definition momentum.H:279
virtual void writeFileHeader(Ostream &os)
Output file header information.
Definition momentum.C:222
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition momentum.C:586
virtual bool execute()
Execute the function-object operations.
Definition momentum.C:478
virtual bool write()
Write the function-object results.
Definition momentum.C:504
coordSystem::cylindrical csys_
Coordinate system for evaluating angular momentum.
Definition momentum.H:309
bool writeVelocity_
Write fields flag.
Definition momentum.H:324
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition readFields.H:146
const ObjectType * cfindObject(const word &fieldName) const
Return const pointer to the object (eg, a field) in the (sub) objectRegistry.
bool foundObject(const word &fieldName) const
Find object (eg, a field) in the (sub) objectRegistry.
const objectRegistry & obr_
Reference to the region objectRegistry.
virtual const objectRegistry & obr() const
The region or sub-region registry being used.
const ObjectType * findObject(const word &fieldName) const
Return const pointer to the object (eg, a field) in the (sub) objectRegistry.
void setResult(const word &entryName, const Type &value)
Add result.
const Time & time_
Reference to the time database.
bool useAllCells() const noexcept
Use all cells, not the cellIDs.
Definition volRegionI.H:24
scalar V() const
Return total volume of the selected region.
Definition volRegionI.H:52
volRegion(const fvMesh &mesh, const dictionary &dict)
Construct from fvMesh and dictionary.
Definition volRegion.C:143
wordRe regionName_
Region name (cellSet, cellZone, ...).
Definition volRegion.H:188
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition volRegion.C:171
static const Enum< regionTypes > regionTypeNames_
Region type names.
Definition volRegion.H:130
bool update()
Update the cached values as required.
Definition volRegion.C:243
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
Definition volRegion.C:261
regionTypes regionType_
Region type.
Definition volRegion.H:183
virtual void updateMesh(const mapPolyMesh &)
Update for changes of mesh.
Definition volRegion.C:255
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
bool writtenHeader_
Flag to identify whether the header has been written.
Definition writeFile.H:157
virtual OFstream & file()
Return access to the file (if only 1).
Definition writeFile.C:270
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
virtual bool writeToFile() const
Flag to allow writing to file.
Definition writeFile.C:286
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Registry of regIOobjects.
bool checkOut(regIOobject *io) const
Remove a regIOobject from registry and free memory if the object is ownedByRegistry....
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Mesh consisting of general polyhedral cells.
Definition polyMesh.H:79
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
bool store()
Register object with its registry and transfer ownership to the registry.
A class for handling words, derived from Foam::string.
Definition word.H:66
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition zero.H:58
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
U
Definition pEqn.H:72
volScalarField & p
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
auto & name
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
const dimensionSet dimPressure
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
void readFields(const typename GeoFieldType::Mesh &mesh, const IOobjectList &objects, const NameMatchPredicate &selectedFields, DynamicList< regIOobject * > &storedObjects)
Read the selected GeometricFields of the templated type and store on the objectRegistry.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
static void writeHeader(Ostream &os, const word &fieldName)
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere).
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
dimensionedScalar log(const dimensionedScalar &ds)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
vector point
Point is a vector.
Definition point.H:37
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
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
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
constexpr char tab
The tab '\t' character(0x09).
Definition Ostream.H:49
static bool initialised_(false)
Info<< "Reading mechanical properties\n"<< endl;IOdictionary mechanicalProperties(IOobject("mechanicalProperties", runTime.constant(), mesh, IOobject::MUST_READ_IF_MODIFIED, IOobject::NO_WRITE));const dictionary &rhoDict(mechanicalProperties.subDict("rho"));word rhoType(rhoDict.get< word >("type"));autoPtr< volScalarField > rhoPtr
dictionary dict
const pointField & pts
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299