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-2024 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
39 namespace Foam
40 {
41 namespace functionObjects
42 {
43  defineTypeNameAndDebug(forces, 0);
44  addToRunTimeSelectionTable(functionObject, forces, dictionary);
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  (
98  IOobject
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  (
126  IOobject
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);
181  }
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 == dimensionedVector(force.dimensions(), Zero);
204  moment == dimensionedVector(moment.dimensions(), Zero);
205  }
206  else
207  {
208  constexpr bool updateAccessTime = false;
209  for (const label patchi : patchIDs_)
210  {
211  force.boundaryFieldRef(updateAccessTime)[patchi] = Zero;
212  moment.boundaryFieldRef(updateAccessTime)[patchi] = Zero;
213  }
214  }
215 }
216 
217 
220 (
221  const tensorField& gradUp,
222  const label patchi
223 ) const
224 {
225  typedef incompressible::turbulenceModel icoTurbModel;
226  typedef compressible::turbulenceModel cmpTurbModel;
227 
228  if (foundObject<icoTurbModel>(icoTurbModel::propertiesName))
229  {
230  const auto& turb =
231  lookupObject<icoTurbModel>(icoTurbModel::propertiesName);
232 
233  return -rho(patchi)*turb.nuEff(patchi)*devTwoSymm(gradUp);
234  }
235  else if (foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
236  {
237  const auto& turb =
238  lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
239 
240  return -turb.muEff(patchi)*devTwoSymm(gradUp);
241  }
242  else if (foundObject<fluidThermo>(fluidThermo::dictName))
243  {
244  const auto& thermo = lookupObject<fluidThermo>(fluidThermo::dictName);
245 
246  return -thermo.mu(patchi)*devTwoSymm(gradUp);
247  }
248  else if (foundObject<transportModel>("transportProperties"))
249  {
250  const auto& laminarT =
251  lookupObject<transportModel>("transportProperties");
252 
253  return -rho(patchi)*laminarT.nu(patchi)*devTwoSymm(gradUp);
254  }
255  else if (foundObject<dictionary>("transportProperties"))
256  {
257  const auto& transportProperties =
258  lookupObject<dictionary>("transportProperties");
259 
261 
262  return -rho(patchi)*nu.value()*devTwoSymm(gradUp);
263  }
264  else
265  {
267  << "No valid model for viscous stress calculation"
269 
270  return volSymmTensorField::null();
271  }
272 }
273 
274 
276 {
277  if (foundObject<fluidThermo>(basicThermo::dictName))
278  {
279  const auto& thermo = lookupObject<fluidThermo>(basicThermo::dictName);
280 
281  return thermo.mu();
282  }
283  else if (foundObject<transportModel>("transportProperties"))
284  {
285  const auto& laminarT =
286  lookupObject<transportModel>("transportProperties");
287 
288  return rho()*laminarT.nu();
289  }
290  else if (foundObject<dictionary>("transportProperties"))
291  {
292  const auto& transportProperties =
293  lookupObject<dictionary>("transportProperties");
294 
296 
297  return rho()*nu;
298  }
299  else
300  {
302  << "No valid model for dynamic viscosity calculation"
304 
305  return volScalarField::null();
306  }
307 }
308 
309 
311 {
312  if (rhoName_ == "rhoInf")
313  {
314  return volScalarField::New
315  (
316  "rho",
318  mesh_,
319  dimensionedScalar(dimDensity, rhoRef_)
320  );
321  }
322 
323  return (lookupObject<volScalarField>(rhoName_));
324 }
325 
326 
328 Foam::functionObjects::forces::rho(const label patchi) const
329 {
330  if (rhoName_ == "rhoInf")
331  {
332  return tmp<scalarField>::New
333  (
334  mesh_.boundary()[patchi].size(),
335  rhoRef_
336  );
337  }
338 
339  const auto& rho = lookupObject<volScalarField>(rhoName_);
340  return rho.boundaryField()[patchi];
341 }
342 
343 
344 Foam::scalar Foam::functionObjects::forces::rho(const volScalarField& p) const
345 {
346  if (p.dimensions() == dimPressure)
347  {
348  return 1;
349  }
350 
351  if (rhoName_ != "rhoInf")
352  {
354  << "Dynamic pressure is expected but kinematic is provided."
355  << exit(FatalError);
356  }
357 
358  return rhoRef_;
359 }
360 
361 
363 (
364  const label patchi,
365  const vectorField& Md,
366  const vectorField& fP,
367  const vectorField& fV
368 )
369 {
370  constexpr bool updateAccessTime = false;
371 
372  sumPatchForcesP_ += sum(fP);
373  sumPatchForcesV_ += sum(fV);
374  force().boundaryFieldRef(updateAccessTime)[patchi] += fP + fV;
375 
376  const vectorField mP(Md^fP);
377  const vectorField mV(Md^fV);
378 
379  sumPatchMomentsP_ += sum(mP);
380  sumPatchMomentsV_ += sum(mV);
381  moment().boundaryFieldRef(updateAccessTime)[patchi] += mP + mV;
382 }
383 
384 
386 (
387  const labelList& cellIDs,
388  const vectorField& Md,
389  const vectorField& f
390 )
391 {
392  auto& force = this->force();
393  auto& moment = this->moment();
394 
395  forAll(cellIDs, i)
396  {
397  const label celli = cellIDs[i];
398 
399  sumInternalForces_ += f[i];
400  force[celli] += f[i];
401 
402  const vector m(Md[i]^f[i]);
403  sumInternalMoments_ += m;
404  moment[celli] = m;
405  }
406 }
407 
408 
410 {
411  if (!forceFilePtr_)
412  {
413  forceFilePtr_ = newFileAtStartTime("force");
414  writeIntegratedDataFileHeader("Force", forceFilePtr_());
415  }
416 
417  if (!momentFilePtr_)
418  {
419  momentFilePtr_ = newFileAtStartTime("moment");
420  writeIntegratedDataFileHeader("Moment", momentFilePtr_());
421  }
422 }
423 
424 
426 (
427  const word& header,
428  OFstream& os
429 ) const
430 {
431  const auto& coordSys = coordSysPtr_();
432  const auto vecDesc = [](const word& root)->string
433  {
434  return root + "_x " + root + "_y " + root + "_z";
435  };
436  writeHeader(os, header);
437  writeHeaderValue(os, "CofR", coordSys.origin());
438  writeHeader(os, "");
439  writeCommented(os, "Time");
440  writeTabbed(os, vecDesc("total"));
441  writeTabbed(os, vecDesc("pressure"));
442  writeTabbed(os, vecDesc("viscous"));
443 
444  if (porosity_)
445  {
446  writeTabbed(os, vecDesc("porous"));
447  }
448 
449  os << endl;
450 }
451 
452 
454 {
455  const auto& coordSys = coordSysPtr_();
456 
457  writeIntegratedDataFile
458  (
459  coordSys.localVector(sumPatchForcesP_),
460  coordSys.localVector(sumPatchForcesV_),
461  coordSys.localVector(sumInternalForces_),
462  forceFilePtr_()
463  );
464 
465  writeIntegratedDataFile
466  (
467  coordSys.localVector(sumPatchMomentsP_),
468  coordSys.localVector(sumPatchMomentsV_),
469  coordSys.localVector(sumInternalMoments_),
470  momentFilePtr_()
471  );
472 }
473 
474 
476 (
477  const vector& pres,
478  const vector& vis,
479  const vector& internal,
480  OFstream& os
481 ) const
482 {
483  writeCurrentTime(os);
484 
485  writeValue(os, pres + vis + internal);
486  writeValue(os, pres);
487  writeValue(os, vis);
488 
489  if (porosity_)
490  {
491  writeValue(os, internal);
492  }
493 
494  os << endl;
495 }
496 
497 
499 (
500  const string& descriptor,
501  const vector& pres,
502  const vector& vis,
503  const vector& internal
504 ) const
505 {
506  if (!log)
507  {
508  return;
509  }
510 
511  Log << " Sum of " << descriptor.c_str() << nl
512  << " Total : " << (pres + vis + internal) << nl
513  << " Pressure : " << pres << nl
514  << " Viscous : " << vis << nl;
515 
516  if (porosity_)
517  {
518  Log << " Porous : " << internal << nl;
519  }
520 }
521 
522 
523 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
524 
526 (
527  const word& name,
528  const Time& runTime,
529  const dictionary& dict,
530  bool readFields
531 )
532 :
533  fvMeshFunctionObject(name, runTime, dict),
534  writeFile(mesh_, name),
535  sumPatchForcesP_(Zero),
536  sumPatchForcesV_(Zero),
537  sumPatchMomentsP_(Zero),
538  sumPatchMomentsV_(Zero),
539  sumInternalForces_(Zero),
540  sumInternalMoments_(Zero),
541  forceFilePtr_(),
542  momentFilePtr_(),
543  coordSysPtr_(nullptr),
544  rhoRef_(VGREAT),
545  pRef_(0),
546  pName_("p"),
547  UName_("U"),
548  rhoName_("rho"),
549  fDName_("fD"),
550  directForceDensity_(false),
551  porosity_(false),
552  writeFields_(false),
553  initialised_(false)
554 {
555  if (readFields)
556  {
557  read(dict);
559  Log << endl;
560  }
561 }
562 
563 
565 (
566  const word& name,
567  const objectRegistry& obr,
568  const dictionary& dict,
569  bool readFields
570 )
571 :
572  fvMeshFunctionObject(name, obr, dict),
573  writeFile(mesh_, name),
574  sumPatchForcesP_(Zero),
575  sumPatchForcesV_(Zero),
576  sumPatchMomentsP_(Zero),
577  sumPatchMomentsV_(Zero),
578  sumInternalForces_(Zero),
579  sumInternalMoments_(Zero),
580  forceFilePtr_(),
581  momentFilePtr_(),
582  coordSysPtr_(nullptr),
583  rhoRef_(VGREAT),
584  pRef_(0),
585  pName_("p"),
586  UName_("U"),
587  rhoName_("rho"),
588  fDName_("fD"),
589  directForceDensity_(false),
590  porosity_(false),
591  writeFields_(false),
592  initialised_(false)
593 {
594  if (readFields)
595  {
596  read(dict);
598  Log << endl;
599  }
600 }
601 
602 
603 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
604 
606 {
607  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
608 
610  {
611  return false;
612  }
613 
614  initialised_ = false;
615 
616  Info<< type() << ' ' << name() << ':' << endl;
617 
618  // Can also use pbm.indices(), but no warnings...
619  patchIDs_ = pbm.patchSet(dict.get<wordRes>("patches")).sortedToc();
620 
621  dict.readIfPresent("directForceDensity", directForceDensity_);
622  if (directForceDensity_)
623  {
624  // Optional name entry for fD
625  if (dict.readIfPresent<word>("fD", fDName_))
626  {
627  Info<< " fD: " << fDName_ << endl;
628  }
629  }
630  else
631  {
632  // Optional field name entries
633  if (dict.readIfPresent<word>("p", pName_))
634  {
635  Info<< " p: " << pName_ << endl;
636  }
637  if (dict.readIfPresent<word>("U", UName_))
638  {
639  Info<< " U: " << UName_ << endl;
640  }
641  if (dict.readIfPresent<word>("rho", rhoName_))
642  {
643  Info<< " rho: " << rhoName_ << endl;
644  }
645 
646  // Reference density needed for incompressible calculations
647  if (rhoName_ == "rhoInf")
648  {
649  rhoRef_ = dict.getCheck<scalar>("rhoInf", scalarMinMax::ge(SMALL));
650  Info<< " Freestream density (rhoInf) set to " << rhoRef_ << endl;
651  }
652 
653  // Reference pressure, 0 by default
654  if (dict.readIfPresent<scalar>("pRef", pRef_))
655  {
656  Info<< " Reference pressure (pRef) set to " << pRef_ << endl;
657  }
658  }
659 
660  dict.readIfPresent("porosity", porosity_);
661  if (porosity_)
662  {
663  Info<< " Including porosity effects" << endl;
664  }
665  else
666  {
667  Info<< " Not including porosity effects" << endl;
668  }
669 
670  writeFields_ = dict.getOrDefault("writeFields", false);
671  if (writeFields_)
672  {
673  Info<< " Fields will be written" << endl;
674  }
675 
676 
677  return true;
678 }
679 
680 
682 {
683  initialise();
684 
685  reset();
686 
687  const point& origin = coordSysPtr_->origin();
688 
689  if (directForceDensity_)
690  {
691  const auto& fD = lookupObject<volVectorField>(fDName_);
692  const auto& fDb = fD.boundaryField();
693 
694  const auto& Sfb = mesh_.Sf().boundaryField();
695  const auto& magSfb = mesh_.magSf().boundaryField();
696  const auto& Cb = mesh_.C().boundaryField();
697 
698  for (const label patchi : patchIDs_)
699  {
700  const vectorField Md(Cb[patchi] - origin);
701 
702  // Pressure force = surfaceUnitNormal*(surfaceNormal & forceDensity)
703  const vectorField fP
704  (
705  Sfb[patchi]/magSfb[patchi]
706  *(
707  Sfb[patchi] & fDb[patchi]
708  )
709  );
710 
711  // Viscous force (total force minus pressure fP)
712  const vectorField fV(magSfb[patchi]*fDb[patchi] - fP);
713 
714  addToPatchFields(patchi, Md, fP, fV);
715  }
716  }
717  else
718  {
719  const auto& p = lookupObject<volScalarField>(pName_);
720  const auto& pb = p.boundaryField();
721 
722  const auto& Sfb = mesh_.Sf().boundaryField();
723  const auto& Cb = mesh_.C().boundaryField();
724 
725  const auto& U = lookupObject<volVectorField>(UName_);
726  tmp<volTensorField> tgradU = fvc::grad(U);
727  const volTensorField& gradU = tgradU();
728  const auto& gradUb = gradU.boundaryField();
729 
730  // Scale pRef by density for incompressible simulations
731  const scalar rhoRef = rho(p);
732  const scalar pRef = pRef_/rhoRef;
733 
734  for (const label patchi : patchIDs_)
735  {
736  const vectorField Md(Cb[patchi] - origin);
737 
738  const vectorField fP(rhoRef*Sfb[patchi]*(pb[patchi] - pRef));
739 
740  const vectorField fV
741  (
742  Sfb[patchi] & devRhoReff(gradUb[patchi], patchi)
743  );
744 
745  addToPatchFields(patchi, Md, fP, fV);
746  }
747  }
748 
749  if (porosity_)
750  {
751  const auto& U = lookupObject<volVectorField>(UName_);
752  const volScalarField rho(this->rho());
753  const volScalarField mu(this->mu());
754 
755  const UPtrList<const porosityModel> models
756  (
757  obr_.csorted<porosityModel>()
758  );
759 
760  if (models.empty())
761  {
763  << "Porosity effects requested, "
764  << "but no porosity models found in the database"
765  << endl;
766  }
767 
768  for (const porosityModel& mdl : models)
769  {
770  // Non-const access required if mesh is changing
771  auto& pm = const_cast<porosityModel&>(mdl);
772 
773  const vectorField fPTot(pm.force(U, rho, mu));
774 
775  const labelList& cellZoneIDs = pm.cellZoneIDs();
776 
777  for (const label zonei : cellZoneIDs)
778  {
779  const cellZone& cZone = mesh_.cellZones()[zonei];
780 
781  const vectorField d(mesh_.C(), cZone);
782  const vectorField fP(fPTot, cZone);
783  const vectorField Md(d - origin);
784 
785  addToInternalField(cZone, Md, fP);
786  }
787  }
788  }
789 
790  reduce(sumPatchForcesP_, sumOp<vector>());
791  reduce(sumPatchForcesV_, sumOp<vector>());
792  reduce(sumPatchMomentsP_, sumOp<vector>());
793  reduce(sumPatchMomentsV_, sumOp<vector>());
794  reduce(sumInternalForces_, sumOp<vector>());
795  reduce(sumInternalMoments_, sumOp<vector>());
796 }
797 
800 {
801  return sumPatchForcesP_ + sumPatchForcesV_ + sumInternalForces_;
802 }
803 
806 {
807  return sumPatchMomentsP_ + sumPatchMomentsV_ + sumInternalMoments_;
808 }
809 
810 
812 {
813  calcForcesMoments();
814 
815  Log << type() << " " << name() << " write:" << nl;
816 
817  const auto& coordSys = coordSysPtr_();
818 
819  const auto localFp(coordSys.localVector(sumPatchForcesP_));
820  const auto localFv(coordSys.localVector(sumPatchForcesV_));
821  const auto localFi(coordSys.localVector(sumInternalForces_));
822 
823  logIntegratedData("forces", localFp, localFv, localFi);
824 
825  const auto localMp(coordSys.localVector(sumPatchMomentsP_));
826  const auto localMv(coordSys.localVector(sumPatchMomentsV_));
827  const auto localMi(coordSys.localVector(sumInternalMoments_));
828 
829  logIntegratedData("moments", localMp, localMv, localMi);
830 
831  setResult("pressureForce", localFp);
832  setResult("viscousForce", localFv);
833  setResult("internalForce", localFi);
834  setResult("pressureMoment", localMp);
835  setResult("viscousMoment", localMv);
836  setResult("internalMoment", localMi);
837 
838  return true;
839 }
840 
841 
843 {
844  if (writeToFile())
845  {
846  Log << " writing force and moment files." << endl;
847 
848  createIntegratedDataFiles();
849  writeIntegratedDataFiles();
850  }
851 
852  if (writeFields_)
853  {
854  Log << " writing force and moment fields." << endl;
855 
856  force().write();
857  moment().write();
858  }
859 
860  Log << endl;
861 
862  return true;
863 }
864 
865 
866 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
virtual vector momentEff() const
Return the total moment.
Definition: forces.C:798
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:53
const polyBoundaryMesh & pbm
dictionary dict
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
defineTypeNameAndDebug(ObukhovLength, 0)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
static bool initialised_(false)
static void writeHeader(Ostream &os, const word &fieldName)
dimensionedScalar log(const dimensionedScalar &ds)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:88
Output to file stream as an OSstream, normally using std::ofstream for the actual output...
Definition: OFstream.H:71
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
const dimensionSet dimViscosity
compressible::turbulenceModel & turb
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
engineTime & runTime
ThermalDiffusivity< CompressibleTurbulenceModel< fluidThermo > > turbulenceModel
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
volVectorField & force()
Return access to the force field.
Definition: forces.C:83
An indexed form of CGAL::Triangulation_cell_base_3<K> used to keep track of the Delaunay cells (tets)...
Definition: indexedCell.H:62
bool store()
Register object with its registry and transfer ownership to the registry.
Definition: regIOobjectI.H:36
virtual vector forceEff() const
Return the total force.
Definition: forces.C:792
IncompressibleTurbulenceModel< transportModel > turbulenceModel
Ignore writing from objectRegistry::writeObject()
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
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:356
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
virtual bool read(const dictionary &dict)
Read the dictionary.
Definition: forces.C:598
void writeIntegratedDataFiles()
Write integrated data to files.
Definition: forces.C:446
Macros for easy insertion into run-time selection tables.
tmp< symmTensorField > devRhoReff(const tensorField &gradUp, const label patchi) const
Return the effective stress (viscous + turbulent) for patch.
Definition: forces.C:213
void writeIntegratedDataFile(const vector &pres, const vector &vis, const vector &internal, OFstream &os) const
Write integrated data to a file.
Definition: forces.C:469
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
A Cartesian coordinate system.
Definition: cartesianCS.H:65
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
psiReactionThermo & thermo
Definition: createFields.H:28
virtual bool execute()
Execute the function object.
Definition: forces.C:804
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
virtual void calcForcesMoments()
Calculate forces and moments.
Definition: forces.C:674
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
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.
Calculate the gradient of the given field.
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< volScalarField > mu() const
Return dynamic viscosity field.
Definition: forces.C:268
static MinMax< scalar > ge(const scalar &minVal)
A semi-infinite range from minVal to the type max.
Definition: MinMaxI.H:24
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
static const GeometricField< symmTensor, fvPatchField, volMesh > & null() noexcept
Return a null GeometricField (reference to a nullObject).
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Vector< scalar > vector
Definition: vector.H:57
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:379
tmp< volScalarField > rho() const
Return rho if specified otherwise rhoRef.
Definition: forces.C:303
const dimensionSet dimPressure
void createIntegratedDataFiles()
Create the integrated-data files.
Definition: forces.C:402
volVectorField & moment()
Return access to the moment field.
Definition: forces.C:111
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
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())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
const dimensionSet dimForce
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
void reset()
Reset containers and fields.
Definition: forces.C:181
void initialise()
Initialise containers and fields.
Definition: forces.C:139
OBJstream os(runTime.globalPath()/outputName)
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:241
labelList f(nPoints)
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
const dimensionSet dimDensity
const dimensionedScalar mu
Atomic mass unit.
Field< tensor > tensorField
Specialisation of Field<T> for tensor.
autoPtr< coordinateSystem > coordSysPtr_
Coordinate system used when evaluating forces and moments.
Definition: forces.H:377
U
Definition: pEqn.H:72
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
forces(const word &name, const Time &runTime, const dictionary &dict, const bool readFields=true)
Construct from Time and dictionary.
Definition: forces.C:519
const objectRegistry & obr_
Reference to the region objectRegistry.
Nothing to be read.
#define Log
Definition: PDRblock.C:28
messageStream Info
Information stream (stdout output on master, null elsewhere)
virtual bool read(const dictionary &dict)
Read optional controls.
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
void writeIntegratedDataFileHeader(const word &header, OFstream &os) const
Write header for an integrated-data file.
Definition: forces.C:419
Registry of regIOobjects.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
void logIntegratedData(const string &descriptor, const vector &pres, const vector &vis, const vector &internal) const
Write integrated data to stream.
Definition: forces.C:492
Request registration (bool: true)
labelList cellIDs
volScalarField & nu
Do not request registration (bool: false)
void reduce(T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) using linear/tree communication schedule.
virtual bool write()
Write to data files/fields and to streams.
Definition: forces.C:835
Namespace for OpenFOAM.
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:491
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127