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-2018 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  {
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 tmp<volScalarField>
335  (
336  new volScalarField
337  (
338  IOobject
339  (
340  "rho",
341  mesh_.time().timeName(),
342  mesh_
343  ),
344  mesh_,
345  dimensionedScalar("rho", dimDensity, rhoRef_)
346  )
347  );
348  }
349  else
350  {
351  return(obr_.lookupObject<volScalarField>(rhoName_));
352  }
353 }
354 
355 
357 {
358  if (p.dimensions() == dimPressure)
359  {
360  return 1.0;
361  }
362  else
363  {
364  if (rhoName_ != "rhoInf")
365  {
367  << "Dynamic pressure is expected but kinematic is provided."
368  << exit(FatalError);
369  }
370 
371  return rhoRef_;
372  }
373 }
374 
375 
377 (
378  const vectorField& Md,
379  const vectorField& fN,
380  const vectorField& fT,
381  const vectorField& fP,
382  const vectorField& d
383 )
384 {
385  if (nBin_ == 1)
386  {
387  force_[0][0] += sum(fN);
388  force_[1][0] += sum(fT);
389  force_[2][0] += sum(fP);
390  moment_[0][0] += sum(Md^fN);
391  moment_[1][0] += sum(Md^fT);
392  moment_[2][0] += sum(Md^fP);
393  }
394  else
395  {
396  scalarField dd((d & binDir_) - binMin_);
397 
398  forAll(dd, i)
399  {
400  label bini = min(max(floor(dd[i]/binDx_), 0), force_[0].size() - 1);
401 
402  force_[0][bini] += fN[i];
403  force_[1][bini] += fT[i];
404  force_[2][bini] += fP[i];
405  moment_[0][bini] += Md[i]^fN[i];
406  moment_[1][bini] += Md[i]^fT[i];
407  moment_[2][bini] += Md[i]^fP[i];
408  }
409  }
410 }
411 
412 
414 {
415  Log << type() << " " << name() << " write:" << nl
416  << " sum of forces:" << nl
417  << " pressure : " << sum(force_[0]) << nl
418  << " viscous : " << sum(force_[1]) << nl
419  << " porous : " << sum(force_[2]) << nl
420  << " sum of moments:" << nl
421  << " pressure : " << sum(moment_[0]) << nl
422  << " viscous : " << sum(moment_[1]) << nl
423  << " porous : " << sum(moment_[2])
424  << endl;
425 
426  writeTime(file(MAIN_FILE));
427 
428  file(MAIN_FILE) << tab << setw(1) << '('
429  << sum(force_[0]) << setw(1) << ' '
430  << sum(force_[1]) << setw(1) << ' '
431  << sum(force_[2]) << setw(3) << ") ("
432  << sum(moment_[0]) << setw(1) << ' '
433  << sum(moment_[1]) << setw(1) << ' '
434  << sum(moment_[2]) << setw(1) << ')';
435 
436  if (localSystem_)
437  {
438  vectorField localForceN(coordSys_.localVector(force_[0]));
439  vectorField localForceT(coordSys_.localVector(force_[1]));
440  vectorField localForceP(coordSys_.localVector(force_[2]));
441  vectorField localMomentN(coordSys_.localVector(moment_[0]));
442  vectorField localMomentT(coordSys_.localVector(moment_[1]));
443  vectorField localMomentP(coordSys_.localVector(moment_[2]));
444 
445  file(MAIN_FILE) << tab << setw(1) << '('
446  << sum(localForceN) << setw(1) << ' '
447  << sum(localForceT) << setw(1) << ' '
448  << sum(localForceP) << setw(3) << ") ("
449  << sum(localMomentN) << setw(1) << ' '
450  << sum(localMomentT) << setw(1) << ' '
451  << sum(localMomentP) << setw(1) << ')';
452  }
453 
454  file(MAIN_FILE) << endl;
455 }
456 
457 
459 {
460  if (nBin_ == 1)
461  {
462  return;
463  }
464 
465  List<Field<vector>> f(force_);
466  List<Field<vector>> m(moment_);
467 
468  if (binCumulative_)
469  {
470  for (label i = 1; i < f[0].size(); i++)
471  {
472  f[0][i] += f[0][i-1];
473  f[1][i] += f[1][i-1];
474  f[2][i] += f[2][i-1];
475 
476  m[0][i] += m[0][i-1];
477  m[1][i] += m[1][i-1];
478  m[2][i] += m[2][i-1];
479  }
480  }
481 
482  writeTime(file(BINS_FILE));
483 
484  forAll(f[0], i)
485  {
486  file(BINS_FILE)
487  << tab << setw(1) << '('
488  << f[0][i] << setw(1) << ' '
489  << f[1][i] << setw(1) << ' '
490  << f[2][i] << setw(3) << ") ("
491  << m[0][i] << setw(1) << ' '
492  << m[1][i] << setw(1) << ' '
493  << m[2][i] << setw(1) << ')';
494  }
495 
496  if (localSystem_)
497  {
498  List<Field<vector>> lf(3);
499  List<Field<vector>> lm(3);
500  lf[0] = coordSys_.localVector(force_[0]);
501  lf[1] = coordSys_.localVector(force_[1]);
502  lf[2] = coordSys_.localVector(force_[2]);
503  lm[0] = coordSys_.localVector(moment_[0]);
504  lm[1] = coordSys_.localVector(moment_[1]);
505  lm[2] = coordSys_.localVector(moment_[2]);
506 
507  if (binCumulative_)
508  {
509  for (label i = 1; i < lf[0].size(); i++)
510  {
511  lf[0][i] += lf[0][i-1];
512  lf[1][i] += lf[1][i-1];
513  lf[2][i] += lf[2][i-1];
514  lm[0][i] += lm[0][i-1];
515  lm[1][i] += lm[1][i-1];
516  lm[2][i] += lm[2][i-1];
517  }
518  }
519 
520  forAll(lf[0], i)
521  {
522  file(BINS_FILE)
523  << tab << setw(1) << '('
524  << lf[0][i] << setw(1) << ' '
525  << lf[1][i] << setw(1) << ' '
526  << lf[2][i] << setw(3) << ") ("
527  << lm[0][i] << setw(1) << ' '
528  << lm[1][i] << setw(1) << ' '
529  << lm[2][i] << setw(1) << ')';
530  }
531  }
532 
533  file(BINS_FILE) << endl;
534 }
535 
536 
537 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
538 
540 (
541  const word& name,
542  const Time& runTime,
543  const dictionary& dict
544 )
545 :
546  fvMeshFunctionObject(name, runTime, dict),
547  logFiles(obr_, name),
548  force_(3),
549  moment_(3),
550  patchSet_(),
551  pName_(word::null),
552  UName_(word::null),
553  rhoName_(word::null),
554  directForceDensity_(false),
555  fDName_(""),
556  rhoRef_(vGreat),
557  pRef_(0),
558  coordSys_(),
559  localSystem_(false),
560  porosity_(false),
561  nBin_(1),
562  binDir_(Zero),
563  binDx_(0.0),
564  binMin_(great),
565  binPoints_(),
566  binCumulative_(true),
567  initialised_(false)
568 {
569  read(dict);
570  resetNames(createFileNames(dict));
571 }
572 
573 
575 (
576  const word& name,
577  const objectRegistry& obr,
578  const dictionary& dict
579 )
580 :
581  fvMeshFunctionObject(name, obr, dict),
582  logFiles(obr_, name),
583  force_(3),
584  moment_(3),
585  patchSet_(),
586  pName_(word::null),
587  UName_(word::null),
588  rhoName_(word::null),
589  directForceDensity_(false),
590  fDName_(""),
591  rhoRef_(vGreat),
592  pRef_(0),
593  coordSys_(),
594  localSystem_(false),
595  porosity_(false),
596  nBin_(1),
597  binDir_(Zero),
598  binDx_(0.0),
599  binMin_(great),
600  binPoints_(),
601  binCumulative_(true),
602  initialised_(false)
603 {
604  read(dict);
605  resetNames(createFileNames(dict));
606 }
607 
608 
609 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
610 
612 {}
613 
614 
615 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
616 
618 {
620 
621  initialised_ = false;
622 
623  Log << type() << " " << name() << ":" << nl;
624 
625  directForceDensity_ = dict.lookupOrDefault("directForceDensity", false);
626 
627  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
628 
629  patchSet_ = pbm.patchSet(wordReList(dict.lookup("patches")));
630 
631  if (directForceDensity_)
632  {
633  // Optional entry for fDName
634  fDName_ = dict.lookupOrDefault<word>("fD", "fD");
635  }
636  else
637  {
638  // Optional entries U and p
639  pName_ = dict.lookupOrDefault<word>("p", "p");
640  UName_ = dict.lookupOrDefault<word>("U", "U");
641  rhoName_ = dict.lookupOrDefault<word>("rho", "rho");
642 
643  // Reference density needed for incompressible calculations
644  if (rhoName_ == "rhoInf")
645  {
646  dict.lookup("rhoInf") >> rhoRef_;
647  }
648 
649  // Reference pressure, 0 by default
650  pRef_ = dict.lookupOrDefault<scalar>("pRef", 0.0);
651  }
652 
653  coordSys_.clear();
654 
655  // Centre of rotation for moment calculations
656  // specified directly, from coordinate system, or implicitly (0 0 0)
657  if (!dict.readIfPresent<point>("CofR", coordSys_.origin()))
658  {
659  coordSys_ = coordinateSystem(obr_, dict);
660  localSystem_ = true;
661  }
662 
663  dict.readIfPresent("porosity", porosity_);
664  if (porosity_)
665  {
666  Log << " Including porosity effects" << endl;
667  }
668  else
669  {
670  Log << " Not including porosity effects" << endl;
671  }
672 
673  if (dict.found("binData"))
674  {
675  const dictionary& binDict(dict.subDict("binData"));
676  binDict.lookup("nBin") >> nBin_;
677 
678  if (nBin_ < 0)
679  {
681  << "Number of bins (nBin) must be zero or greater"
682  << exit(FatalIOError);
683  }
684  else if ((nBin_ == 0) || (nBin_ == 1))
685  {
686  nBin_ = 1;
687  forAll(force_, i)
688  {
689  force_[i].setSize(1);
690  moment_[i].setSize(1);
691  }
692  }
693 
694  if (nBin_ > 1)
695  {
696  binDict.lookup("direction") >> binDir_;
697  binDir_ /= mag(binDir_);
698 
699  binMin_ = great;
700  scalar binMax = -great;
701  forAllConstIter(labelHashSet, patchSet_, iter)
702  {
703  label patchi = iter.key();
704  const polyPatch& pp = pbm[patchi];
705  scalarField d(pp.faceCentres() & binDir_);
706  binMin_ = min(min(d), binMin_);
707  binMax = max(max(d), binMax);
708  }
709  reduce(binMin_, minOp<scalar>());
710  reduce(binMax, maxOp<scalar>());
711 
712  // slightly boost binMax so that region of interest is fully
713  // within bounds
714  binMax = 1.0001*(binMax - binMin_) + binMin_;
715 
716  binDx_ = (binMax - binMin_)/scalar(nBin_);
717 
718  // create the bin points used for writing
719  binPoints_.setSize(nBin_);
720  forAll(binPoints_, i)
721  {
722  binPoints_[i] = (i + 0.5)*binDir_*binDx_;
723  }
724 
725  binDict.lookup("cumulative") >> binCumulative_;
726 
727  // allocate storage for forces and moments
728  forAll(force_, i)
729  {
730  force_[i].setSize(nBin_);
731  moment_[i].setSize(nBin_);
732  }
733  }
734  }
735 
736  if (nBin_ == 1)
737  {
738  // allocate storage for forces and moments
739  force_[0].setSize(1);
740  force_[1].setSize(1);
741  force_[2].setSize(1);
742  moment_[0].setSize(1);
743  moment_[1].setSize(1);
744  moment_[2].setSize(1);
745  }
746 
747  return true;
748 }
749 
750 
752 {
753  initialise();
754 
755  force_[0] = Zero;
756  force_[1] = Zero;
757  force_[2] = Zero;
758 
759  moment_[0] = Zero;
760  moment_[1] = Zero;
761  moment_[2] = Zero;
762 
763  if (directForceDensity_)
764  {
765  const volVectorField& fD = obr_.lookupObject<volVectorField>(fDName_);
766 
767  const surfaceVectorField::Boundary& Sfb =
768  mesh_.Sf().boundaryField();
769 
770  forAllConstIter(labelHashSet, patchSet_, iter)
771  {
772  label patchi = iter.key();
773 
774  vectorField Md
775  (
776  mesh_.C().boundaryField()[patchi] - coordSys_.origin()
777  );
778 
779  scalarField sA(mag(Sfb[patchi]));
780 
781  // Normal force = surfaceUnitNormal*(surfaceNormal & forceDensity)
782  vectorField fN
783  (
784  Sfb[patchi]/sA
785  *(
786  Sfb[patchi] & fD.boundaryField()[patchi]
787  )
788  );
789 
790  // Tangential force (total force minus normal fN)
791  vectorField fT(sA*fD.boundaryField()[patchi] - fN);
792 
793  //- Porous force
794  vectorField fP(Md.size(), Zero);
795 
796  applyBins(Md, fN, fT, fP, mesh_.C().boundaryField()[patchi]);
797  }
798  }
799  else
800  {
801  const volScalarField& p = obr_.lookupObject<volScalarField>(pName_);
802 
803  const surfaceVectorField::Boundary& Sfb =
804  mesh_.Sf().boundaryField();
805 
806  tmp<volSymmTensorField> tdevRhoReff = devRhoReff();
807  const volSymmTensorField::Boundary& devRhoReffb
808  = tdevRhoReff().boundaryField();
809 
810  // Scale pRef by density for incompressible simulations
811  scalar pRef = pRef_/rho(p);
812 
813  forAllConstIter(labelHashSet, patchSet_, iter)
814  {
815  label patchi = iter.key();
816 
817  vectorField Md
818  (
819  mesh_.C().boundaryField()[patchi] - coordSys_.origin()
820  );
821 
822  vectorField fN
823  (
824  rho(p)*Sfb[patchi]*(p.boundaryField()[patchi] - pRef)
825  );
826 
827  vectorField fT(Sfb[patchi] & devRhoReffb[patchi]);
828 
829  vectorField fP(Md.size(), Zero);
830 
831  applyBins(Md, fN, fT, fP, mesh_.C().boundaryField()[patchi]);
832  }
833  }
834 
835  if (porosity_)
836  {
837  const volVectorField& U = obr_.lookupObject<volVectorField>(UName_);
838  const volScalarField rho(this->rho());
839  const volScalarField mu(this->mu());
840 
841  const HashTable<const porosityModel*> models =
842  obr_.lookupClass<porosityModel>();
843 
844  if (models.empty())
845  {
847  << "Porosity effects requested, but no porosity models found "
848  << "in the database"
849  << endl;
850  }
851 
853  {
854  // non-const access required if mesh is changing
855  porosityModel& pm = const_cast<porosityModel&>(*iter());
856 
857  vectorField fPTot(pm.force(U, rho, mu));
858 
859  const labelList& cellZoneIDs = pm.cellZoneIDs();
860 
861  forAll(cellZoneIDs, i)
862  {
863  label zoneI = cellZoneIDs[i];
864  const cellZone& cZone = mesh_.cellZones()[zoneI];
865 
866  const vectorField d(mesh_.C(), cZone);
867  const vectorField fP(fPTot, cZone);
868  const vectorField Md(d - coordSys_.origin());
869 
870  const vectorField fDummy(Md.size(), Zero);
871 
872  applyBins(Md, fDummy, fDummy, fP, d);
873  }
874  }
875  }
876 
881 }
882 
883 
885 {
886  return sum(force_[0]) + sum(force_[1]) + sum(force_[2]);
887 }
888 
889 
891 {
892  return sum(moment_[0]) + sum(moment_[1]) + sum(moment_[2]);
893 }
894 
895 
897 {
898  return true;
899 }
900 
901 
903 {
904  calcForcesMoment();
905 
906  if (Pstream::master())
907  {
908  logFiles::write();
909 
910  writeForces();
911 
912  writeBins();
913 
914  Log << endl;
915  }
916 
917  return true;
918 }
919 
920 
921 // ************************************************************************* //
virtual vector momentEff() const
Return the total moment.
Definition: forces.C:890
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].
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: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.
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:60
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
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:421
virtual vector forceEff() const
Return the total force.
Definition: forces.C:884
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:617
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:458
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:896
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:292
static const word null
An empty word.
Definition: word.H:77
virtual void calcForcesMoment()
Calculate the forces and moments.
Definition: forces.C:751
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:377
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
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: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)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
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:413
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. ...
#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:611
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:206
volScalarField & nu
Top level model for porosity models.
Definition: porosityModel.H:55
virtual bool write()
Write the forces.
Definition: forces.C:902
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