forces.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-2019 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 "forces.H"
27 #include "fvcGrad.H"
28 #include "porosityModel.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace functionObjects
38 {
39  defineTypeNameAndDebug(forces, 0);
40 
41  addToRunTimeSelectionTable(functionObject, forces, dictionary);
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
47 
49 (
50  const dictionary& dict
51 ) const
52 {
53  DynamicList<word> names(1);
54 
55  const word forceType(dict.lookup("type"));
56 
57  // Name for file(fileID::mainFile=0)
58  names.append(forceType);
59 
60  if (dict.found("binData"))
61  {
62  const dictionary& binDict(dict.subDict("binData"));
63  label nb = readLabel(binDict.lookup("nBin"));
64  if (nb > 0)
65  {
66  // Name for file(fileID::binsFile=1)
67  names.append(forceType + "_bins");
68  }
69  }
70 
71  return move(names);
72 }
73 
74 
76 {
77  switch (fileID(i))
78  {
79  case fileID::mainFile:
80  {
81  // force data
82 
83  writeHeader(file(i), "Forces");
84  writeHeaderValue(file(i), "CofR", coordSys_.origin());
85  writeCommented(file(i), "Time");
86 
87  const word forceTypes("(pressure viscous porous)");
88  file(i)
89  << "forces" << forceTypes << tab
90  << "moments" << forceTypes;
91 
92  if (localSystem_)
93  {
94  file(i)
95  << tab
96  << "localForces" << forceTypes << tab
97  << "localMoments" << forceTypes;
98  }
99 
100  break;
101  }
102  case fileID::binsFile:
103  {
104  // bin data
105 
106  writeHeader(file(i), "Force bins");
107  writeHeaderValue(file(i), "bins", nBin_);
108  writeHeaderValue(file(i), "start", binMin_);
109  writeHeaderValue(file(i), "delta", binDx_);
110  writeHeaderValue(file(i), "direction", binDir_);
111 
112  vectorField binPoints(nBin_);
113  writeCommented(file(i), "x co-ords :");
114  forAll(binPoints, pointi)
115  {
116  binPoints[pointi] = (binMin_ + (pointi + 1)*binDx_)*binDir_;
117  file(i) << tab << binPoints[pointi].x();
118  }
119  file(i) << nl;
120 
121  writeCommented(file(i), "y co-ords :");
122  forAll(binPoints, pointi)
123  {
124  file(i) << tab << binPoints[pointi].y();
125  }
126  file(i) << nl;
127 
128  writeCommented(file(i), "z co-ords :");
129  forAll(binPoints, pointi)
130  {
131  file(i) << tab << binPoints[pointi].z();
132  }
133  file(i) << nl;
134 
135  writeCommented(file(i), "Time");
136 
137  const word binForceTypes("[pressure,viscous,porous]");
138  for (label j = 0; j < nBin_; j++)
139  {
140  const word jn('(' + Foam::name(j) + ')');
141  const word f("forces" + jn + binForceTypes);
142  const word m("moments" + jn + binForceTypes);
143 
144  file(i)<< tab << f << tab << m;
145  }
146  if (localSystem_)
147  {
148  for (label j = 0; j < nBin_; j++)
149  {
150  const word jn('(' + Foam::name(j) + ')');
151  const word f("localForces" + jn + binForceTypes);
152  const word m("localMoments" + jn + binForceTypes);
153 
154  file(i)<< tab << f << tab << m;
155  }
156  }
157 
158  break;
159  }
160  default:
161  {
163  << "Unhandled file index: " << i
164  << abort(FatalError);
165  }
166  }
167 
168  file(i)<< endl;
169 }
170 
171 
173 {
174  if (initialised_)
175  {
176  return;
177  }
178 
179  if (directForceDensity_)
180  {
181  if (!obr_.foundObject<volVectorField>(fDName_))
182  {
184  << "Could not find " << fDName_ << " in database."
185  << exit(FatalError);
186  }
187  }
188  else
189  {
190  if
191  (
192  !obr_.foundObject<volVectorField>(UName_)
193  || !obr_.foundObject<volScalarField>(pName_)
194 
195  )
196  {
198  << "Could not find " << UName_ << ", " << pName_
199  << exit(FatalError);
200  }
201 
202  if
203  (
204  rhoName_ != "rhoInf"
205  && !obr_.foundObject<volScalarField>(rhoName_)
206  )
207  {
209  << "Could not find " << rhoName_
210  << exit(FatalError);
211  }
212  }
213 
214  initialised_ = true;
215 }
216 
217 
220 {
221  typedef compressible::turbulenceModel cmpTurbModel;
222  typedef incompressible::turbulenceModel icoTurbModel;
223 
224  if (obr_.foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
225  {
226  const cmpTurbModel& turb =
227  obr_.lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
228 
229  return turb.devRhoReff();
230  }
231  else if (obr_.foundObject<icoTurbModel>(icoTurbModel::propertiesName))
232  {
233  const incompressible::turbulenceModel& turb =
234  obr_.lookupObject<icoTurbModel>(icoTurbModel::propertiesName);
235 
236  return rho()*turb.devReff();
237  }
238  else if (obr_.foundObject<fluidThermo>(fluidThermo::dictName))
239  {
240  const fluidThermo& thermo =
241  obr_.lookupObject<fluidThermo>(fluidThermo::dictName);
242 
243  const volVectorField& U = obr_.lookupObject<volVectorField>(UName_);
244 
245  return -thermo.mu()*dev(twoSymm(fvc::grad(U)));
246  }
247  else if
248  (
249  obr_.foundObject<transportModel>("transportProperties")
250  )
251  {
252  const transportModel& laminarT =
253  obr_.lookupObject<transportModel>("transportProperties");
254 
255  const volVectorField& U = obr_.lookupObject<volVectorField>(UName_);
256 
257  return -rho()*laminarT.nu()*dev(twoSymm(fvc::grad(U)));
258  }
259  else if (obr_.foundObject<dictionary>("transportProperties"))
260  {
261  const dictionary& transportProperties =
262  obr_.lookupObject<dictionary>("transportProperties");
263 
265  (
266  "nu",
267  dimViscosity,
268  transportProperties.lookup("nu")
269  );
270 
271  const volVectorField& U = obr_.lookupObject<volVectorField>(UName_);
272 
273  return -rho()*nu*dev(twoSymm(fvc::grad(U)));
274  }
275  else
276  {
278  << "No valid model for viscous stress calculation"
279  << exit(FatalError);
280 
281  return volSymmTensorField::null();
282  }
283 }
284 
285 
287 {
288  if (obr_.foundObject<fluidThermo>(basicThermo::dictName))
289  {
290  const fluidThermo& thermo =
291  obr_.lookupObject<fluidThermo>(basicThermo::dictName);
292 
293  return thermo.mu();
294  }
295  else if
296  (
297  obr_.foundObject<transportModel>("transportProperties")
298  )
299  {
300  const transportModel& laminarT =
301  obr_.lookupObject<transportModel>("transportProperties");
302 
303  return rho()*laminarT.nu();
304  }
305  else if (obr_.foundObject<dictionary>("transportProperties"))
306  {
307  const dictionary& transportProperties =
308  obr_.lookupObject<dictionary>("transportProperties");
309 
311  (
312  "nu",
313  dimViscosity,
314  transportProperties.lookup("nu")
315  );
316 
317  return rho()*nu;
318  }
319  else
320  {
322  << "No valid model for dynamic viscosity calculation"
323  << exit(FatalError);
324 
325  return volScalarField::null();
326  }
327 }
328 
329 
331 {
332  if (rhoName_ == "rhoInf")
333  {
334  return volScalarField::New
335  (
336  "rho",
337  mesh_,
338  dimensionedScalar(dimDensity, rhoRef_)
339  );
340  }
341  else
342  {
343  return(obr_.lookupObject<volScalarField>(rhoName_));
344  }
345 }
346 
347 
349 {
350  if (p.dimensions() == dimPressure)
351  {
352  return 1.0;
353  }
354  else
355  {
356  if (rhoName_ != "rhoInf")
357  {
359  << "Dynamic pressure is expected but kinematic is provided."
360  << exit(FatalError);
361  }
362 
363  return rhoRef_;
364  }
365 }
366 
367 
369 (
370  const vectorField& Md,
371  const vectorField& fN,
372  const vectorField& fT,
373  const vectorField& fP,
374  const vectorField& d
375 )
376 {
377  if (nBin_ == 1)
378  {
379  force_[0][0] += sum(fN);
380  force_[1][0] += sum(fT);
381  force_[2][0] += sum(fP);
382  moment_[0][0] += sum(Md^fN);
383  moment_[1][0] += sum(Md^fT);
384  moment_[2][0] += sum(Md^fP);
385  }
386  else
387  {
388  scalarField dd((d & binDir_) - binMin_);
389 
390  forAll(dd, i)
391  {
392  label bini = min(max(floor(dd[i]/binDx_), 0), force_[0].size() - 1);
393 
394  force_[0][bini] += fN[i];
395  force_[1][bini] += fT[i];
396  force_[2][bini] += fP[i];
397  moment_[0][bini] += Md[i]^fN[i];
398  moment_[1][bini] += Md[i]^fT[i];
399  moment_[2][bini] += Md[i]^fP[i];
400  }
401  }
402 }
403 
404 
406 {
407  Log << type() << " " << name() << " write:" << nl
408  << " sum of forces:" << nl
409  << " pressure : " << sum(force_[0]) << nl
410  << " viscous : " << sum(force_[1]) << nl
411  << " porous : " << sum(force_[2]) << nl
412  << " sum of moments:" << nl
413  << " pressure : " << sum(moment_[0]) << nl
414  << " viscous : " << sum(moment_[1]) << nl
415  << " porous : " << sum(moment_[2])
416  << endl;
417 
418  writeTime(file(fileID::mainFile));
419 
420  file(fileID::mainFile) << tab << setw(1) << '('
421  << sum(force_[0]) << setw(1) << ' '
422  << sum(force_[1]) << setw(1) << ' '
423  << sum(force_[2]) << setw(3) << ") ("
424  << sum(moment_[0]) << setw(1) << ' '
425  << sum(moment_[1]) << setw(1) << ' '
426  << sum(moment_[2]) << setw(1) << ')';
427 
428  if (localSystem_)
429  {
430  vectorField localForceN(coordSys_.localVector(force_[0]));
431  vectorField localForceT(coordSys_.localVector(force_[1]));
432  vectorField localForceP(coordSys_.localVector(force_[2]));
433  vectorField localMomentN(coordSys_.localVector(moment_[0]));
434  vectorField localMomentT(coordSys_.localVector(moment_[1]));
435  vectorField localMomentP(coordSys_.localVector(moment_[2]));
436 
437  file(fileID::mainFile) << tab << setw(1) << '('
438  << sum(localForceN) << setw(1) << ' '
439  << sum(localForceT) << setw(1) << ' '
440  << sum(localForceP) << setw(3) << ") ("
441  << sum(localMomentN) << setw(1) << ' '
442  << sum(localMomentT) << setw(1) << ' '
443  << sum(localMomentP) << setw(1) << ')';
444  }
445 
446  file(fileID::mainFile) << endl;
447 }
448 
449 
451 {
452  if (nBin_ == 1)
453  {
454  return;
455  }
456 
457  List<Field<vector>> f(force_);
458  List<Field<vector>> m(moment_);
459 
460  if (binCumulative_)
461  {
462  for (label i = 1; i < f[0].size(); i++)
463  {
464  f[0][i] += f[0][i-1];
465  f[1][i] += f[1][i-1];
466  f[2][i] += f[2][i-1];
467 
468  m[0][i] += m[0][i-1];
469  m[1][i] += m[1][i-1];
470  m[2][i] += m[2][i-1];
471  }
472  }
473 
474  writeTime(file(fileID::binsFile));
475 
476  forAll(f[0], i)
477  {
478  file(fileID::binsFile)
479  << tab << setw(1) << '('
480  << f[0][i] << setw(1) << ' '
481  << f[1][i] << setw(1) << ' '
482  << f[2][i] << setw(3) << ") ("
483  << m[0][i] << setw(1) << ' '
484  << m[1][i] << setw(1) << ' '
485  << m[2][i] << setw(1) << ')';
486  }
487 
488  if (localSystem_)
489  {
490  List<Field<vector>> lf(3);
491  List<Field<vector>> lm(3);
492  lf[0] = coordSys_.localVector(force_[0]);
493  lf[1] = coordSys_.localVector(force_[1]);
494  lf[2] = coordSys_.localVector(force_[2]);
495  lm[0] = coordSys_.localVector(moment_[0]);
496  lm[1] = coordSys_.localVector(moment_[1]);
497  lm[2] = coordSys_.localVector(moment_[2]);
498 
499  if (binCumulative_)
500  {
501  for (label i = 1; i < lf[0].size(); i++)
502  {
503  lf[0][i] += lf[0][i-1];
504  lf[1][i] += lf[1][i-1];
505  lf[2][i] += lf[2][i-1];
506  lm[0][i] += lm[0][i-1];
507  lm[1][i] += lm[1][i-1];
508  lm[2][i] += lm[2][i-1];
509  }
510  }
511 
512  forAll(lf[0], i)
513  {
514  file(fileID::binsFile)
515  << tab << setw(1) << '('
516  << lf[0][i] << setw(1) << ' '
517  << lf[1][i] << setw(1) << ' '
518  << lf[2][i] << setw(3) << ") ("
519  << lm[0][i] << setw(1) << ' '
520  << lm[1][i] << setw(1) << ' '
521  << lm[2][i] << setw(1) << ')';
522  }
523  }
524 
525  file(fileID::binsFile) << endl;
526 }
527 
528 
529 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
530 
532 (
533  const word& name,
534  const Time& runTime,
535  const dictionary& dict
536 )
537 :
538  fvMeshFunctionObject(name, runTime, dict),
539  logFiles(obr_, name),
540  force_(3),
541  moment_(3),
542  patchSet_(),
543  pName_(word::null),
544  UName_(word::null),
545  rhoName_(word::null),
546  directForceDensity_(false),
547  fDName_(""),
548  rhoRef_(vGreat),
549  pRef_(0),
550  coordSys_(),
551  localSystem_(false),
552  porosity_(false),
553  nBin_(1),
554  binDir_(Zero),
555  binDx_(0.0),
556  binMin_(great),
557  binPoints_(),
558  binCumulative_(true),
559  initialised_(false)
560 {
561  read(dict);
562  resetNames(createFileNames(dict));
563 }
564 
565 
567 (
568  const word& name,
569  const objectRegistry& obr,
570  const dictionary& dict
571 )
572 :
573  fvMeshFunctionObject(name, obr, dict),
574  logFiles(obr_, name),
575  force_(3),
576  moment_(3),
577  patchSet_(),
578  pName_(word::null),
579  UName_(word::null),
580  rhoName_(word::null),
581  directForceDensity_(false),
582  fDName_(""),
583  rhoRef_(vGreat),
584  pRef_(0),
585  coordSys_(),
586  localSystem_(false),
587  porosity_(false),
588  nBin_(1),
589  binDir_(Zero),
590  binDx_(0.0),
591  binMin_(great),
592  binPoints_(),
593  binCumulative_(true),
594  initialised_(false)
595 {
596  read(dict);
597  resetNames(createFileNames(dict));
598 }
599 
600 
601 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
602 
604 {}
605 
606 
607 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
608 
610 {
612 
613  initialised_ = false;
614 
615  Log << type() << " " << name() << ":" << nl;
616 
617  directForceDensity_ = dict.lookupOrDefault("directForceDensity", false);
618 
619  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
620 
621  patchSet_ = pbm.patchSet(wordReList(dict.lookup("patches")));
622 
623  if (directForceDensity_)
624  {
625  // Optional entry for fDName
626  fDName_ = dict.lookupOrDefault<word>("fD", "fD");
627  }
628  else
629  {
630  // Optional entries U and p
631  pName_ = dict.lookupOrDefault<word>("p", "p");
632  UName_ = dict.lookupOrDefault<word>("U", "U");
633  rhoName_ = dict.lookupOrDefault<word>("rho", "rho");
634 
635  // Reference density needed for incompressible calculations
636  if (rhoName_ == "rhoInf")
637  {
638  dict.lookup("rhoInf") >> rhoRef_;
639  }
640 
641  // Reference pressure, 0 by default
642  pRef_ = dict.lookupOrDefault<scalar>("pRef", 0.0);
643  }
644 
645  coordSys_.clear();
646 
647  // Centre of rotation for moment calculations
648  // specified directly, from coordinate system, or implicitly (0 0 0)
649  if (!dict.readIfPresent<point>("CofR", coordSys_.origin()))
650  {
651  coordSys_ = coordinateSystem(obr_, dict);
652  localSystem_ = true;
653  }
654 
655  dict.readIfPresent("porosity", porosity_);
656  if (porosity_)
657  {
658  Log << " Including porosity effects" << endl;
659  }
660  else
661  {
662  Log << " Not including porosity effects" << endl;
663  }
664 
665  if (dict.found("binData"))
666  {
667  const dictionary& binDict(dict.subDict("binData"));
668  binDict.lookup("nBin") >> nBin_;
669 
670  if (nBin_ < 0)
671  {
673  << "Number of bins (nBin) must be zero or greater"
674  << exit(FatalIOError);
675  }
676  else if ((nBin_ == 0) || (nBin_ == 1))
677  {
678  nBin_ = 1;
679  forAll(force_, i)
680  {
681  force_[i].setSize(1);
682  moment_[i].setSize(1);
683  }
684  }
685 
686  if (nBin_ > 1)
687  {
688  binDict.lookup("direction") >> binDir_;
689  binDir_ /= mag(binDir_);
690 
691  binMin_ = great;
692  scalar binMax = -great;
693  forAllConstIter(labelHashSet, patchSet_, iter)
694  {
695  label patchi = iter.key();
696  const polyPatch& pp = pbm[patchi];
697  scalarField d(pp.faceCentres() & binDir_);
698  binMin_ = min(min(d), binMin_);
699  binMax = max(max(d), binMax);
700  }
701  reduce(binMin_, minOp<scalar>());
702  reduce(binMax, maxOp<scalar>());
703 
704  // slightly boost binMax so that region of interest is fully
705  // within bounds
706  binMax = 1.0001*(binMax - binMin_) + binMin_;
707 
708  binDx_ = (binMax - binMin_)/scalar(nBin_);
709 
710  // create the bin points used for writing
711  binPoints_.setSize(nBin_);
712  forAll(binPoints_, i)
713  {
714  binPoints_[i] = (i + 0.5)*binDir_*binDx_;
715  }
716 
717  binDict.lookup("cumulative") >> binCumulative_;
718 
719  // allocate storage for forces and moments
720  forAll(force_, i)
721  {
722  force_[i].setSize(nBin_);
723  moment_[i].setSize(nBin_);
724  }
725  }
726  }
727 
728  if (nBin_ == 1)
729  {
730  // allocate storage for forces and moments
731  force_[0].setSize(1);
732  force_[1].setSize(1);
733  force_[2].setSize(1);
734  moment_[0].setSize(1);
735  moment_[1].setSize(1);
736  moment_[2].setSize(1);
737  }
738 
739  return true;
740 }
741 
742 
744 {
745  initialise();
746 
747  force_[0] = Zero;
748  force_[1] = Zero;
749  force_[2] = Zero;
750 
751  moment_[0] = Zero;
752  moment_[1] = Zero;
753  moment_[2] = Zero;
754 
755  if (directForceDensity_)
756  {
757  const volVectorField& fD = obr_.lookupObject<volVectorField>(fDName_);
758 
759  const surfaceVectorField::Boundary& Sfb =
760  mesh_.Sf().boundaryField();
761 
762  forAllConstIter(labelHashSet, patchSet_, iter)
763  {
764  label patchi = iter.key();
765 
766  vectorField Md
767  (
768  mesh_.C().boundaryField()[patchi] - coordSys_.origin()
769  );
770 
771  scalarField sA(mag(Sfb[patchi]));
772 
773  // Normal force = surfaceUnitNormal*(surfaceNormal & forceDensity)
774  vectorField fN
775  (
776  Sfb[patchi]/sA
777  *(
778  Sfb[patchi] & fD.boundaryField()[patchi]
779  )
780  );
781 
782  // Tangential force (total force minus normal fN)
783  vectorField fT(sA*fD.boundaryField()[patchi] - fN);
784 
785  //- Porous force
786  vectorField fP(Md.size(), Zero);
787 
788  applyBins(Md, fN, fT, fP, mesh_.C().boundaryField()[patchi]);
789  }
790  }
791  else
792  {
793  const volScalarField& p = obr_.lookupObject<volScalarField>(pName_);
794 
795  const surfaceVectorField::Boundary& Sfb =
796  mesh_.Sf().boundaryField();
797 
798  tmp<volSymmTensorField> tdevRhoReff = devRhoReff();
799  const volSymmTensorField::Boundary& devRhoReffb
800  = tdevRhoReff().boundaryField();
801 
802  // Scale pRef by density for incompressible simulations
803  scalar pRef = pRef_/rho(p);
804 
805  forAllConstIter(labelHashSet, patchSet_, iter)
806  {
807  label patchi = iter.key();
808 
809  vectorField Md
810  (
811  mesh_.C().boundaryField()[patchi] - coordSys_.origin()
812  );
813 
814  vectorField fN
815  (
816  rho(p)*Sfb[patchi]*(p.boundaryField()[patchi] - pRef)
817  );
818 
819  vectorField fT(Sfb[patchi] & devRhoReffb[patchi]);
820 
821  vectorField fP(Md.size(), Zero);
822 
823  applyBins(Md, fN, fT, fP, mesh_.C().boundaryField()[patchi]);
824  }
825  }
826 
827  if (porosity_)
828  {
829  const volVectorField& U = obr_.lookupObject<volVectorField>(UName_);
830  const volScalarField rho(this->rho());
831  const volScalarField mu(this->mu());
832 
833  const HashTable<const porosityModel*> models =
834  obr_.lookupClass<porosityModel>();
835 
836  if (models.empty())
837  {
839  << "Porosity effects requested, but no porosity models found "
840  << "in the database"
841  << endl;
842  }
843 
845  {
846  // non-const access required if mesh is changing
847  porosityModel& pm = const_cast<porosityModel&>(*iter());
848 
849  vectorField fPTot(pm.force(U, rho, mu));
850 
851  const labelList& cellZoneIDs = pm.cellZoneIDs();
852 
853  forAll(cellZoneIDs, i)
854  {
855  label zoneI = cellZoneIDs[i];
856  const cellZone& cZone = mesh_.cellZones()[zoneI];
857 
858  const vectorField d(mesh_.C(), cZone);
859  const vectorField fP(fPTot, cZone);
860  const vectorField Md(d - coordSys_.origin());
861 
862  const vectorField fDummy(Md.size(), Zero);
863 
864  applyBins(Md, fDummy, fDummy, fP, d);
865  }
866  }
867  }
868 
873 }
874 
875 
877 {
878  return sum(force_[0]) + sum(force_[1]) + sum(force_[2]);
879 }
880 
881 
883 {
884  return sum(moment_[0]) + sum(moment_[1]) + sum(moment_[2]);
885 }
886 
887 
889 {
890  return true;
891 }
892 
893 
895 {
896  calcForcesMoment();
897 
898  if (Pstream::master())
899  {
900  logFiles::write();
901 
902  writeForces();
903 
904  writeBins();
905 
906  Log << endl;
907  }
908 
909  return true;
910 }
911 
912 
913 // ************************************************************************* //
virtual vector momentEff() const
Return the total moment.
Definition: forces.C:882
Base class for other coordinate system specifications.
virtual void writeFileHeader(const label i)
Output file header information.
Definition: forces.C:75
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:438
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual tmp< volScalarField > mu() const =0
Dynamic viscosity of mixture [kg/m/s].
virtual bool write()
Write function.
Definition: logFiles.C:180
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const char tab
Definition: Ostream.H:264
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
forces(const word &name, const Time &runTime, const dictionary &dict)
Construct from Time and dictionary.
Definition: forces.C:532
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
addToRunTimeSelectionTable(functionObject, Qdot, dictionary)
virtual tmp< vectorField > force(const volVectorField &U, const volScalarField &rho, const volScalarField &mu)
Return the force over the cell zone(s)
bool empty() const
Return true if the hash table is empty.
Definition: HashTableI.H:72
const dimensionSet dimViscosity
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
virtual vector forceEff() const
Return the total force.
Definition: forces.C:876
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
rhoReactionThermo & thermo
Definition: createFields.H:28
const labelList & cellZoneIDs() const
Return const access to the cell zone IDs.
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool usePatchGroups=true) const
Return the set of patch IDs corresponding to the given names.
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
virtual bool read(const dictionary &)
Read the forces data.
Definition: forces.C:609
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
Macros for easy insertion into run-time selection tables.
void writeBins()
Helper function to write bin data.
Definition: forces.C:450
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const word dictName() const
Return the local dictionary name (final part of scoped name)
Definition: dictionary.H:123
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, const Mesh &, const dimensionSet &, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return a temporary field constructed from name, mesh, dimensionSet.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:699
Templated abstract base class for single-phase incompressible turbulence models.
virtual bool execute()
Execute, currently does nothing.
Definition: forces.C:888
const dimensionSet & dimensions() const
Return dimensions.
virtual bool read(const dictionary &)
Read optional controls.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Calculate the gradient of the given field.
A class for handling words, derived from string.
Definition: word.H:59
tmp< volScalarField > mu() const
Dynamic viscosity field.
Definition: forces.C:286
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:49
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
static const word null
An empty word.
Definition: word.H:77
virtual void calcForcesMoment()
Calculate the forces and moments.
Definition: forces.C:743
static const zero Zero
Definition: zero.H:97
tmp< volScalarField > rho() const
Return rho if specified otherwise rhoRef.
Definition: forces.C:330
void applyBins(const vectorField &Md, const vectorField &fN, const vectorField &fT, const vectorField &fP, const vectorField &d)
Accumulate bin data.
Definition: forces.C:369
An STL-conforming hash table.
Definition: HashTable.H:61
const dimensionSet dimPressure
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
errorManip< error > abort(error &err)
Definition: errorManip.H:131
label readLabel(Istream &is)
Definition: label.H:64
Foam::polyBoundaryMesh.
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor.
static const char nl
Definition: Ostream.H:265
tmp< volSymmTensorField > devRhoReff() const
Return the effective viscous stress (laminar + turbulent).
Definition: forces.C:219
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
void initialise()
Initialise the fields.
Definition: forces.C:172
scalar pRef
Definition: createFields.H:19
labelList f(nPoints)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
A subset of mesh cells.
Definition: cellZone.H:61
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
defineTypeNameAndDebug(Qdot, 0)
const dimensionSet dimDensity
static const GeometricField< symmTensor, fvPatchField, volMesh > & null()
Return a null geometric field.
void writeForces()
Helper function to write force data.
Definition: forces.C:405
const dimensionedScalar mu
Atomic mass unit.
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
label patchi
U
Definition: pEqn.H:72
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Base-class for all transport models used by the incompressible turbulence models. ...
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
#define Log
Report write to Foam::Info if the local log switch is true.
dimensioned< scalar > mag(const dimensioned< Type > &)
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Omanip< int > setw(const int i)
Definition: IOmanip.H:199
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
volScalarField & p
A class for managing temporary objects.
Definition: PtrList.H:53
Registry of regIOobjects.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
virtual ~forces()
Destructor.
Definition: forces.C:603
fileID
Enumeration for ensuring the right file is accessed.
Definition: forces.H:206
volScalarField & nu
Top level model for porosity models.
Definition: porosityModel.H:55
virtual bool write()
Write the forces.
Definition: forces.C:894
wordList createFileNames(const dictionary &dict) const
Create file names for forces and bins.
Definition: forces.C:49
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:583
functionObject base class for creating, maintaining and writing log files e.g. integrated of averaged...
Definition: logFiles.H:57
IOerror FatalIOError
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:284