surfaceFieldValue.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-2017 OpenFOAM Foundation
9  Copyright (C) 2017-2021 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 "surfaceFieldValue.H"
30 #include "fvMesh.H"
31 #include "emptyPolyPatch.H"
32 #include "coupledPolyPatch.H"
33 #include "sampledSurface.H"
34 #include "mergePoints.H"
36 #include "PatchTools.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43 namespace functionObjects
44 {
45 namespace fieldValues
46 {
47  defineTypeNameAndDebug(surfaceFieldValue, 0);
48  addToRunTimeSelectionTable(fieldValue, surfaceFieldValue, runTime);
49  addToRunTimeSelectionTable(functionObject, surfaceFieldValue, dictionary);
50 }
51 }
52 }
53 
54 
55 const Foam::Enum
56 <
58 >
60 ({
61  { regionTypes::stFaceZone, "faceZone" },
62  { regionTypes::stPatch, "patch" },
63  { regionTypes::stObject, "functionObjectSurface" },
64  { regionTypes::stSampled, "sampledSurface" },
65 });
66 
67 
68 const Foam::Enum
69 <
71 >
73 ({
74  // Normal operations
75  { operationType::opNone, "none" },
76  { operationType::opMin, "min" },
77  { operationType::opMax, "max" },
78  { operationType::opSum, "sum" },
79  { operationType::opSumMag, "sumMag" },
80  { operationType::opSumDirection, "sumDirection" },
81  { operationType::opSumDirectionBalance, "sumDirectionBalance" },
82  { operationType::opAverage, "average" },
83  { operationType::opAreaAverage, "areaAverage" },
84  { operationType::opAreaIntegrate, "areaIntegrate" },
85  { operationType::opCoV, "CoV" },
86  { operationType::opAreaNormalAverage, "areaNormalAverage" },
87  { operationType::opAreaNormalIntegrate, "areaNormalIntegrate" },
88  { operationType::opUniformity, "uniformity" },
89 
90  // Using weighting
91  { operationType::opWeightedSum, "weightedSum" },
92  { operationType::opWeightedAverage, "weightedAverage" },
93  { operationType::opWeightedAreaAverage, "weightedAreaAverage" },
94  { operationType::opWeightedAreaIntegrate, "weightedAreaIntegrate" },
95  { operationType::opWeightedUniformity, "weightedUniformity" },
96 
97  // Using absolute weighting
98  { operationType::opAbsWeightedSum, "absWeightedSum" },
99  { operationType::opAbsWeightedAverage, "absWeightedAverage" },
100  { operationType::opAbsWeightedAreaAverage, "absWeightedAreaAverage" },
101  { operationType::opAbsWeightedAreaIntegrate, "absWeightedAreaIntegrate" },
102  { operationType::opAbsWeightedUniformity, "absWeightedUniformity" },
103 });
104 
105 const Foam::Enum
106 <
108 >
110 ({
111  { postOperationType::postOpNone, "none" },
112  { postOperationType::postOpMag, "mag" },
113  { postOperationType::postOpSqrt, "sqrt" },
114 });
115 
116 
117 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
118 
121 {
122  if (stObject == regionType_)
123  {
125  }
126 
127  return mesh_;
128 }
129 
130 
131 void Foam::functionObjects::fieldValues::surfaceFieldValue::setFaceZoneFaces()
132 {
133  // Indices for all matches, already sorted
134  const labelList zoneIds
135  (
136  mesh_.faceZones().indices(selectionNames_)
137  );
138 
139  // Total number of faces selected
140  label numFaces = 0;
141  for (const label zoneId : zoneIds)
142  {
143  numFaces += mesh_.faceZones()[zoneId].size();
144  }
145 
146  if (zoneIds.empty())
147  {
149  << type() << ' ' << name() << ": "
150  << regionTypeNames_[regionType_] << '(' << regionName_ << "):" << nl
151  << " No matching face zone(s): "
152  << flatOutput(selectionNames_) << nl
153  << " Known face zones: "
154  << flatOutput(mesh_.faceZones().names()) << nl
155  << exit(FatalError);
156  }
157 
158  // Could also check this
159  #if 0
160  if (!returnReduce(bool(numFaces), orOp<bool>()))
161  {
163  << type() << ' ' << name() << ": "
164  << regionTypeNames_[regionType_] << '(' << regionName_ << "):" << nl
165  << " The faceZone specification: "
166  << flatOutput(selectionNames_) << nl
167  << " resulted in 0 faces" << nl
168  << exit(FatalError);
169  }
170  #endif
171 
172  faceId_.resize(numFaces);
173  facePatchId_.resize(numFaces);
174  faceFlip_.resize(numFaces);
175 
176  numFaces = 0;
177 
178  for (const label zoneId : zoneIds)
179  {
180  const faceZone& fZone = mesh_.faceZones()[zoneId];
181 
182  forAll(fZone, i)
183  {
184  const label meshFacei = fZone[i];
185  const bool isFlip = fZone.flipMap()[i];
186 
187  // Internal faces
188  label faceId = meshFacei;
189  label facePatchId = -1;
190 
191  // Boundary faces
192  if (!mesh_.isInternalFace(meshFacei))
193  {
194  facePatchId = mesh_.boundaryMesh().whichPatch(meshFacei);
195  const polyPatch& pp = mesh_.boundaryMesh()[facePatchId];
196  const auto* cpp = isA<coupledPolyPatch>(pp);
197 
198  if (cpp)
199  {
200  faceId = (cpp->owner() ? pp.whichFace(meshFacei) : -1);
201  }
202  else if (!isA<emptyPolyPatch>(pp))
203  {
204  faceId = pp.whichFace(meshFacei);
205  }
206  else
207  {
208  faceId = -1;
209  facePatchId = -1;
210  }
211  }
212 
213  if (faceId >= 0)
214  {
215  faceId_[numFaces] = faceId;
216  facePatchId_[numFaces] = facePatchId;
217  faceFlip_[numFaces] = isFlip;
218 
219  ++numFaces;
220  }
221  }
222  }
223 
224  // Shrink to size used
225  faceId_.resize(numFaces);
226  facePatchId_.resize(numFaces);
227  faceFlip_.resize(numFaces);
228  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
229 }
230 
231 
232 void Foam::functionObjects::fieldValues::surfaceFieldValue::setPatchFaces()
233 {
234  // Patch indices for all matches
236 
237  // Total number of faces selected
238  label numFaces = 0;
239 
240  labelList selected
241  (
242  mesh_.boundaryMesh().patchSet
243  (
244  selectionNames_,
245  false // warnNotFound - we do that ourselves
246  ).sortedToc()
247  );
248 
249  DynamicList<label> bad;
250  for (const label patchi : selected)
251  {
252  const polyPatch& pp = mesh_.boundaryMesh()[patchi];
253 
254  if (isA<emptyPolyPatch>(pp))
255  {
256  bad.append(patchi);
257  }
258  else
259  {
260  numFaces += pp.size();
261  }
262  }
263 
264  if (bad.size())
265  {
266  label nGood = (selected.size() - bad.size());
267 
268  auto& os = (nGood > 0 ? WarningInFunction : FatalErrorInFunction);
269 
270  os << "Cannot sample an empty patch" << nl;
271 
272  for (const label patchi : bad)
273  {
274  os << " "
275  << mesh_.boundaryMesh()[patchi].name() << nl;
276  }
277 
278  if (nGood)
279  {
280  os << "No non-empty patches selected" << endl
281  << exit(FatalError);
282  }
283  else
284  {
285  os << "Selected " << nGood << " non-empty patches" << nl;
286  }
287 
288  patchIds.resize(nGood);
289  nGood = 0;
290  for (const label patchi : selected)
291  {
292  if (!bad.found(patchi))
293  {
294  patchIds[nGood] = patchi;
295  ++nGood;
296  }
297  }
298  }
299  else
300  {
301  patchIds = std::move(selected);
302  }
303 
304  if (patchIds.empty())
305  {
307  << type() << ' ' << name() << ": "
308  << regionTypeNames_[regionType_] << '(' << regionName_ << "):" << nl
309  << " No matching patch name(s): "
310  << flatOutput(selectionNames_) << nl
311  << " Known patch names:" << nl
312  << mesh_.boundaryMesh().names() << nl
313  << exit(FatalError);
314  }
315 
316  // Could also check this
317  #if 0
318  if (!returnReduce(bool(numFaces), orOp<bool>()))
319  {
321  << type() << ' ' << name() << ": "
322  << regionTypeNames_[regionType_] << '(' << regionName_ << "):" << nl
323  << " The patch specification: "
324  << flatOutput(selectionNames_) << nl
325  << " resulted in 0 faces" << nl
326  << exit(FatalError);
327  }
328  #endif
329 
330  faceId_.resize(numFaces);
331  facePatchId_.resize(numFaces);
332  faceFlip_.resize(numFaces);
333  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
334 
335  numFaces = 0;
336  for (const label patchi : patchIds)
337  {
338  const polyPatch& pp = mesh_.boundaryMesh()[patchi];
339  const label len = pp.size();
340 
341  SubList<label>(faceId_, len, numFaces) = identity(len);
342  SubList<label>(facePatchId_, len, numFaces) = patchi;
343  SubList<bool>(faceFlip_, len, numFaces) = false;
344 
345  numFaces += len;
346  }
347 }
348 
349 
350 void Foam::functionObjects::fieldValues::surfaceFieldValue::combineMeshGeometry
351 (
352  faceList& faces,
354 ) const
355 {
356  List<faceList> allFaces(Pstream::nProcs());
357  List<pointField> allPoints(Pstream::nProcs());
358 
359  {
360  IndirectList<face> selectedFaces(mesh_.faces(), labelList(faceId_));
361  labelList& meshFaceIds = selectedFaces.addressing();
362 
363  forAll(meshFaceIds, i)
364  {
365  const label patchi = facePatchId_[i];
366  if (patchi != -1)
367  {
368  meshFaceIds[i] += mesh_.boundaryMesh()[patchi].start();
369  }
370  }
371 
372  // Add local faces and points to the all* lists
373  uindirectPrimitivePatch pp(selectedFaces, mesh_.points());
374 
375  allFaces[Pstream::myProcNo()] = pp.localFaces();
376  allPoints[Pstream::myProcNo()] = pp.localPoints();
377  }
378 
379  Pstream::gatherList(allFaces);
381 
382  // Renumber and flatten
383  label nFaces = 0;
384  label nPoints = 0;
385  forAll(allFaces, proci)
386  {
387  nFaces += allFaces[proci].size();
388  nPoints += allPoints[proci].size();
389  }
390 
391  faces.resize(nFaces);
392  points.resize(nPoints);
393 
394  nFaces = 0;
395  nPoints = 0;
396 
397  // My data first
398  {
399  for (const face& f : allFaces[Pstream::myProcNo()])
400  {
401  faces[nFaces++] = offsetOp<face>()(f, nPoints);
402  }
403 
404  for (const point& p : allPoints[Pstream::myProcNo()])
405  {
406  points[nPoints++] = p;
407  }
408  }
409 
410  // Other proc data follows
411  forAll(allFaces, proci)
412  {
413  if (proci == Pstream::myProcNo())
414  {
415  continue;
416  }
417 
418  for (const face& f : allFaces[proci])
419  {
420  faces[nFaces++] = offsetOp<face>()(f, nPoints);
421  }
422 
423  for (const point& p : allPoints[proci])
424  {
425  points[nPoints++] = p;
426  }
427  }
428 
429  // Merge
430  labelList oldToNew;
431  pointField newPoints;
432  bool hasMerged = mergePoints
433  (
434  points,
435  SMALL,
436  false,
437  oldToNew,
438  newPoints
439  );
440 
441  if (hasMerged)
442  {
443  if (debug)
444  {
445  Pout<< "Merged from " << points.size()
446  << " down to " << newPoints.size() << " points" << endl;
447  }
448 
449  points.transfer(newPoints);
450  for (face& f : faces)
451  {
452  inplaceRenumber(oldToNew, f);
453  }
454  }
455 }
456 
457 
458 void Foam::functionObjects::fieldValues::surfaceFieldValue::
459 combineSurfaceGeometry
460 (
461  faceList& faces,
463 ) const
464 {
465  if (stObject == regionType_)
466  {
467  const polySurface& s = dynamicCast<const polySurface>(obr());
468 
469  if (Pstream::parRun())
470  {
471  // Dimension as fraction of surface
472  const scalar mergeDim = 1e-10*boundBox(s.points(), true).mag();
473 
474  labelList pointsMap;
475 
477  (
478  mergeDim,
479  primitivePatch(SubList<face>(s.faces()), s.points()),
480  points,
481  faces,
482  pointsMap
483  );
484  }
485  else
486  {
487  faces = s.faces();
488  points = s.points();
489  }
490  }
491  else if (sampledPtr_)
492  {
493  const sampledSurface& s = *sampledPtr_;
494 
495  if (Pstream::parRun())
496  {
497  // Dimension as fraction of mesh bounding box
498  const scalar mergeDim = 1e-10*mesh_.bounds().mag();
499 
500  labelList pointsMap;
501 
503  (
504  mergeDim,
505  primitivePatch(SubList<face>(s.faces()), s.points()),
506  points,
507  faces,
508  pointsMap
509  );
510  }
511  else
512  {
513  faces = s.faces();
514  points = s.points();
515  }
516  }
517 }
518 
519 
520 Foam::scalar
521 Foam::functionObjects::fieldValues::surfaceFieldValue::totalArea() const
522 {
523  scalar totalArea = 0;
524 
525  if (stObject == regionType_)
526  {
527  const polySurface& s = dynamicCast<const polySurface>(obr());
528 
529  totalArea = gSum(s.magSf());
530  }
531  else if (sampledPtr_)
532  {
533  totalArea = gSum(sampledPtr_->magSf());
534  }
535  else
536  {
537  totalArea = gSum(filterField(mesh_.magSf()));
538  }
539 
540  return totalArea;
541 }
542 
543 
544 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
545 
547 const noexcept
548 {
549  // Few operations do not require the Sf field
550  switch (operation_)
551  {
552  case opNone:
553  case opMin:
554  case opMax:
555  case opSum:
556  case opSumMag:
557  case opAverage:
558  return false;
559 
560  default:
561  return true;
562  }
563 }
564 
565 
567 {
568  if (sampledPtr_)
569  {
570  sampledPtr_->update();
571  }
572 
573  if (!needsUpdate_)
574  {
575  return false;
576  }
577 
578  switch (regionType_)
579  {
580  case stFaceZone:
581  {
582  setFaceZoneFaces();
583  break;
584  }
585  case stPatch:
586  {
587  setPatchFaces();
588  break;
589  }
590  case stObject:
591  {
592  const polySurface& s = dynamicCast<const polySurface>(obr());
593  nFaces_ = returnReduce(s.size(), sumOp<label>());
594  break;
595  }
596  case stSampled:
597  {
598  nFaces_ = returnReduce(sampledPtr_->faces().size(), sumOp<label>());
599  break;
600  }
601 
602  // Compiler warning if we forgot an enumeration
603  }
604 
605  if (nFaces_ == 0)
606  {
608  << type() << ' ' << name() << ": "
609  << regionTypeNames_[regionType_] << '(' << regionName_ << "):" << nl
610  << " Region has no faces" << exit(FatalError);
611  }
612 
613  totalArea_ = totalArea();
614 
615  Log << " total faces = " << nFaces_ << nl
616  << " total area = " << totalArea_ << nl
617  << endl;
618 
619  writeFileHeader(file());
620 
621  needsUpdate_ = false;
622  return true;
623 }
624 
625 
627 (
628  Ostream& os
629 )
630 {
631  if (canWriteHeader() && (operation_ != opNone))
632  {
633  writeCommented(os, "Region type : ");
634  os << regionTypeNames_[regionType_] << ' ' << regionName_ << nl;
635 
636  writeHeaderValue(os, "Faces", nFaces_);
637  writeHeaderValue(os, "Area", totalArea_);
638  writeHeaderValue(os, "Scale factor", scaleFactor_);
639 
640  if (weightFieldNames_.size())
641  {
642  writeHeaderValue
643  (
644  os,
645  "Weight field",
646  flatOutput(weightFieldNames_, FlatOutput::BareComma{})
647  );
648  }
649 
650  writeCommented(os, "Time");
651  if (writeArea_)
652  {
653  os << tab << "Area";
654  }
655 
656  // TBD: add in postOperation information?
657 
658  for (const word& fieldName : fields_)
659  {
660  os << tab << operationTypeNames_[operation_]
661  << '(' << fieldName << ')';
662  }
663 
664  os << endl;
665  }
666 
667  writtenHeader_ = true;
668 }
669 
670 
671 template<>
672 Foam::scalar
674 (
675  const Field<scalar>& values,
676  const vectorField& Sf,
677  const scalarField& weightField
678 ) const
679 {
680  switch (operation_)
681  {
682  case opSumDirection:
683  {
684  const vector n(dict_.get<vector>("direction"));
685  return gSum(pos0(values*(Sf & n))*mag(values));
686  }
687  case opSumDirectionBalance:
688  {
689  const vector n(dict_.get<vector>("direction"));
690  const scalarField nv(values*(Sf & n));
691 
692  return gSum(pos0(nv)*mag(values) - neg(nv)*mag(values));
693  }
694 
695  case opUniformity:
696  case opWeightedUniformity:
697  case opAbsWeightedUniformity:
698  {
699  const scalar areaTotal = gSum(mag(Sf));
700  tmp<scalarField> areaVal(values * mag(Sf));
701 
702  scalar mean, numer;
703 
704  if (is_weightedOp() && canWeight(weightField))
705  {
706  // Weighted quantity = (Weight * phi * dA)
707 
708  tmp<scalarField> weight
709  (
710  weightingFactor(weightField, is_magOp())
711  );
712 
713  // Mean weighted value (area-averaged)
714  mean = gSum(weight()*areaVal()) / areaTotal;
715 
716  // Abs. deviation from weighted mean value
717  numer = gSum(mag(weight*areaVal - (mean * mag(Sf))));
718  }
719  else
720  {
721  // Unweighted quantity = (1 * phi * dA)
722 
723  // Mean value (area-averaged)
724  mean = gSum(areaVal()) / areaTotal;
725 
726  // Abs. deviation from mean value
727  numer = gSum(mag(areaVal - (mean * mag(Sf))));
728  }
729 
730  // Uniformity index
731  const scalar ui = 1 - numer/(2*mag(mean*areaTotal) + ROOTVSMALL);
732 
733  return min(max(ui, 0), 1);
734  }
735 
736  default:
737  {
738  // Fall through to other operations
739  return processSameTypeValues(values, Sf, weightField);
740  }
741  }
742 }
743 
744 
745 template<>
748 (
749  const Field<vector>& values,
750  const vectorField& Sf,
751  const scalarField& weightField
752 ) const
753 {
754  switch (operation_)
755  {
756  case opSumDirection:
757  {
758  const vector n(dict_.get<vector>("direction").normalise());
759 
760  const scalarField nv(n & values);
761  return gSum(pos0(nv)*n*(nv));
762  }
763  case opSumDirectionBalance:
764  {
765  const vector n(dict_.get<vector>("direction").normalise());
766 
767  const scalarField nv(n & values);
768  return gSum(pos0(nv)*n*(nv));
769  }
770  case opAreaNormalAverage:
771  {
772  const scalar val = gSum(values & Sf)/gSum(mag(Sf));
773  return vector(val, 0, 0);
774  }
775  case opAreaNormalIntegrate:
776  {
777  const scalar val = gSum(values & Sf);
778  return vector(val, 0, 0);
779  }
780 
781  case opUniformity:
782  case opWeightedUniformity:
783  case opAbsWeightedUniformity:
784  {
785  const scalar areaTotal = gSum(mag(Sf));
786  tmp<scalarField> areaVal(values & Sf);
787 
788  scalar mean, numer;
789 
790  if (is_weightedOp() && canWeight(weightField))
791  {
792  // Weighted quantity = (Weight * phi . dA)
793 
794  tmp<scalarField> weight
795  (
796  weightingFactor(weightField, is_magOp())
797  );
798 
799  // Mean weighted value (area-averaged)
800  mean = gSum(weight()*areaVal()) / areaTotal;
801 
802  // Abs. deviation from weighted mean value
803  numer = gSum(mag(weight*areaVal - (mean * mag(Sf))));
804  }
805  else
806  {
807  // Unweighted quantity = (1 * phi . dA)
808 
809  // Mean value (area-averaged)
810  mean = gSum(areaVal()) / areaTotal;
811 
812  // Abs. deviation from mean value
813  numer = gSum(mag(areaVal - (mean * mag(Sf))));
814  }
815 
816  // Uniformity index
817  const scalar ui = 1 - numer/(2*mag(mean*areaTotal) + ROOTVSMALL);
818 
819  return vector(min(max(ui, 0), 1), 0, 0);
820  }
821 
822  default:
823  {
824  // Fall through to other operations
825  return processSameTypeValues(values, Sf, weightField);
826  }
827  }
828 }
829 
830 
831 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
832 
833 template<>
836 (
837  const Field<scalar>& weightField,
838  const bool useMag
839 )
840 {
841  if (useMag)
842  {
843  return mag(weightField);
844  }
845 
846  // pass through
847  return weightField;
848 }
849 
850 
851 template<>
854 (
855  const Field<scalar>& weightField,
856  const vectorField& Sf,
857  const bool useMag
858 )
859 {
860  // scalar * unit-normal
861 
862  // Can skip this check - already used canWeight()
868 
869  if (useMag)
870  {
871  return mag(weightField);
872  }
873 
874  // pass through
875  return weightField;
876 }
877 
878 
879 template<>
882 (
883  const Field<scalar>& weightField,
884  const vectorField& Sf,
885  const bool useMag
886 )
887 {
888  // scalar * Area
889 
890  // Can skip this check - already used canWeight()
896 
897  if (useMag)
898  {
899  return mag(weightField * mag(Sf));
900  }
901 
902  return (weightField * mag(Sf));
903 }
904 
905 
906 template<>
909 (
910  const Field<vector>& weightField,
911  const vectorField& Sf,
912  const bool useMag
913 )
914 {
915  // vector (dot) unit-normal
916 
917  // Can skip this check - already used canWeight()
923 
924  const label len = weightField.size();
925 
926  auto tresult = tmp<scalarField>::New(weightField.size());
927  auto& result = tresult.ref();
928 
929  for (label facei=0; facei < len; ++facei)
930  {
931  const vector unitNormal(normalised(Sf[facei]));
932  result[facei] = (weightField[facei] & unitNormal);
933  }
934 
935  if (useMag)
936  {
937  for (scalar& val : result)
938  {
939  val = mag(val);
940  }
941  }
942 
943  return tresult;
944 }
945 
946 
947 template<>
950 (
951  const Field<vector>& weightField,
952  const vectorField& Sf,
953  const bool useMag
954 )
955 {
956  // vector (dot) Area
957 
958  // Can skip this check - already used canWeight()
964 
965  if (useMag)
966  {
967  return mag(weightField & Sf);
968  }
969 
970  return (weightField & Sf);
971 }
972 
973 
974 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
975 
977 (
978  const word& name,
979  const Time& runTime,
980  const dictionary& dict
981 )
982 :
983  fieldValue(name, runTime, dict, typeName),
984  regionType_(regionTypeNames_.get("regionType", dict)),
985  operation_(operationTypeNames_.get("operation", dict)),
986  postOperation_
987  (
988  postOperationTypeNames_.getOrDefault
989  (
990  "postOperation",
991  dict,
992  postOperationType::postOpNone,
993  true // Failsafe behaviour
994  )
995  ),
996  needsUpdate_(true),
997  writeArea_(false),
998  selectionNames_(),
999  weightFieldNames_(),
1000  totalArea_(0),
1001  nFaces_(0),
1002  faceId_(),
1003  facePatchId_(),
1004  faceFlip_()
1005 {
1006  read(dict);
1007 }
1008 
1009 
1012  const word& name,
1013  const objectRegistry& obr,
1014  const dictionary& dict
1015 )
1016 :
1017  fieldValue(name, obr, dict, typeName),
1018  regionType_(regionTypeNames_.get("regionType", dict)),
1019  operation_(operationTypeNames_.get("operation", dict)),
1020  postOperation_
1021  (
1022  postOperationTypeNames_.getOrDefault
1023  (
1024  "postOperation",
1025  dict,
1026  postOperationType::postOpNone,
1027  true // Failsafe behaviour
1028  )
1029  ),
1030  needsUpdate_(true),
1031  writeArea_(false),
1032  selectionNames_(),
1033  weightFieldNames_(),
1034  totalArea_(0),
1035  nFaces_(0),
1036  faceId_(),
1037  facePatchId_(),
1038  faceFlip_()
1039 {
1040  read(dict);
1041 }
1042 
1043 
1044 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
1045 
1046 // Needs completed sampledSurface, surfaceWriter
1048 {}
1049 
1050 
1051 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1052 
1055  const dictionary& dict
1056 )
1057 {
1059 
1060  needsUpdate_ = true;
1061  writeArea_ = dict.getOrDefault("writeArea", false);
1062  weightFieldNames_.clear();
1063  // future?
1064  // sampleFaceScheme_ = dict.getOrDefault<word>("sampleScheme", "cell");
1065 
1066  totalArea_ = 0;
1067  nFaces_ = 0;
1068  faceId_.clear();
1069  facePatchId_.clear();
1070  faceFlip_.clear();
1071  sampledPtr_.reset(nullptr);
1072  surfaceWriterPtr_.reset(nullptr);
1073 
1074  // Can have "name" (word) and/or "names" (wordRes)
1075  //
1076  // If "names" exists AND contains a literal (non-regex) that can be used
1077  // as a suitable value for "name", the "name" entry becomes optional.
1078 
1079  regionName_.clear();
1080  selectionNames_.clear();
1081 
1082  {
1083  dict.readIfPresent("names", selectionNames_);
1084 
1085  for (const auto& item : selectionNames_)
1086  {
1087  if (item.isLiteral())
1088  {
1089  regionName_ = item;
1090  break;
1091  }
1092  }
1093 
1094  // Mandatory if we didn't pick up a value from selectionNames_
1095  dict.readEntry
1096  (
1097  "name",
1098  regionName_,
1100  regionName_.empty()
1101  );
1102 
1103  // Ensure there is always content for selectionNames_
1104  if (selectionNames_.empty())
1105  {
1106  selectionNames_.resize(1);
1107  selectionNames_.first() = regionName_;
1108  }
1109  }
1110 
1111 
1112  // Create sampled surface, but leave 'expired' (ie, no update) since it
1113  // may depend on fields or data that do not yet exist
1114  if (stSampled == regionType_)
1115  {
1116  sampledPtr_ = sampledSurface::New
1117  (
1118  name(),
1119  mesh_,
1120  dict.subDict("sampledSurfaceDict")
1121  );
1122 
1123  // Internal consistency. Want face values, never point values!
1124  sampledPtr_->isPointData(false);
1125  }
1126 
1127  Info<< type() << ' ' << name() << ':' << nl
1128  << " operation = ";
1129 
1130  if (postOperation_ != postOpNone)
1131  {
1132  Info<< postOperationTypeNames_[postOperation_] << '('
1133  << operationTypeNames_[operation_] << ')' << nl;
1134  }
1135  else
1136  {
1137  Info<< operationTypeNames_[operation_] << nl;
1138  }
1139 
1140  if (is_weightedOp())
1141  {
1142  // Can have "weightFields" or "weightField"
1143 
1144  bool missing = true;
1145  if (dict.readIfPresent("weightFields", weightFieldNames_))
1146  {
1147  missing = false;
1148  }
1149  else
1150  {
1151  weightFieldNames_.resize(1);
1152 
1153  if (dict.readIfPresent("weightField", weightFieldNames_.first()))
1154  {
1155  missing = false;
1156  if ("none" == weightFieldNames_.first())
1157  {
1158  // "none" == no weighting
1159  weightFieldNames_.clear();
1160  }
1161  }
1162  }
1163 
1164  if (missing)
1165  {
1166  // Suggest possible alternative unweighted operation?
1168  << "The '" << operationTypeNames_[operation_]
1169  << "' operation is missing a weightField." << nl
1170  << "Either provide the weightField, "
1171  << "use weightField 'none' to suppress weighting," << nl
1172  << "or use a different operation."
1173  << exit(FatalIOError);
1174  }
1175 
1176  Info<< " weight field = ";
1177  if (weightFieldNames_.empty())
1178  {
1179  Info<< "none" << nl;
1180  }
1181  else
1182  {
1183  Info<< flatOutput(weightFieldNames_) << nl;
1184  }
1185  }
1186 
1187  if (stSampled == regionType_ && sampledPtr_)
1188  {
1189  Info<< " sampled surface: ";
1190  sampledPtr_->print(Info, 0);
1191  Info<< nl;
1192  }
1193 
1194  if (writeFields_)
1195  {
1196  const word formatName(dict.get<word>("surfaceFormat"));
1197 
1198  surfaceWriterPtr_.reset
1199  (
1201  (
1202  formatName,
1203  dict.subOrEmptyDict("formatOptions").subOrEmptyDict(formatName)
1204  )
1205  );
1206 
1207  // Propagate field counts (per surface)
1208  surfaceWriterPtr_->nFields(fields_.size());
1209 
1210  if (debug)
1211  {
1212  surfaceWriterPtr_->verbose(true);
1213  }
1214 
1215  if (surfaceWriterPtr_->enabled())
1216  {
1217  Info<< " surfaceFormat = " << formatName << nl;
1218  }
1219  else
1220  {
1221  surfaceWriterPtr_->clear();
1222  }
1223  }
1224 
1225  Info<< nl << endl;
1226 
1227  return true;
1228 }
1229 
1230 
1232 {
1233  if (needsUpdate_ || operation_ != opNone)
1234  {
1236  }
1237 
1238  update();
1239 
1240  if (operation_ != opNone)
1241  {
1242  writeCurrentTime(file());
1243  }
1244 
1245  if (writeArea_)
1246  {
1247  totalArea_ = totalArea();
1248  Log << " total area = " << totalArea_ << endl;
1249 
1250  if (operation_ != opNone && Pstream::master())
1251  {
1252  file() << tab << totalArea_;
1253  }
1254  }
1255 
1256  // Many operations use the Sf field
1257  vectorField Sf;
1258  if (usesSf())
1259  {
1260  if (stObject == regionType_)
1261  {
1262  const polySurface& s = dynamicCast<const polySurface>(obr());
1263  Sf = s.Sf();
1264  }
1265  else if (sampledPtr_)
1266  {
1267  Sf = sampledPtr_->Sf();
1268  }
1269  else
1270  {
1271  Sf = filterField(mesh_.Sf());
1272  }
1273  }
1274 
1275  // Faces and points for surface format (if specified)
1276  faceList faces;
1278 
1279  if (surfaceWriterPtr_)
1280  {
1281  if (withTopologicalMerge())
1282  {
1283  combineMeshGeometry(faces, points);
1284  }
1285  else
1286  {
1287  combineSurfaceGeometry(faces, points);
1288  }
1289  }
1290 
1291 
1292  // Check availability and type of weight field
1293  // Only support a few weight types:
1294  // scalar: 0-N fields
1295  // vector: 0-1 fields
1296 
1297  // Default is a zero-size scalar weight field (ie, weight = 1)
1298  scalarField scalarWeights;
1299  vectorField vectorWeights;
1300 
1301  for (const word& weightName : weightFieldNames_)
1302  {
1303  if (validField<scalar>(weightName))
1304  {
1305  tmp<scalarField> tfld = getFieldValues<scalar>(weightName, true);
1306 
1307  if (scalarWeights.empty())
1308  {
1309  scalarWeights = tfld;
1310  }
1311  else
1312  {
1313  scalarWeights *= tfld;
1314  }
1315  }
1316  else if (validField<vector>(weightName))
1317  {
1318  tmp<vectorField> tfld = getFieldValues<vector>(weightName, true);
1319 
1320  if (vectorWeights.empty())
1321  {
1322  vectorWeights = tfld;
1323  }
1324  else
1325  {
1327  << "weightField " << weightName
1328  << " - only one vector weight field allowed. " << nl
1329  << "weights: " << flatOutput(weightFieldNames_) << nl
1330  << abort(FatalError);
1331  }
1332  }
1333  else if (weightName != "none")
1334  {
1335  // Silently ignore "none", flag everything else as an error
1336 
1337  // TBD: treat missing "rho" like incompressible with rho = 1
1338  // and/or provided rhoRef value
1339 
1341  << "weightField " << weightName
1342  << " not found or an unsupported type" << nl
1343  << abort(FatalError);
1344  }
1345  }
1346 
1347 
1348  // Process the fields
1349  if (returnReduce(!vectorWeights.empty(), orOp<bool>()))
1350  {
1351  if (scalarWeights.size())
1352  {
1353  vectorWeights *= scalarWeights;
1354  }
1355 
1356  writeAll(Sf, vectorWeights, points, faces);
1357  }
1358  else
1359  {
1360  writeAll(Sf, scalarWeights, points, faces);
1361  }
1362 
1363 
1364  if (operation_ != opNone)
1365  {
1366  file() << endl;
1367  Log << endl;
1368  }
1369 
1370 
1371  return true;
1372 }
1373 
1374 
1377  const mapPolyMesh& mpm
1378 )
1379 {
1380  needsUpdate_ = true;
1381 }
1382 
1383 
1386  const polyMesh& mesh
1387 )
1388 {
1389  needsUpdate_ = true;
1390 }
1391 
1392 
1393 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
runTime
engineTime & runTime
Definition: createEngineTime.H:13
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::sampledSurface::New
static autoPtr< sampledSurface > New(const word &name, const polyMesh &mesh, const dictionary &dict)
Return a reference to the selected surface.
Definition: sampledSurface.C:62
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:57
Log
#define Log
Definition: PDRblock.C:35
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::functionObjects::fieldValue::read
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: fieldValue.C:89
Foam::functionObjects::fieldValues::defineTypeNameAndDebug
defineTypeNameAndDebug(multiFieldValue, 0)
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:94
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
uindirectPrimitivePatch.H
update
mesh update()
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::functionObjects::fieldValues::surfaceFieldValue::operationType
operationType
Operation type enumeration.
Definition: surfaceFieldValue.H:482
Foam::functionObjects::fieldValues::surfaceFieldValue::read
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: surfaceFieldValue.C:1054
Foam::functionObjects::fieldValues::surfaceFieldValue::postOperationType
postOperationType
Post-operation type enumeration.
Definition: surfaceFieldValue.H:554
PatchTools.H
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
Foam::functionObjects::timeFunctionObject::storedObjects
objectRegistry & storedObjects()
Definition: timeFunctionObject.C:63
Foam::PatchTools::gatherAndMerge
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< FaceList, PointField > &p, Field< typename PrimitivePatch< FaceList, PointField >::point_type > &mergedPoints, List< typename PrimitivePatch< FaceList, PointField >::face_type > &mergedFaces, labelList &pointMergeMap)
Gather points and faces onto master and merge into single patch.
Definition: PatchToolsGatherAndMerge.C:38
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::UPstream::master
static bool master(const label communicator=worldComm)
Am I the master process.
Definition: UPstream.H:457
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::pos0
dimensionedScalar pos0(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:188
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:594
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::Vector::normalise
Vector< Cmpt > & normalise()
Normalise the vector by its magnitude.
Definition: VectorI.H:123
Foam::functionObjects::fieldValues::surfaceFieldValue::regionTypeNames_
static const Enum< regionTypes > regionTypeNames_
Region type names.
Definition: surfaceFieldValue.H:470
Foam::functionObjects::fieldValues::surfaceFieldValue::regionType_
regionTypes regionType_
Region type.
Definition: surfaceFieldValue.H:598
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::inplaceRenumber
void inplaceRenumber(const labelUList &oldToNew, IntListType &input)
Inplace renumber the values (not the indices) of a list.
Definition: ListOpsTemplates.C:61
Foam::functionObjects::fieldValues::surfaceFieldValue::~surfaceFieldValue
virtual ~surfaceFieldValue()
Destructor.
Definition: surfaceFieldValue.C:1047
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::sumOp
Definition: ops.H:213
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::functionObjects::fieldValues::surfaceFieldValue::surfaceFieldValue
surfaceFieldValue(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
Definition: surfaceFieldValue.C:977
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::functionObjects::fieldValues::surfaceFieldValue::stObject
Calculate with function object surface.
Definition: surfaceFieldValue.H:465
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::polySurface
A surface mesh consisting of general polygon faces and capable of holding fields.
Definition: polySurface.H:67
Foam::functionObjects::fieldValues::surfaceFieldValue::update
bool update()
Update the surface and surface information as required.
Definition: surfaceFieldValue.C:566
coupledPolyPatch.H
patchIds
labelList patchIds
Definition: convertProcessorPatches.H:67
Foam::functionObjects::fieldValue::regionName_
word regionName_
Name of region (patch, zone, etc.)
Definition: fieldValue.H:132
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
faceId
label faceId(-1)
sampledSurface.H
Foam::functionObjects::fieldValues::surfaceFieldValue::updateMesh
virtual void updateMesh(const mapPolyMesh &mpm)
Update for changes of mesh.
Definition: surfaceFieldValue.C:1376
Foam::functionObjects::fieldValues::surfaceFieldValue::regionTypes
regionTypes
Region type enumeration.
Definition: surfaceFieldValue.H:461
Foam::functionObjects::fieldValues::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, multiFieldValue, dictionary)
Foam::objectRegistry::lookupObject
const Type & lookupObject(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:434
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
dict
dictionary dict
Definition: searchingEngine.H:14
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
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::flatOutput
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:216
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:144
Foam::functionObjects::fieldValues::surfaceFieldValue::obr
const objectRegistry & obr() const
The volume mesh or surface registry being used.
Definition: surfaceFieldValue.C:120
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::functionObjects::fieldValues::surfaceFieldValue::weightingFactor
static tmp< scalarField > weightingFactor(const Field< WeightType > &weightField, const bool useMag)
Weighting factor.
emptyPolyPatch.H
Foam::functionObjects::fieldValues::surfaceFieldValue::movePoints
virtual void movePoints(const polyMesh &mesh)
Update for changes of mesh.
Definition: surfaceFieldValue.C:1385
Foam::mergePoints
label mergePoints(const PointList &points, const scalar mergeTol, const bool verbose, labelList &pointMap, typename PointList::const_reference origin=PointList::value_type::zero)
Sorts and merges points. All points closer than/equal mergeTol get merged.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::normalised
VectorSpace< Form, Cmpt, Ncmpts > normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:487
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::Pstream::gatherList
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
Definition: gatherScatterList.C:52
Foam::tab
constexpr char tab
Definition: Ostream.H:403
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=worldComm)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:463
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::functionObjects::fieldValues::surfaceFieldValue::operationTypeNames_
static const Enum< operationType > operationTypeNames_
Operation type names.
Definition: surfaceFieldValue.H:550
f
labelList f(nPoints)
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::surfaceWriter::New
static autoPtr< surfaceWriter > New(const word &writeType)
Return a reference to the selected surfaceWriter.
Definition: surfaceWriter.C:64
Foam::FlatOutput::BareComma
Surround with '\0' and '\0' separate with ','.
Definition: FlatOutput.H:84
Foam::Vector< scalar >
Foam::UPstream::parRun
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
Foam::List< label >
Foam::functionObjects::fieldValue
Intermediate class for handling field value-based function objects.
Definition: fieldValue.H:119
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::functionObjects::fvMeshFunctionObject::mesh_
const fvMesh & mesh_
Reference to the fvMesh.
Definition: fvMeshFunctionObject.H:73
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::functionObjects::fieldValues::surfaceFieldValue::writeFileHeader
virtual void writeFileHeader(Ostream &os)
Output file header information.
Definition: surfaceFieldValue.C:627
Foam::identity
labelList identity(const label len, label start=0)
Create identity map of the given length with (map[i] == i)
Definition: labelList.C:38
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
mergePoints.H
Merge points. See below.
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:161
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::primitivePatch
PrimitivePatch< SubList< face >, const pointField & > primitivePatch
A PrimitivePatch with a SubList addressing for the faces, const reference for the point field.
Definition: primitivePatch.H:51
surfaceFieldValue.H
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Foam::DelaunayMeshTools::allPoints
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Foam::uindirectPrimitivePatch
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
Definition: uindirectPrimitivePatch.H:49
Foam::keyType::LITERAL
String literal.
Definition: keyType.H:81
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::functionObjects::fieldValues::surfaceFieldValue::areaWeightingFactor
static tmp< scalarField > areaWeightingFactor(const Field< WeightType > &weightField, const vectorField &Sf, const bool useMag)
Weighting factor, weight field with area factor.
Foam::functionObjects::fieldValue::write
virtual bool write()
Write.
Definition: fieldValue.C:115
Foam::point
vector point
Point is a vector.
Definition: point.H:43
Foam::functionObjects::fieldValues::surfaceFieldValue::write
virtual bool write()
Calculate and write.
Definition: surfaceFieldValue.C:1231
Foam::orOp
Definition: ops.H:234
Foam::functionObjects::fieldValues::surfaceFieldValue::processValues
Type processValues(const Field< Type > &values, const vectorField &Sf, const Field< WeightType > &weightField) const
Apply the 'operation' to the values. Wrapper around.
Definition: surfaceFieldValueTemplates.C:318
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:199
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:328
Foam::functionObjects::fieldValues::surfaceFieldValue::postOperationTypeNames_
static const Enum< postOperationType > postOperationTypeNames_
Operation type names.
Definition: surfaceFieldValue.H:562
Foam::UPstream::nProcs
static label nProcs(const label communicator=worldComm)
Number of processes in parallel run, and 1 for serial run.
Definition: UPstream.H:445
Foam::functionObjects::fieldValues::surfaceFieldValue::usesSf
bool usesSf() const noexcept
True if the operation needs a surface Sf.
Definition: surfaceFieldValue.C:546