forces.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) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2015-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "forces.H"
30#include "fvcGrad.H"
31#include "porosityModel.H"
34#include "cartesianCS.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
41namespace functionObjects
42{
45}
46}
47
48
49// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
50
52(
53 const dictionary& dict,
54 const word& e3Name,
55 const word& e1Name
56)
57{
58 coordSysPtr_.clear();
59
60 point origin(Zero);
61 if (dict.readIfPresent<point>("CofR", origin))
62 {
63 const vector e3 = e3Name == word::null ?
64 vector(0, 0, 1) : dict.get<vector>(e3Name);
65 const vector e1 = e1Name == word::null ?
66 vector(1, 0, 0) : dict.get<vector>(e1Name);
67
68 coordSysPtr_.reset(new coordSystem::cartesian(origin, e3, e1));
69 }
70 else
71 {
72 // The 'coordinateSystem' sub-dictionary is optional,
73 // but enforce use of a cartesian system if not found.
74
75 if (dict.found(coordinateSystem::typeName_()))
76 {
77 // New() for access to indirect (global) coordinate system
80 (
81 obr_,
82 dict,
83 coordinateSystem::typeName_()
84 );
85 }
86 else
87 {
89 }
90 }
91}
92
93
95{
96 auto* forcePtr = mesh_.getObjectPtr<volVectorField>(scopedName("force"));
97
98 if (!forcePtr)
99 {
100 forcePtr = new volVectorField
101 (
103 (
104 scopedName("force"),
105 time_.timeName(),
106 mesh_,
109 ),
110 mesh_,
112 );
113
114 mesh_.objectRegistry::store(forcePtr);
115 }
116
117 return *forcePtr;
118}
119
120
122{
123 auto* momentPtr = mesh_.getObjectPtr<volVectorField>(scopedName("moment"));
124
125 if (!momentPtr)
126 {
127 momentPtr = new volVectorField
128 (
130 (
131 scopedName("moment"),
132 time_.timeName(),
133 mesh_,
136 ),
137 mesh_,
139 );
140
141 mesh_.objectRegistry::store(momentPtr);
142 }
143
144 return *momentPtr;
145}
146
147
149{
150 if (initialised_)
151 {
152 return;
153 }
154
155 if (directForceDensity_)
156 {
157 if (!foundObject<volVectorField>(fDName_))
158 {
160 << "Could not find " << fDName_ << " in database"
161 << exit(FatalError);
162 }
163 }
164 else
165 {
166 if
167 (
168 !foundObject<volVectorField>(UName_)
169 || !foundObject<volScalarField>(pName_)
170 )
171 {
173 << "Could not find U: " << UName_
174 << " or p:" << pName_ << " in database"
175 << exit(FatalError);
176 }
177
178 if (rhoName_ != "rhoInf" && !foundObject<volScalarField>(rhoName_))
179 {
181 << "Could not find rho:" << rhoName_ << " in database"
182 << exit(FatalError);
183 }
184 }
185
186 initialised_ = true;
187}
188
189
191{
192 sumPatchForcesP_ = Zero;
193 sumPatchForcesV_ = Zero;
194 sumPatchMomentsP_ = Zero;
195 sumPatchMomentsV_ = Zero;
196
197 sumInternalForces_ = Zero;
198 sumInternalMoments_ = Zero;
199
200 auto& force = this->force();
201 auto& moment = this->moment();
202 force == dimensionedVector(force.dimensions(), Zero);
203 moment == dimensionedVector(moment.dimensions(), Zero);
204}
205
206
209{
210 typedef compressible::turbulenceModel cmpTurbModel;
211 typedef incompressible::turbulenceModel icoTurbModel;
212
213 if (foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
214 {
215 const auto& turb =
216 lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
217
218 return turb.devRhoReff();
219 }
220 else if (foundObject<icoTurbModel>(icoTurbModel::propertiesName))
221 {
222 const auto& turb =
223 lookupObject<icoTurbModel>(icoTurbModel::propertiesName);
224
225 return rho()*turb.devReff();
226 }
227 else if (foundObject<fluidThermo>(fluidThermo::dictName))
228 {
229 const auto& thermo = lookupObject<fluidThermo>(fluidThermo::dictName);
230
231 const auto& U = lookupObject<volVectorField>(UName_);
232
233 return -thermo.mu()*dev(twoSymm(fvc::grad(U)));
234 }
235 else if (foundObject<transportModel>("transportProperties"))
236 {
237 const auto& laminarT =
238 lookupObject<transportModel>("transportProperties");
239
240 const auto& U = lookupObject<volVectorField>(UName_);
241
242 return -rho()*laminarT.nu()*dev(twoSymm(fvc::grad(U)));
243 }
244 else if (foundObject<dictionary>("transportProperties"))
245 {
246 const auto& transportProperties =
247 lookupObject<dictionary>("transportProperties");
248
250
251 const auto& U = lookupObject<volVectorField>(UName_);
252
253 return -rho()*nu*dev(twoSymm(fvc::grad(U)));
254 }
255 else
256 {
258 << "No valid model for viscous stress calculation"
259 << exit(FatalError);
260
262 }
263}
264
265
267{
268 if (foundObject<fluidThermo>(basicThermo::dictName))
269 {
270 const auto& thermo = lookupObject<fluidThermo>(basicThermo::dictName);
271
272 return thermo.mu();
273 }
274 else if (foundObject<transportModel>("transportProperties"))
275 {
276 const auto& laminarT =
277 lookupObject<transportModel>("transportProperties");
278
279 return rho()*laminarT.nu();
280 }
281 else if (foundObject<dictionary>("transportProperties"))
282 {
283 const auto& transportProperties =
284 lookupObject<dictionary>("transportProperties");
285
287
288 return rho()*nu;
289 }
290 else
291 {
293 << "No valid model for dynamic viscosity calculation"
294 << exit(FatalError);
295
296 return volScalarField::null();
297 }
298}
299
300
302{
303 if (rhoName_ == "rhoInf")
304 {
306 (
308 (
309 "rho",
310 mesh_.time().timeName(),
311 mesh_
312 ),
313 mesh_,
315 );
316 }
317
318 return (lookupObject<volScalarField>(rhoName_));
319}
320
321
323{
324 if (p.dimensions() == dimPressure)
325 {
326 return 1;
327 }
328
329 if (rhoName_ != "rhoInf")
330 {
332 << "Dynamic pressure is expected but kinematic is provided."
333 << exit(FatalError);
334 }
335
336 return rhoRef_;
337}
338
339
341(
342 const label patchi,
343 const vectorField& Md,
344 const vectorField& fP,
345 const vectorField& fV
346)
347{
348 sumPatchForcesP_ += sum(fP);
349 sumPatchForcesV_ += sum(fV);
350 force().boundaryFieldRef()[patchi] += fP + fV;
351
352 const vectorField mP(Md^fP);
353 const vectorField mV(Md^fV);
354
355 sumPatchMomentsP_ += sum(mP);
356 sumPatchMomentsV_ += sum(mV);
357 moment().boundaryFieldRef()[patchi] += mP + mV;
358}
359
360
362(
363 const labelList& cellIDs,
364 const vectorField& Md,
365 const vectorField& f
366)
367{
368 auto& force = this->force();
369 auto& moment = this->moment();
370
371 forAll(cellIDs, i)
372 {
373 const label celli = cellIDs[i];
374
375 sumInternalForces_ += f[i];
376 force[celli] += f[i];
377
378 const vector m(Md[i]^f[i]);
379 sumInternalMoments_ += m;
380 moment[celli] = m;
381 }
382}
383
384
386{
387 if (!forceFilePtr_.valid())
388 {
389 forceFilePtr_ = createFile("force");
390 writeIntegratedDataFileHeader("Force", forceFilePtr_());
391 }
392
393 if (!momentFilePtr_.valid())
394 {
395 momentFilePtr_ = createFile("moment");
396 writeIntegratedDataFileHeader("Moment", momentFilePtr_());
397 }
398}
399
400
402(
403 const word& header,
404 OFstream& os
405) const
406{
407 const auto& coordSys = coordSysPtr_();
408 const auto vecDesc = [](const word& root)->string
409 {
410 return root + "_x " + root + "_y " + root + "_z";
411 };
412 writeHeader(os, header);
413 writeHeaderValue(os, "CofR", coordSys.origin());
414 writeHeader(os, "");
415 writeCommented(os, "Time");
416 writeTabbed(os, vecDesc("total"));
417 writeTabbed(os, vecDesc("pressure"));
418 writeTabbed(os, vecDesc("viscous"));
419
420 if (porosity_)
421 {
422 writeTabbed(os, vecDesc("porous"));
423 }
424
425 os << endl;
426}
427
428
430{
431 const auto& coordSys = coordSysPtr_();
432
433 writeIntegratedDataFile
434 (
435 coordSys.localVector(sumPatchForcesP_),
436 coordSys.localVector(sumPatchForcesV_),
437 coordSys.localVector(sumInternalForces_),
438 forceFilePtr_()
439 );
440
441 writeIntegratedDataFile
442 (
443 coordSys.localVector(sumPatchMomentsP_),
444 coordSys.localVector(sumPatchMomentsV_),
445 coordSys.localVector(sumInternalMoments_),
446 momentFilePtr_()
447 );
448}
449
450
452(
453 const vector& pres,
454 const vector& vis,
455 const vector& internal,
456 OFstream& os
457) const
458{
459 writeCurrentTime(os);
460
461 writeValue(os, pres + vis + internal);
462 writeValue(os, pres);
463 writeValue(os, vis);
464
465 if (porosity_)
466 {
467 writeValue(os, internal);
468 }
469
470 os << endl;
471}
472
473
475(
476 const string& descriptor,
477 const vector& pres,
478 const vector& vis,
479 const vector& internal
480) const
481{
482 if (!log)
483 {
484 return;
485 }
486
487 Log << " Sum of " << descriptor.c_str() << nl
488 << " Total : " << (pres + vis + internal) << nl
489 << " Pressure : " << pres << nl
490 << " Viscous : " << vis << nl;
491
492 if (porosity_)
493 {
494 Log << " Porous : " << internal << nl;
495 }
496}
497
498
499// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
500
502(
503 const word& name,
504 const Time& runTime,
505 const dictionary& dict,
506 bool readFields
507)
508:
510 writeFile(mesh_, name),
511 sumPatchForcesP_(Zero),
512 sumPatchForcesV_(Zero),
513 sumPatchMomentsP_(Zero),
514 sumPatchMomentsV_(Zero),
515 sumInternalForces_(Zero),
516 sumInternalMoments_(Zero),
517 forceFilePtr_(),
518 momentFilePtr_(),
519 coordSysPtr_(nullptr),
520 patchSet_(),
521 rhoRef_(VGREAT),
522 pRef_(0),
523 pName_("p"),
524 UName_("U"),
525 rhoName_("rho"),
526 fDName_("fD"),
527 directForceDensity_(false),
528 porosity_(false),
529 writeFields_(false),
530 initialised_(false)
531{
532 if (readFields)
533 {
534 read(dict);
536 Log << endl;
537 }
538}
539
540
542(
543 const word& name,
544 const objectRegistry& obr,
545 const dictionary& dict,
546 bool readFields
547)
548:
550 writeFile(mesh_, name),
551 sumPatchForcesP_(Zero),
552 sumPatchForcesV_(Zero),
553 sumPatchMomentsP_(Zero),
554 sumPatchMomentsV_(Zero),
555 sumInternalForces_(Zero),
556 sumInternalMoments_(Zero),
557 forceFilePtr_(),
558 momentFilePtr_(),
559 coordSysPtr_(nullptr),
560 patchSet_(),
561 rhoRef_(VGREAT),
562 pRef_(0),
563 pName_("p"),
564 UName_("U"),
565 rhoName_("rho"),
566 fDName_("fD"),
567 directForceDensity_(false),
568 porosity_(false),
569 writeFields_(false),
570 initialised_(false)
571{
572 if (readFields)
573 {
574 read(dict);
576 Log << endl;
577 }
578}
579
580
581// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
582
584{
586 {
587 return false;
588 }
589
590 initialised_ = false;
591
592 Info<< type() << " " << name() << ":" << endl;
593
594 patchSet_ =
595 mesh_.boundaryMesh().patchSet
596 (
597 dict.get<wordRes>("patches")
598 );
599
600 dict.readIfPresent("directForceDensity", directForceDensity_);
601 if (directForceDensity_)
602 {
603 // Optional name entry for fD
604 if (dict.readIfPresent<word>("fD", fDName_))
605 {
606 Info<< " fD: " << fDName_ << endl;
607 }
608 }
609 else
610 {
611 // Optional field name entries
612 if (dict.readIfPresent<word>("p", pName_))
613 {
614 Info<< " p: " << pName_ << endl;
615 }
616 if (dict.readIfPresent<word>("U", UName_))
617 {
618 Info<< " U: " << UName_ << endl;
619 }
620 if (dict.readIfPresent<word>("rho", rhoName_))
621 {
622 Info<< " rho: " << rhoName_ << endl;
623 }
624
625 // Reference density needed for incompressible calculations
626 if (rhoName_ == "rhoInf")
627 {
628 rhoRef_ = dict.getCheck<scalar>("rhoInf", scalarMinMax::ge(SMALL));
629 Info<< " Freestream density (rhoInf) set to " << rhoRef_ << endl;
630 }
631
632 // Reference pressure, 0 by default
633 if (dict.readIfPresent<scalar>("pRef", pRef_))
634 {
635 Info<< " Reference pressure (pRef) set to " << pRef_ << endl;
636 }
637 }
638
639 dict.readIfPresent("porosity", porosity_);
640 if (porosity_)
641 {
642 Info<< " Including porosity effects" << endl;
643 }
644 else
645 {
646 Info<< " Not including porosity effects" << endl;
647 }
648
649 writeFields_ = dict.getOrDefault("writeFields", false);
650 if (writeFields_)
651 {
652 Info<< " Fields will be written" << endl;
653 }
654
655
656 return true;
657}
658
659
661{
662 initialise();
663
664 reset();
665
666 const point& origin = coordSysPtr_->origin();
667
668 if (directForceDensity_)
669 {
670 const auto& fD = lookupObject<volVectorField>(fDName_);
671
672 const auto& Sfb = mesh_.Sf().boundaryField();
673
674 for (const label patchi : patchSet_)
675 {
676 const vectorField& d = mesh_.C().boundaryField()[patchi];
677
678 const vectorField Md(d - origin);
679
680 const scalarField sA(mag(Sfb[patchi]));
681
682 // Pressure force = surfaceUnitNormal*(surfaceNormal & forceDensity)
683 const vectorField fP
684 (
685 Sfb[patchi]/sA
686 *(
687 Sfb[patchi] & fD.boundaryField()[patchi]
688 )
689 );
690
691 // Viscous force (total force minus pressure fP)
692 const vectorField fV(sA*fD.boundaryField()[patchi] - fP);
693
694 addToPatchFields(patchi, Md, fP, fV);
695 }
696 }
697 else
698 {
699 const auto& p = lookupObject<volScalarField>(pName_);
700
701 const auto& Sfb = mesh_.Sf().boundaryField();
702
703 tmp<volSymmTensorField> tdevRhoReff = devRhoReff();
704 const auto& devRhoReffb = tdevRhoReff().boundaryField();
705
706 // Scale pRef by density for incompressible simulations
707 const scalar rhoRef = rho(p);
708 const scalar pRef = pRef_/rhoRef;
709
710 for (const label patchi : patchSet_)
711 {
712 const vectorField& d = mesh_.C().boundaryField()[patchi];
713
714 const vectorField Md(d - origin);
715
716 const vectorField fP
717 (
718 rhoRef*Sfb[patchi]*(p.boundaryField()[patchi] - pRef)
719 );
720
721 const vectorField fV(Sfb[patchi] & devRhoReffb[patchi]);
722
723 addToPatchFields(patchi, Md, fP, fV);
724 }
725 }
726
727 if (porosity_)
728 {
729 const auto& U = lookupObject<volVectorField>(UName_);
730 const volScalarField rho(this->rho());
731 const volScalarField mu(this->mu());
732
733 const auto models = obr_.lookupClass<porosityModel>();
734
735 if (models.empty())
736 {
738 << "Porosity effects requested, "
739 << "but no porosity models found in the database"
740 << endl;
741 }
742
743 forAllConstIters(models, iter)
744 {
745 // Non-const access required if mesh is changing
746 auto& pm = const_cast<porosityModel&>(*iter());
747
748 const vectorField fPTot(pm.force(U, rho, mu));
749
750 const labelList& cellZoneIDs = pm.cellZoneIDs();
751
752 for (const label zonei : cellZoneIDs)
753 {
754 const cellZone& cZone = mesh_.cellZones()[zonei];
755
756 const vectorField d(mesh_.C(), cZone);
757 const vectorField fP(fPTot, cZone);
758 const vectorField Md(d - origin);
759
760 addToInternalField(cZone, Md, fP);
761 }
762 }
763 }
764
765 reduce(sumPatchForcesP_, sumOp<vector>());
766 reduce(sumPatchForcesV_, sumOp<vector>());
767 reduce(sumPatchMomentsP_, sumOp<vector>());
768 reduce(sumPatchMomentsV_, sumOp<vector>());
769 reduce(sumInternalForces_, sumOp<vector>());
770 reduce(sumInternalMoments_, sumOp<vector>());
771}
772
773
775{
776 return sumPatchForcesP_ + sumPatchForcesV_ + sumInternalForces_;
777}
778
779
781{
782 return sumPatchMomentsP_ + sumPatchMomentsV_ + sumInternalMoments_;
783}
784
785
787{
788 calcForcesMoments();
789
790 Log << type() << " " << name() << " write:" << nl;
791
792 const auto& coordSys = coordSysPtr_();
793
794 const auto localFp(coordSys.localVector(sumPatchForcesP_));
795 const auto localFv(coordSys.localVector(sumPatchForcesV_));
796 const auto localFi(coordSys.localVector(sumInternalForces_));
797
798 logIntegratedData("forces", localFp, localFv, localFi);
799
800 const auto localMp(coordSys.localVector(sumPatchMomentsP_));
801 const auto localMv(coordSys.localVector(sumPatchMomentsV_));
802 const auto localMi(coordSys.localVector(sumInternalMoments_));
803
804 logIntegratedData("moments", localMp, localMv, localMi);
805
806 setResult("pressureForce", localFp);
807 setResult("viscousForce", localFv);
808 setResult("internalForce", localFi);
809 setResult("pressureMoment", localMp);
810 setResult("viscousMoment", localMv);
811 setResult("internalMoment", localMi);
812
813 return true;
814}
815
816
818{
819 if (writeToFile())
820 {
821 Log << " writing force and moment files." << endl;
822
823 createIntegratedDataFiles();
824 writeIntegratedDataFiles();
825 }
826
827 if (writeFields_)
828 {
829 Log << " writing force and moment fields." << endl;
830
831 force().write();
832 moment().write();
833 }
834
835 Log << endl;
836
837 return true;
838}
839
840
841// ************************************************************************* //
#define Log
Definition: PDRblock.C:35
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
labelList cellIDs
compressible::turbulenceModel & turb
static const GeometricField< symmTensor, fvPatchField, volMesh > & null()
Return a null geometric field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
Templated abstract base class for single-phase incompressible turbulence models.
static MinMax< scalar > ge(const scalar &minVal)
A semi-infinite range from minVal to the type max.
Definition: MinMaxI.H:31
Output to file stream, using an OSstream.
Definition: OFstream.H:57
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
virtual bool read()
Re-read model coefficients if they have changed.
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
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
static const word dictName
Definition: basicThermo.H:256
A subset of mesh cells.
Definition: cellZone.H:65
A Cartesian coordinate system.
Definition: cartesianCS.H:72
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Abstract base-class for Time/database function objects.
void writeIntegratedDataFile()
Write integrated coefficients to the integrated-coefficient file.
Definition: forceCoeffs.C:250
Computes forces and moments over a given list of patches by integrating pressure and viscous forces a...
Definition: forces.H:325
virtual void calcForcesMoments()
Calculate forces and moments.
Definition: forces.C:660
void writeIntegratedDataFileHeader(const word &header, OFstream &os) const
Write header for an integrated-data file.
Definition: forces.C:402
void initialise()
Initialise containers and fields.
Definition: forces.C:148
void addToPatchFields(const label patchi, const vectorField &Md, const vectorField &fP, const vectorField &fV)
Add patch contributions to force and moment fields.
Definition: forces.C:341
void writeIntegratedDataFiles()
Write integrated data to files.
Definition: forces.C:429
volVectorField & moment()
Return access to the moment field.
Definition: forces.C:121
void createIntegratedDataFiles()
Create the integrated-data files.
Definition: forces.C:385
tmp< volScalarField > mu() const
Return dynamic viscosity field.
Definition: forces.C:266
virtual bool read(const dictionary &dict)
Read the dictionary.
Definition: forces.C:583
virtual vector forceEff() const
Return the total force.
Definition: forces.C:774
void logIntegratedData(const string &descriptor, const vector &pres, const vector &vis, const vector &internal) const
Write integrated data to stream.
Definition: forces.C:475
autoPtr< coordinateSystem > coordSysPtr_
Coordinate system used when evaluating forces and moments.
Definition: forces.H:363
virtual vector momentEff() const
Return the total moment.
Definition: forces.C:780
void addToInternalField(const labelList &cellIDs, const vectorField &Md, const vectorField &f)
Definition: forces.C:362
tmp< volScalarField > rho() const
Return rho if specified otherwise rhoRef.
Definition: forces.C:301
void reset()
Reset containers and fields.
Definition: forces.C:190
virtual bool execute()
Execute the function object.
Definition: forces.C:786
volVectorField & force()
Return access to the force field.
Definition: forces.C:94
virtual bool write()
Write to data files/fields and to streams.
Definition: forces.C:817
tmp< volSymmTensorField > devRhoReff() const
Return the effective stress (viscous + turbulent)
Definition: forces.C:208
void setCoordinateSystem(const dictionary &dict, const word &e3Name=word::null, const word &e1Name=word::null)
Set the co-ordinate system from dictionary and axes names.
Definition: forces.C:52
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Computes the natural logarithm of an input volScalarField.
Definition: log.H:230
Computes the magnitude of an input field.
Definition: mag.H:153
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:158
const objectRegistry & obr_
Reference to the region objectRegistry.
Base class for writing single files from the function objects.
Definition: writeFile.H:120
Registry of regIOobjects.
Top level model for porosity models.
Definition: porosityModel.H:60
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A class for managing temporary objects.
Definition: tmp.H:65
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
U
Definition: pEqn.H:72
volScalarField & p
const volScalarField & mu
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
Calculate the gradient of the given field.
#define WarningInFunction
Report a warning using Foam::Warning.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Namespace for OpenFOAM.
const dimensionSet dimPressure
const dimensionSet dimViscosity
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:83
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static void writeHeader(Ostream &os, const word &fieldName)
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: MSwindows.C:598
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
const dimensionSet dimForce
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Vector< scalar > vector
Definition: vector.H:61
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
volScalarField & nu
dictionary dict
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition: stdFoam.H:278