faceSource.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 | Copyright (C) 2011-2015 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 "faceSource.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  template<>
43  {
44  "faceZone",
45  "patch",
46  "sampledSurface"
47  };
48 
49 
50  template<>
52  {
53  "none",
54  "sum",
55  "sumMag",
56  "sumDirection",
57  "sumDirectionBalance",
58  "average",
59  "weightedAverage",
60  "areaAverage",
61  "weightedAreaAverage",
62  "areaIntegrate",
63  "min",
64  "max",
65  "CoV",
66  "areaNormalAverage",
67  "areaNormalIntegrate"
68  };
69 
70  namespace fieldValues
71  {
72  defineTypeNameAndDebug(faceSource, 0);
74  }
75 }
76 
77 
80 
83 
84 
85 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
86 
87 void Foam::fieldValues::faceSource::setFaceZoneFaces()
88 {
89  label zoneId = mesh().faceZones().findZoneID(sourceName_);
90 
91  if (zoneId < 0)
92  {
93  FatalErrorIn("faceSource::faceSource::setFaceZoneFaces()")
94  << type() << " " << name_ << ": "
95  << sourceTypeNames_[source_] << "(" << sourceName_ << "):" << nl
96  << " Unknown face zone name: " << sourceName_
97  << ". Valid face zones are: " << mesh().faceZones().names()
98  << nl << exit(FatalError);
99  }
100 
101  const faceZone& fZone = mesh().faceZones()[zoneId];
102 
103  DynamicList<label> faceIds(fZone.size());
104  DynamicList<label> facePatchIds(fZone.size());
105  DynamicList<label> faceSigns(fZone.size());
106 
107  forAll(fZone, i)
108  {
109  label faceI = fZone[i];
110 
111  label faceId = -1;
112  label facePatchId = -1;
113  if (mesh().isInternalFace(faceI))
114  {
115  faceId = faceI;
116  facePatchId = -1;
117  }
118  else
119  {
120  facePatchId = mesh().boundaryMesh().whichPatch(faceI);
121  const polyPatch& pp = mesh().boundaryMesh()[facePatchId];
122  if (isA<coupledPolyPatch>(pp))
123  {
124  if (refCast<const coupledPolyPatch>(pp).owner())
125  {
126  faceId = pp.whichFace(faceI);
127  }
128  else
129  {
130  faceId = -1;
131  }
132  }
133  else if (!isA<emptyPolyPatch>(pp))
134  {
135  faceId = faceI - pp.start();
136  }
137  else
138  {
139  faceId = -1;
140  facePatchId = -1;
141  }
142  }
143 
144  if (faceId >= 0)
145  {
146  if (fZone.flipMap()[i])
147  {
148  faceSigns.append(-1);
149  }
150  else
151  {
152  faceSigns.append(1);
153  }
154  faceIds.append(faceId);
155  facePatchIds.append(facePatchId);
156  }
157  }
158 
159  faceId_.transfer(faceIds);
160  facePatchId_.transfer(facePatchIds);
161  faceSign_.transfer(faceSigns);
162  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
163 
164  if (debug)
165  {
166  Pout<< "Original face zone size = " << fZone.size()
167  << ", new size = " << faceId_.size() << endl;
168  }
169 }
170 
171 
172 void Foam::fieldValues::faceSource::setPatchFaces()
173 {
174  const label patchId = mesh().boundaryMesh().findPatchID(sourceName_);
175 
176  if (patchId < 0)
177  {
178  FatalErrorIn("faceSource::constructFaceAddressing()")
179  << type() << " " << name_ << ": "
180  << sourceTypeNames_[source_] << "(" << sourceName_ << "):" << nl
181  << " Unknown patch name: " << sourceName_
182  << ". Valid patch names are: "
183  << mesh().boundaryMesh().names() << nl
184  << exit(FatalError);
185  }
186 
187  const polyPatch& pp = mesh().boundaryMesh()[patchId];
188 
189  label nFaces = pp.size();
190  if (isA<emptyPolyPatch>(pp))
191  {
192  nFaces = 0;
193  }
194 
195  faceId_.setSize(nFaces);
196  facePatchId_.setSize(nFaces);
197  faceSign_.setSize(nFaces);
198  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
199 
200  forAll(faceId_, faceI)
201  {
202  faceId_[faceI] = faceI;
203  facePatchId_[faceI] = patchId;
204  faceSign_[faceI] = 1;
205  }
206 }
207 
208 
209 void Foam::fieldValues::faceSource::sampledSurfaceFaces(const dictionary& dict)
210 {
211  surfacePtr_ = sampledSurface::New
212  (
213  name_,
214  mesh(),
215  dict.subDict("sampledSurfaceDict")
216  );
217  surfacePtr_().update();
218  nFaces_ = returnReduce(surfacePtr_().faces().size(), sumOp<label>());
219 }
220 
221 
222 void Foam::fieldValues::faceSource::combineMeshGeometry
223 (
224  faceList& faces,
225  pointField& points
226 ) const
227 {
228  List<faceList> allFaces(Pstream::nProcs());
230 
231  labelList globalFacesIs(faceId_);
232  forAll(globalFacesIs, i)
233  {
234  if (facePatchId_[i] != -1)
235  {
236  label patchI = facePatchId_[i];
237  globalFacesIs[i] += mesh().boundaryMesh()[patchI].start();
238  }
239  }
240 
241  // Add local faces and points to the all* lists
243  (
244  IndirectList<face>(mesh().faces(), globalFacesIs),
245  mesh().points()
246  );
247  allFaces[Pstream::myProcNo()] = pp.localFaces();
248  allPoints[Pstream::myProcNo()] = pp.localPoints();
249 
250  Pstream::gatherList(allFaces);
251  Pstream::gatherList(allPoints);
252 
253  // Renumber and flatten
254  label nFaces = 0;
255  label nPoints = 0;
256  forAll(allFaces, procI)
257  {
258  nFaces += allFaces[procI].size();
259  nPoints += allPoints[procI].size();
260  }
261 
262  faces.setSize(nFaces);
263  points.setSize(nPoints);
264 
265  nFaces = 0;
266  nPoints = 0;
267 
268  // My own data first
269  {
270  const faceList& fcs = allFaces[Pstream::myProcNo()];
271  forAll(fcs, i)
272  {
273  const face& f = fcs[i];
274  face& newF = faces[nFaces++];
275  newF.setSize(f.size());
276  forAll(f, fp)
277  {
278  newF[fp] = f[fp] + nPoints;
279  }
280  }
281 
282  const pointField& pts = allPoints[Pstream::myProcNo()];
283  forAll(pts, i)
284  {
285  points[nPoints++] = pts[i];
286  }
287  }
288 
289  // Other proc data follows
290  forAll(allFaces, procI)
291  {
292  if (procI != Pstream::myProcNo())
293  {
294  const faceList& fcs = allFaces[procI];
295  forAll(fcs, i)
296  {
297  const face& f = fcs[i];
298  face& newF = faces[nFaces++];
299  newF.setSize(f.size());
300  forAll(f, fp)
301  {
302  newF[fp] = f[fp] + nPoints;
303  }
304  }
305 
306  const pointField& pts = allPoints[procI];
307  forAll(pts, i)
308  {
309  points[nPoints++] = pts[i];
310  }
311  }
312  }
313 
314  // Merge
315  labelList oldToNew;
316  pointField newPoints;
317  bool hasMerged = mergePoints
318  (
319  points,
320  SMALL,
321  false,
322  oldToNew,
323  newPoints
324  );
325 
326  if (hasMerged)
327  {
328  if (debug)
329  {
330  Pout<< "Merged from " << points.size()
331  << " down to " << newPoints.size() << " points" << endl;
332  }
333 
334  points.transfer(newPoints);
335  forAll(faces, i)
336  {
337  inplaceRenumber(oldToNew, faces[i]);
338  }
339  }
340 }
341 
342 
343 void Foam::fieldValues::faceSource::combineSurfaceGeometry
344 (
345  faceList& faces,
346  pointField& points
347 ) const
348 {
349  if (surfacePtr_.valid())
350  {
351  const sampledSurface& s = surfacePtr_();
352 
353  if (Pstream::parRun())
354  {
355  // dimension as fraction of mesh bounding box
356  scalar mergeDim = 1e-10*mesh().bounds().mag();
357 
358  labelList pointsMap;
359 
361  (
362  mergeDim,
364  (
365  SubList<face>(s.faces(), s.faces().size()),
366  s.points()
367  ),
368  points,
369  faces,
370  pointsMap
371  );
372  }
373  else
374  {
375  faces = s.faces();
376  points = s.points();
377  }
378  }
379 }
380 
381 
382 Foam::scalar Foam::fieldValues::faceSource::totalArea() const
383 {
384  scalar totalArea;
385 
386  if (surfacePtr_.valid())
387  {
388  totalArea = gSum(surfacePtr_().magSf());
389  }
390  else
391  {
392  totalArea = gSum(filterField(mesh().magSf(), false));
393  }
394 
395  return totalArea;
396 }
397 
398 
399 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
400 
402 {
403  dict.lookup("sourceName") >> sourceName_;
404 
405  switch (source_)
406  {
407  case stFaceZone:
408  {
409  setFaceZoneFaces();
410  break;
411  }
412  case stPatch:
413  {
414  setPatchFaces();
415  break;
416  }
417  case stSampledSurface:
418  {
419  sampledSurfaceFaces(dict);
420  break;
421  }
422  default:
423  {
424  FatalErrorIn("faceSource::initialise()")
425  << type() << " " << name_ << ": "
426  << sourceTypeNames_[source_] << "(" << sourceName_ << "):"
427  << nl << " Unknown source type. Valid source types are:"
428  << sourceTypeNames_.sortedToc() << nl << exit(FatalError);
429  }
430  }
431 
432  if (nFaces_ == 0)
433  {
434  WarningIn
435  (
436  "Foam::fieldValues::faceSource::initialise(const dictionary&)"
437  )
438  << type() << " " << name_ << ": "
439  << sourceTypeNames_[source_] << "(" << sourceName_ << "):" << nl
440  << " Source has no faces - deactivating" << endl;
441 
442  active_ = false;
443  return;
444  }
445 
446  if (surfacePtr_.valid())
447  {
448  surfacePtr_().update();
449  }
450 
451  totalArea_ = totalArea();
452 
453  Info<< type() << " " << name_ << ":" << nl
454  << " total faces = " << nFaces_
455  << nl
456  << " total area = " << totalArea_
457  << nl;
458 
459  if (dict.readIfPresent("weightField", weightFieldName_))
460  {
461  Info<< " weight field = " << weightFieldName_ << nl;
462 
463  if (source_ == stSampledSurface)
464  {
466  (
467  "void Foam::fieldValues::faceSource::initialise"
468  "("
469  "const dictionary&"
470  ")",
471  dict
472  )
473  << "Cannot use weightField for a sampledSurface"
474  << exit(FatalIOError);
475  }
476  }
477 
478  if (dict.found("orientedWeightField"))
479  {
480  if (weightFieldName_ == "none")
481  {
482  dict.lookup("orientedWeightField") >> weightFieldName_;
483  Info<< " weight field = " << weightFieldName_ << nl;
484  orientWeightField_ = true;
485  }
486  else
487  {
489  (
490  "void Foam::fieldValues::faceSource::initialise"
491  "("
492  "const dictionary&"
493  ")",
494  dict
495  )
496  << "Either weightField or orientedWeightField can be supplied, "
497  << "but not both"
498  << exit(FatalIOError);
499  }
500  }
501 
502  List<word> orientedFields;
503  if (dict.readIfPresent("orientedFields", orientedFields))
504  {
505  orientedFieldsStart_ = fields_.size();
506  fields_.append(orientedFields);
507  }
508 
509  if (dict.readIfPresent("scaleFactor", scaleFactor_))
510  {
511  Info<< " scale factor = " << scaleFactor_ << nl;
512  }
513 
514  Info<< nl << endl;
515 
516  if (valueOutput_)
517  {
518  const word surfaceFormat(dict.lookup("surfaceFormat"));
519 
520  surfaceWriterPtr_.reset
521  (
523  (
524  surfaceFormat,
525  dict.subOrEmptyDict("formatOptions").
526  subOrEmptyDict(surfaceFormat)
527  ).ptr()
528  );
529  }
530 }
531 
532 
534 {
535  writeCommented(file(), "Source : ");
536  file() << sourceTypeNames_[source_] << " " << sourceName_ << endl;
537  writeCommented(file(), "Faces : ");
538  file() << nFaces_ << endl;
539  writeCommented(file(), "Area : ");
540  file() << totalArea_ << endl;
541 
542  writeCommented(file(), "Time");
543  if (writeArea_)
544  {
545  file() << tab << "Area";
546  }
547 
548  forAll(fields_, i)
549  {
550  file()
551  << tab << operationTypeNames_[operation_]
552  << "(" << fields_[i] << ")";
553  }
554 
555  file() << endl;
556 }
557 
558 
559 template<>
561 (
562  const Field<scalar>& values,
563  const vectorField& Sf,
564  const scalarField& weightField
565 ) const
566 {
567  switch (operation_)
568  {
569  case opSumDirection:
570  {
571  vector n(dict_.lookup("direction"));
572  return sum(pos(values*(Sf & n))*mag(values));
573  }
574  case opSumDirectionBalance:
575  {
576  vector n(dict_.lookup("direction"));
577  const scalarField nv(values*(Sf & n));
578 
579  return sum(pos(nv)*mag(values) - neg(nv)*mag(values));
580  }
581  default:
582  {
583  // Fall through to other operations
584  return processSameTypeValues(values, Sf, weightField);
585  }
586  }
587 }
588 
589 
590 template<>
592 (
593  const Field<vector>& values,
594  const vectorField& Sf,
595  const scalarField& weightField
596 ) const
597 {
598  switch (operation_)
599  {
600  case opSumDirection:
601  {
602  vector n(dict_.lookup("direction"));
603  n /= mag(n) + ROOTVSMALL;
604  const scalarField nv(n & values);
605 
606  return sum(pos(nv)*n*(nv));
607  }
608  case opSumDirectionBalance:
609  {
610  vector n(dict_.lookup("direction"));
611  n /= mag(n) + ROOTVSMALL;
612  const scalarField nv(n & values);
613 
614  return sum(pos(nv)*n*(nv));
615  }
616  case opAreaNormalAverage:
617  {
618  scalar result = sum(values & Sf)/sum(mag(Sf));
619  return vector(result, 0.0, 0.0);
620  }
621  case opAreaNormalIntegrate:
622  {
623  scalar result = sum(values & Sf);
624  return vector(result, 0.0, 0.0);
625  }
626  default:
627  {
628  // Fall through to other operations
629  return processSameTypeValues(values, Sf, weightField);
630  }
631  }
632 }
633 
634 
635 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
636 
638 (
639  const word& name,
640  const objectRegistry& obr,
641  const dictionary& dict,
642  const bool loadFromFiles
643 )
644 :
645  fieldValue(name, obr, dict, typeName, loadFromFiles),
646  surfaceWriterPtr_(NULL),
647  source_(sourceTypeNames_.read(dict.lookup("source"))),
648  operation_(operationTypeNames_.read(dict.lookup("operation"))),
649  weightFieldName_("none"),
650  orientWeightField_(false),
651  orientedFieldsStart_(labelMax),
652  scaleFactor_(1.0),
653  writeArea_(dict.lookupOrDefault("writeArea", false)),
654  nFaces_(0),
655  faceId_(),
656  facePatchId_(),
657  faceSign_()
658 {
659  read(dict);
660 }
661 
662 
663 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
664 
666 {}
667 
668 
669 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
670 
672 {
673  fieldValue::read(dict);
674 
675  if (active_)
676  {
677  initialise(dict);
678  }
679 }
680 
681 
683 {
685 
686  if (active_)
687  {
688  if (surfacePtr_.valid())
689  {
690  surfacePtr_().update();
691  }
692 
693  if (Pstream::master())
694  {
695  file() << obr_.time().value();
696  }
697 
698  if (writeArea_)
699  {
700  totalArea_ = totalArea();
701  if (Pstream::master())
702  {
703  file() << tab << totalArea_;
704  }
705  if (log_) Info<< " total area = " << totalArea_ << endl;
706  }
707 
708  // construct weight field. Note: zero size means weight = 1
709  scalarField weightField;
710  if (weightFieldName_ != "none")
711  {
712  weightField =
713  getFieldValues<scalar>
714  (
715  weightFieldName_,
716  true,
717  orientWeightField_
718  );
719  }
720 
721  // Combine onto master
722  combineFields(weightField);
723 
724  // process the fields
725  forAll(fields_, i)
726  {
727  const word& fieldName = fields_[i];
728  bool ok = false;
729 
730  bool orient = i >= orientedFieldsStart_;
731  ok = ok || writeValues<scalar>(fieldName, weightField, orient);
732  ok = ok || writeValues<vector>(fieldName, weightField, orient);
733  ok = ok
734  || writeValues<sphericalTensor>(fieldName, weightField, orient);
735  ok = ok || writeValues<symmTensor>(fieldName, weightField, orient);
736  ok = ok || writeValues<tensor>(fieldName, weightField, orient);
737 
738  if (!ok)
739  {
740  WarningIn("void Foam::fieldValues::faceSource::write()")
741  << "Requested field " << fieldName
742  << " not found in database and not processed"
743  << endl;
744  }
745  }
746 
747  if (Pstream::master())
748  {
749  file()<< endl;
750  }
751 
752  if (log_) Info<< endl;
753  }
754 }
755 
756 
757 // ************************************************************************* //
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:380
A List with indirect addressing.
Definition: IndirectList.H:102
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
const pointField & points
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
static const NamedEnum< sourceType, 3 > sourceTypeNames_
Source type names.
Definition: faceSource.H:321
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensionedScalar neg(const dimensionedScalar &ds)
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 ))
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:390
dimensioned< scalar > mag(const dimensioned< Type > &)
static const NamedEnum< operationType, 15 > operationTypeNames_
Operation type names.
Definition: faceSource.H:345
An abstract class for surfaces with sampling.
labelList f(nPoints)
Merge points. See below.
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
A class for handling words, derived from string.
Definition: word.H:59
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
virtual ~faceSource()
Destructor.
Definition: faceSource.C:665
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:76
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Base class for field value -based function objects.
Definition: fieldValue.H:62
messageStream Info
dynamicFvMesh & mesh
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
virtual void write()
Write to screen/file.
Definition: fieldValue.C:55
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:638
void initialise(const dictionary &dict)
Initialise, e.g. face addressing.
Definition: faceSource.C:401
Namespace for OpenFOAM.
virtual void read(const dictionary &dict)
Read from dictionary.
Definition: fieldValue.C:42
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Type gSum(const FieldField< Field, Type > &f)
label n
virtual void writeFileHeader(const label i)
Output file header information.
Definition: faceSource.C:533
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
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.
const Field< PointType > & localPoints() const
Return pointField of points in patch.
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:386
static const char nl
Definition: Ostream.H:260
const double e
Elementary charge.
Definition: doubleFloat.H:78
void setSize(const label)
Reset size of List.
Definition: List.C:318
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
IOerror FatalIOError
#define WarningIn(functionName)
Report a warning using Foam::Warning.
label faceId(-1)
#define forAll(list, i)
Definition: UList.H:421
faceSource(const word &name, const objectRegistry &obr, const dictionary &dict, const bool loadFromFiles=false)
Construct from components.
Definition: faceSource.C:638
Macros for easy insertion into run-time selection tables.
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:404
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
A list of faces which address into the list of points.
addToRunTimeSelectionTable(fieldValue, cellSource, dictionary)
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
label patchId(-1)
virtual void read(const dictionary &)
Read from dictionary.
Definition: faceSource.C:671
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:300
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.
static const char tab
Definition: Ostream.H:259
static autoPtr< sampledSurface > New(const word &name, const polyMesh &, const dictionary &)
Return a reference to the selected surface.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:452
error FatalError
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Registry of regIOobjects.
defineTypeNameAndDebug(cellSource, 0)
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.
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:398
bool read(const char *, int32_t &)
Definition: int32IO.C:87
virtual const pointField & points() const =0
Points of surface.
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
A List obtained as a section of another List.
Definition: SubList.H:53
virtual const faceList & faces() const =0
Faces of surface.
dimensionedScalar pos(const dimensionedScalar &ds)
label nPoints
#define FatalIOErrorIn(functionName, ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:325
virtual void write()
Calculate and write.
Definition: faceSource.C:682
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:377
static autoPtr< surfaceWriter > New(const word &writeType)
Return a reference to the selected surfaceWriter.
Definition: surfaceWriter.C:55
dictionary subOrEmptyDict(const word &, const bool mustRead=false) const
Find and return a sub-dictionary as a copy, or.
Definition: dictionary.C:675
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
static const label labelMax
Definition: label.H:62