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 | Copyright (C) 2011-2016 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(MAIN_FILE=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(BINS_FILE=1)
67  names.append(forceType + "_bins");
68  }
69  }
70 
71  return names;
72 }
73 
74 
76 {
77  switch (fileID(i))
78  {
79  case MAIN_FILE:
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 BINS_FILE:
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  {
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 
264  dimensionedScalar nu(transportProperties.lookup("nu"));
265 
266  const volVectorField& U = obr_.lookupObject<volVectorField>(UName_);
267 
268  return -rho()*nu*dev(twoSymm(fvc::grad(U)));
269  }
270  else
271  {
273  << "No valid model for viscous stress calculation"
274  << exit(FatalError);
275 
276  return volSymmTensorField::null();
277  }
278 }
279 
280 
282 {
283  if (obr_.foundObject<fluidThermo>(basicThermo::dictName))
284  {
285  const fluidThermo& thermo =
286  obr_.lookupObject<fluidThermo>(basicThermo::dictName);
287 
288  return thermo.mu();
289  }
290  else if
291  (
292  obr_.foundObject<transportModel>("transportProperties")
293  )
294  {
295  const transportModel& laminarT =
296  obr_.lookupObject<transportModel>("transportProperties");
297 
298  return rho()*laminarT.nu();
299  }
300  else if (obr_.foundObject<dictionary>("transportProperties"))
301  {
302  const dictionary& transportProperties =
303  obr_.lookupObject<dictionary>("transportProperties");
304 
306  (
307  "nu",
308  dimViscosity,
309  transportProperties.lookup("nu")
310  );
311 
312  return rho()*nu;
313  }
314  else
315  {
317  << "No valid model for dynamic viscosity calculation"
318  << exit(FatalError);
319 
320  return volScalarField::null();
321  }
322 }
323 
324 
326 {
327  if (rhoName_ == "rhoInf")
328  {
329  return tmp<volScalarField>
330  (
331  new volScalarField
332  (
333  IOobject
334  (
335  "rho",
336  mesh_.time().timeName(),
337  mesh_
338  ),
339  mesh_,
340  dimensionedScalar("rho", dimDensity, rhoRef_)
341  )
342  );
343  }
344  else
345  {
346  return(obr_.lookupObject<volScalarField>(rhoName_));
347  }
348 }
349 
350 
352 {
353  if (p.dimensions() == dimPressure)
354  {
355  return 1.0;
356  }
357  else
358  {
359  if (rhoName_ != "rhoInf")
360  {
362  << "Dynamic pressure is expected but kinematic is provided."
363  << exit(FatalError);
364  }
365 
366  return rhoRef_;
367  }
368 }
369 
370 
372 (
373  const vectorField& Md,
374  const vectorField& fN,
375  const vectorField& fT,
376  const vectorField& fP,
377  const vectorField& d
378 )
379 {
380  if (nBin_ == 1)
381  {
382  force_[0][0] += sum(fN);
383  force_[1][0] += sum(fT);
384  force_[2][0] += sum(fP);
385  moment_[0][0] += sum(Md^fN);
386  moment_[1][0] += sum(Md^fT);
387  moment_[2][0] += sum(Md^fP);
388  }
389  else
390  {
391  scalarField dd((d & binDir_) - binMin_);
392 
393  forAll(dd, i)
394  {
395  label bini = min(max(floor(dd[i]/binDx_), 0), force_[0].size() - 1);
396 
397  force_[0][bini] += fN[i];
398  force_[1][bini] += fT[i];
399  force_[2][bini] += fP[i];
400  moment_[0][bini] += Md[i]^fN[i];
401  moment_[1][bini] += Md[i]^fT[i];
402  moment_[2][bini] += Md[i]^fP[i];
403  }
404  }
405 }
406 
407 
409 {
410  Log << type() << " " << name() << " write:" << nl
411  << " sum of forces:" << nl
412  << " pressure : " << sum(force_[0]) << nl
413  << " viscous : " << sum(force_[1]) << nl
414  << " porous : " << sum(force_[2]) << nl
415  << " sum of moments:" << nl
416  << " pressure : " << sum(moment_[0]) << nl
417  << " viscous : " << sum(moment_[1]) << nl
418  << " porous : " << sum(moment_[2])
419  << endl;
420 
421  writeTime(file(MAIN_FILE));
422  file(MAIN_FILE) << tab << setw(1) << '('
423  << sum(force_[0]) << setw(1) << ' '
424  << sum(force_[1]) << setw(1) << ' '
425  << sum(force_[2]) << setw(3) << ") ("
426  << sum(moment_[0]) << setw(1) << ' '
427  << sum(moment_[1]) << setw(1) << ' '
428  << sum(moment_[2]) << setw(1) << ')'
429  << endl;
430 
431  if (localSystem_)
432  {
433  vectorField localForceN(coordSys_.localVector(force_[0]));
434  vectorField localForceT(coordSys_.localVector(force_[1]));
435  vectorField localForceP(coordSys_.localVector(force_[2]));
436  vectorField localMomentN(coordSys_.localVector(moment_[0]));
437  vectorField localMomentT(coordSys_.localVector(moment_[1]));
438  vectorField localMomentP(coordSys_.localVector(moment_[2]));
439 
440  writeTime(file(MAIN_FILE));
441  file(MAIN_FILE) << tab << setw(1) << '('
442  << sum(localForceN) << setw(1) << ' '
443  << sum(localForceT) << setw(1) << ' '
444  << sum(localForceP) << setw(3) << ") ("
445  << sum(localMomentN) << setw(1) << ' '
446  << sum(localMomentT) << setw(1) << ' '
447  << sum(localMomentP) << setw(1) << ')'
448  << endl;
449  }
450 }
451 
452 
454 {
455  if (nBin_ == 1)
456  {
457  return;
458  }
459 
460  List<Field<vector>> f(force_);
461  List<Field<vector>> m(moment_);
462 
463  if (binCumulative_)
464  {
465  for (label i = 1; i < f[0].size(); i++)
466  {
467  f[0][i] += f[0][i-1];
468  f[1][i] += f[1][i-1];
469  f[2][i] += f[2][i-1];
470 
471  m[0][i] += m[0][i-1];
472  m[1][i] += m[1][i-1];
473  m[2][i] += m[2][i-1];
474  }
475  }
476 
477  writeTime(file(BINS_FILE));
478 
479  forAll(f[0], i)
480  {
481  file(BINS_FILE)
482  << tab << setw(1) << '('
483  << f[0][i] << setw(1) << ' '
484  << f[1][i] << setw(1) << ' '
485  << f[2][i] << setw(3) << ") ("
486  << m[0][i] << setw(1) << ' '
487  << m[1][i] << setw(1) << ' '
488  << m[2][i] << setw(1) << ')';
489  }
490 
491  if (localSystem_)
492  {
493  List<Field<vector>> lf(3);
494  List<Field<vector>> lm(3);
495  lf[0] = coordSys_.localVector(force_[0]);
496  lf[1] = coordSys_.localVector(force_[1]);
497  lf[2] = coordSys_.localVector(force_[2]);
498  lm[0] = coordSys_.localVector(moment_[0]);
499  lm[1] = coordSys_.localVector(moment_[1]);
500  lm[2] = coordSys_.localVector(moment_[2]);
501 
502  if (binCumulative_)
503  {
504  for (label i = 1; i < lf[0].size(); i++)
505  {
506  lf[0][i] += lf[0][i-1];
507  lf[1][i] += lf[1][i-1];
508  lf[2][i] += lf[2][i-1];
509  lm[0][i] += lm[0][i-1];
510  lm[1][i] += lm[1][i-1];
511  lm[2][i] += lm[2][i-1];
512  }
513  }
514 
515  forAll(lf[0], i)
516  {
517  file(BINS_FILE)
518  << tab << setw(1) << '('
519  << lf[0][i] << setw(1) << ' '
520  << lf[1][i] << setw(1) << ' '
521  << lf[2][i] << setw(3) << ") ("
522  << lm[0][i] << setw(1) << ' '
523  << lm[1][i] << setw(1) << ' '
524  << lm[2][i] << setw(1) << ')';
525  }
526  }
527 
528  file(BINS_FILE) << endl;
529 }
530 
531 
532 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
533 
535 (
536  const word& name,
537  const Time& runTime,
538  const dictionary& dict
539 )
540 :
541  fvMeshFunctionObject(name, runTime, dict),
542  logFiles(obr_, name),
543  force_(3),
544  moment_(3),
545  patchSet_(),
546  pName_(word::null),
547  UName_(word::null),
548  rhoName_(word::null),
549  directForceDensity_(false),
550  fDName_(""),
551  rhoRef_(VGREAT),
552  pRef_(0),
553  coordSys_(),
554  localSystem_(false),
555  porosity_(false),
556  nBin_(1),
557  binDir_(Zero),
558  binDx_(0.0),
559  binMin_(GREAT),
560  binPoints_(),
561  binCumulative_(true),
562  initialised_(false)
563 {
564  read(dict);
565  resetNames(createFileNames(dict));
566 }
567 
568 
570 (
571  const word& name,
572  const objectRegistry& obr,
573  const dictionary& dict
574 )
575 :
576  fvMeshFunctionObject(name, obr, dict),
577  logFiles(obr_, name),
578  force_(3),
579  moment_(3),
580  patchSet_(),
581  pName_(word::null),
582  UName_(word::null),
583  rhoName_(word::null),
584  directForceDensity_(false),
585  fDName_(""),
586  rhoRef_(VGREAT),
587  pRef_(0),
588  coordSys_(),
589  localSystem_(false),
590  porosity_(false),
591  nBin_(1),
592  binDir_(Zero),
593  binDx_(0.0),
594  binMin_(GREAT),
595  binPoints_(),
596  binCumulative_(true),
597  initialised_(false)
598 {
599  read(dict);
600  resetNames(createFileNames(dict));
601 }
602 
603 
604 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
605 
607 {}
608 
609 
610 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
611 
613 {
615 
616  initialised_ = false;
617 
618  Log << type() << " " << name() << ":" << nl;
619 
620  directForceDensity_ = dict.lookupOrDefault("directForceDensity", false);
621 
622  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
623 
624  patchSet_ = pbm.patchSet(wordReList(dict.lookup("patches")));
625 
626  if (directForceDensity_)
627  {
628  // Optional entry for fDName
629  fDName_ = dict.lookupOrDefault<word>("fD", "fD");
630  }
631  else
632  {
633  // Optional entries U and p
634  pName_ = dict.lookupOrDefault<word>("p", "p");
635  UName_ = dict.lookupOrDefault<word>("U", "U");
636  rhoName_ = dict.lookupOrDefault<word>("rho", "rho");
637 
638  // Reference density needed for incompressible calculations
639  if (rhoName_ == "rhoInf")
640  {
641  dict.lookup("rhoInf") >> rhoRef_;
642  }
643 
644  // Reference pressure, 0 by default
645  pRef_ = dict.lookupOrDefault<scalar>("pRef", 0.0);
646  }
647 
648  coordSys_.clear();
649 
650  // Centre of rotation for moment calculations
651  // specified directly, from coordinate system, or implicitly (0 0 0)
652  if (!dict.readIfPresent<point>("CofR", coordSys_.origin()))
653  {
654  coordSys_ = coordinateSystem(obr_, dict);
655  localSystem_ = true;
656  }
657 
658  dict.readIfPresent("porosity", porosity_);
659  if (porosity_)
660  {
661  Log << " Including porosity effects" << endl;
662  }
663  else
664  {
665  Log << " Not including porosity effects" << endl;
666  }
667 
668  if (dict.found("binData"))
669  {
670  const dictionary& binDict(dict.subDict("binData"));
671  binDict.lookup("nBin") >> nBin_;
672 
673  if (nBin_ < 0)
674  {
676  << "Number of bins (nBin) must be zero or greater"
677  << exit(FatalIOError);
678  }
679  else if ((nBin_ == 0) || (nBin_ == 1))
680  {
681  nBin_ = 1;
682  forAll(force_, i)
683  {
684  force_[i].setSize(1);
685  moment_[i].setSize(1);
686  }
687  }
688 
689  if (nBin_ > 1)
690  {
691  binDict.lookup("direction") >> binDir_;
692  binDir_ /= mag(binDir_);
693 
694  binMin_ = GREAT;
695  scalar binMax = -GREAT;
696  forAllConstIter(labelHashSet, patchSet_, iter)
697  {
698  label patchi = iter.key();
699  const polyPatch& pp = pbm[patchi];
700  scalarField d(pp.faceCentres() & binDir_);
701  binMin_ = min(min(d), binMin_);
702  binMax = max(max(d), binMax);
703  }
704  reduce(binMin_, minOp<scalar>());
705  reduce(binMax, maxOp<scalar>());
706 
707  // slightly boost binMax so that region of interest is fully
708  // within bounds
709  binMax = 1.0001*(binMax - binMin_) + binMin_;
710 
711  binDx_ = (binMax - binMin_)/scalar(nBin_);
712 
713  // create the bin points used for writing
714  binPoints_.setSize(nBin_);
715  forAll(binPoints_, i)
716  {
717  binPoints_[i] = (i + 0.5)*binDir_*binDx_;
718  }
719 
720  binDict.lookup("cumulative") >> binCumulative_;
721 
722  // allocate storage for forces and moments
723  forAll(force_, i)
724  {
725  force_[i].setSize(nBin_);
726  moment_[i].setSize(nBin_);
727  }
728  }
729  }
730 
731  if (nBin_ == 1)
732  {
733  // allocate storage for forces and moments
734  force_[0].setSize(1);
735  force_[1].setSize(1);
736  force_[2].setSize(1);
737  moment_[0].setSize(1);
738  moment_[1].setSize(1);
739  moment_[2].setSize(1);
740  }
741 
742  return true;
743 }
744 
745 
747 {
748  initialise();
749 
750  force_[0] = Zero;
751  force_[1] = Zero;
752  force_[2] = Zero;
753 
754  moment_[0] = Zero;
755  moment_[1] = Zero;
756  moment_[2] = Zero;
757 
758  if (directForceDensity_)
759  {
760  const volVectorField& fD = obr_.lookupObject<volVectorField>(fDName_);
761 
762  const surfaceVectorField::Boundary& Sfb =
763  mesh_.Sf().boundaryField();
764 
765  forAllConstIter(labelHashSet, patchSet_, iter)
766  {
767  label patchi = iter.key();
768 
769  vectorField Md
770  (
771  mesh_.C().boundaryField()[patchi] - coordSys_.origin()
772  );
773 
774  scalarField sA(mag(Sfb[patchi]));
775 
776  // Normal force = surfaceUnitNormal*(surfaceNormal & forceDensity)
777  vectorField fN
778  (
779  Sfb[patchi]/sA
780  *(
781  Sfb[patchi] & fD.boundaryField()[patchi]
782  )
783  );
784 
785  // Tangential force (total force minus normal fN)
786  vectorField fT(sA*fD.boundaryField()[patchi] - fN);
787 
788  //- Porous force
789  vectorField fP(Md.size(), Zero);
790 
791  applyBins(Md, fN, fT, fP, mesh_.C().boundaryField()[patchi]);
792  }
793  }
794  else
795  {
796  const volScalarField& p = obr_.lookupObject<volScalarField>(pName_);
797 
798  const surfaceVectorField::Boundary& Sfb =
799  mesh_.Sf().boundaryField();
800 
801  tmp<volSymmTensorField> tdevRhoReff = devRhoReff();
802  const volSymmTensorField::Boundary& devRhoReffb
803  = tdevRhoReff().boundaryField();
804 
805  // Scale pRef by density for incompressible simulations
806  scalar pRef = pRef_/rho(p);
807 
808  forAllConstIter(labelHashSet, patchSet_, iter)
809  {
810  label patchi = iter.key();
811 
812  vectorField Md
813  (
814  mesh_.C().boundaryField()[patchi] - coordSys_.origin()
815  );
816 
817  vectorField fN
818  (
819  rho(p)*Sfb[patchi]*(p.boundaryField()[patchi] - pRef)
820  );
821 
822  vectorField fT(Sfb[patchi] & devRhoReffb[patchi]);
823 
824  vectorField fP(Md.size(), Zero);
825 
826  applyBins(Md, fN, fT, fP, mesh_.C().boundaryField()[patchi]);
827  }
828  }
829 
830  if (porosity_)
831  {
832  const volVectorField& U = obr_.lookupObject<volVectorField>(UName_);
833  const volScalarField rho(this->rho());
834  const volScalarField mu(this->mu());
835 
836  const HashTable<const porosityModel*> models =
837  obr_.lookupClass<porosityModel>();
838 
839  if (models.empty())
840  {
842  << "Porosity effects requested, but no porosity models found "
843  << "in the database"
844  << endl;
845  }
846 
848  {
849  // non-const access required if mesh is changing
850  porosityModel& pm = const_cast<porosityModel&>(*iter());
851 
852  vectorField fPTot(pm.force(U, rho, mu));
853 
854  const labelList& cellZoneIDs = pm.cellZoneIDs();
855 
856  forAll(cellZoneIDs, i)
857  {
858  label zoneI = cellZoneIDs[i];
859  const cellZone& cZone = mesh_.cellZones()[zoneI];
860 
861  const vectorField d(mesh_.C(), cZone);
862  const vectorField fP(fPTot, cZone);
863  const vectorField Md(d - coordSys_.origin());
864 
865  const vectorField fDummy(Md.size(), Zero);
866 
867  applyBins(Md, fDummy, fDummy, fP, d);
868  }
869  }
870  }
871 
876 }
877 
878 
880 {
881  return sum(force_[0]) + sum(force_[1]) + sum(force_[2]);
882 }
883 
884 
886 {
887  return sum(moment_[0]) + sum(moment_[1]) + sum(moment_[2]);
888 }
889 
890 
892 {
893  return true;
894 }
895 
896 
898 {
899  calcForcesMoment();
900 
901  if (Pstream::master())
902  {
903  logFiles::write();
904 
905  writeForces();
906 
907  writeBins();
908 
909  Log << endl;
910  }
911 
912  return true;
913 }
914 
915 
916 // ************************************************************************* //
virtual vector momentEff() const
Return the total moment.
Definition: forces.C:885
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:431
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:428
virtual tmp< volScalarField > mu() const =0
Dynamic viscosity of mixture [kg/m/s].
U
Definition: pEqn.H:83
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:261
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
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
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
compressible::turbulenceModel & turb
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:412
virtual vector forceEff() const
Return the total force.
Definition: forces.C:879
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
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:612
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:453
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:115
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:692
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Templated abstract base class for single-phase incompressible turbulence models.
virtual bool execute()
Execute, currently does nothing.
Definition: forces.C:891
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:281
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:292
static const word null
An empty word.
Definition: word.H:77
virtual void calcForcesMoment()
Calculate the forces and moments.
Definition: forces.C:746
static const zero Zero
Definition: zero.H:91
tmp< volScalarField > rho() const
Return rho if specified otherwise rhoRef.
Definition: forces.C:325
void applyBins(const vectorField &Md, const vectorField &fN, const vectorField &fT, const vectorField &fP, const vectorField &d)
Accumulate bin data.
Definition: forces.C:372
An STL-conforming hash table.
Definition: HashTable.H:62
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
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::polyBoundaryMesh.
forces(const forces &)
Disallow default bitwise copy construct.
virtual tmp< volSymmTensorField > devReff() const
Return the effective stress tensor.
static const char nl
Definition: Ostream.H:262
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)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:485
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
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:408
const dimensionedScalar mu
Atomic mass unit.
label patchi
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. ...
#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
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
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:606
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
fileID
Enumeration for ensuring the right file is accessed.
Definition: forces.H:209
volScalarField & nu
Top level model for porosity models.
Definition: porosityModel.H:55
virtual bool write()
Write the forces.
Definition: forces.C:897
wordList createFileNames(const dictionary &dict) const
Create file names for forces and bins.
Definition: forces.C:49
addToRunTimeSelectionTable(functionObject, add, dictionary)
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:576
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