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 
88 {
90 
91  if (zoneId < 0)
92  {
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);
163 
164  if (debug)
165  {
166  Pout<< "Original face zone size = " << fZone.size()
167  << ", new size = " << faceId_.size() << endl;
168  }
169 }
170 
171 
173 {
174  const label patchId = mesh().boundaryMesh().findPatchID(sourceName_);
175 
176  if (patchId < 0)
177  {
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 
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 
223 (
224  faceList& faces,
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();
249 
250  Pstream::gatherList(allFaces);
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 
344 (
345  faceList& faces,
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 
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  switch (source_)
404  {
405  case stFaceZone:
406  {
407  setFaceZoneFaces();
408  break;
409  }
410  case stPatch:
411  {
412  setPatchFaces();
413  break;
414  }
415  case stSampledSurface:
416  {
417  sampledSurfaceFaces(dict);
418  break;
419  }
420  default:
421  {
423  << type() << " " << name_ << ": "
424  << sourceTypeNames_[source_] << "(" << sourceName_ << "):"
425  << nl << " Unknown source type. Valid source types are:"
426  << sourceTypeNames_.sortedToc() << nl << exit(FatalError);
427  }
428  }
429 
430  if (nFaces_ == 0)
431  {
433  << type() << " " << name_ << ": "
434  << sourceTypeNames_[source_] << "(" << sourceName_ << "):" << nl
435  << " Source has no faces - deactivating" << endl;
436 
437  active_ = false;
438  return;
439  }
440 
441  if (surfacePtr_.valid())
442  {
443  surfacePtr_().update();
444  }
445 
446  totalArea_ = totalArea();
447 
448  if (log_) Info
449  << type() << " " << name_ << ":" << nl
450  << " total faces = " << nFaces_ << nl
451  << " total area = " << totalArea_ << nl;
452 
453  if (dict.readIfPresent("weightField", weightFieldName_))
454  {
455  if (log_) Info << " weight field = " << weightFieldName_ << nl;
456 
457  if (source_ == stSampledSurface)
458  {
460  << "Cannot use weightField for a sampledSurface"
461  << exit(FatalIOError);
462  }
463  }
464 
465  if (dict.found("orientedWeightField"))
466  {
467  if (weightFieldName_ == "none")
468  {
469  dict.lookup("orientedWeightField") >> weightFieldName_;
470  if (log_) Info << " weight field = " << weightFieldName_ << nl;
471  orientWeightField_ = true;
472  }
473  else
474  {
476  << "Either weightField or orientedWeightField can be supplied, "
477  << "but not both"
478  << exit(FatalIOError);
479  }
480  }
481 
482  List<word> orientedFields;
483  if (dict.readIfPresent("orientedFields", orientedFields))
484  {
485  orientedFieldsStart_ = fields_.size();
486  fields_.append(orientedFields);
487  }
488 
489  if (log_) Info << nl << endl;
490 
491  if (valueOutput_)
492  {
493  const word surfaceFormat(dict.lookup("surfaceFormat"));
494 
495  surfaceWriterPtr_.reset
496  (
498  (
499  surfaceFormat,
500  dict.subOrEmptyDict("formatOptions").
501  subOrEmptyDict(surfaceFormat)
502  ).ptr()
503  );
504  }
505 }
506 
507 
509 {
510  writeHeaderValue(os, "Source", sourceTypeNames_[source_]);
511  writeHeaderValue(os, "Name", sourceName_);
512  writeHeaderValue(os, "Faces", nFaces_);
513  writeHeaderValue(os, "Total area", totalArea_);
514  writeHeaderValue(os, "Scale factor", scaleFactor_);
515 
516  writeCommented(os, "Time");
517  if (writeArea_)
518  {
519  os << tab << "Area";
520  }
521 
522  forAll(fields_, i)
523  {
524  os << tab << operationTypeNames_[operation_]
525  << "(" << fields_[i] << ")";
526  }
527 
528  os << endl;
529 }
530 
531 
532 template<>
534 (
535  const Field<scalar>& values,
536  const vectorField& Sf,
537  const scalarField& weightField
538 ) const
539 {
540  switch (operation_)
541  {
542  case opSumDirection:
543  {
544  vector n(dict_.lookup("direction"));
545  return gSum(pos(values*(Sf & n))*mag(values));
546  }
547  case opSumDirectionBalance:
548  {
549  vector n(dict_.lookup("direction"));
550  const scalarField nv(values*(Sf & n));
551 
552  return gSum(pos(nv)*mag(values) - neg(nv)*mag(values));
553  }
554  default:
555  {
556  // Fall through to other operations
557  return processSameTypeValues(values, Sf, weightField);
558  }
559  }
560 }
561 
562 
563 template<>
565 (
566  const Field<vector>& values,
567  const vectorField& Sf,
568  const scalarField& weightField
569 ) const
570 {
571  switch (operation_)
572  {
573  case opSumDirection:
574  {
575  vector n(dict_.lookup("direction"));
576  n /= mag(n) + ROOTVSMALL;
577  const scalarField nv(n & values);
578 
579  return gSum(pos(nv)*n*(nv));
580  }
581  case opSumDirectionBalance:
582  {
583  vector n(dict_.lookup("direction"));
584  n /= mag(n) + ROOTVSMALL;
585  const scalarField nv(n & values);
586 
587  return gSum(pos(nv)*n*(nv));
588  }
589  case opAreaNormalAverage:
590  {
591  scalar result = gSum(values & Sf)/gSum(mag(Sf));
592  return vector(result, 0.0, 0.0);
593  }
594  case opAreaNormalIntegrate:
595  {
596  scalar result = gSum(values & Sf);
597  return vector(result, 0.0, 0.0);
598  }
599  default:
600  {
601  // Fall through to other operations
602  return processSameTypeValues(values, Sf, weightField);
603  }
604  }
605 }
606 
607 
608 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
609 
611 (
612  const word& name,
613  const objectRegistry& obr,
614  const dictionary& dict,
615  const bool loadFromFiles
616 )
617 :
618  fieldValue(name, obr, dict, typeName, loadFromFiles),
619  surfaceWriterPtr_(NULL),
620  source_(sourceTypeNames_.read(dict.lookup("source"))),
621  operation_(operationTypeNames_.read(dict.lookup("operation"))),
622  weightFieldName_("none"),
623  orientWeightField_(false),
624  orientedFieldsStart_(labelMax),
625  writeArea_(dict.lookupOrDefault("writeArea", false)),
626  nFaces_(0),
627  faceId_(),
628  facePatchId_(),
629  faceSign_()
630 {
631  if (active_)
632  {
633  read(dict);
634  writeFileHeader(file());
635  }
636 }
637 
638 
639 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
640 
642 {}
643 
644 
645 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
646 
648 {
649  if (active_)
650  {
652 
653  // No additional info to read
654  initialise(dict);
655  }
656 }
657 
658 
660 {
662 
663  if (active_)
664  {
665  if (surfacePtr_.valid())
666  {
667  surfacePtr_().update();
668  }
669 
670  writeTime(file());
671 
672  if (writeArea_)
673  {
674  totalArea_ = totalArea();
675  file() << tab << totalArea_;
676  if (log_) Info<< " total area = " << totalArea_ << endl;
677  }
678 
679  // Construct weight field. Note: zero size indicates unweighted
680  scalarField weightField;
681  if (weightFieldName_ != "none")
682  {
683  weightField =
684  setFieldValues<scalar>
685  (
686  weightFieldName_,
687  true,
688  orientWeightField_
689  );
690  }
691 
692  // Process the fields
693  forAll(fields_, i)
694  {
695  const word& fieldName = fields_[i];
696  bool ok = false;
697 
698  bool orient = i >= orientedFieldsStart_;
699  ok = ok || writeValues<scalar>(fieldName, weightField, orient);
700  ok = ok || writeValues<vector>(fieldName, weightField, orient);
701  ok = ok ||
702  writeValues<sphericalTensor>(fieldName, weightField, orient);
703  ok = ok || writeValues<symmTensor>(fieldName, weightField, orient);
704  ok = ok || writeValues<tensor>(fieldName, weightField, orient);
705 
706  if (!ok)
707  {
709  << "Requested field " << fieldName
710  << " not found in database and not processed"
711  << endl;
712  }
713  }
714 
715  file()<< endl;
716 
717  if (log_) Info<< endl;
718  }
719 }
720 
721 
722 // ************************************************************************* //
Foam::fieldValues::faceSource::write
virtual void write()
Calculate and write.
Definition: faceSource.C:659
Foam::fieldValues::faceSource::nFaces_
label nFaces_
Global number of faces.
Definition: faceSource.H:408
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::fieldValues::faceSource::setPatchFaces
void setPatchFaces()
Set faces to evaluate based on a patch.
Definition: faceSource.C:172
cyclicPolyPatch.H
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::DynamicList< label >
Foam::PrimitivePatch::localPoints
const Field< PointType > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatchTemplate.C:432
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::labelMax
static const label labelMax
Definition: label.H:62
Foam::UPstream::nProcs
static label nProcs(const label communicator=0)
Number of processes in parallel run.
Definition: UPstream.H:387
Foam::fieldValues::faceSource::sampledSurfaceFaces
void sampledSurfaceFaces(const dictionary &)
Set faces according to sampledSurface.
Definition: faceSource.C:209
Foam::UPstream::parRun
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:377
Foam::sampledSurface::New
static autoPtr< sampledSurface > New(const word &name, const polyMesh &, const dictionary &)
Return a reference to the selected surface.
Definition: sampledSurface.C:117
Foam::fieldValue::sourceName_
word sourceName_
Name of source object.
Definition: fieldValue.H:83
Foam::List::transfer
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Foam::fieldValues::faceSource::faceSign_
labelList faceSign_
List of +1/-1 representing face flip map.
Definition: faceSource.H:421
Foam::FatalIOError
IOerror FatalIOError
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::gSum
Type gSum(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:564
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::fieldValues::faceSource::faceId
const labelList & faceId() const
Return the local list of face IDs.
Definition: faceSourceI.H:38
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::inplaceRenumber
void inplaceRenumber(const labelUList &oldToNew, ListType &)
Inplace renumber the values of a list.
Definition: ListOpsTemplates.C:57
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:50
Foam::fieldValue::read
virtual void read(const dictionary &dict)
Read from dictionary.
Definition: fieldValue.C:42
n
label n
Definition: TABSMDCalcMethod2.H:31
coupledPolyPatch.H
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:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::fieldValues::faceSource::processValues
Type processValues(const Field< Type > &values, const vectorField &Sf, const scalarField &weightField) const
Apply the 'operation' to the values. Wrapper around.
Definition: faceSourceTemplates.C:261
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::fieldValues::faceSource::combineSurfaceGeometry
void combineSurfaceGeometry(faceList &faces, pointField &points) const
Combine surface faces and points from multiple processors.
Definition: faceSource.C:344
sampledSurface.H
Foam::fieldValues::faceSource::writeFileHeader
virtual void writeFileHeader(Ostream &os) const
Output file header information.
Definition: faceSource.C:508
Foam::IndirectList
A List with indirect addressing.
Definition: IndirectList.H:102
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:703
Foam::sampledSurface
An abstract class for surfaces with sampling.
Definition: sampledSurface.H:77
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::fieldValues::faceSource::initialise
void initialise(const dictionary &dict)
Initialise, e.g. face addressing.
Definition: faceSource.C:401
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::ZoneMesh::findZoneID
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:348
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::vector
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
s
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))
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
emptyPolyPatch.H
Foam::ZoneMesh::names
wordList names() const
Return a list of zone names.
Definition: ZoneMesh.C:263
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::UPstream::myProcNo
static int myProcNo(const label communicator=0)
Number of this process (starting from masterNo() = 0)
Definition: UPstream.H:405
Foam::fieldValues::faceSource::~faceSource
virtual ~faceSource()
Destructor.
Definition: faceSource.C:641
Foam::fieldValues::faceSource::read
virtual void read(const dictionary &)
Read from dictionary.
Definition: faceSource.C:647
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::fieldValues::faceSource::source_
sourceType source_
Source type.
Definition: faceSource.H:387
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
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:49
Foam::polyPatch::whichFace
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:389
Foam::sumOp
Definition: ops.H:162
Foam::tab
static const char tab
Definition: Ostream.H:259
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::fieldValue
Base class for field value-based function objects.
Definition: fieldValue.H:63
f
labelList f(nPoints)
Foam::functionObjectState::name_
const word name_
Name of model.
Definition: functionObjectState.H:72
Foam::surfaceWriter::New
static autoPtr< surfaceWriter > New(const word &writeType)
Return a reference to the selected surfaceWriter.
Definition: surfaceWriter.C:55
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
points
const pointField & points
Definition: gmvOutputHeader.H:1
mergePoints.H
Merge points. See below.
Foam::PrimitivePatch::localFaces
const List< Face > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatchTemplate.C:372
Foam::fieldValues::faceSource::operationTypeNames_
static const NamedEnum< operationType, 15 > operationTypeNames_
Operation type names.
Definition: faceSource.H:345
Foam::fieldValues::faceSource::faceSource
faceSource(const word &name, const objectRegistry &obr, const dictionary &dict, const bool loadFromFiles=false)
Construct from components.
Definition: faceSource.C:611
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:330
patchId
label patchId(-1)
Foam::DelaunayMeshTools::allPoints
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::fieldValue::write
virtual void write()
Write to screen/file.
Definition: fieldValue.C:58
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::fieldValues::addToRunTimeSelectionTable
addToRunTimeSelectionTable(fieldValue, cellSource, dictionary)
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:43
Foam::fieldValue::mesh
const fvMesh & mesh() const
Helper function to return the reference to the mesh.
Definition: fieldValueI.H:79
Foam::faceZone::flipMap
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
Foam::fieldValues::faceSource::facePatchId_
labelList facePatchId_
Local list of patch ID per face.
Definition: faceSource.H:417
Foam::neg
dimensionedScalar neg(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:201
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::mergePoints
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.
Foam::fieldValues::faceSource::totalArea
scalar totalArea() const
Calculate and return total area of the faceSource: sum(magSf)
Definition: faceSource.C:382
Foam::fieldValues::faceSource::combineMeshGeometry
void combineMeshGeometry(faceList &faces, pointField &points) const
Combine mesh faces and points from multiple processors.
Definition: faceSource.C:223
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::fieldValues::faceSource::faceId_
labelList faceId_
Local list of face IDs.
Definition: faceSource.H:414
Foam::fieldValues::defineTypeNameAndDebug
defineTypeNameAndDebug(cellSource, 0)
Foam::fieldValues::faceSource::setFaceZoneFaces
void setFaceZoneFaces()
Set faces to evaluate based on a face zone.
Definition: faceSource.C:87
faceSource.H
Foam::NamedEnum
Initialise the NamedEnum HashTable from the static list of names.
Definition: NamedEnum.H:52
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88
Foam::fieldValues::faceSource::sourceTypeNames_
static const NamedEnum< sourceType, 3 > sourceTypeNames_
Source type names.
Definition: faceSource.H:321
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190