Loading...
Searching...
No Matches
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-2025 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"
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 point origin(Zero);
59
60 // With objectRegistry for access to indirect (global) coordinate systems
61 coordSysPtr_ = coordinateSystem::NewIfPresent(obr_, dict);
62
63 if (coordSysPtr_)
64 {
65 // Report ...
66 }
67 else if (dict.readIfPresent("CofR", origin))
68 {
69 const vector e3
70 (
71 e3Name.empty() ? vector(0, 0, 1) : dict.get<vector>(e3Name)
72 );
73 const vector e1
74 (
75 e1Name.empty() ? vector(1, 0, 0) : dict.get<vector>(e1Name)
76 );
77
78 coordSysPtr_.reset(new coordSystem::cartesian(origin, e3, e1));
79 }
80 else
81 {
82 // No 'coordinateSystem' or 'CofR'
83 // - enforce a cartesian system
84
86 }
87}
88
89
91{
92 auto* ptr = mesh_.getObjectPtr<volVectorField>(scopedName("force"));
93
94 if (!ptr)
95 {
96 ptr = new volVectorField
97 (
99 (
100 scopedName("force"),
101 time_.timeName(),
102 mesh_.thisDb(),
106 ),
107 mesh_,
109 );
110
112 }
113
114 return *ptr;
115}
116
117
119{
120 auto* ptr = mesh_.getObjectPtr<volVectorField>(scopedName("moment"));
121
122 if (!ptr)
123 {
124 ptr = new volVectorField
125 (
127 (
128 scopedName("moment"),
129 time_.timeName(),
130 mesh_.thisDb(),
134 ),
135 mesh_,
137 );
138
140 }
141
142 return *ptr;
143}
144
145
147{
148 if (initialised_)
149 {
150 return;
151 }
152
153 if (directForceDensity_)
154 {
155 if (!foundObject<volVectorField>(fDName_))
156 {
158 << "Could not find " << fDName_ << " in database"
159 << exit(FatalError);
160 }
161 }
162 else
163 {
164 if
165 (
166 !foundObject<volVectorField>(UName_)
167 || !foundObject<volScalarField>(pName_)
168 )
169 {
171 << "Could not find U: " << UName_
172 << " or p:" << pName_ << " in database"
173 << exit(FatalError);
174 }
175
176 if (rhoName_ != "rhoInf" && !foundObject<volScalarField>(rhoName_))
177 {
179 << "Could not find rho:" << rhoName_ << " in database"
180 << exit(FatalError);
182 }
183
184 initialised_ = true;
185}
186
187
189{
190 sumPatchForcesP_ = Zero;
191 sumPatchForcesV_ = Zero;
192 sumPatchMomentsP_ = Zero;
193 sumPatchMomentsV_ = Zero;
194
195 sumInternalForces_ = Zero;
196 sumInternalMoments_ = Zero;
197
198 auto& force = this->force();
199 auto& moment = this->moment();
200
201 if (porosity_)
202 {
203 force == Zero;
204 moment == Zero;
205 }
206 else
207 {
208 for (const label patchi : patchIDs_)
209 {
210 force.boundaryField().constCast()[patchi] = Zero;
211 moment.boundaryField().constCast()[patchi] = Zero;
212 }
213 }
214}
215
216
219(
220 const tensorField& gradUp,
221 const label patchi
222) const
223{
224 typedef incompressible::turbulenceModel icoModel;
225 typedef compressible::turbulenceModel cmpModel;
226
227 if (const auto* turbp = cfindObject<icoModel>(icoModel::propertiesName))
228 {
229 const auto& turb = *turbp;
230
231 return -rho(patchi)*turb.nuEff(patchi)*devTwoSymm(gradUp);
232 }
233
234 if (const auto* turbp = cfindObject<cmpModel>(cmpModel::propertiesName))
235 {
236 const auto& turb = *turbp;
237
238 return -turb.muEff(patchi)*devTwoSymm(gradUp);
239 }
240
241 if (const auto* thermop = cfindObject<fluidThermo>(fluidThermo::dictName))
242 {
243 const auto& thermo = *thermop;
244
245 return -thermo.mu(patchi)*devTwoSymm(gradUp);
246 }
247
248 if (const auto* props = cfindObject<transportModel>("transportProperties"))
249 {
250 const auto& laminarT = *props;
251
252 return -rho(patchi)*laminarT.nu(patchi)*devTwoSymm(gradUp);
253 }
254
255 if (const auto* props = cfindObject<dictionary>("transportProperties"))
256 {
257 const dimensionedScalar nu("nu", dimViscosity, *props);
258
259 return -rho(patchi)*nu.value()*devTwoSymm(gradUp);
260 }
261
262 // Failed to resolve any model
264 << "No valid model for viscous stress calculation"
265 << exit(FatalError);
266
267 return nullptr;
268}
269
270
272{
273 if (const auto* thermop = cfindObject<fluidThermo>(basicThermo::dictName))
274 {
275 const auto& thermo = *thermop;
276
277 return thermo.mu();
278 }
279
280 if (const auto* props = cfindObject<transportModel>("transportProperties"))
281 {
282 const auto& laminarT = *props;
283
284 return rho()*laminarT.nu();
285 }
286
287 if (const auto* props = cfindObject<dictionary>("transportProperties"))
288 {
289 const dimensionedScalar nu("nu", dimViscosity, *props);
290
291 return rho()*nu;
292 }
293
294 // Failed to resolve any model
296 << "No valid model for dynamic viscosity calculation"
297 << exit(FatalError);
298
299 return nullptr;
300}
301
302
304{
305 if (rhoName_ == "rhoInf")
306 {
308 (
309 "rho",
311 mesh_,
313 );
315
316 return lookupObject<volScalarField>(rhoName_);
317}
318
319
321Foam::functionObjects::forces::rho(const label patchi) const
322{
323 if (rhoName_ == "rhoInf")
324 {
326 (
327 mesh_.boundary()[patchi].size(),
328 rhoRef_
329 );
331
332 const auto& rho = lookupObject<volScalarField>(rhoName_);
333 return rho.boundaryField()[patchi];
334}
335
336
337Foam::scalar Foam::functionObjects::forces::rho(const volScalarField& p) const
338{
339 if (p.dimensions() == dimPressure)
340 {
341 return 1;
342 }
343
344 if (rhoName_ != "rhoInf")
345 {
347 << "Dynamic pressure is expected but kinematic is provided."
349 }
350
351 return rhoRef_;
352}
353
354
356(
357 const label patchi,
358 const vectorField& Md,
359 const vectorField& fP,
360 const vectorField& fV
361)
362{
363 sumPatchForcesP_ += sum(fP);
364 sumPatchForcesV_ += sum(fV);
365 force().boundaryField().constCast()[patchi] += fP + fV;
366
367 const vectorField mP(Md^fP);
368 const vectorField mV(Md^fV);
370 sumPatchMomentsP_ += sum(mP);
371 sumPatchMomentsV_ += sum(mV);
372 moment().boundaryField().constCast()[patchi] += mP + mV;
373}
374
375
377(
378 const labelList& cellIDs,
379 const vectorField& Md,
380 const vectorField& f
381)
382{
383 auto& force = this->force();
384 auto& moment = this->moment();
385
386 forAll(cellIDs, i)
387 {
388 const label celli = cellIDs[i];
389
390 sumInternalForces_ += f[i];
391 force[celli] += f[i];
392
393 const vector m(Md[i]^f[i]);
395 moment[celli] = m;
396 }
397}
398
399
401{
402 if (!forceFilePtr_)
403 {
404 forceFilePtr_ = newFileAtStartTime("force");
405 writeIntegratedDataFileHeader("Force", forceFilePtr_());
406 }
407
408 if (!momentFilePtr_)
412 }
413}
414
415
417(
418 const word& header,
419 OFstream& os
420) const
421{
422 const auto& coordSys = coordSysPtr_();
423 const auto vecDesc = [](const word& root)->string
424 {
425 return root + "_x " + root + "_y " + root + "_z";
426 };
427 writeHeader(os, header);
428 writeHeaderValue(os, "CofR", coordSys.origin());
429 writeHeader(os, "");
430 writeCommented(os, "Time");
431 writeTabbed(os, vecDesc("total"));
432 writeTabbed(os, vecDesc("pressure"));
433 writeTabbed(os, vecDesc("viscous"));
434
435 if (porosity_)
436 {
437 writeTabbed(os, vecDesc("porous"));
438 }
439
440 os << endl;
441}
442
443
445{
446 const auto& coordSys = coordSysPtr_();
447
448 writeIntegratedDataFile
449 (
450 coordSys.localVector(sumPatchForcesP_),
451 coordSys.localVector(sumPatchForcesV_),
452 coordSys.localVector(sumInternalForces_),
453 forceFilePtr_()
454 );
455
456 writeIntegratedDataFile
457 (
458 coordSys.localVector(sumPatchMomentsP_),
459 coordSys.localVector(sumPatchMomentsV_),
460 coordSys.localVector(sumInternalMoments_),
462 );
463}
464
465
467(
468 const vector& pres,
469 const vector& vis,
470 const vector& internal,
471 OFstream& os
472) const
473{
474 writeCurrentTime(os);
475
476 writeValue(os, pres + vis + internal);
477 writeValue(os, pres);
478 writeValue(os, vis);
479
480 if (porosity_)
481 {
482 writeValue(os, internal);
483 }
484
485 os << endl;
486}
487
488
490(
491 const string& descriptor,
492 const vector& pres,
493 const vector& vis,
494 const vector& internal
495) const
496{
497 if (!log)
498 {
499 return;
500 }
501
502 Log << " Sum of " << descriptor.c_str() << nl
503 << " Total : " << (pres + vis + internal) << nl
504 << " Pressure : " << pres << nl
505 << " Viscous : " << vis << nl;
506
507 if (porosity_)
508 {
509 Log << " Porous : " << internal << nl;
510 }
511}
512
513
514// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
515
517(
518 const word& name,
519 const Time& runTime,
520 const dictionary& dict,
521 bool readFields
522)
523:
524 fvMeshFunctionObject(name, runTime, dict),
525 writeFile(mesh_, name),
526 sumPatchForcesP_(Zero),
527 sumPatchForcesV_(Zero),
528 sumPatchMomentsP_(Zero),
529 sumPatchMomentsV_(Zero),
530 sumInternalForces_(Zero),
531 sumInternalMoments_(Zero),
532 forceFilePtr_(),
533 momentFilePtr_(),
534 coordSysPtr_(nullptr),
535 rhoRef_(VGREAT),
536 pRef_(0),
537 pName_("p"),
538 UName_("U"),
539 rhoName_("rho"),
540 fDName_("fD"),
541 directForceDensity_(false),
542 porosity_(false),
543 writeFields_(false),
544 initialised_(false)
545{
546 if (readFields)
547 {
550 Log << endl;
551 }
552}
553
554
556(
557 const word& name,
558 const objectRegistry& obr,
559 const dictionary& dict,
560 bool readFields
561)
562:
563 fvMeshFunctionObject(name, obr, dict),
564 writeFile(mesh_, name),
565 sumPatchForcesP_(Zero),
566 sumPatchForcesV_(Zero),
567 sumPatchMomentsP_(Zero),
568 sumPatchMomentsV_(Zero),
569 sumInternalForces_(Zero),
570 sumInternalMoments_(Zero),
571 forceFilePtr_(),
572 momentFilePtr_(),
573 coordSysPtr_(nullptr),
574 rhoRef_(VGREAT),
575 pRef_(0),
576 pName_("p"),
577 UName_("U"),
578 rhoName_("rho"),
579 fDName_("fD"),
580 directForceDensity_(false),
581 porosity_(false),
582 writeFields_(false),
583 initialised_(false)
584{
585 if (readFields)
586 {
587 read(dict);
590 }
591}
592
593
594// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
595
597{
598 const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
599
601 {
602 return false;
603 }
604
605 initialised_ = false;
606
607 Info<< type() << ' ' << name() << ':' << endl;
608
609 // Can also use pbm.indices(), but no warnings...
610 patchIDs_ = pbm.patchSet(dict.get<wordRes>("patches")).sortedToc();
611
612 dict.readIfPresent("directForceDensity", directForceDensity_);
613 if (directForceDensity_)
614 {
615 // Optional name entry for fD
616 if (dict.readIfPresent<word>("fD", fDName_))
617 {
618 Info<< " fD: " << fDName_ << endl;
619 }
620 }
621 else
622 {
623 // Optional field name entries
624 if (dict.readIfPresent<word>("p", pName_))
625 {
626 Info<< " p: " << pName_ << endl;
627 }
628 if (dict.readIfPresent<word>("U", UName_))
629 {
630 Info<< " U: " << UName_ << endl;
631 }
632 if (dict.readIfPresent<word>("rho", rhoName_))
633 {
634 Info<< " rho: " << rhoName_ << endl;
635 }
636
637 // Reference density needed for incompressible calculations
638 if (rhoName_ == "rhoInf")
639 {
640 rhoRef_ = dict.getCheck<scalar>("rhoInf", scalarMinMax::ge(SMALL));
641 Info<< " Freestream density (rhoInf) set to " << rhoRef_ << endl;
642 }
643
644 // Reference pressure, 0 by default
645 if (dict.readIfPresent<scalar>("pRef", pRef_))
646 {
647 Info<< " Reference pressure (pRef) set to " << pRef_ << endl;
648 }
649 }
650
651 dict.readIfPresent("porosity", porosity_);
652 if (porosity_)
653 {
654 Info<< " Including porosity effects" << endl;
655 }
656 else
657 {
658 Info<< " Not including porosity effects" << endl;
659 }
660
661 writeFields_ = dict.getOrDefault("writeFields", false);
662 if (writeFields_)
663 {
664 Info<< " Fields will be written" << endl;
666
667
668 return true;
669}
670
671
673{
674 initialise();
675
676 reset();
677
678 const point& origin = coordSysPtr_->origin();
679
680 if (directForceDensity_)
681 {
682 const auto& fD = lookupObject<volVectorField>(fDName_);
683 const auto& fDb = fD.boundaryField();
684
685 const auto& Sfb = mesh_.Sf().boundaryField();
686 const auto& magSfb = mesh_.magSf().boundaryField();
687 const auto& Cb = mesh_.C().boundaryField();
688
689 for (const label patchi : patchIDs_)
690 {
691 const vectorField Md(Cb[patchi] - origin);
692
693 // Pressure force = surfaceUnitNormal*(surfaceNormal & forceDensity)
694 const vectorField fP
695 (
696 Sfb[patchi]/magSfb[patchi]
697 *(
698 Sfb[patchi] & fDb[patchi]
699 )
700 );
701
702 // Viscous force (total force minus pressure fP)
703 const vectorField fV(magSfb[patchi]*fDb[patchi] - fP);
704
705 addToPatchFields(patchi, Md, fP, fV);
706 }
707 }
708 else
709 {
710 const auto& p = lookupObject<volScalarField>(pName_);
711 const auto& pb = p.boundaryField();
712
713 const auto& Sfb = mesh_.Sf().boundaryField();
714 const auto& Cb = mesh_.C().boundaryField();
715
716 const auto& U = lookupObject<volVectorField>(UName_);
717 tmp<volTensorField> tgradU = fvc::grad(U);
718 const volTensorField& gradU = tgradU();
719 const auto& gradUb = gradU.boundaryField();
720
721 // Scale pRef by density for incompressible simulations
722 const scalar rhoRef = rho(p);
723 const scalar pRef = pRef_/rhoRef;
724
725 for (const label patchi : patchIDs_)
726 {
727 const vectorField Md(Cb[patchi] - origin);
728
729 const vectorField fP(rhoRef*Sfb[patchi]*(pb[patchi] - pRef));
730
731 const vectorField fV
732 (
733 Sfb[patchi] & devRhoReff(gradUb[patchi], patchi)
734 );
735
736 addToPatchFields(patchi, Md, fP, fV);
737 }
738 }
739
740 if (porosity_)
741 {
742 const auto& U = lookupObject<volVectorField>(UName_);
743 const volScalarField rho(this->rho());
744 const volScalarField mu(this->mu());
745
746 const UPtrList<const porosityModel> models
747 (
748 obr_.csorted<porosityModel>()
749 );
750
751 if (models.empty())
752 {
754 << "Porosity effects requested, "
755 << "but no porosity models found in the database"
756 << endl;
757 }
758
759 for (const porosityModel& mdl : models)
760 {
761 // Non-const access required if mesh is changing
762 auto& pm = const_cast<porosityModel&>(mdl);
763
764 const vectorField fPTot(pm.force(U, rho, mu));
765
766 const labelList& cellZoneIDs = pm.cellZoneIDs();
767
768 for (const label zonei : cellZoneIDs)
769 {
770 const cellZone& cZone = mesh_.cellZones()[zonei];
771
772 const vectorField d(mesh_.C(), cZone);
773 const vectorField fP(fPTot, cZone);
774 const vectorField Md(d - origin);
775
776 addToInternalField(cZone, Md, fP);
777 }
778 }
779 }
780
781 reduce(sumPatchForcesP_, sumOp<vector>());
782 reduce(sumPatchForcesV_, sumOp<vector>());
787}
788
793}
794
797{
799}
800
801
803{
804 calcForcesMoments();
805
806 Log << type() << " " << name() << " write:" << nl;
807
808 const auto& coordSys = coordSysPtr_();
809
810 const auto localFp(coordSys.localVector(sumPatchForcesP_));
811 const auto localFv(coordSys.localVector(sumPatchForcesV_));
812 const auto localFi(coordSys.localVector(sumInternalForces_));
813
814 logIntegratedData("forces", localFp, localFv, localFi);
815
816 const auto localMp(coordSys.localVector(sumPatchMomentsP_));
817 const auto localMv(coordSys.localVector(sumPatchMomentsV_));
818 const auto localMi(coordSys.localVector(sumInternalMoments_));
819
820 logIntegratedData("moments", localMp, localMv, localMi);
821
822 setResult("pressureForce", localFp);
823 setResult("viscousForce", localFv);
824 setResult("internalForce", localFi);
825 setResult("pressureMoment", localMp);
826 setResult("viscousMoment", localMv);
827 setResult("internalMoment", localMi);
828
829 return true;
830}
831
832
834{
835 if (writeToFile())
836 {
837 Log << " writing force and moment files." << endl;
838
839 createIntegratedDataFiles();
840 writeIntegratedDataFiles();
841 }
842
843 if (writeFields_)
844 {
845 Log << " writing force and moment fields." << endl;
846
847 force().write();
848 moment().write();
849 }
850
851 Log << endl;
852
853 return true;
854}
855
856
857// ************************************************************************* //
CGAL::indexedCell< K > Cb
#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
compressible::turbulenceModel & turb
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ 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
static MinMax< scalar > ge(const scalar &minVal)
Output to file stream as an OSstream, normally using std::ofstream for the actual output.
Definition OFstream.H:75
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition UPtrList.H:101
bool empty() const noexcept
True if the list is empty (ie, size() is zero).
Definition UPtrListI.H:99
static const word dictName
The dictionary name ("thermophysicalProperties").
A subset of mesh cells.
Definition cellZone.H:61
A Cartesian coordinate system.
Definition cartesianCS.H:68
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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 forces and moments over a given list of patches by integrating pressure and viscous forces a...
Definition forces.H:321
virtual void calcForcesMoments()
Calculate forces and moments.
Definition forces.C:665
void writeIntegratedDataFileHeader(const word &header, OFstream &os) const
Write header for an integrated-data file.
Definition forces.C:410
word pName_
Name of pressure field.
Definition forces.H:397
scalar rhoRef_
Reference density needed for incompressible calculations.
Definition forces.H:387
bool porosity_
Flag to include porosity effects.
Definition forces.H:422
void initialise()
Initialise containers and fields.
Definition forces.C:139
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:349
bool directForceDensity_
Flag to directly supply force density.
Definition forces.H:417
word UName_
Name of velocity field.
Definition forces.H:402
vector sumPatchForcesV_
Sum of patch viscous forces.
Definition forces.H:336
void writeIntegratedDataFiles()
Write integrated data to files.
Definition forces.C:437
vector sumPatchForcesP_
Sum of patch pressure forces.
Definition forces.H:331
vector sumInternalMoments_
Sum of internal moments.
Definition forces.H:356
volVectorField & moment()
Return access to the moment field.
Definition forces.C:111
void createIntegratedDataFiles()
Create the integrated-data files.
Definition forces.C:393
bool initialised_
Flag of initialisation (internal).
Definition forces.H:432
autoPtr< OFstream > forceFilePtr_
File stream for forces.
Definition forces.H:364
vector sumPatchMomentsV_
Sum of patch viscous moments.
Definition forces.H:346
tmp< volScalarField > mu() const
Return dynamic viscosity field.
Definition forces.C:264
vector sumInternalForces_
Sum of internal forces.
Definition forces.H:351
forces(const word &name, const Time &runTime, const dictionary &dict, const bool readFields=true)
Construct from name, Time and dictionary.
Definition forces.C:510
virtual vector forceEff() const
Return the total force.
Definition forces.C:783
void logIntegratedData(const string &descriptor, const vector &pres, const vector &vis, const vector &internal) const
Write integrated data to stream.
Definition forces.C:483
autoPtr< coordinateSystem > coordSysPtr_
Coordinate system used when evaluating forces and moments.
Definition forces.H:377
virtual vector momentEff() const
Return the total moment.
Definition forces.C:789
void addToInternalField(const labelList &cellIDs, const vectorField &Md, const vectorField &f)
Add cell contributions to force and moment fields, and include porosity effects.
Definition forces.C:370
autoPtr< OFstream > momentFilePtr_
File stream for moments.
Definition forces.H:369
vector sumPatchMomentsP_
Sum of patch pressure moments.
Definition forces.H:341
word rhoName_
Name of density field.
Definition forces.H:407
void writeIntegratedDataFile(const vector &pres, const vector &vis, const vector &internal, OFstream &os) const
Write integrated data to a file.
Definition forces.C:460
tmp< symmTensorField > devRhoReff(const tensorField &gradUp, const label patchi) const
Return the effective stress (viscous + turbulent) for patch.
Definition forces.C:212
bool writeFields_
Flag to write force and moment fields.
Definition forces.H:427
tmp< volScalarField > rho() const
Return rho if specified otherwise rhoRef.
Definition forces.C:296
void reset()
Reset containers and fields.
Definition forces.C:181
virtual bool execute()
Execute the function-object operations.
Definition forces.C:795
volVectorField & force()
Return access to the force field.
Definition forces.C:83
virtual bool write()
Write the function-object results.
Definition forces.C:826
scalar pRef_
Reference pressure.
Definition forces.H:392
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:45
word fDName_
Name of force density field.
Definition forces.H:412
labelList patchIDs_
Selected operand patches.
Definition forces.H:382
virtual bool read(const dictionary &)
Read the function-object dictionary.
Definition forces.C:589
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
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition readFields.H:146
const ObjectType & lookupObject(const word &fieldName) const
Lookup and return object (eg, a field) from the (sub) objectRegistry.
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.
void setResult(const word &entryName, const Type &value)
Add result.
const Time & time_
Reference to the time database.
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
void writeValue(Ostream &os, const Type &val) const
Write a given value to stream with the space delimiter.
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
Registry of regIOobjects.
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
Top level model for porosity models.
bool store()
Register object with its registry and transfer ownership to the registry.
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:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
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
U
Definition pEqn.H:72
volScalarField & p
limits reset(1/(limits.max()+VSMALL), 1/(limits.min()+VSMALL))
engineTime & runTime
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
Calculate the gradient of the given field.
auto & name
#define WarningInFunction
Report a warning using Foam::Warning.
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
IncompressibleTurbulenceModel< transportModel > turbulenceModel
Namespace for OpenFOAM.
const dimensionSet dimPressure
const dimensionSet dimViscosity
GeometricField< vector, fvPatchField, volMesh > volVectorField
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.
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
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 POSIX.C:801
dimensionedScalar log(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
const dimensionSet dimForce
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
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.
Vector< scalar > vector
Definition vector.H:57
const dimensionSet dimDensity
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
static bool initialised_(false)
labelList f(nPoints)
volScalarField & nu
dictionary dict
psiReactionThermo & thermo
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299