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-2024 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"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace functionObjects
41 {
43 }
44 }
45 
46 
47 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
48 
50 (
51  const dictionary& dict
52 ) const
53 {
55 
56  const word forceType(dict.lookup("type"));
57 
58  // Name for file(fileID::mainFile=0)
59  names.append(forceType);
60 
61  if (dict.found("binData"))
62  {
63  const dictionary& binDict(dict.subDict("binData"));
64  const label nb = binDict.lookup<label>("nBin");
65  if (nb > 0)
66  {
67  // Name for file(fileID::binsFile=1)
68  names.append(forceType + "_bins");
69  }
70  }
71 
72  return names;
73 }
74 
75 
77 {
78  const word forceTypes
79  (
80  porosity_
81  ? "(pressure viscous porous)"
82  : "(pressure viscous)"
83  );
84 
85  switch (fileID(i))
86  {
87  case fileID::mainFile:
88  {
89  // force data
90 
91  writeHeader(file(i), "Forces");
92  writeCoRValueHeader(file(i));
93  writeCommented(file(i), "Time");
94 
95  writeCoRHeader(file(i));
96 
97  file(i)
98  << "forces" << forceTypes << tab
99  << "moments" << forceTypes;
100 
101  break;
102  }
103  case fileID::binsFile:
104  {
105  // bin data
106 
107  writeHeader(file(i), "Force bins");
108  writeHeaderValue(file(i), "bins", nBin_);
109  writeHeaderValue(file(i), "start", binMin_);
110  writeHeaderValue(file(i), "delta", binDx_);
111  writeHeaderValue(file(i), "direction", binDir_);
112 
113  vectorField binPoints(nBin_);
114  writeCommented(file(i), "x co-ords :");
115  forAll(binPoints, pointi)
116  {
117  binPoints[pointi] = (binMin_ + (pointi + 1)*binDx_)*binDir_;
118  file(i) << tab << binPoints[pointi].x();
119  }
120  file(i) << nl;
121 
122  writeCommented(file(i), "y co-ords :");
123  forAll(binPoints, pointi)
124  {
125  file(i) << tab << binPoints[pointi].y();
126  }
127  file(i) << nl;
128 
129  writeCommented(file(i), "z co-ords :");
130  forAll(binPoints, pointi)
131  {
132  file(i) << tab << binPoints[pointi].z();
133  }
134  file(i) << nl;
135 
136  writeCommented(file(i), "Time");
137 
138  for (label j = 0; j < nBin_; j++)
139  {
140  const word jn('(' + Foam::name(j) + ')');
141  const word f("forces" + jn + forceTypes);
142  const word m("moments" + jn + forceTypes);
143 
144  file(i)<< tab << f << tab << m;
145  }
146 
147  break;
148  }
149  default:
150  {
152  << "Unhandled file index: " << i
153  << abort(FatalError);
154  }
155  }
156 
157  file(i)<< endl;
158 }
159 
160 
162 {
163  if (initialised_)
164  {
165  return;
166  }
167 
168  if (directForceDensity_)
169  {
170  if (!obr_.foundObject<volVectorField>(fDName_))
171  {
173  << "Could not find " << fDName_ << " in database."
174  << exit(FatalError);
175  }
176  }
177  else
178  {
179  if
180  (
181  !obr_.foundObject<volVectorField>(UName_)
182  || !obr_.foundObject<volScalarField>(pName_)
183 
184  )
185  {
187  << "Could not find " << UName_ << ", " << pName_
188  << exit(FatalError);
189  }
190 
191  if
192  (
193  rhoName_ != "rhoInf"
194  && !obr_.foundObject<volScalarField>(rhoName_)
195  )
196  {
198  << "Could not find " << rhoName_
199  << exit(FatalError);
200  }
201  }
202 
203  initialised_ = true;
204 }
205 
206 
209 {
211  typedef compressible::momentumTransportModel cmpModel;
212  typedef phaseIncompressible::momentumTransportModel phaseIcoModel;
213  typedef phaseCompressible::momentumTransportModel phaseCmpModel;
214 
215  const word& modelName = momentumTransportModel::typeName;
216  const word phaseModelName =
217  phaseName_ == word::null
218  ? word::null
219  : IOobject::groupName(momentumTransportModel::typeName, phaseName_);
220 
221  if (obr_.foundObject<icoModel>(modelName))
222  {
224  obr_.lookupObject<icoModel>(modelName);
225 
226  return alpha()*rho()*model.devSigma();
227  }
228  else if (obr_.foundObject<cmpModel>(modelName))
229  {
230  const cmpModel& model =
231  obr_.lookupObject<cmpModel>(modelName);
232 
233  return alpha()*model.devTau();
234  }
235  else if (obr_.foundObject<phaseIcoModel>(phaseModelName))
236  {
237  const phaseIcoModel& model =
238  obr_.lookupObject<phaseIcoModel>(phaseModelName);
239 
240  return rho()*model.devSigma();
241  }
242  else if (obr_.foundObject<phaseCmpModel>(phaseModelName))
243  {
244  const phaseCmpModel& model =
245  obr_.lookupObject<phaseCmpModel>(phaseModelName);
246 
247  return model.devTau();
248  }
249  else if (obr_.foundObject<dictionary>("physicalProperties"))
250  {
251  // Legacy support for icoFoam
252 
254  obr_.lookupObject<dictionary>("physicalProperties");
255 
256  const dimensionedScalar nu
257  (
258  "nu",
261  );
262 
263  const volVectorField& U = obr_.lookupObject<volVectorField>(UName_);
264 
265  return -rho()*nu*dev(twoSymm(fvc::grad(U)));
266  }
267  else
268  {
270  << "No valid model for viscous stress calculation"
271  << exit(FatalError);
272 
273  return volSymmTensorField::null();
274  }
275 }
276 
277 
279 {
281  typedef compressible::momentumTransportModel cmpModel;
282  typedef phaseIncompressible::momentumTransportModel phaseIcoModel;
283  typedef phaseCompressible::momentumTransportModel phaseCmpModel;
284 
285  const word& modelName = momentumTransportModel::typeName;
286  const word phaseModelName =
287  phaseName_ == word::null
288  ? word::null
289  : IOobject::groupName(momentumTransportModel::typeName, phaseName_);
290 
291  if (obr_.foundObject<icoModel>(modelName))
292  {
294  obr_.lookupObject<icoModel>(modelName);
295 
296  return rho()*model.nu();
297  }
298  else if (obr_.foundObject<cmpModel>(modelName))
299  {
300  const cmpModel& model =
301  obr_.lookupObject<cmpModel>(modelName);
302 
303  return model.rho()*model.nu();
304  }
305  else if (obr_.foundObject<phaseIcoModel>(phaseModelName))
306  {
307  const phaseIcoModel& model =
308  obr_.lookupObject<phaseIcoModel>(phaseModelName);
309 
310  return rho()*model.nu();
311  }
312  else if (obr_.foundObject<phaseCmpModel>(phaseModelName))
313  {
314  const phaseCmpModel& model =
315  obr_.lookupObject<phaseCmpModel>(phaseModelName);
316 
317  return model.rho()*model.nu();
318  }
319  else if (obr_.foundObject<dictionary>("physicalProperties"))
320  {
321  // Legacy support for icoFoam
322 
324  obr_.lookupObject<dictionary>("physicalProperties");
325 
326  const dimensionedScalar nu
327  (
328  "nu",
331  );
332 
333  return rho()*nu;
334  }
335  else
336  {
338  << "No valid model for dynamic viscosity calculation"
339  << exit(FatalError);
340 
341  return volScalarField::null();
342  }
343 }
344 
345 
347 {
348  if (rhoName_ == "rhoInf")
349  {
350  return volScalarField::New
351  (
352  "rho",
353  mesh_,
354  dimensionedScalar(dimDensity, rhoRef_)
355  );
356  }
357  else
358  {
359  return(obr_.lookupObject<volScalarField>(rhoName_));
360  }
361 }
362 
363 
365 (
366  const volScalarField& p
367 ) const
368 {
369  if (p.dimensions() == dimPressure)
370  {
371  return 1.0;
372  }
373  else
374  {
375  if (rhoName_ != "rhoInf")
376  {
378  << "Dynamic pressure is expected but kinematic is provided."
379  << exit(FatalError);
380  }
381 
382  return rhoRef_;
383  }
384 }
385 
386 
388 {
389  if (phaseName_ == word::null)
390  {
391  return volScalarField::New
392  (
393  "alpha",
394  mesh_,
396  );
397  }
398  else
399  {
400  return obr_.lookupObject<volScalarField>
401  (
402  IOobject::groupName("alpha", phaseName_)
403  );
404  }
405 }
406 
407 
409 (
410  const label patchi
411 ) const
412 {
413  if (phaseName_ == word::null)
414  {
415  return tmp<scalarField>
416  (
417  new scalarField(mesh_.boundary()[patchi].size(), 1)
418  );
419  }
420  else
421  {
422  return obr_.lookupObject<volScalarField>
423  (
424  IOobject::groupName("alpha", phaseName_)
425  ).boundaryField()[patchi];
426  }
427 }
428 
429 
431 (
432  const vectorField& Md,
433  const vectorField& fN,
434  const vectorField& fT,
435  const vectorField& fP,
436  const vectorField& d
437 )
438 {
439  if (nBin_ == 1)
440  {
441  force_[0][0] += sum(fN);
442  force_[1][0] += sum(fT);
443  force_[2][0] += sum(fP);
444  moment_[0][0] += sum(Md^fN);
445  moment_[1][0] += sum(Md^fT);
446  moment_[2][0] += sum(Md^fP);
447  }
448  else
449  {
450  scalarField dd((d & binDir_) - binMin_);
451 
452  forAll(dd, i)
453  {
454  label bini = min(max(floor(dd[i]/binDx_), 0), force_[0].size() - 1);
455 
456  force_[0][bini] += fN[i];
457  force_[1][bini] += fT[i];
458  force_[2][bini] += fP[i];
459  moment_[0][bini] += Md[i]^fN[i];
460  moment_[1][bini] += Md[i]^fT[i];
461  moment_[2][bini] += Md[i]^fP[i];
462  }
463  }
464 }
465 
466 
468 {}
469 
470 
472 {
473  file << "CofR" << tab;
474 }
475 
476 
478 {
479  file << CofR();
480 }
481 
482 
484 {
485  Log << type() << " " << name() << " write:" << nl
486  << " sum of forces:" << nl
487  << " pressure : " << sum(force_[0]) << nl
488  << " viscous : " << sum(force_[1]) << nl
489  << " porous : " << sum(force_[2]) << nl
490  << " sum of moments:" << nl
491  << " pressure : " << sum(moment_[0]) << nl
492  << " viscous : " << sum(moment_[1]) << nl
493  << " porous : " << sum(moment_[2])
494  << endl;
495 
496  writeTime(file(fileID::mainFile));
497  writeCofR(file(fileID::mainFile));
498 
499  if (porosity_)
500  {
501  file(fileID::mainFile) << tab << setw(1) << '('
502  << sum(force_[0]) << setw(1) << ' '
503  << sum(force_[1]) << setw(1) << ' '
504  << sum(force_[2]) << setw(3) << ") ("
505  << sum(moment_[0]) << setw(1) << ' '
506  << sum(moment_[1]) << setw(1) << ' '
507  << sum(moment_[2]) << setw(1) << ')';
508  }
509  else
510  {
511  file(fileID::mainFile) << tab << setw(1) << '('
512  << sum(force_[0]) << setw(1) << ' '
513  << sum(force_[1]) << setw(3) << ") ("
514  << sum(moment_[0]) << setw(1) << ' '
515  << sum(moment_[1]) << setw(1) << ')';
516  }
517 
518  file(fileID::mainFile) << endl;
519 }
520 
521 
523 {
524  if (nBin_ == 1)
525  {
526  return;
527  }
528 
529  List<Field<vector>> f(force_);
530  List<Field<vector>> m(moment_);
531 
532  if (binCumulative_)
533  {
534  for (label i = 1; i < f[0].size(); i++)
535  {
536  f[0][i] += f[0][i-1];
537  f[1][i] += f[1][i-1];
538  f[2][i] += f[2][i-1];
539 
540  m[0][i] += m[0][i-1];
541  m[1][i] += m[1][i-1];
542  m[2][i] += m[2][i-1];
543  }
544  }
545 
546  writeTime(file(fileID::binsFile));
547 
548 
549  forAll(f[0], i)
550  {
551  if (porosity_)
552  {
553  file(fileID::binsFile)
554  << tab << setw(1) << '('
555  << f[0][i] << setw(1) << ' '
556  << f[1][i] << setw(1) << ' '
557  << f[2][i] << setw(3) << ") ("
558  << m[0][i] << setw(1) << ' '
559  << m[1][i] << setw(1) << ' '
560  << m[2][i] << setw(1) << ')';
561  }
562  else
563  {
564  file(fileID::binsFile)
565  << tab << setw(1) << '('
566  << f[0][i] << setw(1) << ' '
567  << f[1][i] << setw(3) << ") ("
568  << m[0][i] << setw(1) << ' '
569  << m[1][i] << setw(1) << ')';
570  }
571  }
572 
573  file(fileID::binsFile) << endl;
574 }
575 
576 
577 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
578 
580 (
581  const word& name,
582  const Time& runTime,
583  const dictionary& dict
584 )
585 :
586  fvMeshFunctionObject(name, runTime, dict),
587  logFiles(obr_, name),
588  force_(3),
589  moment_(3),
590  patchSet_(),
591  pName_(word::null),
592  UName_(word::null),
593  rhoName_(word::null),
594  phaseName_(word::null),
595  directForceDensity_(false),
596  fDName_(""),
597  rhoRef_(vGreat),
598  pRef_(0),
599  porosity_(false),
600  nBin_(1),
601  binDir_(Zero),
602  binDx_(0.0),
603  binMin_(great),
604  binPoints_(),
605  binCumulative_(true),
606  initialised_(false)
607 {
608  read(dict);
609 }
610 
611 
613 (
614  const word& name,
615  const objectRegistry& obr,
616  const dictionary& dict
617 )
618 :
620  logFiles(obr_, name),
621  force_(3),
622  moment_(3),
623  patchSet_(),
624  pName_(word::null),
625  UName_(word::null),
626  rhoName_(word::null),
627  phaseName_(word::null),
628  directForceDensity_(false),
629  fDName_(""),
630  rhoRef_(vGreat),
631  pRef_(0),
632  porosity_(false),
633  nBin_(1),
634  binDir_(Zero),
635  binDx_(0.0),
636  binMin_(great),
637  binPoints_(),
638  binCumulative_(true),
639  initialised_(false)
640 {
641  read(dict);
642 }
643 
644 
645 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
646 
648 {}
649 
650 
651 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
652 
654 {
656 
657  initialised_ = false;
658 
659  Log << type() << " " << name() << ":" << nl;
660 
661  directForceDensity_ = dict.lookupOrDefault("directForceDensity", false);
662 
663  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
664 
665  patchSet_ = patchSet(dict);
666 
667  if (directForceDensity_)
668  {
669  // Optional entry for fDName
670  fDName_ = dict.lookupOrDefault<word>("fD", "fD");
671  }
672  else
673  {
674  // Optional phase entry
675  phaseName_ = dict.lookupOrDefault<word>("phase", word::null);
676 
677  // Optional U, p and rho entries
678  pName_ =
679  dict.lookupOrDefault<word>
680  (
681  "p",
682  IOobject::groupName("p", phaseName_)
683  );
684  UName_ =
685  dict.lookupOrDefault<word>
686  (
687  "U",
688  IOobject::groupName("U", phaseName_)
689  );
690  rhoName_ =
691  dict.lookupOrDefault<word>
692  (
693  "rho",
694  IOobject::groupName("rho", phaseName_)
695  );
696 
697  // Reference density needed for incompressible calculations
698  if (rhoName_ == "rhoInf")
699  {
700  dict.lookup("rhoInf") >> rhoRef_;
701  }
702 
703  // Reference pressure, 0 by default
704  pRef_ = dict.lookupOrDefault<scalar>("pRef", 0.0);
705  }
706 
707  dict.readIfPresent("porosity", porosity_);
708  if (porosity_)
709  {
710  Log << " Including porosity effects" << endl;
711  }
712  else
713  {
714  Log << " Not including porosity effects" << endl;
715  }
716 
717  if (dict.found("binData"))
718  {
719  const dictionary& binDict(dict.subDict("binData"));
720  binDict.lookup("nBin") >> nBin_;
721 
722  if (nBin_ < 0)
723  {
725  << "Number of bins (nBin) must be zero or greater"
726  << exit(FatalIOError);
727  }
728  else if ((nBin_ == 0) || (nBin_ == 1))
729  {
730  nBin_ = 1;
731  forAll(force_, i)
732  {
733  force_[i].setSize(1);
734  moment_[i].setSize(1);
735  }
736  }
737 
738  if (nBin_ > 1)
739  {
740  binDict.lookup("direction") >> binDir_;
741  binDir_ /= mag(binDir_);
742 
743  binMin_ = great;
744  scalar binMax = -great;
745  forAllConstIter(labelHashSet, patchSet_, iter)
746  {
747  const label patchi = iter.key();
748  const polyPatch& pp = pbm[patchi];
749  const scalarField d(pp.faceCentres() & binDir_);
750  binMin_ = min(min(d), binMin_);
751  binMax = max(max(d), binMax);
752  }
753  reduce(binMin_, minOp<scalar>());
754  reduce(binMax, maxOp<scalar>());
755 
756  // slightly boost binMax so that region of interest is fully
757  // within bounds
758  binMax = 1.0001*(binMax - binMin_) + binMin_;
759 
760  binDx_ = (binMax - binMin_)/scalar(nBin_);
761 
762  // create the bin points used for writing
763  binPoints_.setSize(nBin_);
764  forAll(binPoints_, i)
765  {
766  binPoints_[i] = (i + 0.5)*binDir_*binDx_;
767  }
768 
769  binDict.lookup("cumulative") >> binCumulative_;
770 
771  // allocate storage for forces and moments
772  forAll(force_, i)
773  {
774  force_[i].setSize(nBin_);
775  moment_[i].setSize(nBin_);
776  }
777  }
778  }
779 
780  if (nBin_ == 1)
781  {
782  // allocate storage for forces and moments
783  force_[0].setSize(1);
784  force_[1].setSize(1);
785  force_[2].setSize(1);
786  moment_[0].setSize(1);
787  moment_[1].setSize(1);
788  moment_[2].setSize(1);
789  }
790 
791  resetNames(createFileNames(dict));
792 
793  return true;
794 }
795 
796 
798 {
799  initialise();
800 
801  force_[0] = Zero;
802  force_[1] = Zero;
803  force_[2] = Zero;
804 
805  moment_[0] = Zero;
806  moment_[1] = Zero;
807  moment_[2] = Zero;
808 
809  if (directForceDensity_)
810  {
811  const volVectorField& fD = obr_.lookupObject<volVectorField>(fDName_);
812 
813  const surfaceVectorField::Boundary& Sfb =
814  mesh_.Sf().boundaryField();
815 
816  forAllConstIter(labelHashSet, patchSet_, iter)
817  {
818  const label patchi = iter.key();
819 
820  const vectorField Md
821  (
822  mesh_.C().boundaryField()[patchi] - CofR
823  );
824 
825  const scalarField sA(mag(Sfb[patchi]));
826 
827  // Normal force = surfaceUnitNormal*(surfaceNormal & forceDensity)
828  const vectorField fN
829  (
830  Sfb[patchi]/sA
831  *(
832  Sfb[patchi] & fD.boundaryField()[patchi]
833  )
834  );
835 
836  // Tangential force (total force minus normal fN)
837  const vectorField fT(sA*fD.boundaryField()[patchi] - fN);
838 
839  //- Porous force
840  const vectorField fP(Md.size(), Zero);
841 
842  applyBins(Md, fN, fT, fP, mesh_.C().boundaryField()[patchi]);
843  }
844  }
845  else
846  {
847  const volScalarField& p = obr_.lookupObject<volScalarField>(pName_);
848 
849  const surfaceVectorField::Boundary& Sfb =
850  mesh_.Sf().boundaryField();
851 
852  tmp<volSymmTensorField> tdevTau = devTau();
853  const volSymmTensorField::Boundary& devTaub =
854  tdevTau().boundaryField();
855 
856  // Scale pRef by density for incompressible simulations
857  const scalar pRef = pRef_/rho(p);
858 
859  forAllConstIter(labelHashSet, patchSet_, iter)
860  {
861  const label patchi = iter.key();
862 
863  const vectorField Md
864  (
865  mesh_.C().boundaryField()[patchi] - CofR
866  );
867 
868  const vectorField fN
869  (
870  alpha(patchi)
871  *rho(p)
872  *Sfb[patchi]
873  *(p.boundaryField()[patchi] - pRef)
874  );
875 
876  const vectorField fT(Sfb[patchi] & devTaub[patchi]);
877 
878  const vectorField fP(Md.size(), Zero);
879 
880  applyBins(Md, fN, fT, fP, mesh_.C().boundaryField()[patchi]);
881  }
882  }
883 
884  if (porosity_)
885  {
886  const volVectorField& U = obr_.lookupObject<volVectorField>(UName_);
887  const volScalarField rho(this->rho());
888  const volScalarField mu(this->mu());
889 
890  const HashTable<const porosityModel*> models =
891  obr_.lookupClass<porosityModel>();
892 
893  if (models.empty())
894  {
896  << "Porosity effects requested, but no porosity models found "
897  << "in the database"
898  << endl;
899  }
900 
902  {
903  const porosityModel& pm = *iter();
904 
905  const vectorField fPTot(pm.force(U, rho, mu));
906 
907  const cellZone& cZone = mesh_.cellZones()[pm.zoneName()];
908  const vectorField d(mesh_.C(), cZone);
909  const vectorField fP(fPTot, cZone);
910  const vectorField Md(d - CofR);
911 
912  const vectorField fDummy(Md.size(), Zero);
913 
914  applyBins(Md, fDummy, fDummy, fP, d);
915  }
916  }
917 
922 }
923 
924 
926 {
928 }
929 
930 
932 {
933  return sum(force_[0]) + sum(force_[1]) + sum(force_[2]);
934 }
935 
936 
938 {
939  return sum(moment_[0]) + sum(moment_[1]) + sum(moment_[2]);
940 }
941 
942 
944 {
945  return true;
946 }
947 
948 
950 {
951  calcForcesMoments(CofR());
952 
953  if (Pstream::master())
954  {
955  logFiles::write();
956 
957  writeForces();
958 
959  writeBins();
960 
961  Log << endl;
962  }
963 
964  return true;
965 }
966 
967 
968 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
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, PatchField, GeoMesh > & null()
Return a null geometric field.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, 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
A subset of mesh cells.
Definition: cellZone.H:58
Base class for single-phase compressible turbulence models.
A list of keyword definitions, which are a keyword followed by any number of values (e....
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:710
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:925
tmp< volScalarField > alpha() const
Get the volume fraction field.
Definition: forcesBase.C:387
void initialise()
Initialise the fields.
Definition: forcesBase.C:161
virtual void writeCoRHeader(Ostream &file)
Write the time varying centre of rotation column header.
Definition: forcesBase.C:471
tmp< volSymmTensorField > devTau() const
Return the effective viscous stress (laminar + turbulent).
Definition: forcesBase.C:208
void applyBins(const vectorField &Md, const vectorField &fN, const vectorField &fT, const vectorField &fP, const vectorField &d)
Accumulate bin data.
Definition: forcesBase.C:431
tmp< volScalarField > mu() const
Dynamic viscosity field.
Definition: forcesBase.C:278
virtual vector forceEff() const
Return the total force.
Definition: forcesBase.C:931
forcesBase(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: forcesBase.C:580
virtual vector momentEff() const
Return the total moment.
Definition: forcesBase.C:937
virtual void writeCoRValueHeader(Ostream &file)
Write the constant centre of rotation value in the header.
Definition: forcesBase.C:467
void writeBins()
Helper function to write bin data.
Definition: forcesBase.C:522
virtual ~forcesBase()
Destructor.
Definition: forcesBase.C:647
virtual void writeFileHeader(const label i)
Output file header information.
Definition: forcesBase.C:76
void writeForces()
Helper function to write force data.
Definition: forcesBase.C:483
tmp< volScalarField > rho() const
Return rho if specified otherwise rhoRef.
Definition: forcesBase.C:346
fileID
Enumeration for ensuring the right file is accessed.
Definition: forcesBase.H:72
virtual bool execute()
Execute, currently does nothing.
Definition: forcesBase.C:943
virtual bool write()
Write the forces.
Definition: forcesBase.C:949
virtual void writeCofR(Ostream &file)
Write the time varying centre of rotation.
Definition: forcesBase.C:477
virtual bool read(const dictionary &)
Read the forces data.
Definition: forcesBase.C:653
wordList createFileNames(const dictionary &dict) const
Create file names for forces and bins.
Definition: forcesBase.C:50
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< volSymmTensorField > devSigma() const
Return the effective stress tensor.
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.
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)
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
void writeHeader(std::ostream &, const bool isBinary, const std::string &title)
Write header.
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
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:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimKinematicViscosity
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
static const char tab
Definition: Ostream.H:265
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const dimensionSet dimless
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
dimensioned< scalar > mag(const dimensioned< Type > &)
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOerror FatalIOError
error FatalError
static const char nl
Definition: Ostream.H:266
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