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-2020 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  // Note: Only possible to create bin files after bins have been initialised
54 
55  if (writeToFile() && !forceFilePtr_)
56  {
57  forceFilePtr_ = createFile("force");
59  momentFilePtr_ = createFile("moment");
61 
62  if (nBin_ > 1)
63  {
64  forceBinFilePtr_ = createFile("forceBin");
65  writeBinHeader("Force", forceBinFilePtr_());
66  momentBinFilePtr_ = createFile("momentBin");
67  writeBinHeader("Moment", momentBinFilePtr_());
68  }
69  }
70 }
71 
72 
74 (
75  const word& header,
76  Ostream& os
77 ) const
78 {
79  writeHeader(os, header);
80  writeHeaderValue(os, "CofR", coordSysPtr_->origin());
81  writeHeader(os, "");
82  writeCommented(os, "Time");
83  writeTabbed(os, "(total_x total_y total_z)");
84  writeTabbed(os, "(pressure_x pressure_y pressure_z)");
85  writeTabbed(os, "(viscous_x viscous_y viscous_z)");
86 
87  if (porosity_)
88  {
89  writeTabbed(os, "(porous_x porous_y porous_z)");
90  }
91 
92  os << endl;
93 }
94 
95 
97 (
98  const word& header,
99  Ostream& os
100 ) const
101 {
102  writeHeader(os, header + " bins");
103  writeHeaderValue(os, "bins", nBin_);
104  writeHeaderValue(os, "start", binMin_);
105  writeHeaderValue(os, "end", binMax_);
106  writeHeaderValue(os, "delta", binDx_);
107  writeHeaderValue(os, "direction", binDir_);
108 
109  vectorField binPoints(nBin_);
110  writeCommented(os, "x co-ords :");
111  forAll(binPoints, pointi)
112  {
113  binPoints[pointi] = (binMin_ + (pointi + 1)*binDx_)*binDir_;
114  os << tab << binPoints[pointi].x();
115  }
116  os << nl;
117 
118  writeCommented(os, "y co-ords :");
119  forAll(binPoints, pointi)
120  {
121  os << tab << binPoints[pointi].y();
122  }
123  os << nl;
124 
125  writeCommented(os, "z co-ords :");
126  forAll(binPoints, pointi)
127  {
128  os << tab << binPoints[pointi].z();
129  }
130  os << nl;
131 
132  writeHeader(os, "");
133  writeCommented(os, "Time");
134 
135  for (label j = 0; j < nBin_; j++)
136  {
137  const word jn(Foam::name(j) + ':');
138  os << tab << jn << "(total_x total_y total_z)"
139  << tab << jn << "(pressure_x pressure_y pressure_z)"
140  << tab << jn << "(viscous_x viscous_y viscous_z)";
141 
142  if (porosity_)
143  {
144  os << tab << jn << "(porous_x porous_y porous_z)";
145  }
146  }
147 
148  os << endl;
149 }
150 
151 
153 (
154  const dictionary& dict,
155  const word& e3Name,
156  const word& e1Name
157 )
158 {
159  coordSysPtr_.clear();
160 
161  point origin(Zero);
162  if (dict.readIfPresent<point>("CofR", origin))
163  {
164  const vector e3 = e3Name == word::null ?
165  vector(0, 0, 1) : dict.get<vector>(e3Name);
166  const vector e1 = e1Name == word::null ?
167  vector(1, 0, 0) : dict.get<vector>(e1Name);
168 
169  coordSysPtr_.reset(new coordSystem::cartesian(origin, e3, e1));
170  }
171  else
172  {
173  // The 'coordinateSystem' sub-dictionary is optional,
174  // but enforce use of a cartesian system if not found.
175 
176  if (dict.found(coordinateSystem::typeName_()))
177  {
178  // New() for access to indirect (global) coordinate system
179  coordSysPtr_ =
181  (
182  obr_,
183  dict,
184  coordinateSystem::typeName_()
185  );
186  }
187  else
188  {
189  coordSysPtr_.reset(new coordSystem::cartesian(dict));
190  }
191  }
192 
193 }
194 
195 
197 {
198  if (initialised_)
199  {
200  return;
201  }
202 
203  if (directForceDensity_)
204  {
205  if (!foundObject<volVectorField>(fDName_))
206  {
208  << "Could not find " << fDName_ << " in database"
209  << exit(FatalError);
210  }
211  }
212  else
213  {
214  if
215  (
216  !foundObject<volVectorField>(UName_)
217  || !foundObject<volScalarField>(pName_)
218 
219  )
220  {
222  << "Could not find U: " << UName_ << " or p:" << pName_
223  << " in database"
224  << exit(FatalError);
225  }
226 
227  if (rhoName_ != "rhoInf" && !foundObject<volScalarField>(rhoName_))
228  {
230  << "Could not find rho:" << rhoName_
231  << exit(FatalError);
232  }
233  }
234 
235  initialiseBins();
236 
237  initialised_ = true;
238 }
239 
240 
242 {
243  if (nBin_ > 1)
244  {
245  const polyBoundaryMesh& pbm = mesh_.boundaryMesh();
246 
247  // Determine extents of patches
248  scalar geomMin = GREAT;
249  scalar geomMax = -GREAT;
250  for (const label patchi : patchSet_)
251  {
252  const polyPatch& pp = pbm[patchi];
253  scalarField d(pp.faceCentres() & binDir_);
254  geomMin = min(min(d), geomMin);
255  geomMax = max(max(d), geomMax);
256  }
257 
258  // Include porosity
259  if (porosity_)
260  {
261  const HashTable<const porosityModel*> models =
262  obr_.lookupClass<porosityModel>();
263 
264  const scalarField dd(mesh_.C() & binDir_);
265 
266  forAllConstIters(models, iter)
267  {
268  const porosityModel& pm = *iter();
269  const labelList& cellZoneIDs = pm.cellZoneIDs();
270 
271  for (const label zonei : cellZoneIDs)
272  {
273  const cellZone& cZone = mesh_.cellZones()[zonei];
274  const scalarField d(dd, cZone);
275  geomMin = min(min(d), geomMin);
276  geomMax = max(max(d), geomMax);
277  }
278  }
279  }
280 
281  reduce(geomMin, minOp<scalar>());
282  reduce(geomMax, maxOp<scalar>());
283 
284  // Slightly boost max so that region of interest is fully within bounds
285  geomMax = 1.0001*(geomMax - geomMin) + geomMin;
286 
287  // Use geometry limits if not specified by the user
288  if (binMin_ == GREAT)
289  {
290  binMin_ = geomMin;
291  }
292  if (binMax_ == GREAT)
293  {
294  binMax_ = geomMax;
295  }
296 
297  binDx_ = (binMax_ - binMin_)/scalar(nBin_);
298 
299  // Create the bin mid-points used for writing
300  binPoints_.setSize(nBin_);
301  forAll(binPoints_, i)
302  {
303  binPoints_[i] = (i + 0.5)*binDir_*binDx_;
304  }
305  }
306 
307  // Allocate storage for forces and moments
308  forAll(force_, i)
309  {
310  force_[i].setSize(nBin_, vector::zero);
311  moment_[i].setSize(nBin_, vector::zero);
312  }
313 }
314 
315 
317 {
318  force_[0] = Zero;
319  force_[1] = Zero;
320  force_[2] = Zero;
321 
322  moment_[0] = Zero;
323  moment_[1] = Zero;
324  moment_[2] = Zero;
325 
326  if (writeFields_)
327  {
328  volVectorField& force =
329  lookupObjectRef<volVectorField>(scopedName("force"));
330 
331  force == dimensionedVector(force.dimensions(), Zero);
332 
333  volVectorField& moment =
334  lookupObjectRef<volVectorField>(scopedName("moment"));
335 
336  moment == dimensionedVector(moment.dimensions(), Zero);
337  }
338 }
339 
340 
343 {
344  typedef compressible::turbulenceModel cmpTurbModel;
345  typedef incompressible::turbulenceModel icoTurbModel;
346 
347  const auto& U = lookupObject<volVectorField>(UName_);
348 
349  if (foundObject<cmpTurbModel>(cmpTurbModel::propertiesName))
350  {
351  const cmpTurbModel& turb =
352  lookupObject<cmpTurbModel>(cmpTurbModel::propertiesName);
353 
354  return turb.devRhoReff(U);
355  }
356  else if (foundObject<icoTurbModel>(icoTurbModel::propertiesName))
357  {
359  lookupObject<icoTurbModel>(icoTurbModel::propertiesName);
360 
361  return rho()*turb.devReff(U);
362  }
363  else if (foundObject<fluidThermo>(fluidThermo::dictName))
364  {
365  const fluidThermo& thermo =
366  lookupObject<fluidThermo>(fluidThermo::dictName);
367 
368  return -thermo.mu()*dev(twoSymm(fvc::grad(U)));
369  }
370  else if (foundObject<transportModel>("transportProperties"))
371  {
372  const transportModel& laminarT =
373  lookupObject<transportModel>("transportProperties");
374 
375  return -rho()*laminarT.nu()*dev(twoSymm(fvc::grad(U)));
376  }
377  else if (foundObject<dictionary>("transportProperties"))
378  {
380  lookupObject<dictionary>("transportProperties");
381 
383 
384  return -rho()*nu*dev(twoSymm(fvc::grad(U)));
385  }
386  else
387  {
389  << "No valid model for viscous stress calculation"
390  << exit(FatalError);
391 
392  return volSymmTensorField::null();
393  }
394 }
395 
396 
398 {
399  if (foundObject<fluidThermo>(basicThermo::dictName))
400  {
401  const fluidThermo& thermo =
402  lookupObject<fluidThermo>(basicThermo::dictName);
403 
404  return thermo.mu();
405  }
406  else if (foundObject<transportModel>("transportProperties"))
407  {
408  const transportModel& laminarT =
409  lookupObject<transportModel>("transportProperties");
410 
411  return rho()*laminarT.nu();
412  }
413  else if (foundObject<dictionary>("transportProperties"))
414  {
416  lookupObject<dictionary>("transportProperties");
417 
419 
420  return rho()*nu;
421  }
422  else
423  {
425  << "No valid model for dynamic viscosity calculation"
426  << exit(FatalError);
427 
428  return volScalarField::null();
429  }
430 }
431 
432 
434 {
435  if (rhoName_ == "rhoInf")
436  {
438  (
439  IOobject
440  (
441  "rho",
442  mesh_.time().timeName(),
443  mesh_
444  ),
445  mesh_,
446  dimensionedScalar("rho", dimDensity, rhoRef_)
447  );
448  }
449 
450  return(lookupObject<volScalarField>(rhoName_));
451 }
452 
453 
455 {
456  if (p.dimensions() == dimPressure)
457  {
458  return 1.0;
459  }
460 
461  if (rhoName_ != "rhoInf")
462  {
464  << "Dynamic pressure is expected but kinematic is provided."
465  << exit(FatalError);
466  }
467 
468  return rhoRef_;
469 }
470 
471 
473 (
474  const vectorField& Md,
475  const vectorField& fN,
476  const vectorField& fT,
477  const vectorField& fP,
478  const vectorField& d
479 )
480 {
481  if (nBin_ == 1)
482  {
483  force_[0][0] += sum(fN);
484  force_[1][0] += sum(fT);
485  force_[2][0] += sum(fP);
486  moment_[0][0] += sum(Md^fN);
487  moment_[1][0] += sum(Md^fT);
488  moment_[2][0] += sum(Md^fP);
489  }
490  else
491  {
492  scalarField dd((d & binDir_) - binMin_);
493 
494  forAll(dd, i)
495  {
496  label bini = min(max(floor(dd[i]/binDx_), 0), force_[0].size() - 1);
497 
498  force_[0][bini] += fN[i];
499  force_[1][bini] += fT[i];
500  force_[2][bini] += fP[i];
501  moment_[0][bini] += Md[i]^fN[i];
502  moment_[1][bini] += Md[i]^fT[i];
503  moment_[2][bini] += Md[i]^fP[i];
504  }
505  }
506 }
507 
508 
510 (
511  const label patchi,
512  const vectorField& Md,
513  const vectorField& fN,
514  const vectorField& fT,
515  const vectorField& fP
516 )
517 {
518  if (!writeFields_)
519  {
520  return;
521  }
522 
523  auto& force = lookupObjectRef<volVectorField>(scopedName("force"));
524  vectorField& pf = force.boundaryFieldRef()[patchi];
525  pf += fN + fT + fP;
526 
527  auto& moment = lookupObjectRef<volVectorField>(scopedName("moment"));
528  vectorField& pm = moment.boundaryFieldRef()[patchi];
529  pm = Md^pf;
530 }
531 
532 
534 (
535  const labelList& cellIDs,
536  const vectorField& Md,
537  const vectorField& fN,
538  const vectorField& fT,
539  const vectorField& fP
540 )
541 {
542  if (!writeFields_)
543  {
544  return;
545  }
546 
547  auto& force = lookupObjectRef<volVectorField>(scopedName("force"));
548  auto& moment = lookupObjectRef<volVectorField>(scopedName("moment"));
549 
550  forAll(cellIDs, i)
551  {
552  label celli = cellIDs[i];
553  force[celli] += fN[i] + fT[i] + fP[i];
554  moment[celli] = Md[i]^force[celli];
555  }
556 }
557 
558 
560 (
561  const string& descriptor,
562  const vectorField& fm0,
563  const vectorField& fm1,
564  const vectorField& fm2,
565  autoPtr<OFstream>& osPtr
566 ) const
567 {
568  vector pressure = sum(fm0);
569  vector viscous = sum(fm1);
570  vector porous = sum(fm2);
571  vector total = pressure + viscous + porous;
572 
573  Log << " Sum of " << descriptor.c_str() << nl
574  << " Total : " << total << nl
575  << " Pressure : " << pressure << nl
576  << " Viscous : " << viscous << nl;
577 
578  if (porosity_)
579  {
580  Log << " Porous : " << porous << nl;
581  }
582 
583  if (writeToFile())
584  {
585  Ostream& os = osPtr();
586 
587  writeCurrentTime(os);
588 
589  os << tab << total
590  << tab << pressure
591  << tab << viscous;
592 
593  if (porosity_)
594  {
595  os << tab << porous;
596  }
597 
598  os << endl;
599  }
600 }
601 
602 
604 {
605  Log << type() << " " << name() << " write:" << nl;
606 
607  const auto& coordSys = coordSysPtr_();
608 
609  writeIntegratedForceMoment
610  (
611  "forces",
612  coordSys.localVector(force_[0]),
613  coordSys.localVector(force_[1]),
614  coordSys.localVector(force_[2]),
615  forceFilePtr_
616  );
617 
618  writeIntegratedForceMoment
619  (
620  "moments",
621  coordSys.localVector(moment_[0]),
622  coordSys.localVector(moment_[1]),
623  coordSys.localVector(moment_[2]),
624  momentFilePtr_
625  );
626 
627  Log << endl;
628 }
629 
630 
632 (
633  const List<Field<vector>>& fm,
634  autoPtr<OFstream>& osPtr
635 ) const
636 {
637  if ((nBin_ == 1) || !writeToFile())
638  {
639  return;
640  }
641 
642  List<Field<vector>> f(fm);
643 
644  if (binCumulative_)
645  {
646  for (label i = 1; i < f[0].size(); i++)
647  {
648  f[0][i] += f[0][i-1];
649  f[1][i] += f[1][i-1];
650  f[2][i] += f[2][i-1];
651  }
652  }
653 
654  Ostream& os = osPtr();
655 
656  writeCurrentTime(os);
657 
658  forAll(f[0], i)
659  {
660  vector total = f[0][i] + f[1][i] + f[2][i];
661 
662  os << tab << total
663  << tab << f[0][i]
664  << tab << f[1][i];
665 
666  if (porosity_)
667  {
668  os << tab << f[2][i];
669  }
670  }
671 
672  os << nl;
673 }
674 
675 
677 {
678  const auto& coordSys = coordSysPtr_();
679 
680  List<Field<vector>> lf(3);
681  List<Field<vector>> lm(3);
682  lf[0] = coordSys.localVector(force_[0]);
683  lf[1] = coordSys.localVector(force_[1]);
684  lf[2] = coordSys.localVector(force_[2]);
685  lm[0] = coordSys.localVector(moment_[0]);
686  lm[1] = coordSys.localVector(moment_[1]);
687  lm[2] = coordSys.localVector(moment_[2]);
688 
689  writeBinnedForceMoment(lf, forceBinFilePtr_);
690  writeBinnedForceMoment(lm, momentBinFilePtr_);
691 }
692 
693 
694 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
695 
697 (
698  const word& name,
699  const Time& runTime,
700  const dictionary& dict,
701  bool readFields
702 )
703 :
705  writeFile(mesh_, name),
706  force_(3),
707  moment_(3),
708  forceFilePtr_(),
709  momentFilePtr_(),
710  forceBinFilePtr_(),
711  momentBinFilePtr_(),
712  patchSet_(),
713  pName_("p"),
714  UName_("U"),
715  rhoName_("rho"),
716  directForceDensity_(false),
717  fDName_("fD"),
718  rhoRef_(VGREAT),
719  pRef_(0),
720  coordSysPtr_(nullptr),
721  porosity_(false),
722  nBin_(1),
723  binDir_(Zero),
724  binDx_(0),
725  binMin_(GREAT),
726  binMax_(GREAT),
727  binPoints_(),
728  binCumulative_(true),
729  writeFields_(false),
730  initialised_(false)
731 {
732  if (readFields)
733  {
734  read(dict);
735  setCoordinateSystem(dict);
736  Log << endl;
737  }
738 }
739 
740 
742 (
743  const word& name,
744  const objectRegistry& obr,
745  const dictionary& dict,
746  bool readFields
747 )
748 :
750  writeFile(mesh_, name),
751  force_(3),
752  moment_(3),
753  forceFilePtr_(),
754  momentFilePtr_(),
755  forceBinFilePtr_(),
756  momentBinFilePtr_(),
757  patchSet_(),
758  pName_("p"),
759  UName_("U"),
760  rhoName_("rho"),
761  directForceDensity_(false),
762  fDName_("fD"),
763  rhoRef_(VGREAT),
764  pRef_(0),
765  coordSysPtr_(nullptr),
766  porosity_(false),
767  nBin_(1),
768  binDir_(Zero),
769  binDx_(0),
770  binMin_(GREAT),
771  binMax_(GREAT),
772  binPoints_(),
773  binCumulative_(true),
774  writeFields_(false),
775  initialised_(false)
776 {
777  if (readFields)
778  {
779  read(dict);
780  setCoordinateSystem(dict);
781  Log << endl;
782  }
783 
784 /*
785  // Turn off writing to file
786  writeToFile_ = false;
787 
788  forAll(force_, i)
789  {
790  force_[i].setSize(nBin_, vector::zero);
791  moment_[i].setSize(nBin_, vector::zero);
792  }
793 */
794 }
795 
796 
797 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
798 
800 {
803 
804  initialised_ = false;
805 
806  Info<< type() << " " << name() << ":" << nl;
807 
808  directForceDensity_ = dict.getOrDefault("directForceDensity", false);
809 
810  patchSet_ =
811  mesh_.boundaryMesh().patchSet
812  (
813  dict.get<wordRes>("patches")
814  );
815 
816  if (directForceDensity_)
817  {
818  // Optional entry for fDName
819  if (dict.readIfPresent<word>("fD", fDName_))
820  {
821  Info<< " fD: " << fDName_ << endl;
822  }
823  }
824  else
825  {
826  // Optional field name entries
827  if (dict.readIfPresent<word>("p", pName_))
828  {
829  Info<< " p: " << pName_ << endl;
830  }
831  if (dict.readIfPresent<word>("U", UName_))
832  {
833  Info<< " U: " << UName_ << endl;
834  }
835  if (dict.readIfPresent<word>("rho", rhoName_))
836  {
837  Info<< " rho: " << rhoName_ << endl;
838  }
839 
840  // Reference density needed for incompressible calculations
841  if (rhoName_ == "rhoInf")
842  {
843  rhoRef_ = dict.get<scalar>("rhoInf");
844  Info<< " Freestream density (rhoInf) set to " << rhoRef_ << endl;
845  }
846 
847  // Reference pressure, 0 by default
848  if (dict.readIfPresent<scalar>("pRef", pRef_))
849  {
850  Info<< " Reference pressure (pRef) set to " << pRef_ << endl;
851  }
852  }
853 
854  dict.readIfPresent("porosity", porosity_);
855  if (porosity_)
856  {
857  Info<< " Including porosity effects" << endl;
858  }
859  else
860  {
861  Info<< " Not including porosity effects" << endl;
862  }
863 
864  if (dict.found("binData"))
865  {
866  Info<< " Activated data bins" << endl;
867  const dictionary& binDict(dict.subDict("binData"));
868  nBin_ = binDict.get<label>("nBin");
869 
870  if (nBin_ < 0)
871  {
873  << "Number of bins (nBin) must be zero or greater"
874  << exit(FatalIOError);
875  }
876  else if (nBin_ == 0)
877  {
878  // Case of no bins equates to a single bin to collect all data
879  nBin_ = 1;
880  }
881  else
882  {
883  Info<< " Employing " << nBin_ << " bins" << endl;
884  if (binDict.readIfPresent("min", binMin_))
885  {
886  Info<< " - min : " << binMin_ << endl;
887  }
888  if (binDict.readIfPresent("max", binMax_))
889  {
890  Info<< " - max : " << binMax_ << endl;
891  }
892 
893  binCumulative_ = binDict.get<bool>("cumulative");
894  Info<< " - cumuluative : " << binCumulative_ << endl;
895 
896  binDir_ = binDict.get<vector>("direction");
897  binDir_.normalise();
898  Info<< " - direction : " << binDir_ << endl;
899  }
900  }
901 
902  writeFields_ = dict.getOrDefault("writeFields", false);
903 
904  if (writeFields_)
905  {
906  Info<< " Fields will be written" << endl;
907 
908  volVectorField* forcePtr
909  (
910  new volVectorField
911  (
912  IOobject
913  (
914  scopedName("force"),
915  time_.timeName(),
916  mesh_,
919  ),
920  mesh_,
922  )
923  );
924 
925  mesh_.objectRegistry::store(forcePtr);
926 
927  volVectorField* momentPtr
928  (
929  new volVectorField
930  (
931  IOobject
932  (
933  scopedName("moment"),
934  time_.timeName(),
935  mesh_,
938  ),
939  mesh_,
941  )
942  );
943 
944  mesh_.objectRegistry::store(momentPtr);
945  }
946 
947  return true;
948 }
949 
950 
952 {
953  initialise();
954 
955  resetFields();
956 
957  const point& origin = coordSysPtr_->origin();
958 
959  if (directForceDensity_)
960  {
961  const volVectorField& fD = lookupObject<volVectorField>(fDName_);
962 
963  const surfaceVectorField::Boundary& Sfb = mesh_.Sf().boundaryField();
964 
965  for (const label patchi : patchSet_)
966  {
967  vectorField Md(mesh_.C().boundaryField()[patchi] - origin);
968 
969  scalarField sA(mag(Sfb[patchi]));
970 
971  // Normal force = surfaceUnitNormal*(surfaceNormal & forceDensity)
972  vectorField fN
973  (
974  Sfb[patchi]/sA
975  *(
976  Sfb[patchi] & fD.boundaryField()[patchi]
977  )
978  );
979 
980  // Tangential force (total force minus normal fN)
981  vectorField fT(sA*fD.boundaryField()[patchi] - fN);
982 
983  // Porous force
984  vectorField fP(Md.size(), Zero);
985 
986  addToFields(patchi, Md, fN, fT, fP);
987 
988  applyBins(Md, fN, fT, fP, mesh_.C().boundaryField()[patchi]);
989  }
990  }
991  else
992  {
993  const volScalarField& p = lookupObject<volScalarField>(pName_);
994 
995  const surfaceVectorField::Boundary& Sfb = mesh_.Sf().boundaryField();
996 
997  tmp<volSymmTensorField> tdevRhoReff = devRhoReff();
998  const volSymmTensorField::Boundary& devRhoReffb
999  = tdevRhoReff().boundaryField();
1000 
1001  // Scale pRef by density for incompressible simulations
1002  scalar pRef = pRef_/rho(p);
1003 
1004  for (const label patchi : patchSet_)
1005  {
1006  vectorField Md(mesh_.C().boundaryField()[patchi] - origin);
1007 
1008  vectorField fN
1009  (
1010  rho(p)*Sfb[patchi]*(p.boundaryField()[patchi] - pRef)
1011  );
1012 
1013  vectorField fT(Sfb[patchi] & devRhoReffb[patchi]);
1014 
1015  vectorField fP(Md.size(), Zero);
1016 
1017  addToFields(patchi, Md, fN, fT, fP);
1018 
1019  applyBins(Md, fN, fT, fP, mesh_.C().boundaryField()[patchi]);
1020  }
1021  }
1022 
1023  if (porosity_)
1024  {
1025  const volVectorField& U = lookupObject<volVectorField>(UName_);
1026  const volScalarField rho(this->rho());
1027  const volScalarField mu(this->mu());
1028 
1029  const HashTable<const porosityModel*> models =
1030  obr_.lookupClass<porosityModel>();
1031 
1032  if (models.empty())
1033  {
1035  << "Porosity effects requested, but no porosity models found "
1036  << "in the database"
1037  << endl;
1038  }
1039 
1040  forAllConstIters(models, iter)
1041  {
1042  // Non-const access required if mesh is changing
1043  porosityModel& pm = const_cast<porosityModel&>(*iter());
1044 
1045  vectorField fPTot(pm.force(U, rho, mu));
1046 
1047  const labelList& cellZoneIDs = pm.cellZoneIDs();
1048 
1049  for (const label zonei : cellZoneIDs)
1050  {
1051  const cellZone& cZone = mesh_.cellZones()[zonei];
1052 
1053  const vectorField d(mesh_.C(), cZone);
1054  const vectorField fP(fPTot, cZone);
1055  const vectorField Md(d - origin);
1056 
1057  const vectorField fDummy(Md.size(), Zero);
1058 
1059  addToFields(cZone, Md, fDummy, fDummy, fP);
1060 
1061  applyBins(Md, fDummy, fDummy, fP, d);
1062  }
1063  }
1064  }
1065 
1069  Pstream::listCombineScatter(moment_);
1070 }
1071 
1072 
1074 {
1075  return sum(force_[0]) + sum(force_[1]) + sum(force_[2]);
1076 }
1077 
1078 
1080 {
1081  return sum(moment_[0]) + sum(moment_[1]) + sum(moment_[2]);
1082 }
1083 
1084 
1086 {
1087  calcForcesMoment();
1088 
1089  if (Pstream::master())
1090  {
1091  createFiles();
1092 
1093  writeForces();
1094 
1095  writeBins();
1096 
1097  Log << endl;
1098  }
1099 
1100  // Write state/results information
1101  setResult("normalForce", sum(force_[0]));
1102  setResult("tangentialForce", sum(force_[1]));
1103  setResult("porousForce", sum(force_[2]));
1104 
1105  setResult("normalMoment", sum(moment_[0]));
1106  setResult("tangentialMoment", sum(moment_[1]));
1107  setResult("porousMoment", sum(moment_[2]));
1108 
1109  return true;
1110 }
1111 
1112 
1114 {
1115  if (writeFields_)
1116  {
1117  lookupObject<volVectorField>(scopedName("force")).write();
1118  lookupObject<volVectorField>(scopedName("moment")).write();
1119  }
1120 
1121  return true;
1122 }
1123 
1124 
1125 // ************************************************************************* //
Foam::functionObjects::forces::forceEff
virtual vector forceEff() const
Return the total force.
Definition: forces.C:1073
Foam::IOobject::NO_WRITE
Definition: IOobject.H:195
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:169
p
volScalarField & p
Definition: createFieldRefs.H:8
Log
#define Log
Definition: PDRblock.C:35
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:397
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
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:603
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
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:52
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimDensity
const dimensionSet dimDensity
Foam::functionObjects::forces::write
virtual bool write()
Write the forces.
Definition: forces.C:1113
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
Foam::functionObjects::forces::momentBinFilePtr_
autoPtr< OFstream > momentBinFilePtr_
Moment bins.
Definition: forces.H:264
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:97
turbulentTransportModel.H
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::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::dictionary::get
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:107
porosityModel.H
Foam::writeHeader
static void writeHeader(Ostream &os, const word &fieldName)
Definition: rawSurfaceWriterImpl.C:66
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:560
Foam::functionObjects::forces::writeBins
void writeBins()
Write binned data.
Definition: forces.C:676
Foam::dimForce
const dimensionSet dimForce
Foam::Vector::normalise
Vector< Cmpt > & normalise()
Normalise the vector by its magnitude.
Definition: VectorI.H:123
Foam::functionObjects::pressure
Provides several methods to convert an input pressure field into derived forms, including:
Definition: pressure.H:349
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:88
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:153
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
Foam::functionObjects::forces::forceFilePtr_
autoPtr< OFstream > forceFilePtr_
Forces.
Definition: forces.H:255
Foam::functionObjects::forces::momentFilePtr_
autoPtr< OFstream > momentFilePtr_
Moments.
Definition: forces.H:258
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
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:951
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:241
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:433
Foam::functionObjects::writeFile::read
virtual bool read(const dictionary &dict)
Read.
Definition: writeFile.C:213
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:1079
Foam::functionObjects::forces::resetFields
void resetFields()
Reset the fields prior to accumulation of force/moments.
Definition: forces.C:316
Foam::Field< vector >
Foam::functionObjects::forces::nBin_
label nBin_
Number of bins.
Definition: forces.H:303
Foam::functionObjects::forces::initialise
void initialise()
Initialise the fields.
Definition: forces.C:196
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
Foam::functionObjects::forces::devRhoReff
tmp< volSymmTensorField > devRhoReff() const
Return the effective viscous stress (laminar + turbulent).
Definition: forces.C:342
Foam::coordinateSystem::New
static autoPtr< coordinateSystem > New(word modelType, const objectRegistry &obr, const dictionary &dict)
Definition: coordinateSystemNew.C:84
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::functionObjects::forces::read
virtual bool read(const dictionary &)
Read the forces data.
Definition: forces.C:799
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::functionObjects::forces::forceBinFilePtr_
autoPtr< OFstream > forceBinFilePtr_
Force bins.
Definition: forces.H:261
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:510
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:123
Foam::functionObjects::forces::execute
virtual bool execute()
Execute, currently does nothing.
Definition: forces.C:1085
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:173
os
OBJstream os(runTime.globalPath()/outputName)
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:74
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:632
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::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::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:453
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:473
Foam::tab
constexpr char tab
Definition: Ostream.H:403
Foam::nl
constexpr char nl
Definition: Ostream.H:404
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::polyPatch::faceCentres
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:321
Foam::functionObjects::writeFile::createFile
virtual autoPtr< OFstream > createFile(const word &name, scalar timeValue) const
Return autoPtr to a new file for a given time.
Definition: writeFile.C:83
f
labelList f(nPoints)
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
Foam::coordSystem::cartesian
A Cartesian coordinate system.
Definition: cartesianCS.H:69
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:51
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::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::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(ObukhovLength, 0)
Foam::word::null
static const word null
An empty word.
Definition: word.H:80
Foam::functionObjects::readFields
Reads fields from the time directories and adds them to the mesh database for further post-processing...
Definition: readFields.H:155
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::VectorSpace< Vector< scalar >, scalar, 3 >::zero
static const Vector< scalar > zero
Definition: VectorSpace.H:115
Foam::functionObjects::writeFile
Base class for writing single files from the function objects.
Definition: writeFile.H:119
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:473
Foam::HashTable::empty
bool empty() const noexcept
True if the hash table is empty.
Definition: HashTableI.H:59
Foam::functionObjects::writeFile::writeToFile
virtual bool writeToFile() const
Flag to allow writing to file.
Definition: writeFile.C:253
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:188
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
Foam::functionObjects::forces::forces
forces(const forces &)=delete
No copy construct.
turbulentFluidThermoModel.H
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
Foam::dictionary::dictName
word dictName() const
The local dictionary name (final part of scoped name)
Definition: dictionaryI.H:60