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