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-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 "surfaceFieldValue.H"
30 #include "fvMesh.H"
31 #include "emptyPolyPatch.H"
32 #include "coupledPolyPatch.H"
33 #include "sampledSurface.H"
34 #include "mergePoints.H"
35 #include "indirectPrimitivePatch.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, dictionary);
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::postOpSqrt, "sqrt" },
113 });
114 
115 
116 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
117 
120 {
121  if (stObject == regionType_)
122  {
124  }
125 
126  return mesh_;
127 }
128 
129 
130 void Foam::functionObjects::fieldValues::surfaceFieldValue::setFaceZoneFaces()
131 {
132  const label zoneId = mesh_.faceZones().findZoneID(regionName_);
133 
134  if (zoneId < 0)
135  {
137  << type() << " " << name() << ": "
138  << regionTypeNames_[regionType_] << '(' << regionName_ << "):" << nl
139  << " Unknown face zone name: " << regionName_
140  << ". Valid face zones are: " << mesh_.faceZones().names()
141  << nl << exit(FatalError);
142  }
143 
144  const faceZone& fZone = mesh_.faceZones()[zoneId];
145 
146  DynamicList<label> faceIds(fZone.size());
147  DynamicList<label> facePatchIds(fZone.size());
148  DynamicList<bool> faceFlip(fZone.size());
149 
150  forAll(fZone, i)
151  {
152  const label facei = fZone[i];
153 
154  label faceId = -1;
155  label facePatchId = -1;
156  if (mesh_.isInternalFace(facei))
157  {
158  faceId = facei;
159  facePatchId = -1;
160  }
161  else
162  {
163  facePatchId = mesh_.boundaryMesh().whichPatch(facei);
164  const polyPatch& pp = mesh_.boundaryMesh()[facePatchId];
165  if (isA<coupledPolyPatch>(pp))
166  {
167  if (refCast<const coupledPolyPatch>(pp).owner())
168  {
169  faceId = pp.whichFace(facei);
170  }
171  else
172  {
173  faceId = -1;
174  }
175  }
176  else if (!isA<emptyPolyPatch>(pp))
177  {
178  faceId = facei - pp.start();
179  }
180  else
181  {
182  faceId = -1;
183  facePatchId = -1;
184  }
185  }
186 
187  if (faceId >= 0)
188  {
189  faceIds.append(faceId);
190  facePatchIds.append(facePatchId);
191  faceFlip.append(fZone.flipMap()[i] ? true : false);
192  }
193  }
194 
195  faceId_.transfer(faceIds);
196  facePatchId_.transfer(facePatchIds);
197  faceFlip_.transfer(faceFlip);
198  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
199 
200  if (debug)
201  {
202  Pout<< "Original face zone size = " << fZone.size()
203  << ", new size = " << faceId_.size() << endl;
204  }
205 }
206 
207 
208 void Foam::functionObjects::fieldValues::surfaceFieldValue::setPatchFaces()
209 {
210  const label patchid = mesh_.boundaryMesh().findPatchID(regionName_);
211 
212  if (patchid < 0)
213  {
215  << type() << " " << name() << ": "
216  << regionTypeNames_[regionType_] << '(' << regionName_ << "):" << nl
217  << " Unknown patch name: " << regionName_
218  << ". Valid patch names are: "
219  << mesh_.boundaryMesh().names() << nl
220  << exit(FatalError);
221  }
222 
223  const polyPatch& pp = mesh_.boundaryMesh()[patchid];
224 
225  label nFaces = pp.size();
226  if (isA<emptyPolyPatch>(pp))
227  {
228  nFaces = 0;
229  }
230 
231  faceId_.setSize(nFaces);
232  facePatchId_.setSize(nFaces);
233  faceFlip_.setSize(nFaces);
234  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
235 
236  forAll(faceId_, facei)
237  {
238  faceId_[facei] = facei;
239  facePatchId_[facei] = patchid;
240  faceFlip_[facei] = false;
241  }
242 }
243 
244 
245 void Foam::functionObjects::fieldValues::surfaceFieldValue::combineMeshGeometry
246 (
247  faceList& faces,
249 ) const
250 {
251  List<faceList> allFaces(Pstream::nProcs());
252  List<pointField> allPoints(Pstream::nProcs());
253 
254  labelList globalFacesIs(faceId_);
255  forAll(globalFacesIs, i)
256  {
257  if (facePatchId_[i] != -1)
258  {
259  const label patchi = facePatchId_[i];
260  globalFacesIs[i] += mesh_.boundaryMesh()[patchi].start();
261  }
262  }
263 
264  // Add local faces and points to the all* lists
266  (
267  IndirectList<face>(mesh_.faces(), globalFacesIs),
268  mesh_.points()
269  );
270  allFaces[Pstream::myProcNo()] = pp.localFaces();
271  allPoints[Pstream::myProcNo()] = pp.localPoints();
272 
273  Pstream::gatherList(allFaces);
275 
276  // Renumber and flatten
277  label nFaces = 0;
278  label nPoints = 0;
279  forAll(allFaces, proci)
280  {
281  nFaces += allFaces[proci].size();
282  nPoints += allPoints[proci].size();
283  }
284 
285  faces.setSize(nFaces);
286  points.setSize(nPoints);
287 
288  nFaces = 0;
289  nPoints = 0;
290 
291  // My own data first
292  {
293  const faceList& fcs = allFaces[Pstream::myProcNo()];
294  for (const face& f : fcs)
295  {
296  face& newF = faces[nFaces++];
297  newF.setSize(f.size());
298  forAll(f, fp)
299  {
300  newF[fp] = f[fp] + nPoints;
301  }
302  }
303 
304  const pointField& pts = allPoints[Pstream::myProcNo()];
305  for (const point& pt : pts)
306  {
307  points[nPoints++] = pt;
308  }
309  }
310 
311  // Other proc data follows
312  forAll(allFaces, proci)
313  {
314  if (proci != Pstream::myProcNo())
315  {
316  const faceList& fcs = allFaces[proci];
317  for (const face& f : fcs)
318  {
319  face& newF = faces[nFaces++];
320  newF.setSize(f.size());
321  forAll(f, fp)
322  {
323  newF[fp] = f[fp] + nPoints;
324  }
325  }
326 
327  const pointField& pts = allPoints[proci];
328  for (const point& pt : pts)
329  {
330  points[nPoints++] = pt;
331  }
332  }
333  }
334 
335  // Merge
336  labelList oldToNew;
337  pointField newPoints;
338  bool hasMerged = mergePoints
339  (
340  points,
341  SMALL,
342  false,
343  oldToNew,
344  newPoints
345  );
346 
347  if (hasMerged)
348  {
349  if (debug)
350  {
351  Pout<< "Merged from " << points.size()
352  << " down to " << newPoints.size() << " points" << endl;
353  }
354 
355  points.transfer(newPoints);
356  for (face& f : faces)
357  {
358  inplaceRenumber(oldToNew, f);
359  }
360  }
361 }
362 
363 
364 void Foam::functionObjects::fieldValues::surfaceFieldValue::
365 combineSurfaceGeometry
366 (
367  faceList& faces,
369 ) const
370 {
371  if (stObject == regionType_)
372  {
373  const polySurface& s = dynamicCast<const polySurface>(obr());
374 
375  if (Pstream::parRun())
376  {
377  // Dimension as fraction of surface
378  const scalar mergeDim = 1e-10*boundBox(s.points(), true).mag();
379 
380  labelList pointsMap;
381 
383  (
384  mergeDim,
386  (
387  SubList<face>(s.faces(), s.faces().size()),
388  s.points()
389  ),
390  points,
391  faces,
392  pointsMap
393  );
394  }
395  else
396  {
397  faces = s.faces();
398  points = s.points();
399  }
400  }
401  else if (sampledPtr_.valid())
402  {
403  const sampledSurface& s = sampledPtr_();
404 
405  if (Pstream::parRun())
406  {
407  // Dimension as fraction of mesh bounding box
408  const scalar mergeDim = 1e-10*mesh_.bounds().mag();
409 
410  labelList pointsMap;
411 
413  (
414  mergeDim,
416  (
417  SubList<face>(s.faces(), s.faces().size()),
418  s.points()
419  ),
420  points,
421  faces,
422  pointsMap
423  );
424  }
425  else
426  {
427  faces = s.faces();
428  points = s.points();
429  }
430  }
431 }
432 
433 
434 Foam::scalar
435 Foam::functionObjects::fieldValues::surfaceFieldValue::totalArea() const
436 {
437  scalar totalArea = 0;
438 
439  if (stObject == regionType_)
440  {
441  const polySurface& s = dynamicCast<const polySurface>(obr());
442 
443  totalArea = gSum(s.magSf());
444  }
445  else if (sampledPtr_.valid())
446  {
447  totalArea = gSum(sampledPtr_().magSf());
448  }
449  else
450  {
451  totalArea = gSum(filterField(mesh_.magSf()));
452  }
453 
454  return totalArea;
455 }
456 
457 
458 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
459 
461 {
462  // Only a few operations do not require the Sf field
463  switch (operation_)
464  {
465  case opNone:
466  case opMin:
467  case opMax:
468  case opSum:
469  case opSumMag:
470  case opAverage:
471  {
472  return false;
473  }
474  default:
475  {
476  return true;
477  }
478  }
479 }
480 
481 
483 {
484  if (sampledPtr_.valid())
485  {
486  sampledPtr_->update();
487  }
488 
489  if (!needsUpdate_)
490  {
491  return false;
492  }
493 
494  switch (regionType_)
495  {
496  case stFaceZone:
497  {
498  setFaceZoneFaces();
499  break;
500  }
501  case stPatch:
502  {
503  setPatchFaces();
504  break;
505  }
506  case stObject:
507  {
508  const polySurface& s = dynamicCast<const polySurface>(obr());
509  nFaces_ = returnReduce(s.size(), sumOp<label>());
510  break;
511  }
512  case stSampled:
513  {
514  nFaces_ = returnReduce(sampledPtr_->faces().size(), sumOp<label>());
515  break;
516  }
517 
518  // Compiler warning if we forgot an enumeration
519  }
520 
521  if (nFaces_ == 0)
522  {
524  << type() << " " << name() << ": "
525  << regionTypeNames_[regionType_] << '(' << regionName_ << "):" << nl
526  << " Region has no faces" << exit(FatalError);
527  }
528 
529  totalArea_ = totalArea();
530 
531  Log << " total faces = " << nFaces_ << nl
532  << " total area = " << totalArea_ << nl
533  << endl;
534 
535  writeFileHeader(file());
536 
537  needsUpdate_ = false;
538  return true;
539 }
540 
541 
543 (
544  Ostream& os
545 ) const
546 {
547  if (operation_ != opNone)
548  {
549  writeCommented(os, "Region type : ");
550  os << regionTypeNames_[regionType_] << " " << regionName_ << endl;
551 
552  writeHeaderValue(os, "Faces", nFaces_);
553  writeHeaderValue(os, "Area", totalArea_);
554  writeHeaderValue(os, "Scale factor", scaleFactor_);
555 
556  if (weightFieldName_ != "none")
557  {
558  writeHeaderValue(os, "Weight field", weightFieldName_);
559  }
560 
561  writeCommented(os, "Time");
562  if (writeArea_)
563  {
564  os << tab << "Area";
565  }
566 
567  for (const word& fieldName : fields_)
568  {
569  os << tab << operationTypeNames_[operation_]
570  << "(" << fieldName << ")";
571  }
572 
573  os << endl;
574  }
575 }
576 
577 
578 template<>
579 Foam::scalar
581 (
582  const Field<scalar>& values,
583  const vectorField& Sf,
584  const scalarField& weightField
585 ) const
586 {
587  switch (operation_)
588  {
589  case opSumDirection:
590  {
591  const vector n(dict_.get<vector>("direction"));
592  return gSum(pos0(values*(Sf & n))*mag(values));
593  }
594  case opSumDirectionBalance:
595  {
596  const vector n(dict_.get<vector>("direction"));
597  const scalarField nv(values*(Sf & n));
598 
599  return gSum(pos0(nv)*mag(values) - neg(nv)*mag(values));
600  }
601 
602  case opUniformity:
603  case opWeightedUniformity:
604  case opAbsWeightedUniformity:
605  {
606  const scalar areaTotal = gSum(mag(Sf));
607  tmp<scalarField> areaVal(values * mag(Sf));
608 
609  scalar mean, numer;
610 
611  if (canWeight(weightField))
612  {
613  // Weighted quantity = (Weight * phi * dA)
614 
615  tmp<scalarField> weight(weightingFactor(weightField));
616 
617  // Mean weighted value (area-averaged)
618  mean = gSum(weight()*areaVal()) / areaTotal;
619 
620  // Abs. deviation from weighted mean value
621  numer = gSum(mag(weight*areaVal - (mean * mag(Sf))));
622  }
623  else
624  {
625  // Unweighted quantity = (1 * phi * dA)
626 
627  // Mean value (area-averaged)
628  mean = gSum(areaVal()) / areaTotal;
629 
630  // Abs. deviation from mean value
631  numer = gSum(mag(areaVal - (mean * mag(Sf))));
632  }
633 
634  // Uniformity index
635  const scalar ui = 1 - numer/(2*mag(mean*areaTotal) + ROOTVSMALL);
636 
637  return min(max(ui, 0), 1);
638  }
639 
640  default:
641  {
642  // Fall through to other operations
643  return processSameTypeValues(values, Sf, weightField);
644  }
645  }
646 }
647 
648 
649 template<>
652 (
653  const Field<vector>& values,
654  const vectorField& Sf,
655  const scalarField& weightField
656 ) const
657 {
658  switch (operation_)
659  {
660  case opSumDirection:
661  {
662  const vector n(dict_.get<vector>("direction").normalise());
663 
664  const scalarField nv(n & values);
665  return gSum(pos0(nv)*n*(nv));
666  }
667  case opSumDirectionBalance:
668  {
669  const vector n(dict_.get<vector>("direction").normalise());
670 
671  const scalarField nv(n & values);
672  return gSum(pos0(nv)*n*(nv));
673  }
674  case opAreaNormalAverage:
675  {
676  const scalar val = gSum(values & Sf)/gSum(mag(Sf));
677  return vector(val, 0, 0);
678  }
679  case opAreaNormalIntegrate:
680  {
681  const scalar val = gSum(values & Sf);
682  return vector(val, 0, 0);
683  }
684 
685  case opUniformity:
686  case opWeightedUniformity:
687  case opAbsWeightedUniformity:
688  {
689  const scalar areaTotal = gSum(mag(Sf));
690  tmp<scalarField> areaVal(values & Sf);
691 
692  scalar mean, numer;
693 
694  if (canWeight(weightField))
695  {
696  // Weighted quantity = (Weight * phi . dA)
697 
698  tmp<scalarField> weight(weightingFactor(weightField));
699 
700  // Mean weighted value (area-averaged)
701  mean = gSum(weight()*areaVal()) / areaTotal;
702 
703  // Abs. deviation from weighted mean value
704  numer = gSum(mag(weight*areaVal - (mean * mag(Sf))));
705  }
706  else
707  {
708  // Unweighted quantity = (1 * phi . dA)
709 
710  // Mean value (area-averaged)
711  mean = gSum(areaVal()) / areaTotal;
712 
713  // Abs. deviation from mean value
714  numer = gSum(mag(areaVal - (mean * mag(Sf))));
715  }
716 
717  // Uniformity index
718  const scalar ui = 1 - numer/(2*mag(mean*areaTotal) + ROOTVSMALL);
719 
720  return vector(min(max(ui, 0), 1), 0, 0);
721  }
722 
723  default:
724  {
725  // Fall through to other operations
726  return processSameTypeValues(values, Sf, weightField);
727  }
728  }
729 }
730 
731 
732 template<>
735 (
736  const Field<scalar>& weightField
737 ) const
738 {
739  if (usesMag())
740  {
741  return mag(weightField);
742  }
743 
744  // pass through
745  return weightField;
746 }
747 
748 
749 template<>
752 (
753  const Field<scalar>& weightField,
754  const vectorField& Sf
755 ) const
756 {
757  // scalar * Area
758  if (returnReduce(weightField.empty(), andOp<bool>()))
759  {
760  // No weight field - revert to unweighted form
761  return mag(Sf);
762  }
763  else if (usesMag())
764  {
765  return mag(weightField * mag(Sf));
766  }
767 
768  return (weightField * mag(Sf));
769 }
770 
771 
772 template<>
775 (
776  const Field<vector>& weightField,
777  const vectorField& Sf
778 ) const
779 {
780  // vector (dot) Area
781  if (returnReduce(weightField.empty(), andOp<bool>()))
782  {
783  // No weight field - revert to unweighted form
784  return mag(Sf);
785  }
786  else if (usesMag())
787  {
788  return mag(weightField & Sf);
789  }
790 
791  return (weightField & Sf);
792 }
793 
794 
795 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
796 
798 (
799  const word& name,
800  const Time& runTime,
801  const dictionary& dict
802 )
803 :
804  fieldValue(name, runTime, dict, typeName),
805  regionType_(regionTypeNames_.get("regionType", dict)),
806  operation_(operationTypeNames_.get("operation", dict)),
807  postOperation_
808  (
809  postOperationTypeNames_.lookupOrDefault
810  (
811  "postOperation",
812  dict,
813  postOperationType::postOpNone,
814  true // Failsafe behaviour
815  )
816  ),
817  weightFieldName_("none"),
818  needsUpdate_(true),
819  writeArea_(false),
820  totalArea_(0),
821  nFaces_(0),
822  faceId_(),
823  facePatchId_(),
824  faceFlip_()
825 {
826  read(dict);
827 }
828 
829 
831 (
832  const word& name,
833  const objectRegistry& obr,
834  const dictionary& dict
835 )
836 :
837  fieldValue(name, obr, dict, typeName),
838  regionType_(regionTypeNames_.get("regionType", dict)),
839  operation_(operationTypeNames_.get("operation", dict)),
840  postOperation_
841  (
842  postOperationTypeNames_.lookupOrDefault
843  (
844  "postOperation",
845  dict,
846  postOperationType::postOpNone,
847  true // Failsafe behaviour
848  )
849  ),
850  weightFieldName_("none"),
851  needsUpdate_(true),
852  writeArea_(false),
853  totalArea_(0),
854  nFaces_(0),
855  faceId_(),
856  facePatchId_(),
857  faceFlip_()
858 {
859  read(dict);
860 }
861 
862 
863 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
864 
866 (
867  const dictionary& dict
868 )
869 {
871 
872  weightFieldName_ = "none";
873  needsUpdate_ = true;
874  writeArea_ = dict.lookupOrDefault("writeArea", false);
875  totalArea_ = 0;
876  nFaces_ = 0;
877  faceId_.clear();
878  facePatchId_.clear();
879  faceFlip_.clear();
880  sampledPtr_.clear();
881  surfaceWriterPtr_.clear();
882 
883  dict.readEntry("name", regionName_);
884 
885  // Create sampled surface, but leave 'expired' (ie, no update) since it
886  // may depend on fields or data that do not yet exist
887  if (stSampled == regionType_)
888  {
889  sampledPtr_ = sampledSurface::New
890  (
891  name(),
892  mesh_,
893  dict.subDict("sampledSurfaceDict")
894  );
895  }
896 
897  Info<< type() << " " << name() << ":" << nl
898  << " operation = ";
899 
900  if (postOperation_ != postOpNone)
901  {
902  Info<< postOperationTypeNames_[postOperation_] << '('
903  << operationTypeNames_[operation_] << ')' << nl;
904  }
905  else
906  {
907  Info<< operationTypeNames_[operation_] << nl;
908  }
909 
910  if (usesWeight())
911  {
912  if (stSampled == regionType_)
913  {
915  << "Cannot use weighted operation '"
916  << operationTypeNames_[operation_]
917  << "' for sampledSurface"
918  << exit(FatalIOError);
919  }
920 
921  if (dict.readIfPresent("weightField", weightFieldName_))
922  {
923  Info<< " weight field = " << weightFieldName_ << nl;
924  }
925  else
926  {
927  // Suggest possible alternative unweighted operation?
929  << "The '" << operationTypeNames_[operation_]
930  << "' operation is missing a weightField." << nl
931  << "Either provide the weightField, "
932  << "use weightField 'none' to suppress weighting," << nl
933  << "or use a different operation."
934  << exit(FatalIOError);
935  }
936  }
937 
938  // Backwards compatibility for v1612 and older
939  List<word> orientedFields;
940  if (dict.readIfPresent("orientedFields", orientedFields))
941  {
942  fields_.append(orientedFields);
943 
945  << "The 'orientedFields' option is deprecated. These fields can "
946  << "and have been added to the standard 'fields' list."
947  << endl;
948  }
949 
950  if (writeFields_)
951  {
952  const word formatName(dict.get<word>("surfaceFormat"));
953 
954  surfaceWriterPtr_.reset
955  (
957  (
958  formatName,
959  dict.subOrEmptyDict("formatOptions").subOrEmptyDict(formatName)
960  )
961  );
962 
963  if (debug)
964  {
965  surfaceWriterPtr_->verbose() = true;
966  }
967 
968  if (surfaceWriterPtr_->enabled())
969  {
970  Info<< " surfaceFormat = " << formatName << nl;
971  }
972  else
973  {
974  surfaceWriterPtr_->clear();
975  }
976  }
977 
978  Info<< nl << endl;
979 
980  return true;
981 }
982 
983 
985 {
986  if (needsUpdate_ || operation_ != opNone)
987  {
989  }
990 
991  update();
992 
993  if (operation_ != opNone)
994  {
995  writeCurrentTime(file());
996  }
997 
998  if (writeArea_)
999  {
1000  totalArea_ = totalArea();
1001  Log << " total area = " << totalArea_ << endl;
1002 
1003  if (operation_ != opNone && Pstream::master())
1004  {
1005  file() << tab << totalArea_;
1006  }
1007  }
1008 
1009  // Many operations use the Sf field
1010  vectorField Sf;
1011  if (usesSf())
1012  {
1013  if (stObject == regionType_)
1014  {
1015  const polySurface& s = dynamicCast<const polySurface>(obr());
1016  Sf = s.Sf();
1017  }
1018  else if (sampledPtr_.valid())
1019  {
1020  Sf = sampledPtr_().Sf();
1021  }
1022  else
1023  {
1024  Sf = filterField(mesh_.Sf());
1025  }
1026  }
1027 
1028  // Faces and points for surface format (if specified)
1029  faceList faces;
1031 
1032  if (surfaceWriterPtr_.valid())
1033  {
1034  if (withTopologicalMerge())
1035  {
1036  combineMeshGeometry(faces, points);
1037  }
1038  else
1039  {
1040  combineSurfaceGeometry(faces, points);
1041  }
1042  }
1043 
1044  // Only a few weight types (scalar, vector)
1045  if (weightFieldName_ != "none")
1046  {
1047  if (validField<scalar>(weightFieldName_))
1048  {
1049  scalarField weightField
1050  (
1051  getFieldValues<scalar>(weightFieldName_, true)
1052  );
1053 
1054  // Process the fields
1055  writeAll(Sf, weightField, points, faces);
1056  }
1057  else if (validField<vector>(weightFieldName_))
1058  {
1059  vectorField weightField
1060  (
1061  getFieldValues<vector>(weightFieldName_, true)
1062  );
1063 
1064  // Process the fields
1065  writeAll(Sf, weightField, points, faces);
1066  }
1067  else
1068  {
1070  << "weightField " << weightFieldName_
1071  << " not found or an unsupported type"
1072  << abort(FatalError);
1073  }
1074  }
1075  else
1076  {
1077  // Default is a zero-size scalar weight field (ie, weight = 1)
1078  scalarField weightField;
1079 
1080  // Process the fields
1081  writeAll(Sf, weightField, points, faces);
1082  }
1083 
1084  if (operation_ != opNone)
1085  {
1086  file() << endl;
1087  Log << endl;
1088  }
1089 
1090 
1091  return true;
1092 }
1093 
1094 
1097  const mapPolyMesh& mpm
1098 )
1099 {
1100  needsUpdate_ = true;
1101 }
1102 
1103 
1106  const polyMesh& mesh
1107 )
1108 {
1109  needsUpdate_ = true;
1110 }
1111 
1112 
1113 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::val
label ListType::const_reference val
Definition: ListOps.H:407
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:64
Foam::Enum
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: IOstreamOption.H:51
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::functionObjects::fieldValues::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, fieldValueDelta, dictionary)
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::functionObjects::fieldValue::read
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: fieldValue.C:91
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
update
mesh update()
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::functionObjects::fieldValues::surfaceFieldValue::operationType
operationType
Operation type enumeration.
Definition: surfaceFieldValue.H:415
Foam::functionObjects::fieldValues::surfaceFieldValue::read
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: surfaceFieldValue.C:866
Foam::functionObjects::fieldValues::surfaceFieldValue::postOperationType
postOperationType
Post-operation type enumeration.
Definition: surfaceFieldValue.H:487
PatchTools.H
Foam::DynamicList< label >
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::UPstream::nProcs
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:426
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:414
Foam::functionObjects::fieldValues::defineTypeNameAndDebug
defineTypeNameAndDebug(fieldValueDelta, 0)
Foam::functionObjects::fieldValues::surfaceFieldValue::usesSf
bool usesSf() const
True if the operation needs a surface Sf.
Definition: surfaceFieldValue.C:460
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:182
Foam::FatalIOError
IOerror FatalIOError
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
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
An Ostream wrapper for parallel output to std::cout.
Foam::Vector::normalise
Vector< Cmpt > & normalise()
Normalise the vector by its magnitude.
Definition: VectorI.H:128
Foam::functionObjects::fieldValues::surfaceFieldValue::regionTypeNames_
static const Enum< regionTypes > regionTypeNames_
Region type names.
Definition: surfaceFieldValue.H:403
Foam::functionObjects::fieldValues::surfaceFieldValue::regionType_
regionTypes regionType_
Region type.
Definition: surfaceFieldValue.H:530
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::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:290
Foam::functionObjects::fieldValues::surfaceFieldValue::surfaceFieldValue
surfaceFieldValue(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
Definition: surfaceFieldValue.C:798
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::functionObjects::fieldValues::surfaceFieldValue::stObject
Calculate with function object surface.
Definition: surfaceFieldValue.H:398
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
Foam::functionObjects::fieldValues::surfaceFieldValue::weightingFactor
tmp< scalarField > weightingFactor(const Field< WeightType > &weightField) const
Weighting factor.
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:482
coupledPolyPatch.H
Foam::functionObjects::fieldValue::regionName_
word regionName_
Name of region (patch, zone, etc.)
Definition: fieldValue.H:82
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::Field< scalar >
Foam::polyPatch::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:305
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:65
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
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:1096
Foam::functionObjects::fieldValues::surfaceFieldValue::regionTypes
regionTypes
Region type enumeration.
Definition: surfaceFieldValue.H:394
Foam::andOp
Definition: ops.H:233
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
indirectPrimitivePatch.H
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:121
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::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::functionObjects::fieldValues::surfaceFieldValue::obr
const objectRegistry & obr() const
The volume mesh or surface registry being used.
Definition: surfaceFieldValue.C:119
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:51
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:311
emptyPolyPatch.H
Foam::functionObjects::fieldValues::surfaceFieldValue::movePoints
virtual void movePoints(const polyMesh &mesh)
Update for changes of mesh.
Definition: surfaceFieldValue.C:1105
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::functionObjects::fieldValues::surfaceFieldValue::writeFileHeader
virtual void writeFileHeader(Ostream &os) const
Output file header information.
Definition: surfaceFieldValue.C:543
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:444
Foam::UPstream::master
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:438
Foam::indirectPrimitivePatch
PrimitivePatch< face, IndirectList, const pointField & > indirectPrimitivePatch
A PrimitivePatch using an IndirectList for the faces.
Definition: indirectPrimitivePatch.H:47
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
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::polyPatch::whichFace
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:398
Foam::tab
constexpr char tab
Definition: Ostream.H:371
Foam::nl
constexpr char nl
Definition: Ostream.H:372
Foam::functionObjects::fieldValues::surfaceFieldValue::operationTypeNames_
static const Enum< operationType > operationTypeNames_
Operation type names.
Definition: surfaceFieldValue.H:483
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::Vector< scalar >
Foam::List< word >
Foam::functionObjects::fieldValue
Base class for field value-based function objects.
Definition: fieldValue.H:65
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
mergePoints.H
Merge points. See below.
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:160
Foam::primitivePatch
PrimitivePatch< face, SubList, const pointField & > primitivePatch
Addressing for a faceList slice.
Definition: primitivePatch.H:47
surfaceFieldValue.H
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:375
Foam::DelaunayMeshTools::allPoints
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
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:984
Foam::PatchTools::gatherAndMerge
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< Face, FaceList, PointField, PointType > &p, Field< PointType > &mergedPoints, List< Face > &mergedFaces, labelList &pointMergeMap)
Gather points and faces onto master and merge into single patch.
Definition: PatchToolsGatherAndMerge.C:44
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::faceZone::flipMap
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:271
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:293
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:199
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::functionObjects::fieldValues::surfaceFieldValue::postOperationTypeNames_
static const Enum< postOperationType > postOperationTypeNames_
Operation type names.
Definition: surfaceFieldValue.H:494
Log
#define Log
Report write to Foam::Info if the local log switch is true.
Definition: messageStream.H:332