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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "forces.H"
30 #include "fvcGrad.H"
31 #include "porosityModel.H"
35 #include "cartesianCS.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace functionObjects
42 {
43  defineTypeNameAndDebug(forces, 0);
44  addToRunTimeSelectionTable(functionObject, forces, dictionary);
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
50 
52 {
53  return this->name() + ":" + name;
54 }
55 
56 
58 {
59  // Note: Only possible to create bin files after bins have been initialised
60 
61  if (writeToFile() && !forceFilePtr_.valid())
62  {
63  forceFilePtr_ = createFile("force");
64  writeIntegratedHeader("Force", forceFilePtr_());
65  momentFilePtr_ = createFile("moment");
66  writeIntegratedHeader("Moment", momentFilePtr_());
67 
68  if (nBin_ > 1)
69  {
70  forceBinFilePtr_ = createFile("forceBin");
71  writeBinHeader("Force", forceBinFilePtr_());
72  momentBinFilePtr_ = createFile("momentBin");
73  writeBinHeader("Moment", momentBinFilePtr_());
74  }
75  }
76 }
77 
78 
80 (
81  const word& header,
82  Ostream& os
83 ) const
84 {
85  writeHeader(os, header);
86  writeHeaderValue(os, "CofR", coordSys_.origin());
87  writeHeader(os, "");
88  writeCommented(os, "Time");
89  writeTabbed(os, "(total_x total_y total_z)");
90  writeTabbed(os, "(pressure_x pressure_y pressure_z)");
91  writeTabbed(os, "(viscous_x viscous_y viscous_z)");
92 
93  if (porosity_)
94  {
95  writeTabbed(os, "(porous_x porous_y porous_z)");
96  }
97 
98  os << endl;
99 }
100 
101 
103 (
104  const word& header,
105  Ostream& os
106 ) const
107 {
108  writeHeader(os, header + " bins");
109  writeHeaderValue(os, "bins", nBin_);
110  writeHeaderValue(os, "start", binMin_);
111  writeHeaderValue(os, "end", binMax_);
112  writeHeaderValue(os, "delta", binDx_);
113  writeHeaderValue(os, "direction", binDir_);
114 
115  vectorField binPoints(nBin_);
116  writeCommented(os, "x co-ords :");
117  forAll(binPoints, pointi)
118  {
119  binPoints[pointi] = (binMin_ + (pointi + 1)*binDx_)*binDir_;
120  os << tab << binPoints[pointi].x();
121  }
122  os << nl;
123 
124  writeCommented(os, "y co-ords :");
125  forAll(binPoints, pointi)
126  {
127  os << tab << binPoints[pointi].y();
128  }
129  os << nl;
130 
131  writeCommented(os, "z co-ords :");
132  forAll(binPoints, pointi)
133  {
134  os << tab << binPoints[pointi].z();
135  }
136  os << nl;
137 
138  writeHeader(os, "");
139  writeCommented(os, "Time");
140 
141  for (label j = 0; j < nBin_; j++)
142  {
143  const word jn(Foam::name(j) + ':');
144  os << tab << jn << "(total_x total_y total_z)"
145  << tab << jn << "(pressure_x pressure_y pressure_z)"
146  << tab << jn << "(viscous_x viscous_y viscous_z)";
147 
148  if (porosity_)
149  {
150  os << tab << jn << "(porous_x porous_y porous_z)";
151  }
152  }
153 
154  os << endl;
155 }
156 
157 
159 (
160  const dictionary& dict,
161  const word& e3Name,
162  const word& e1Name
163 )
164 {
165  coordSys_.clear();
166 
167  if (dict.readIfPresent<point>("CofR", coordSys_.origin()))
168  {
169  const vector e3 = e3Name == word::null ?
170  vector(0, 0, 1) : dict.get<vector>(e3Name);
171  const vector e1 = e1Name == word::null ?
172  vector(1, 0, 0) : dict.get<vector>(e1Name);
173 
174  coordSys_ =
175  coordSystem::cartesian(coordSys_.origin(), e3, e1);
176  }
177  else
178  {
179  // The 'coordinateSystem' sub-dictionary is optional,
180  // but enforce use of a cartesian system if not found.
181 
182  if (dict.found(coordinateSystem::typeName_()))
183  {
184  // New() for access to indirect (global) coordinate system
185  coordSys_ =
187  (
188  obr_,
189  dict,
190  coordinateSystem::typeName_()
191  );
192  }
193  else
194  {
195  coordSys_ = coordSystem::cartesian(dict);
196  }
197  }
198 
199 }
200 
201 
203 {
204  if (initialised_)
205  {
206  return;
207  }
208 
209  if (directForceDensity_)
210  {
211  if (!foundObject<volVectorField>(fDName_))
212  {
214  << "Could not find " << fDName_ << " in database"
215  << exit(FatalError);
216  }
217  }
218  else
219  {
220  if
221  (
222  !foundObject<volVectorField>(UName_)
223  || !foundObject<volScalarField>(pName_)
224 
225  )
226  {
228  << "Could not find U: " << UName_ << " or p:" << pName_
229  << " in database"
230  << exit(FatalError);
231  }
232 
233  if (rhoName_ != "rhoInf" && !foundObject<volScalarField>(rhoName_))
234  {
236  << "Could not find rho:" << rhoName_
237  << exit(FatalError);
238  }
239  }
240 
241  initialiseBins();
242 
243  initialised_ = true;
244 }
245 
246 
248 {
249  if (nBin_ > 1)
250  {
251  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
252 
253  // Determine extents of patches
254  scalar geomMin = GREAT;
255  scalar geomMax = -GREAT;
256  for (const label patchi : patchSet_)
257  {
258  const polyPatch& pp = pbm[patchi];
259  scalarField d(pp.faceCentres() & binDir_);
260  geomMin = min(min(d), geomMin);
261  geomMax = max(max(d), geomMax);
262  }
263 
264  // Include porosity
265  if (porosity_)
266  {
267  const HashTable<const porosityModel*> models =
268  obr_.lookupClass<porosityModel>();
269 
270  const scalarField dd(mesh_.C() & binDir_);
271 
272  forAllConstIters(models, iter)
273  {
274  const porosityModel& pm = *iter();
275  const labelList& cellZoneIDs = pm.cellZoneIDs();
276 
277  for (const label zonei : cellZoneIDs)
278  {
279  const cellZone& cZone = mesh_.cellZones()[zonei];
280  const scalarField d(dd, cZone);
281  geomMin = min(min(d), geomMin);
282  geomMax = max(max(d), geomMax);
283  }
284  }
285  }
286 
287  reduce(geomMin, minOp<scalar>());
288  reduce(geomMax, maxOp<scalar>());
289 
290  // Slightly boost max so that region of interest is fully within bounds
291  geomMax = 1.0001*(geomMax - geomMin) + geomMin;
292 
293  // Use geometry limits if not specified by the user
294  if (binMin_ == GREAT)
295  {
296  binMin_ = geomMin;
297  }
298  if (binMax_ == GREAT)
299  {
300  binMax_ = geomMax;
301  }
302 
303  binDx_ = (binMax_ - binMin_)/scalar(nBin_);
304 
305  // Create the bin mid-points used for writing
306  binPoints_.setSize(nBin_);
307  forAll(binPoints_, i)
308  {
309  binPoints_[i] = (i + 0.5)*binDir_*binDx_;
310  }
311  }
312 
313  // Allocate storage for forces and moments
314  forAll(force_, i)
315  {
316  force_[i].setSize(nBin_, vector::zero);
317  moment_[i].setSize(nBin_, vector::zero);
318  }
319 }
320 
321 
323 {
324  force_[0] = Zero;
325  force_[1] = Zero;
326  force_[2] = Zero;
327 
328  moment_[0] = Zero;
329  moment_[1] = Zero;
330  moment_[2] = Zero;
331 
332  if (writeFields_)
333  {
334  volVectorField& force =
335  lookupObjectRef<volVectorField>(fieldName("force"));
336 
337  force == dimensionedVector(force.dimensions(), Zero);
338 
339  volVectorField& moment =
340  lookupObjectRef<volVectorField>(fieldName("moment"));
341 
342  moment == dimensionedVector(moment.dimensions(), Zero);
343  }
344 }
345 
346 
349 {
350  typedef compressible::turbulenceModel cmpTurbModel;
351  typedef incompressible::turbulenceModel icoTurbModel;
352 
353  if (foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
354  {
355  const cmpTurbModel& turb =
356  lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
357 
358  return turb.devRhoReff();
359  }
360  else if (foundObject<icoTurbModel>(icoTurbModel::propertiesName))
361  {
363  lookupObject<icoTurbModel>(icoTurbModel::propertiesName);
364 
365  return rho()*turb.devReff();
366  }
367  else if (foundObject<fluidThermo>(fluidThermo::dictName))
368  {
369  const fluidThermo& thermo =
370  lookupObject<fluidThermo>(fluidThermo::dictName);
371 
372  const volVectorField& U = lookupObject<volVectorField>(UName_);
373 
374  return -thermo.mu()*dev(twoSymm(fvc::grad(U)));
375  }
376  else if (foundObject<transportModel>("transportProperties"))
377  {
378  const transportModel& laminarT =
379  lookupObject<transportModel>("transportProperties");
380 
381  const volVectorField& U = lookupObject<volVectorField>(UName_);
382 
383  return -rho()*laminarT.nu()*dev(twoSymm(fvc::grad(U)));
384  }
385  else if (foundObject<dictionary>("transportProperties"))
386  {
388  lookupObject<dictionary>("transportProperties");
389 
391 
392  const volVectorField& U = lookupObject<volVectorField>(UName_);
393 
394  return -rho()*nu*dev(twoSymm(fvc::grad(U)));
395  }
396  else
397  {
399  << "No valid model for viscous stress calculation"
400  << exit(FatalError);
401 
402  return volSymmTensorField::null();
403  }
404 }
405 
406 
408 {
409  if (foundObject<fluidThermo>(basicThermo::dictName))
410  {
411  const fluidThermo& thermo =
412  lookupObject<fluidThermo>(basicThermo::dictName);
413 
414  return thermo.mu();
415  }
416  else if (foundObject<transportModel>("transportProperties"))
417  {
418  const transportModel& laminarT =
419  lookupObject<transportModel>("transportProperties");
420 
421  return rho()*laminarT.nu();
422  }
423  else if (foundObject<dictionary>("transportProperties"))
424  {
426  lookupObject<dictionary>("transportProperties");
427 
429 
430  return rho()*nu;
431  }
432  else
433  {
435  << "No valid model for dynamic viscosity calculation"
436  << exit(FatalError);
437 
438  return volScalarField::null();
439  }
440 }
441 
442 
444 {
445  if (rhoName_ == "rhoInf")
446  {
448  (
449  IOobject
450  (
451  "rho",
452  mesh_.time().timeName(),
453  mesh_
454  ),
455  mesh_,
456  dimensionedScalar("rho", dimDensity, rhoRef_)
457  );
458  }
459 
460  return(lookupObject<volScalarField>(rhoName_));
461 }
462 
463 
465 {
466  if (p.dimensions() == dimPressure)
467  {
468  return 1.0;
469  }
470 
471  if (rhoName_ != "rhoInf")
472  {
474  << "Dynamic pressure is expected but kinematic is provided."
475  << exit(FatalError);
476  }
477 
478  return rhoRef_;
479 }
480 
481 
483 (
484  const vectorField& Md,
485  const vectorField& fN,
486  const vectorField& fT,
487  const vectorField& fP,
488  const vectorField& d
489 )
490 {
491  if (nBin_ == 1)
492  {
493  force_[0][0] += sum(fN);
494  force_[1][0] += sum(fT);
495  force_[2][0] += sum(fP);
496  moment_[0][0] += sum(Md^fN);
497  moment_[1][0] += sum(Md^fT);
498  moment_[2][0] += sum(Md^fP);
499  }
500  else
501  {
502  scalarField dd((d & binDir_) - binMin_);
503 
504  forAll(dd, i)
505  {
506  label bini = min(max(floor(dd[i]/binDx_), 0), force_[0].size() - 1);
507 
508  force_[0][bini] += fN[i];
509  force_[1][bini] += fT[i];
510  force_[2][bini] += fP[i];
511  moment_[0][bini] += Md[i]^fN[i];
512  moment_[1][bini] += Md[i]^fT[i];
513  moment_[2][bini] += Md[i]^fP[i];
514  }
515  }
516 }
517 
518 
520 (
521  const label patchi,
522  const vectorField& Md,
523  const vectorField& fN,
524  const vectorField& fT,
525  const vectorField& fP
526 )
527 {
528  if (!writeFields_)
529  {
530  return;
531  }
532 
533  volVectorField& force =
534  lookupObjectRef<volVectorField>(fieldName("force"));
535 
536  vectorField& pf = force.boundaryFieldRef()[patchi];
537  pf += fN + fT + fP;
538 
539  volVectorField& moment =
540  lookupObjectRef<volVectorField>(fieldName("moment"));
541 
542  vectorField& pm = moment.boundaryFieldRef()[patchi];
543  pm += Md;
544 }
545 
546 
548 (
549  const labelList& cellIDs,
550  const vectorField& Md,
551  const vectorField& fN,
552  const vectorField& fT,
553  const vectorField& fP
554 )
555 {
556  if (!writeFields_)
557  {
558  return;
559  }
560 
561  volVectorField& force =
562  lookupObjectRef<volVectorField>(fieldName("force"));
563 
564  volVectorField& moment =
565  lookupObjectRef<volVectorField>(fieldName("moment"));
566 
567  forAll(cellIDs, i)
568  {
569  label celli = cellIDs[i];
570  force[celli] += fN[i] + fT[i] + fP[i];
571  moment[celli] += Md[i];
572  }
573 }
574 
575 
577 (
578  const string& descriptor,
579  const vectorField& fm0,
580  const vectorField& fm1,
581  const vectorField& fm2,
582  autoPtr<OFstream>& osPtr
583 ) const
584 {
585  vector pressure = sum(fm0);
586  vector viscous = sum(fm1);
587  vector porous = sum(fm2);
588  vector total = pressure + viscous + porous;
589 
590  Log << " Sum of " << descriptor.c_str() << nl
591  << " Total : " << total << nl
592  << " Pressure : " << pressure << nl
593  << " Viscous : " << viscous << nl;
594 
595  if (porosity_)
596  {
597  Log << " Porous : " << porous << nl;
598  }
599 
600  if (writeToFile())
601  {
602  Ostream& os = osPtr();
603 
604  writeCurrentTime(os);
605 
606  os << tab << total
607  << tab << pressure
608  << tab << viscous;
609 
610  if (porosity_)
611  {
612  os << tab << porous;
613  }
614 
615  os << endl;
616  }
617 }
618 
619 
621 {
622  Log << type() << " " << name() << " write:" << nl;
623 
624  writeIntegratedForceMoment
625  (
626  "forces",
627  coordSys_.localVector(force_[0]),
628  coordSys_.localVector(force_[1]),
629  coordSys_.localVector(force_[2]),
630  forceFilePtr_
631  );
632 
633  writeIntegratedForceMoment
634  (
635  "moments",
636  coordSys_.localVector(moment_[0]),
637  coordSys_.localVector(moment_[1]),
638  coordSys_.localVector(moment_[2]),
639  momentFilePtr_
640  );
641 
642  Log << endl;
643 }
644 
645 
647 (
648  const List<Field<vector>>& fm,
649  autoPtr<OFstream>& osPtr
650 ) const
651 {
652  if ((nBin_ == 1) || !writeToFile())
653  {
654  return;
655  }
656 
657  List<Field<vector>> f(fm);
658 
659  if (binCumulative_)
660  {
661  for (label i = 1; i < f[0].size(); i++)
662  {
663  f[0][i] += f[0][i-1];
664  f[1][i] += f[1][i-1];
665  f[2][i] += f[2][i-1];
666  }
667  }
668 
669  Ostream& os = osPtr();
670 
671  writeCurrentTime(os);
672 
673  forAll(f[0], i)
674  {
675  vector total = f[0][i] + f[1][i] + f[2][i];
676 
677  os << tab << total
678  << tab << f[0][i]
679  << tab << f[1][i];
680 
681  if (porosity_)
682  {
683  os << tab << f[2][i];
684  }
685  }
686 
687  os << nl;
688 }
689 
690 
692 {
693  List<Field<vector>> lf(3);
694  List<Field<vector>> lm(3);
695  lf[0] = coordSys_.localVector(force_[0]);
696  lf[1] = coordSys_.localVector(force_[1]);
697  lf[2] = coordSys_.localVector(force_[2]);
698  lm[0] = coordSys_.localVector(moment_[0]);
699  lm[1] = coordSys_.localVector(moment_[1]);
700  lm[2] = coordSys_.localVector(moment_[2]);
701 
702  writeBinnedForceMoment(lf, forceBinFilePtr_);
703  writeBinnedForceMoment(lm, momentBinFilePtr_);
704 }
705 
706 
707 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
708 
710 (
711  const word& name,
712  const Time& runTime,
713  const dictionary& dict,
714  bool readFields
715 )
716 :
718  writeFile(mesh_, name),
719  force_(3),
720  moment_(3),
721  forceFilePtr_(),
722  momentFilePtr_(),
723  forceBinFilePtr_(),
724  momentBinFilePtr_(),
725  patchSet_(),
726  pName_("p"),
727  UName_("U"),
728  rhoName_("rho"),
729  directForceDensity_(false),
730  fDName_("fD"),
731  rhoRef_(VGREAT),
732  pRef_(0),
733  coordSys_(),
734  porosity_(false),
735  nBin_(1),
736  binDir_(Zero),
737  binDx_(0),
738  binMin_(GREAT),
739  binMax_(GREAT),
740  binPoints_(),
741  binCumulative_(true),
742  writeFields_(false),
743  initialised_(false)
744 {
745  if (readFields)
746  {
747  read(dict);
748  setCoordinateSystem(dict);
749  Log << endl;
750  }
751 }
752 
753 
755 (
756  const word& name,
757  const objectRegistry& obr,
758  const dictionary& dict,
759  bool readFields
760 )
761 :
763  writeFile(mesh_, name),
764  force_(3),
765  moment_(3),
766  forceFilePtr_(),
767  momentFilePtr_(),
768  forceBinFilePtr_(),
769  momentBinFilePtr_(),
770  patchSet_(),
771  pName_("p"),
772  UName_("U"),
773  rhoName_("rho"),
774  directForceDensity_(false),
775  fDName_("fD"),
776  rhoRef_(VGREAT),
777  pRef_(0),
778  coordSys_(),
779  porosity_(false),
780  nBin_(1),
781  binDir_(Zero),
782  binDx_(0),
783  binMin_(GREAT),
784  binMax_(GREAT),
785  binPoints_(),
786  binCumulative_(true),
787  writeFields_(false),
788  initialised_(false)
789 {
790  if (readFields)
791  {
792  read(dict);
793  setCoordinateSystem(dict);
794  Log << endl;
795  }
796 
797 /*
798  // Turn off writing to file
799  writeToFile_ = false;
800 
801  forAll(force_, i)
802  {
803  force_[i].setSize(nBin_, vector::zero);
804  moment_[i].setSize(nBin_, vector::zero);
805  }
806 */
807 }
808 
809 
810 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
811 
813 {
816 
817  initialised_ = false;
818 
819  Info<< type() << " " << name() << ":" << nl;
820 
821  directForceDensity_ = dict.lookupOrDefault("directForceDensity", false);
822 
823  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
824 
825  patchSet_ = pbm.patchSet(dict.get<wordRes>("patches"));
826 
827  if (directForceDensity_)
828  {
829  // Optional entry for fDName
830  if (dict.readIfPresent<word>("fD", fDName_))
831  {
832  Info<< " fD: " << fDName_ << endl;
833  }
834  }
835  else
836  {
837  // Optional field name entries
838  if (dict.readIfPresent<word>("p", pName_))
839  {
840  Info<< " p: " << pName_ << endl;
841  }
842  if (dict.readIfPresent<word>("U", UName_))
843  {
844  Info<< " U: " << UName_ << endl;
845  }
846  if (dict.readIfPresent<word>("rho", rhoName_))
847  {
848  Info<< " rho: " << rhoName_ << endl;
849  }
850 
851  // Reference density needed for incompressible calculations
852  if (rhoName_ == "rhoInf")
853  {
854  rhoRef_ = dict.get<scalar>("rhoInf");
855  Info<< " Freestream density (rhoInf) set to " << rhoRef_ << endl;
856  }
857 
858  // Reference pressure, 0 by default
859  if (dict.readIfPresent<scalar>("pRef", pRef_))
860  {
861  Info<< " Reference pressure (pRef) set to " << pRef_ << endl;
862  }
863  }
864 
865  dict.readIfPresent("porosity", porosity_);
866  if (porosity_)
867  {
868  Info<< " Including porosity effects" << endl;
869  }
870  else
871  {
872  Info<< " Not including porosity effects" << endl;
873  }
874 
875  if (dict.found("binData"))
876  {
877  Info<< " Activated data bins" << endl;
878  const dictionary& binDict(dict.subDict("binData"));
879  nBin_ = binDict.get<label>("nBin");
880 
881  if (nBin_ < 0)
882  {
884  << "Number of bins (nBin) must be zero or greater"
885  << exit(FatalIOError);
886  }
887  else if (nBin_ == 0)
888  {
889  // Case of no bins equates to a single bin to collect all data
890  nBin_ = 1;
891  }
892  else
893  {
894  Info<< " Employing " << nBin_ << " bins" << endl;
895  if (binDict.readIfPresent("min", binMin_))
896  {
897  Info<< " - min : " << binMin_ << endl;
898  }
899  if (binDict.readIfPresent("max", binMax_))
900  {
901  Info<< " - max : " << binMax_ << endl;
902  }
903 
904  binCumulative_ = binDict.get<bool>("cumulative");
905  Info<< " - cumuluative : " << binCumulative_ << endl;
906 
907  binDir_ = binDict.get<vector>("direction");
908  binDir_.normalise();
909  Info<< " - direction : " << binDir_ << endl;
910  }
911  }
912 
913  writeFields_ = dict.lookupOrDefault("writeFields", false);
914 
915  if (writeFields_)
916  {
917  Info<< " Fields will be written" << endl;
918 
919  volVectorField* forcePtr
920  (
921  new volVectorField
922  (
923  IOobject
924  (
925  fieldName("force"),
926  time_.timeName(),
927  mesh_,
930  ),
931  mesh_,
933  )
934  );
935 
936  mesh_.objectRegistry::store(forcePtr);
937 
938  volVectorField* momentPtr
939  (
940  new volVectorField
941  (
942  IOobject
943  (
944  fieldName("moment"),
945  time_.timeName(),
946  mesh_,
949  ),
950  mesh_,
952  )
953  );
954 
955  mesh_.objectRegistry::store(momentPtr);
956  }
957 
958  return true;
959 }
960 
961 
963 {
964  initialise();
965 
966  resetFields();
967 
968  if (directForceDensity_)
969  {
970  const volVectorField& fD = lookupObject<volVectorField>(fDName_);
971 
972  const surfaceVectorField::Boundary& Sfb = mesh_.Sf().boundaryField();
973 
974  for (const label patchi : patchSet_)
975  {
976  vectorField Md
977  (
978  mesh_.C().boundaryField()[patchi] - coordSys_.origin()
979  );
980 
981  scalarField sA(mag(Sfb[patchi]));
982 
983  // Normal force = surfaceUnitNormal*(surfaceNormal & forceDensity)
984  vectorField fN
985  (
986  Sfb[patchi]/sA
987  *(
988  Sfb[patchi] & fD.boundaryField()[patchi]
989  )
990  );
991 
992  // Tangential force (total force minus normal fN)
993  vectorField fT(sA*fD.boundaryField()[patchi] - fN);
994 
995  // Porous force
996  vectorField fP(Md.size(), Zero);
997 
998  addToFields(patchi, Md, fN, fT, fP);
999 
1000  applyBins(Md, fN, fT, fP, mesh_.C().boundaryField()[patchi]);
1001  }
1002  }
1003  else
1004  {
1005  const volScalarField& p = lookupObject<volScalarField>(pName_);
1006 
1007  const surfaceVectorField::Boundary& Sfb = mesh_.Sf().boundaryField();
1008 
1009  tmp<volSymmTensorField> tdevRhoReff = devRhoReff();
1010  const volSymmTensorField::Boundary& devRhoReffb
1011  = tdevRhoReff().boundaryField();
1012 
1013  // Scale pRef by density for incompressible simulations
1014  scalar pRef = pRef_/rho(p);
1015 
1016  for (const label patchi : patchSet_)
1017  {
1018  vectorField Md
1019  (
1020  mesh_.C().boundaryField()[patchi] - coordSys_.origin()
1021  );
1022 
1023  vectorField fN
1024  (
1025  rho(p)*Sfb[patchi]*(p.boundaryField()[patchi] - pRef)
1026  );
1027 
1028  vectorField fT(Sfb[patchi] & devRhoReffb[patchi]);
1029 
1030  vectorField fP(Md.size(), Zero);
1031 
1032  addToFields(patchi, Md, fN, fT, fP);
1033 
1034  applyBins(Md, fN, fT, fP, mesh_.C().boundaryField()[patchi]);
1035  }
1036  }
1037 
1038  if (porosity_)
1039  {
1040  const volVectorField& U = lookupObject<volVectorField>(UName_);
1041  const volScalarField rho(this->rho());
1042  const volScalarField mu(this->mu());
1043 
1044  const HashTable<const porosityModel*> models =
1045  obr_.lookupClass<porosityModel>();
1046 
1047  if (models.empty())
1048  {
1050  << "Porosity effects requested, but no porosity models found "
1051  << "in the database"
1052  << endl;
1053  }
1054 
1055  forAllConstIters(models, iter)
1056  {
1057  // Non-const access required if mesh is changing
1058  porosityModel& pm = const_cast<porosityModel&>(*iter());
1059 
1060  vectorField fPTot(pm.force(U, rho, mu));
1061 
1062  const labelList& cellZoneIDs = pm.cellZoneIDs();
1063 
1064  for (const label zonei : cellZoneIDs)
1065  {
1066  const cellZone& cZone = mesh_.cellZones()[zonei];
1067 
1068  const vectorField d(mesh_.C(), cZone);
1069  const vectorField fP(fPTot, cZone);
1070  const vectorField Md(d - coordSys_.origin());
1071 
1072  const vectorField fDummy(Md.size(), Zero);
1073 
1074  addToFields(cZone, Md, fDummy, fDummy, fP);
1075 
1076  applyBins(Md, fDummy, fDummy, fP, d);
1077  }
1078  }
1079  }
1080 
1084  Pstream::listCombineScatter(moment_);
1085 }
1086 
1087 
1089 {
1090  return sum(force_[0]) + sum(force_[1]) + sum(force_[2]);
1091 }
1092 
1093 
1095 {
1096  return sum(moment_[0]) + sum(moment_[1]) + sum(moment_[2]);
1097 }
1098 
1099 
1101 {
1102  calcForcesMoment();
1103 
1104  if (Pstream::master())
1105  {
1106  createFiles();
1107 
1108  writeForces();
1109 
1110  writeBins();
1111 
1112  Log << endl;
1113  }
1114 
1115  // Write state/results information
1116  setResult("normalForce", sum(force_[0]));
1117  setResult("tangentialForce", sum(force_[1]));
1118  setResult("porousForce", sum(force_[2]));
1119 
1120  setResult("normalMoment", sum(moment_[0]));
1121  setResult("tangentialMoment", sum(moment_[1]));
1122  setResult("porousMoment", sum(moment_[2]));
1123 
1124  return true;
1125 }
1126 
1127 
1129 {
1130  if (writeFields_)
1131  {
1132  lookupObject<volVectorField>(fieldName("force")).write();
1133  lookupObject<volVectorField>(fieldName("moment")).write();
1134  }
1135 
1136  return true;
1137 }
1138 
1139 
1140 // ************************************************************************* //
Foam::functionObjects::forces::forceEff
virtual vector forceEff() const
Return the total force.
Definition: forces.C:1088
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::dimPressure
const dimensionSet dimPressure
Foam::maxOp
Definition: ops.H:223
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::functionObjects::forces::mu
tmp< volScalarField > mu() const
Dynamic viscosity field.
Definition: forces.C:407
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::porosityModel::cellZoneIDs
const labelList & cellZoneIDs() const
Return const access to the cell zone IDs.
Definition: porosityModelI.H:62
Foam::functionObjects::forces::writeForces
void writeForces()
Write force data.
Definition: forces.C:620
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::transportModel::nu
virtual tmp< volScalarField > nu() const =0
Return the laminar viscosity.
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFieldRefs.H:4
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::dimDensity
const dimensionSet dimDensity
Foam::functionObjects::forces::write
virtual bool write()
Write the forces.
Definition: forces.C:1128
Foam::minOp
Definition: ops.H:224
Foam::jn
dimensionedScalar jn(const int n, const dimensionedScalar &ds)
Definition: dimensionedScalar.C:305
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
transportProperties
IOdictionary transportProperties(IOobject("transportProperties", runTime.constant(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::functionObjects::forces::writeBinHeader
void writeBinHeader(const word &header, Ostream &os) const
Write header for binned data.
Definition: forces.C:103
turbulentTransportModel.H
Foam::functionObjects::forces::fieldName
word fieldName(const word &name) const
Create a field name.
Definition: forces.C:51
Foam::fluidThermo
Fundamental fluid thermodynamic properties.
Definition: fluidThermo.H:52
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:81
porosityModel.H
Foam::writeHeader
static void writeHeader(Ostream &os, const word &fieldName)
Definition: rawSurfaceWriterImpl.C:49
Foam::functionObjects::forces::writeIntegratedForceMoment
void writeIntegratedForceMoment(const string &descriptor, const vectorField &fm0, const vectorField &fm1, const vectorField &fm2, autoPtr< OFstream > &osPtr) const
Helper function to write integrated forces and moments.
Definition: forces.C:577
Foam::functionObjects::forces::writeBins
void writeBins()
Write binned data.
Definition: forces.C:691
Foam::dimForce
const dimensionSet dimForce
Foam::Vector::normalise
Vector< Cmpt > & normalise()
Normalise the vector by its magnitude.
Definition: VectorI.H:128
Foam::functionObjects::pressure
Includes tools to manipulate the pressure into different forms.
Definition: pressure.H:256
cartesianCS.H
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
rho
rho
Definition: readInitialConditions.H:96
Foam::functionObjects::forces::setCoordinateSystem
void setCoordinateSystem(const dictionary &dict, const word &e3Name=word::null, const word &e1Name=word::null)
Set the co-ordinate system from dictionary and axes names.
Definition: forces.C:159
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, add, dictionary)
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::cellZone
A subset of mesh cells.
Definition: cellZone.H:62
Foam::functionObjects::forces::calcForcesMoment
virtual void calcForcesMoment()
Calculate the forces and moments.
Definition: forces.C:962
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::functionObjects::forces::initialiseBins
void initialiseBins()
Initialise the collection bins.
Definition: forces.C:247
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::functionObjects::forces::rho
tmp< volScalarField > rho() const
Return rho if specified otherwise rhoRef.
Definition: forces.C:443
Foam::functionObjects::writeFile::read
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:212
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::porosityModel::force
virtual tmp< vectorField > force(const volVectorField &U, const volScalarField &rho, const volScalarField &mu)
Return the force over the cell zone(s)
Foam::functionObjects::forces::momentEff
virtual vector momentEff() const
Return the total moment.
Definition: forces.C:1094
Foam::label
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:62
Foam::functionObjects::forces::resetFields
void resetFields()
Reset the fields prior to accumulation of force/moments.
Definition: forces.C:322
Foam::Field< vector >
Foam::functionObjects::forces::initialise
void initialise()
Initialise the fields.
Definition: forces.C:202
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::functionObjects::forces::devRhoReff
tmp< volSymmTensorField > devRhoReff() const
Return the effective viscous stress (laminar + turbulent).
Definition: forces.C:348
Foam::coordinateSystem::New
static autoPtr< coordinateSystem > New(word modelType, const objectRegistry &obr, const dictionary &dict)
Definition: coordinateSystemNew.C:82
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::functionObjects::forces::read
virtual bool read(const dictionary &)
Read the forces data.
Definition: forces.C:812
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionary.H:458
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dimViscosity
const dimensionSet dimViscosity
Foam::functionObjects::forces::addToFields
void addToFields(const label patchi, const vectorField &Md, const vectorField &fN, const vectorField &fT, const vectorField &fP)
Add patch contributions to force and moment fields.
Definition: forces.C:520
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::functionObjects::forces::execute
virtual bool execute()
Execute, currently does nothing.
Definition: forces.C:1100
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:166
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
Foam::functionObjects::forces::writeIntegratedHeader
void writeIntegratedHeader(const word &header, Ostream &os) const
Write header for integrated data.
Definition: forces.C:80
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::functionObjects::forces::writeBinnedForceMoment
void writeBinnedForceMoment(const List< Field< vector >> &fm, autoPtr< OFstream > &osPtr) const
Helper function to write binned forces and moments.
Definition: forces.C:647
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::transportModel
Base-class for all transport models used by the incompressible turbulence models.
Definition: transportModel.H:53
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:438
Foam::HashTable
A HashTable similar to std::unordered_map.
Definition: HashTable.H:105
Foam::Pstream::listCombineGather
static void listCombineGather(const List< commsStruct > &comms, List< T > &Value, const CombineOp &cop, const int tag, const label comm)
Definition: combineGatherScatter.C:290
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:53
U
U
Definition: pEqn.H:72
Foam::porosityModel
Top level model for porosity models.
Definition: porosityModel.H:57
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
forces.H
Foam::functionObjects::forces::applyBins
void applyBins(const vectorField &Md, const vectorField &fN, const vectorField &fT, const vectorField &fP, const vectorField &d)
Accumulate bin data.
Definition: forces.C:483
Foam::tab
constexpr char tab
Definition: Ostream.H:371
Foam::nl
constexpr char nl
Definition: Ostream.H:372
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::polyPatch::faceCentres
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:311
Foam::functionObject::name
const word & name() const
Return the name of this functionObject.
Definition: functionObject.C:131
f
labelList f(nPoints)
Foam::coordSystem::cartesian
A Cartesian coordinate system.
Definition: cartesianCS.H:69
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:752
Foam::Vector< scalar >
Foam::ThermalDiffusivity
Templated wrapper class to provide compressible turbulence models thermal diffusivity based thermal t...
Definition: phaseCompressibleTurbulenceModelFwd.H:47
Foam::List< label >
Foam::functionObjects::forces::createFiles
void createFiles()
Create the output files.
Definition: forces.C:57
fvcGrad.H
Calculate the gradient of the given field.
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
Foam::polyBoundaryMesh::patchSet
labelHashSet patchSet(const UList< wordRe > &patchNames, const bool warnNotFound=true, const bool useGroups=true) const
Return the set of patch IDs corresponding to the given names.
Definition: polyBoundaryMesh.C:857
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
Foam::sum
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
Definition: DimensionedFieldFunctions.C:327
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
Foam::functionObjects::readFields
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:108
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
Foam::functionObjects::writeFile
functionObject base class for writing single files
Definition: writeFile.H:59
Foam::IncompressibleTurbulenceModel
Templated abstract base class for single-phase incompressible turbulence models.
Definition: IncompressibleTurbulenceModel.H:55
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:375
Foam::HashTable::empty
bool empty() const noexcept
True if the hash table is empty.
Definition: HashTableI.H:59
Foam::plusEqOp
Definition: ops.H:72
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::GeometricField< symmTensor, fvPatchField, volMesh >::null
static const GeometricField< symmTensor, fvPatchField, volMesh > & null()
Return a null geometric field.
Definition: GeometricFieldI.H:32
Foam::Pstream::listCombineScatter
static void listCombineScatter(const List< commsStruct > &comms, List< T > &Value, const int tag, const label comm)
Scatter data. Reverse of combineGather.
Definition: combineGatherScatter.C:432
Foam::GeometricField< vector, fvPatchField, volMesh >
turb
compressible::turbulenceModel & turb
Definition: setRegionFluidFields.H:10
Foam::IOobject::NO_READ
Definition: IOobject.H:123
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
Foam::functionObjects::forces::forces
forces(const forces &)=delete
No copy construct.
turbulentFluidThermoModel.H
Log
#define Log
Report write to Foam::Info if the local log switch is true.
Definition: messageStream.H:332
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106