forcesBase.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "forcesBase.H"
27 #include "fvcGrad.H"
28 #include "porosityModel.H"
33 #include "fluidThermo.H"
34 #include "surfaceInterpolate.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace functionObjects
42 {
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49 
51 (
52  const dictionary& dict
53 ) const
54 {
56 
57  const word forceType(dict.lookup("type"));
58 
59  // Name for file(fileID::mainFile=0)
60  names.append(forceType);
61 
62  if (dict.found("binData"))
63  {
64  const dictionary& binDict(dict.subDict("binData"));
65  const label nb = binDict.lookup<label>("nBin");
66  if (nb > 0)
67  {
68  // Name for file(fileID::binsFile=1)
69  names.append(forceType + "_bins");
70  }
71  }
72 
73  return names;
74 }
75 
76 
78 {
79  const word forceTypes
80  (
81  porosity_
82  ? "(pressure viscous porous)"
83  : "(pressure viscous)"
84  );
85 
86  switch (fileID(i))
87  {
88  case fileID::mainFile:
89  {
90  // force data
91 
92  writeHeader(file(i), "Forces");
93  writeCoRValueHeader(file(i));
94  writeCommented(file(i), "Time");
95 
96  writeCoRHeader(file(i));
97 
98  file(i)
99  << "forces" << forceTypes << tab
100  << "moments" << forceTypes;
101 
102  break;
103  }
104  case fileID::binsFile:
105  {
106  // bin data
107 
108  writeHeader(file(i), "Force bins");
109  writeHeaderValue(file(i), "bins", nBin_);
110  writeHeaderValue(file(i), "start", binMin_);
111  writeHeaderValue(file(i), "delta", binDx_);
112  writeHeaderValue(file(i), "direction", binDir_);
113 
114  vectorField binPoints(nBin_);
115  writeCommented(file(i), "x co-ords :");
116  forAll(binPoints, pointi)
117  {
118  binPoints[pointi] = (binMin_ + (pointi + 1)*binDx_)*binDir_;
119  file(i) << tab << binPoints[pointi].x();
120  }
121  file(i) << nl;
122 
123  writeCommented(file(i), "y co-ords :");
124  forAll(binPoints, pointi)
125  {
126  file(i) << tab << binPoints[pointi].y();
127  }
128  file(i) << nl;
129 
130  writeCommented(file(i), "z co-ords :");
131  forAll(binPoints, pointi)
132  {
133  file(i) << tab << binPoints[pointi].z();
134  }
135  file(i) << nl;
136 
137  writeCommented(file(i), "Time");
138 
139  for (label j = 0; j < nBin_; j++)
140  {
141  const word jn('(' + Foam::name(j) + ')');
142  const word f("forces" + jn + forceTypes);
143  const word m("moments" + jn + forceTypes);
144 
145  file(i)<< tab << f << tab << m;
146  }
147 
148  break;
149  }
150  default:
151  {
153  << "Unhandled file index: " << i
154  << abort(FatalError);
155  }
156  }
157 
158  file(i)<< endl;
159 }
160 
161 
163 {
164  if (initialised_)
165  {
166  return;
167  }
168 
169  if (directForceDensity_)
170  {
171  if (!obr_.foundObject<volVectorField>(fDName_))
172  {
174  << "Could not find " << fDName_ << " in database."
175  << exit(FatalError);
176  }
177  }
178  else
179  {
180  if
181  (
182  !obr_.foundObject<volVectorField>(UName_)
183  || !obr_.foundObject<volScalarField>(pName_)
184 
185  )
186  {
188  << "Could not find " << UName_ << ", " << pName_
189  << exit(FatalError);
190  }
191 
192  if
193  (
194  rhoName_ != "rhoInf"
195  && !obr_.foundObject<volScalarField>(rhoName_)
196  )
197  {
199  << "Could not find " << rhoName_
200  << exit(FatalError);
201  }
202  }
203 
204  initialised_ = true;
205 }
206 
207 
210 {
212  typedef compressible::momentumTransportModel cmpModel;
213  typedef phaseIncompressible::momentumTransportModel phaseIcoModel;
214  typedef phaseCompressible::momentumTransportModel phaseCmpModel;
215 
216  const word& modelName = momentumTransportModel::typeName;
217  const word phaseModelName =
218  phaseName_ == word::null
219  ? word::null
220  : IOobject::groupName(momentumTransportModel::typeName, phaseName_);
221 
222  if (obr_.foundObject<icoModel>(modelName))
223  {
225  obr_.lookupObject<icoModel>(modelName);
226 
227  return fvc::interpolate(alpha()*rho())*model.devSigma();
228  }
229  else if (obr_.foundObject<cmpModel>(modelName))
230  {
231  const cmpModel& model =
232  obr_.lookupObject<cmpModel>(modelName);
233 
234  return fvc::interpolate(alpha())*model.devTau();
235  }
236  else if (obr_.foundObject<phaseIcoModel>(phaseModelName))
237  {
238  const phaseIcoModel& model =
239  obr_.lookupObject<phaseIcoModel>(phaseModelName);
240 
241  return fvc::interpolate(rho())*model.devSigma();
242  }
243  else if (obr_.foundObject<phaseCmpModel>(phaseModelName))
244  {
245  const phaseCmpModel& model =
246  obr_.lookupObject<phaseCmpModel>(phaseModelName);
247 
248  return model.devTau();
249  }
250  else
251  {
253  << "No valid model for viscous stress calculation"
254  << exit(FatalError);
255 
256  return surfaceVectorField::null();
257  }
258 }
259 
260 
262 {
264  typedef compressible::momentumTransportModel cmpModel;
265  typedef phaseIncompressible::momentumTransportModel phaseIcoModel;
266  typedef phaseCompressible::momentumTransportModel phaseCmpModel;
267 
268  const word& modelName = momentumTransportModel::typeName;
269  const word phaseModelName =
270  phaseName_ == word::null
271  ? word::null
272  : IOobject::groupName(momentumTransportModel::typeName, phaseName_);
273 
274  if (obr_.foundObject<icoModel>(modelName))
275  {
277  obr_.lookupObject<icoModel>(modelName);
278 
279  return rho()*model.nu();
280  }
281  else if (obr_.foundObject<cmpModel>(modelName))
282  {
283  const cmpModel& model =
284  obr_.lookupObject<cmpModel>(modelName);
285 
286  return model.rho()*model.nu();
287  }
288  else if (obr_.foundObject<phaseIcoModel>(phaseModelName))
289  {
290  const phaseIcoModel& model =
291  obr_.lookupObject<phaseIcoModel>(phaseModelName);
292 
293  return rho()*model.nu();
294  }
295  else if (obr_.foundObject<phaseCmpModel>(phaseModelName))
296  {
297  const phaseCmpModel& model =
298  obr_.lookupObject<phaseCmpModel>(phaseModelName);
299 
300  return model.rho()*model.nu();
301  }
302  else if (obr_.foundObject<dictionary>("physicalProperties"))
303  {
304  // Legacy support for icoFoam
305 
307  obr_.lookupObject<dictionary>("physicalProperties");
308 
309  const dimensionedScalar nu
310  (
311  "nu",
314  );
315 
316  return rho()*nu;
317  }
318  else
319  {
321  << "No valid model for dynamic viscosity calculation"
322  << exit(FatalError);
323 
324  return volScalarField::null();
325  }
326 }
327 
328 
330 {
331  if (rhoName_ == "rhoInf")
332  {
333  return volScalarField::New
334  (
335  "rho",
336  mesh_,
337  dimensionedScalar(dimDensity, rhoRef_)
338  );
339  }
340  else
341  {
342  return(obr_.lookupObject<volScalarField>(rhoName_));
343  }
344 }
345 
346 
348 (
349  const volScalarField& p
350 ) const
351 {
352  if (p.dimensions() == dimPressure)
353  {
354  return 1.0;
355  }
356  else
357  {
358  if (rhoName_ != "rhoInf")
359  {
361  << "Dynamic pressure is expected but kinematic is provided."
362  << exit(FatalError);
363  }
364 
365  return rhoRef_;
366  }
367 }
368 
369 
371 {
372  if (phaseName_ == word::null)
373  {
374  return volScalarField::New
375  (
376  "alpha",
377  mesh_,
379  );
380  }
381  else
382  {
383  return obr_.lookupObject<volScalarField>
384  (
385  IOobject::groupName("alpha", phaseName_)
386  );
387  }
388 }
389 
390 
392 (
393  const label patchi
394 ) const
395 {
396  if (phaseName_ == word::null)
397  {
398  return tmp<scalarField>
399  (
400  new scalarField(mesh_.boundary()[patchi].size(), 1)
401  );
402  }
403  else
404  {
405  return obr_.lookupObject<volScalarField>
406  (
407  IOobject::groupName("alpha", phaseName_)
408  ).boundaryField()[patchi];
409  }
410 }
411 
412 
414 (
415  const vectorField& Md,
416  const vectorField& fN,
417  const vectorField& fT,
418  const vectorField& fP,
419  const vectorField& d
420 )
421 {
422  if (nBin_ == 1)
423  {
424  force_[0][0] += sum(fN);
425  force_[1][0] += sum(fT);
426  force_[2][0] += sum(fP);
427  moment_[0][0] += sum(Md^fN);
428  moment_[1][0] += sum(Md^fT);
429  moment_[2][0] += sum(Md^fP);
430  }
431  else
432  {
433  scalarField dd((d & binDir_) - binMin_);
434 
435  forAll(dd, i)
436  {
437  label bini = min(max(floor(dd[i]/binDx_), 0), force_[0].size() - 1);
438 
439  force_[0][bini] += fN[i];
440  force_[1][bini] += fT[i];
441  force_[2][bini] += fP[i];
442  moment_[0][bini] += Md[i]^fN[i];
443  moment_[1][bini] += Md[i]^fT[i];
444  moment_[2][bini] += Md[i]^fP[i];
445  }
446  }
447 }
448 
449 
451 {}
452 
453 
455 {
456  file << "CofR" << tab;
457 }
458 
459 
461 {
462  file << CofR();
463 }
464 
465 
467 {
468  Log << type() << " " << name() << " write:" << nl
469  << " sum of forces:" << nl
470  << " pressure : " << sum(force_[0]) << nl
471  << " viscous : " << sum(force_[1]) << nl
472  << " porous : " << sum(force_[2]) << nl
473  << " sum of moments:" << nl
474  << " pressure : " << sum(moment_[0]) << nl
475  << " viscous : " << sum(moment_[1]) << nl
476  << " porous : " << sum(moment_[2])
477  << endl;
478 
479  writeTime(file(fileID::mainFile));
480  writeCofR(file(fileID::mainFile));
481 
482  if (porosity_)
483  {
484  file(fileID::mainFile) << tab << setw(1) << '('
485  << sum(force_[0]) << setw(1) << ' '
486  << sum(force_[1]) << setw(1) << ' '
487  << sum(force_[2]) << setw(3) << ") ("
488  << sum(moment_[0]) << setw(1) << ' '
489  << sum(moment_[1]) << setw(1) << ' '
490  << sum(moment_[2]) << setw(1) << ')';
491  }
492  else
493  {
494  file(fileID::mainFile) << tab << setw(1) << '('
495  << sum(force_[0]) << setw(1) << ' '
496  << sum(force_[1]) << setw(3) << ") ("
497  << sum(moment_[0]) << setw(1) << ' '
498  << sum(moment_[1]) << setw(1) << ')';
499  }
500 
501  file(fileID::mainFile) << endl;
502 }
503 
504 
506 {
507  if (nBin_ == 1)
508  {
509  return;
510  }
511 
512  List<Field<vector>> f(force_);
513  List<Field<vector>> m(moment_);
514 
515  if (binCumulative_)
516  {
517  for (label i = 1; i < f[0].size(); i++)
518  {
519  f[0][i] += f[0][i-1];
520  f[1][i] += f[1][i-1];
521  f[2][i] += f[2][i-1];
522 
523  m[0][i] += m[0][i-1];
524  m[1][i] += m[1][i-1];
525  m[2][i] += m[2][i-1];
526  }
527  }
528 
529  writeTime(file(fileID::binsFile));
530 
531 
532  forAll(f[0], i)
533  {
534  if (porosity_)
535  {
536  file(fileID::binsFile)
537  << tab << setw(1) << '('
538  << f[0][i] << setw(1) << ' '
539  << f[1][i] << setw(1) << ' '
540  << f[2][i] << setw(3) << ") ("
541  << m[0][i] << setw(1) << ' '
542  << m[1][i] << setw(1) << ' '
543  << m[2][i] << setw(1) << ')';
544  }
545  else
546  {
547  file(fileID::binsFile)
548  << tab << setw(1) << '('
549  << f[0][i] << setw(1) << ' '
550  << f[1][i] << setw(3) << ") ("
551  << m[0][i] << setw(1) << ' '
552  << m[1][i] << setw(1) << ')';
553  }
554  }
555 
556  file(fileID::binsFile) << endl;
557 }
558 
559 
560 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
561 
563 (
564  const word& name,
565  const Time& runTime,
566  const dictionary& dict
567 )
568 :
569  fvMeshFunctionObject(name, runTime, dict),
570  logFiles(obr_, name),
571  force_(3),
572  moment_(3),
573  patchSet_(),
574  pName_(word::null),
575  UName_(word::null),
576  rhoName_(word::null),
577  phaseName_(word::null),
578  directForceDensity_(false),
579  fDName_(""),
580  rhoRef_(vGreat),
581  pRef_(0),
582  porosity_(false),
583  nBin_(1),
584  binDir_(Zero),
585  binDx_(0.0),
586  binMin_(great),
587  binPoints_(),
588  binCumulative_(true),
589  initialised_(false)
590 {
591  read(dict);
592 }
593 
594 
596 (
597  const word& name,
598  const objectRegistry& obr,
599  const dictionary& dict
600 )
601 :
603  logFiles(obr_, name),
604  force_(3),
605  moment_(3),
606  patchSet_(),
607  pName_(word::null),
608  UName_(word::null),
609  rhoName_(word::null),
610  phaseName_(word::null),
611  directForceDensity_(false),
612  fDName_(""),
613  rhoRef_(vGreat),
614  pRef_(0),
615  porosity_(false),
616  nBin_(1),
617  binDir_(Zero),
618  binDx_(0.0),
619  binMin_(great),
620  binPoints_(),
621  binCumulative_(true),
622  initialised_(false)
623 {
624  read(dict);
625 }
626 
627 
628 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
629 
631 {}
632 
633 
634 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
635 
637 {
639 
640  initialised_ = false;
641 
642  Log << type() << " " << name() << ":" << nl;
643 
644  directForceDensity_ = dict.lookupOrDefault("directForceDensity", false);
645 
646  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
647 
648  patchSet_ = mesh_.boundaryMesh().patchSet(dict);
649 
650  if (directForceDensity_)
651  {
652  // Optional entry for fDName
653  fDName_ = dict.lookupOrDefault<word>("fD", "fD");
654  }
655  else
656  {
657  // Optional phase entry
658  phaseName_ = dict.lookupOrDefault<word>("phase", word::null);
659 
660  // Optional U, p and rho entries
661  pName_ =
662  dict.lookupOrDefault<word>
663  (
664  "p",
665  IOobject::groupName("p", phaseName_)
666  );
667  UName_ =
668  dict.lookupOrDefault<word>
669  (
670  "U",
671  IOobject::groupName("U", phaseName_)
672  );
673  rhoName_ =
674  dict.lookupOrDefault<word>
675  (
676  "rho",
677  IOobject::groupName("rho", phaseName_)
678  );
679 
680  // Reference density needed for incompressible calculations
681  if (rhoName_ == "rhoInf")
682  {
683  dict.lookup("rhoInf") >> rhoRef_;
684  }
685 
686  // Reference pressure, 0 by default
687  pRef_ = dict.lookupOrDefault<scalar>("pRef", 0.0);
688  }
689 
690  dict.readIfPresent("porosity", porosity_);
691  if (porosity_)
692  {
693  Log << " Including porosity effects" << endl;
694  }
695  else
696  {
697  Log << " Not including porosity effects" << endl;
698  }
699 
700  if (dict.found("binData"))
701  {
702  const dictionary& binDict(dict.subDict("binData"));
703  binDict.lookup("nBin") >> nBin_;
704 
705  if (nBin_ < 0)
706  {
708  << "Number of bins (nBin) must be zero or greater"
709  << exit(FatalIOError);
710  }
711  else if ((nBin_ == 0) || (nBin_ == 1))
712  {
713  nBin_ = 1;
714  forAll(force_, i)
715  {
716  force_[i].setSize(1);
717  moment_[i].setSize(1);
718  }
719  }
720 
721  if (nBin_ > 1)
722  {
723  binDict.lookup("direction") >> binDir_;
724  binDir_ /= mag(binDir_);
725 
726  binMin_ = great;
727  scalar binMax = -great;
728  forAllConstIter(labelHashSet, patchSet_, iter)
729  {
730  const label patchi = iter.key();
731  const polyPatch& pp = pbm[patchi];
732  const scalarField d(pp.faceCentres() & binDir_);
733  binMin_ = min(min(d), binMin_);
734  binMax = max(max(d), binMax);
735  }
736  reduce(binMin_, minOp<scalar>());
737  reduce(binMax, maxOp<scalar>());
738 
739  // slightly boost binMax so that region of interest is fully
740  // within bounds
741  binMax = 1.0001*(binMax - binMin_) + binMin_;
742 
743  binDx_ = (binMax - binMin_)/scalar(nBin_);
744 
745  // create the bin points used for writing
746  binPoints_.setSize(nBin_);
747  forAll(binPoints_, i)
748  {
749  binPoints_[i] = (i + 0.5)*binDir_*binDx_;
750  }
751 
752  binDict.lookup("cumulative") >> binCumulative_;
753 
754  // allocate storage for forces and moments
755  forAll(force_, i)
756  {
757  force_[i].setSize(nBin_);
758  moment_[i].setSize(nBin_);
759  }
760  }
761  }
762 
763  if (nBin_ == 1)
764  {
765  // allocate storage for forces and moments
766  force_[0].setSize(1);
767  force_[1].setSize(1);
768  force_[2].setSize(1);
769  moment_[0].setSize(1);
770  moment_[1].setSize(1);
771  moment_[2].setSize(1);
772  }
773 
774  resetNames(createFileNames(dict));
775 
776  return true;
777 }
778 
779 
781 {
782  initialise();
783 
784  force_[0] = Zero;
785  force_[1] = Zero;
786  force_[2] = Zero;
787 
788  moment_[0] = Zero;
789  moment_[1] = Zero;
790  moment_[2] = Zero;
791 
792  if (directForceDensity_)
793  {
794  const volVectorField& fD = obr_.lookupObject<volVectorField>(fDName_);
795 
796  const surfaceVectorField::Boundary& Sfb =
797  mesh_.Sf().boundaryField();
798 
799  forAllConstIter(labelHashSet, patchSet_, iter)
800  {
801  const label patchi = iter.key();
802 
803  const vectorField Md
804  (
805  mesh_.C().boundaryField()[patchi] - CofR
806  );
807 
808  const scalarField sA(mag(Sfb[patchi]));
809 
810  // Normal force = surfaceUnitNormal*(surfaceNormal & forceDensity)
811  const vectorField fN
812  (
813  Sfb[patchi]/sA
814  *(
815  Sfb[patchi] & fD.boundaryField()[patchi]
816  )
817  );
818 
819  // Tangential force (total force minus normal fN)
820  const vectorField fT(sA*fD.boundaryField()[patchi] - fN);
821 
822  //- Porous force
823  const vectorField fP(Md.size(), Zero);
824 
825  applyBins(Md, fN, fT, fP, mesh_.C().boundaryField()[patchi]);
826  }
827  }
828  else
829  {
830  const volScalarField& p = obr_.lookupObject<volScalarField>(pName_);
831 
832  const surfaceVectorField::Boundary& Sfb =
833  mesh_.Sf().boundaryField();
834 
835  const surfaceScalarField::Boundary& magSfb =
836  mesh_.magSf().boundaryField();
837 
838  tmp<surfaceVectorField> tdevTau = devTau();
839  const surfaceVectorField::Boundary& devTaub = tdevTau().boundaryField();
840 
841  // Scale pRef by density for incompressible simulations
842  const scalar pRef = pRef_/rho(p);
843 
844  forAllConstIter(labelHashSet, patchSet_, iter)
845  {
846  const label patchi = iter.key();
847 
848  const vectorField Md
849  (
850  mesh_.C().boundaryField()[patchi] - CofR
851  );
852 
853  const vectorField fN
854  (
855  alpha(patchi)
856  *rho(p)
857  *Sfb[patchi]
858  *(p.boundaryField()[patchi] - pRef)
859  );
860 
861  const vectorField fT(magSfb[patchi] * devTaub[patchi]);
862 
863  const vectorField fP(Md.size(), Zero);
864 
865  applyBins(Md, fN, fT, fP, mesh_.C().boundaryField()[patchi]);
866  }
867  }
868 
869  if (porosity_)
870  {
871  const volVectorField& U = obr_.lookupObject<volVectorField>(UName_);
872  const volScalarField rho(this->rho());
873  const volScalarField mu(this->mu());
874 
875  const HashTable<const porosityModel*> models =
876  obr_.lookupClass<porosityModel>();
877 
878  if (models.empty())
879  {
881  << "Porosity effects requested, but no porosity models found "
882  << "in the database"
883  << endl;
884  }
885 
887  {
888  const porosityModel& pm = *iter();
889 
890  const vectorField fPTot(pm.force(U, rho, mu));
891 
892  const cellZone& cZone = mesh_.cellZones()[pm.zoneName()];
893  const vectorField d(mesh_.C(), cZone);
894  const vectorField fP(fPTot, cZone);
895  const vectorField Md(d - CofR);
896 
897  const vectorField fDummy(Md.size(), Zero);
898 
899  applyBins(Md, fDummy, fDummy, fP, d);
900  }
901  }
902 
907 }
908 
909 
911 {
913 }
914 
915 
917 {
918  return sum(force_[0]) + sum(force_[1]) + sum(force_[2]);
919 }
920 
921 
923 {
924  return sum(moment_[0]) + sum(moment_[1]) + sum(moment_[2]);
925 }
926 
927 
929 {
930  return true;
931 }
932 
933 
935 {
936  calcForcesMoments(CofR());
937 
938  if (Pstream::master())
939  {
940  logFiles::write();
941 
942  writeForces();
943 
944  writeBins();
945 
946  Log << endl;
947  }
948 
949  return true;
950 }
951 
952 
953 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
Macros for easy insertion into run-time selection tables.
Generic GeometricBoundaryField class.
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
static const GeometricField< Type, GeoMesh, PrimitiveField > & null()
Return a null geometric field.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
An STL-conforming hash table.
Definition: HashTable.H:127
bool empty() const
Return true if the hash table is empty.
Definition: HashTableI.H:72
static word groupName(Name name, const word &group)
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
Named list of cell indices representing a sub-set of the mesh.
Definition: cellZone.H:61
Base class for single-phase compressible turbulence models.
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
Calculates the forces and moments by integrating the pressure and skin-friction forces over a given l...
Definition: forcesBase.H:64
virtual void calcForcesMoments()
Calculate the forces and moments.
Definition: forcesBase.C:910
tmp< volScalarField > alpha() const
Get the volume fraction field.
Definition: forcesBase.C:370
void initialise()
Initialise the fields.
Definition: forcesBase.C:162
virtual void writeCoRHeader(Ostream &file)
Write the time varying centre of rotation column header.
Definition: forcesBase.C:454
tmp< surfaceVectorField > devTau() const
Return the effective surface stress.
Definition: forcesBase.C:209
void applyBins(const vectorField &Md, const vectorField &fN, const vectorField &fT, const vectorField &fP, const vectorField &d)
Accumulate bin data.
Definition: forcesBase.C:414
tmp< volScalarField > mu() const
Dynamic viscosity field.
Definition: forcesBase.C:261
virtual vector forceEff() const
Return the total force.
Definition: forcesBase.C:916
forcesBase(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: forcesBase.C:563
virtual vector momentEff() const
Return the total moment.
Definition: forcesBase.C:922
virtual void writeCoRValueHeader(Ostream &file)
Write the constant centre of rotation value in the header.
Definition: forcesBase.C:450
void writeBins()
Helper function to write bin data.
Definition: forcesBase.C:505
virtual ~forcesBase()
Destructor.
Definition: forcesBase.C:630
virtual void writeFileHeader(const label i)
Output file header information.
Definition: forcesBase.C:77
void writeForces()
Helper function to write force data.
Definition: forcesBase.C:466
tmp< volScalarField > rho() const
Return rho if specified otherwise rhoRef.
Definition: forcesBase.C:329
fileID
Enumeration for ensuring the right file is accessed.
Definition: forcesBase.H:72
virtual bool execute()
Execute, currently does nothing.
Definition: forcesBase.C:928
virtual bool write()
Write the forces.
Definition: forcesBase.C:934
virtual void writeCofR(Ostream &file)
Write the time varying centre of rotation.
Definition: forcesBase.C:460
virtual bool read(const dictionary &)
Read the forces data.
Definition: forcesBase.C:636
wordList createFileNames(const dictionary &dict) const
Create file names for forces and bins.
Definition: forcesBase.C:51
Specialisation of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:60
const wordList & names() const
Return the list of log file names.
Definition: logFiles.C:107
virtual bool write()
Write function.
Definition: logFiles.C:173
virtual bool read(const dictionary &)
Read optional controls.
Base class for single-phase incompressible turbulence models.
virtual tmp< surfaceVectorField > devSigma() const
Return the effective surface stress.
virtual tmp< volScalarField > nu() const
Return the laminar viscosity.
Registry of regIOobjects.
Templated abstract base class for multiphase compressible turbulence models.
Templated abstract base class for multiphase incompressible turbulence models.
A base class for physical properties.
Foam::polyBoundaryMesh.
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the patch set corresponding to the given names.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:276
Top level model for porosity models.
Definition: porosityModel.H:57
const word & zoneName() const
Return const access to the cell zone name.
virtual tmp< vectorField > force(const volVectorField &U, const volScalarField &rho, const volScalarField &mu) const
Return the force over the cell zone(s)
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
static const word null
An empty word.
Definition: word.H:77
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:346
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
Calculate the gradient of the given field.
label patchi
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
#define WarningInFunction
Report a warning using Foam::Warning.
#define Log
Report write to Foam::Info if the local log switch is true.
const dimensionedScalar mu
Atomic mass unit.
defineTypeNameAndDebug(adjustTimeStepToCombustion, 0)
static tmp< SurfaceField< Type > > interpolate(const VolField< Type > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
void writeHeader(std::ostream &, const bool isBinary, const std::string &title)
Write header.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
const dimensionSet dimPressure
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimKinematicViscosity
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
static const char tab
Definition: Ostream.H:266
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
const dimensionSet dimDensity
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh, PrimitiveField > &df)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
static const char nl
Definition: Ostream.H:267
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
labelList f(nPoints)
dictionary dict
volScalarField & p