surfaceFieldValue.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "surfaceFieldValue.H"
27 #include "fvMesh.H"
28 #include "cyclicPolyPatch.H"
29 #include "emptyPolyPatch.H"
30 #include "coupledPolyPatch.H"
31 #include "sampledSurface.H"
32 #include "mergePoints.H"
33 #include "indirectPrimitivePatch.H"
34 #include "PatchTools.H"
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 namespace Foam
40 {
41 namespace functionObjects
42 {
43 namespace fieldValues
44 {
45  defineTypeNameAndDebug(surfaceFieldValue, 0);
46  addToRunTimeSelectionTable(fieldValue, surfaceFieldValue, dictionary);
47  addToRunTimeSelectionTable(functionObject, surfaceFieldValue, dictionary);
48 }
49 }
50 }
51 
52 template<>
53 const char* Foam::NamedEnum
54 <
56  3
57 >::names[] =
58 {
59  "faceZone",
60  "patch",
61  "sampledSurface"
62 };
63 
64 template<>
65 const char* Foam::NamedEnum
66 <
68  17
69 >::names[] =
70 {
71  "none",
72  "sum",
73  "weightedSum",
74  "sumMag",
75  "sumDirection",
76  "sumDirectionBalance",
77  "average",
78  "weightedAverage",
79  "areaAverage",
80  "weightedAreaAverage",
81  "areaIntegrate",
82  "weightedAreaIntegrate",
83  "min",
84  "max",
85  "CoV",
86  "areaNormalAverage",
87  "areaNormalIntegrate"
88 };
89 
90 const Foam::NamedEnum
91 <
93  3
95 
96 const Foam::NamedEnum
97 <
99  17
101 
102 
103 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
104 
105 void Foam::functionObjects::fieldValues::surfaceFieldValue::setFaceZoneFaces()
106 {
107  label zoneId = mesh_.faceZones().findZoneID(regionName_);
108 
109  if (zoneId < 0)
110  {
112  << type() << " " << name() << ": "
113  << regionTypeNames_[regionType_] << "(" << regionName_ << "):" << nl
114  << " Unknown face zone name: " << regionName_
115  << ". Valid face zones are: " << mesh_.faceZones().names()
116  << nl << exit(FatalError);
117  }
118 
119  const faceZone& fZone = mesh_.faceZones()[zoneId];
120 
121  DynamicList<label> faceIds(fZone.size());
122  DynamicList<label> facePatchIds(fZone.size());
123  DynamicList<label> faceSigns(fZone.size());
124 
125  forAll(fZone, i)
126  {
127  label facei = fZone[i];
128 
129  label faceId = -1;
130  label facePatchId = -1;
131  if (mesh_.isInternalFace(facei))
132  {
133  faceId = facei;
134  facePatchId = -1;
135  }
136  else
137  {
138  facePatchId = mesh_.boundaryMesh().whichPatch(facei);
139  const polyPatch& pp = mesh_.boundaryMesh()[facePatchId];
140  if (isA<coupledPolyPatch>(pp))
141  {
142  if (refCast<const coupledPolyPatch>(pp).owner())
143  {
144  faceId = pp.whichFace(facei);
145  }
146  else
147  {
148  faceId = -1;
149  }
150  }
151  else if (!isA<emptyPolyPatch>(pp))
152  {
153  faceId = facei - pp.start();
154  }
155  else
156  {
157  faceId = -1;
158  facePatchId = -1;
159  }
160  }
161 
162  if (faceId >= 0)
163  {
164  if (fZone.flipMap()[i])
165  {
166  faceSigns.append(-1);
167  }
168  else
169  {
170  faceSigns.append(1);
171  }
172  faceIds.append(faceId);
173  facePatchIds.append(facePatchId);
174  }
175  }
176 
177  faceId_.transfer(faceIds);
178  facePatchId_.transfer(facePatchIds);
179  faceSign_.transfer(faceSigns);
180  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
181 
182  if (debug)
183  {
184  Pout<< "Original face zone size = " << fZone.size()
185  << ", new size = " << faceId_.size() << endl;
186  }
187 }
188 
189 
190 void Foam::functionObjects::fieldValues::surfaceFieldValue::setPatchFaces()
191 {
192  const label patchid = mesh_.boundaryMesh().findPatchID(regionName_);
193 
194  if (patchid < 0)
195  {
197  << type() << " " << name() << ": "
198  << regionTypeNames_[regionType_] << "(" << regionName_ << "):" << nl
199  << " Unknown patch name: " << regionName_
200  << ". Valid patch names are: "
201  << mesh_.boundaryMesh().names() << nl
202  << exit(FatalError);
203  }
204 
205  const polyPatch& pp = mesh_.boundaryMesh()[patchid];
206 
207  label nFaces = pp.size();
208  if (isA<emptyPolyPatch>(pp))
209  {
210  nFaces = 0;
211  }
212 
213  faceId_.setSize(nFaces);
214  facePatchId_.setSize(nFaces);
215  faceSign_.setSize(nFaces);
216  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
217 
218  forAll(faceId_, facei)
219  {
220  faceId_[facei] = facei;
221  facePatchId_[facei] = patchid;
222  faceSign_[facei] = 1;
223  }
224 }
225 
226 
227 void Foam::functionObjects::fieldValues::surfaceFieldValue::sampledSurfaceFaces
228 (
229  const dictionary& dict
230 )
231 {
232  surfacePtr_ = sampledSurface::New
233  (
234  name(),
235  mesh_,
236  dict.subDict("sampledSurfaceDict")
237  );
238  surfacePtr_().update();
239  nFaces_ = returnReduce(surfacePtr_().faces().size(), sumOp<label>());
240 }
241 
242 
243 void Foam::functionObjects::fieldValues::surfaceFieldValue::combineMeshGeometry
244 (
245  faceList& faces,
246  pointField& points
247 ) const
248 {
249  List<faceList> allFaces(Pstream::nProcs());
251 
252  labelList globalFacesIs(faceId_);
253  forAll(globalFacesIs, i)
254  {
255  if (facePatchId_[i] != -1)
256  {
257  label patchi = facePatchId_[i];
258  globalFacesIs[i] += mesh_.boundaryMesh()[patchi].start();
259  }
260  }
261 
262  // Add local faces and points to the all* lists
264  (
265  IndirectList<face>(mesh_.faces(), globalFacesIs),
266  mesh_.points()
267  );
268  allFaces[Pstream::myProcNo()] = pp.localFaces();
269  allPoints[Pstream::myProcNo()] = pp.localPoints();
270 
271  Pstream::gatherList(allFaces);
272  Pstream::gatherList(allPoints);
273 
274  // Renumber and flatten
275  label nFaces = 0;
276  label nPoints = 0;
277  forAll(allFaces, proci)
278  {
279  nFaces += allFaces[proci].size();
280  nPoints += allPoints[proci].size();
281  }
282 
283  faces.setSize(nFaces);
284  points.setSize(nPoints);
285 
286  nFaces = 0;
287  nPoints = 0;
288 
289  // My own data first
290  {
291  const faceList& fcs = allFaces[Pstream::myProcNo()];
292  forAll(fcs, i)
293  {
294  const face& f = fcs[i];
295  face& newF = faces[nFaces++];
296  newF.setSize(f.size());
297  forAll(f, fp)
298  {
299  newF[fp] = f[fp] + nPoints;
300  }
301  }
302 
303  const pointField& pts = allPoints[Pstream::myProcNo()];
304  forAll(pts, i)
305  {
306  points[nPoints++] = pts[i];
307  }
308  }
309 
310  // Other proc data follows
311  forAll(allFaces, proci)
312  {
313  if (proci != Pstream::myProcNo())
314  {
315  const faceList& fcs = allFaces[proci];
316  forAll(fcs, i)
317  {
318  const face& f = fcs[i];
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  forAll(pts, i)
329  {
330  points[nPoints++] = pts[i];
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  forAll(faces, i)
357  {
358  inplaceRenumber(oldToNew, faces[i]);
359  }
360  }
361 }
362 
363 
364 void Foam::functionObjects::fieldValues::surfaceFieldValue::
365 combineSurfaceGeometry
366 (
367  faceList& faces,
368  pointField& points
369 ) const
370 {
371  if (surfacePtr_.valid())
372  {
373  const sampledSurface& s = surfacePtr_();
374 
375  if (Pstream::parRun())
376  {
377  // Dimension as fraction of mesh bounding box
378  scalar mergeDim = 1e-10*mesh_.bounds().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 }
402 
403 
404 Foam::scalar
405 Foam::functionObjects::fieldValues::surfaceFieldValue::totalArea() const
406 {
407  scalar totalArea;
408 
409  if (surfacePtr_.valid())
410  {
411  totalArea = gSum(surfacePtr_().magSf());
412  }
413  else
414  {
415  totalArea = gSum(filterField(mesh_.magSf(), false));
416  }
417 
418  return totalArea;
419 }
420 
421 
422 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
423 
425 (
426  const dictionary& dict
427 )
428 {
429  dict.lookup("name") >> regionName_;
430 
431  switch (regionType_)
432  {
433  case regionTypes::faceZone:
434  {
435  setFaceZoneFaces();
436  break;
437  }
438  case regionTypes::patch:
439  {
440  setPatchFaces();
441  break;
442  }
443  case regionTypes::sampledSurface:
444  {
445  sampledSurfaceFaces(dict);
446  break;
447  }
448  default:
449  {
451  << type() << " " << name() << ": "
452  << regionTypeNames_[regionType_] << "(" << regionName_ << "):"
453  << nl << " Unknown region type. Valid region types are:"
454  << regionTypeNames_.sortedToc() << nl << exit(FatalError);
455  }
456  }
457 
458  if (nFaces_ == 0)
459  {
461  << type() << " " << name() << ": "
462  << regionTypeNames_[regionType_] << "(" << regionName_ << "):" << nl
463  << " Region has no faces" << exit(FatalError);
464  }
465 
466  if (surfacePtr_.valid())
467  {
468  surfacePtr_().update();
469  }
470 
471  totalArea_ = totalArea();
472 
473  Info<< type() << " " << name() << ":" << nl
474  << " total faces = " << nFaces_
475  << nl
476  << " total area = " << totalArea_
477  << nl;
478 
479  if (dict.readIfPresent("weightField", weightFieldName_))
480  {
481  Info<< " weight field = " << weightFieldName_ << nl;
482 
483  if (regionType_ == regionTypes::sampledSurface)
484  {
486  << "Cannot use weightField for a sampledSurface"
487  << exit(FatalIOError);
488  }
489  }
490 
491  if (dict.found("orientedWeightField"))
492  {
493  if (weightFieldName_ == "none")
494  {
495  dict.lookup("orientedWeightField") >> weightFieldName_;
496  Info<< " weight field = " << weightFieldName_ << nl;
497  orientWeightField_ = true;
498  }
499  else
500  {
502  << "Either weightField or orientedWeightField can be supplied, "
503  << "but not both"
504  << exit(FatalIOError);
505  }
506  }
507 
508  List<word> orientedFields;
509  if (dict.readIfPresent("orientedFields", orientedFields))
510  {
511  orientedFieldsStart_ = fields_.size();
512  fields_.append(orientedFields);
513  }
514 
515  if (dict.readIfPresent("scaleFactor", scaleFactor_))
516  {
517  Info<< " scale factor = " << scaleFactor_ << nl;
518  }
519 
520  Info<< nl << endl;
521 
522  if (writeFields_)
523  {
524  const word surfaceFormat(dict.lookup("surfaceFormat"));
525 
526  surfaceWriterPtr_.reset
527  (
528  surfaceWriter::New(surfaceFormat, dict).ptr()
529  );
530  }
531 }
532 
533 
535 (
536  const label i
537 )
538 {
539  if (operation_ != operationType::none)
540  {
541  writeCommented(file(), "Region type : ");
542  file() << regionTypeNames_[regionType_] << " " << regionName_ << endl;
543  writeCommented(file(), "Faces : ");
544  file() << nFaces_ << endl;
545  writeCommented(file(), "Area : ");
546  file() << totalArea_ << endl;
547 
548  writeCommented(file(), "Time");
549  if (writeArea_)
550  {
551  file() << tab << "Area";
552  }
553 
554  forAll(fields_, fieldi)
555  {
556  file()
557  << tab << operationTypeNames_[operation_]
558  << "(" << fields_[fieldi] << ")";
559  }
560 
561  file() << endl;
562  }
563 }
564 
565 
566 template<>
569 (
570  const Field<scalar>& values,
571  const vectorField& Sf,
572  const scalarField& weightField
573 ) const
574 {
575  switch (operation_)
576  {
577  case operationType::sumDirection:
578  {
579  vector n(dict_.lookup("direction"));
580  return sum(pos0(values*(Sf & n))*mag(values));
581  }
582  case operationType::sumDirectionBalance:
583  {
584  vector n(dict_.lookup("direction"));
585  const scalarField nv(values*(Sf & n));
586 
587  return sum(pos0(nv)*mag(values) - neg(nv)*mag(values));
588  }
589  default:
590  {
591  // Fall through to other operations
592  return processSameTypeValues(values, Sf, weightField);
593  }
594  }
595 }
596 
597 
598 template<>
601 (
602  const Field<vector>& values,
603  const vectorField& Sf,
604  const scalarField& weightField
605 ) const
606 {
607  switch (operation_)
608  {
609  case operationType::sumDirection:
610  {
611  vector n(dict_.lookup("direction"));
612  n /= mag(n) + rootVSmall;
613  const scalarField nv(n & values);
614 
615  return sum(pos0(nv)*n*(nv));
616  }
617  case operationType::sumDirectionBalance:
618  {
619  vector n(dict_.lookup("direction"));
620  n /= mag(n) + rootVSmall;
621  const scalarField nv(n & values);
622 
623  return sum(pos0(nv)*n*(nv));
624  }
625  case operationType::areaNormalAverage:
626  {
627  scalar result = sum(values & Sf)/sum(mag(Sf));
628  return vector(result, 0.0, 0.0);
629  }
630  case operationType::areaNormalIntegrate:
631  {
632  scalar result = sum(values & Sf);
633  return vector(result, 0.0, 0.0);
634  }
635  default:
636  {
637  // Fall through to other operations
638  return processSameTypeValues(values, Sf, weightField);
639  }
640  }
641 }
642 
643 
644 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
645 
647 (
648  const word& name,
649  const Time& runTime,
650  const dictionary& dict
651 )
652 :
653  fieldValue(name, runTime, dict, typeName),
654  surfaceWriterPtr_(nullptr),
655  regionType_(regionTypeNames_.read(dict.lookup("regionType"))),
656  operation_(operationTypeNames_.read(dict.lookup("operation"))),
657  weightFieldName_("none"),
658  orientWeightField_(false),
659  orientedFieldsStart_(labelMax),
660  scaleFactor_(1.0),
661  writeArea_(dict.lookupOrDefault("writeArea", false)),
662  nFaces_(0),
663  faceId_(),
664  facePatchId_(),
665  faceSign_()
666 {
667  read(dict);
668 }
669 
671 (
672  const word& name,
673  const objectRegistry& obr,
674  const dictionary& dict
675 )
676 :
677  fieldValue(name, obr, dict, typeName),
678  surfaceWriterPtr_(nullptr),
679  regionType_(regionTypeNames_.read(dict.lookup("regionType"))),
680  operation_(operationTypeNames_.read(dict.lookup("operation"))),
681  weightFieldName_("none"),
682  orientWeightField_(false),
683  orientedFieldsStart_(labelMax),
684  scaleFactor_(1.0),
685  writeArea_(dict.lookupOrDefault("writeArea", false)),
686  nFaces_(0),
687  faceId_(),
688  facePatchId_(),
689  faceSign_()
690 {
691  read(dict);
692 }
693 
694 
695 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
696 
698 {}
699 
700 
701 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
702 
704 (
705  const dictionary& dict
706 )
707 {
708  fieldValue::read(dict);
709  initialise(dict);
710 
711  return true;
712 }
713 
714 
716 {
717  if (operation_ != operationType::none)
718  {
720  }
721 
722  if (surfacePtr_.valid())
723  {
724  surfacePtr_().update();
725  }
726 
727  if (operation_ != operationType::none && Pstream::master())
728  {
729  writeTime(file());
730  }
731 
732  if (writeArea_)
733  {
734  totalArea_ = totalArea();
735  if (operation_ != operationType::none && Pstream::master())
736  {
737  file() << tab << totalArea_;
738  }
739  Log << " total area = " << totalArea_ << endl;
740  }
741 
742  // Write the surface geometry
743  if (surfaceWriterPtr_.valid())
744  {
745  faceList faces;
747 
748  if (surfacePtr_.valid())
749  {
750  combineSurfaceGeometry(faces, points);
751  }
752  else
753  {
754  combineMeshGeometry(faces, points);
755  }
756 
757  if (Pstream::master())
758  {
759  surfaceWriterPtr_->write
760  (
761  outputDir(),
762  regionTypeNames_[regionType_] + ("_" + regionName_),
763  points,
764  faces
765  );
766  }
767  }
768 
769  // Construct weight field. Note: zero size means weight = 1
770  scalarField weightField;
771  if (weightFieldName_ != "none")
772  {
773  weightField =
774  getFieldValues<scalar>
775  (
776  weightFieldName_,
777  true,
778  orientWeightField_
779  );
780  }
781 
782  // Combine onto master
783  combineFields(weightField);
784 
785  // Process the fields
786  forAll(fields_, i)
787  {
788  const word& fieldName = fields_[i];
789  bool ok = false;
790 
791  bool orient = i >= orientedFieldsStart_;
792  ok = ok || writeValues<scalar>(fieldName, weightField, orient);
793  ok = ok || writeValues<vector>(fieldName, weightField, orient);
794  ok = ok
795  || writeValues<sphericalTensor>(fieldName, weightField, orient);
796  ok = ok || writeValues<symmTensor>(fieldName, weightField, orient);
797  ok = ok || writeValues<tensor>(fieldName, weightField, orient);
798 
799  if (!ok)
800  {
802  << "Requested field " << fieldName
803  << " not found in database and not processed"
804  << endl;
805  }
806  }
807 
808  if (operation_ != operationType::none && Pstream::master())
809  {
810  file() << endl;
811  }
812 
813  Log << endl;
814 
815  return true;
816 }
817 
818 
819 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:667
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
label faceId(-1)
surfaceFieldValue(const word &name, const Time &runTime, const dictionary &dict)
Construct from name, Time and dictionary.
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
void initialise(const dictionary &dict)
Initialise, e.g. face addressing.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
static const char tab
Definition: Ostream.H:259
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
addToRunTimeSelectionTable(functionObject, fieldValueDelta, dictionary)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
An abstract class for surfaces with sampling.
static autoPtr< surfaceWriter > New(const word &writeType, const IOstream::streamFormat writeFormat)
Return a reference to the selected surfaceWriter.
Definition: surfaceWriter.C:45
static const NamedEnum< operationType, 17 > operationTypeNames_
Operation type names.
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:278
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:248
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:429
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
virtual bool read(const dictionary &dict)
Read from dictionary.
Definition: fieldValue.C:90
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:51
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:934
virtual const faceList & faces() const =0
Faces of surface.
A list of faces which address into the list of points.
A List obtained as a section of another List.
Definition: SubList.H:53
bool read(const char *, int32_t &)
Definition: int32IO.C:85
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().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))
Type gSum(const FieldField< Field, Type > &f)
const pointField & points
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
defineTypeNameAndDebug(fieldValueDelta, 0)
static void gatherAndMerge(const scalar mergeDist, const PrimitivePatch< FaceList, PointField > &p, Field< typename PrimitivePatch< FaceList, PointField >::PointType > &mergedPoints, List< typename PrimitivePatch< FaceList, PointField >::FaceType > &mergedFaces, labelList &pointMergeMap)
Gather points and faces onto master and merge into single patch.
bool readIfPresent(const word &, T &, bool recursive=false, bool patternMatch=true) const
Find an entry if present, and assign to T.
Type processValues(const Field< Type > &values, const vectorField &Sf, const scalarField &weightField) const
Apply the &#39;operation&#39; to the values. Wrapper around.
static const label labelMax
Definition: label.H:62
static autoPtr< sampledSurface > New(const word &name, const polyMesh &, const dictionary &)
Return a reference to the selected surface.
dimensionedScalar pos0(const dimensionedScalar &ds)
static const char nl
Definition: Ostream.H:260
Merge points. See below.
labelList f(nPoints)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void setSize(const label)
Reset size of List.
Definition: List.C:281
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:411
label patchi
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
#define WarningInFunction
Report a warning using Foam::Warning.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:331
virtual bool read(const dictionary &)
Read from dictionary.
virtual const pointField & points() const =0
Points of surface.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:303
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
label mergePoints(const UList< Type > &points, const scalar mergeTol, const bool verbose, labelList &pointMap, const Type &origin=Type::zero)
Sorts and merges points. All points closer than/equal mergeTol get merged.
#define Log
Report write to Foam::Info if the local log switch is true.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
label n
Base class for field value -based function objects.
Definition: fieldValue.H:62
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
static const NamedEnum< regionTypes, 3 > regionTypeNames_
region type names
Registry of regIOobjects.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
A List with indirect addressing.
Definition: IndirectList.H:101
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
Namespace for OpenFOAM.
virtual void writeFileHeader(const label i)
Output file header information.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:812
virtual bool write()
Write.
Definition: fieldValue.C:112
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:383
IOerror FatalIOError