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-2021 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 "emptyPolyPatch.H"
28 #include "coupledPolyPatch.H"
29 #include "sampledSurface.H"
30 #include "mergePoints.H"
31 #include "indirectPrimitivePatch.H"
32 #include "PatchTools.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace functionObjects
40 {
41 namespace fieldValues
42 {
43  defineTypeNameAndDebug(surfaceFieldValue, 0);
44  addToRunTimeSelectionTable(fieldValue, surfaceFieldValue, dictionary);
45  addToRunTimeSelectionTable(functionObject, surfaceFieldValue, dictionary);
46 }
47 }
48 }
49 
50 template<>
51 const char* Foam::NamedEnum
52 <
54  3
55 >::names[] =
56 {
57  "faceZone",
58  "patch",
59  "sampledSurface"
60 };
61 
62 template<>
63 const char* Foam::NamedEnum
64 <
66  16
67 >::names[] =
68 {
69  "none",
70  "sum",
71  "sumMag",
72  "sumDirection",
73  "sumDirectionBalance",
74  "orientedSum",
75  "average",
76  "areaAverage",
77  "areaIntegrate",
78  "min",
79  "max",
80  "minMag",
81  "maxMag",
82  "CoV",
83  "areaNormalAverage",
84  "areaNormalIntegrate"
85 };
86 
87 const Foam::NamedEnum
88 <
90  3
92 
93 const Foam::NamedEnum
94 <
96  16
98 
99 
100 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
101 
102 void Foam::functionObjects::fieldValues::surfaceFieldValue::setFaceZoneFaces()
103 {
104  label zoneId = mesh_.faceZones().findZoneID(regionName_);
105 
106  if (zoneId < 0)
107  {
109  << type() << " " << name() << ": "
110  << regionTypeNames_[regionType_] << "(" << regionName_ << "):" << nl
111  << " Unknown face zone name: " << regionName_
112  << ". Valid face zones are: " << mesh_.faceZones().names()
113  << nl << exit(FatalError);
114  }
115 
116  const faceZone& fZone = mesh_.faceZones()[zoneId];
117 
118  DynamicList<label> faceIds(fZone.size());
119  DynamicList<label> facePatchIds(fZone.size());
120  DynamicList<label> faceSigns(fZone.size());
121 
122  forAll(fZone, i)
123  {
124  label facei = fZone[i];
125 
126  label faceId = -1;
127  label facePatchId = -1;
128  if (mesh_.isInternalFace(facei))
129  {
130  faceId = facei;
131  facePatchId = -1;
132  }
133  else
134  {
135  facePatchId = mesh_.boundaryMesh().whichPatch(facei);
136  const polyPatch& pp = mesh_.boundaryMesh()[facePatchId];
137  if (isA<coupledPolyPatch>(pp))
138  {
139  if (refCast<const coupledPolyPatch>(pp).owner())
140  {
141  faceId = pp.whichFace(facei);
142  }
143  else
144  {
145  faceId = -1;
146  }
147  }
148  else if (!isA<emptyPolyPatch>(pp))
149  {
150  faceId = facei - pp.start();
151  }
152  else
153  {
154  faceId = -1;
155  facePatchId = -1;
156  }
157  }
158 
159  if (faceId >= 0)
160  {
161  if (fZone.flipMap()[i])
162  {
163  faceSigns.append(-1);
164  }
165  else
166  {
167  faceSigns.append(1);
168  }
169  faceIds.append(faceId);
170  facePatchIds.append(facePatchId);
171  }
172  }
173 
174  faceId_.transfer(faceIds);
175  facePatchId_.transfer(facePatchIds);
176  faceSign_.transfer(faceSigns);
177  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
178 
179  if (debug)
180  {
181  Pout<< "Original face zone size = " << fZone.size()
182  << ", new size = " << faceId_.size() << endl;
183  }
184 }
185 
186 
187 void Foam::functionObjects::fieldValues::surfaceFieldValue::setPatchFaces()
188 {
189  const label patchid = mesh_.boundaryMesh().findPatchID(regionName_);
190 
191  if (patchid < 0)
192  {
194  << type() << " " << name() << ": "
195  << regionTypeNames_[regionType_] << "(" << regionName_ << "):" << nl
196  << " Unknown patch name: " << regionName_
197  << ". Valid patch names are: "
198  << mesh_.boundaryMesh().names() << nl
199  << exit(FatalError);
200  }
201 
202  const polyPatch& pp = mesh_.boundaryMesh()[patchid];
203 
204  label nFaces = pp.size();
205  if (isA<emptyPolyPatch>(pp))
206  {
207  nFaces = 0;
208  }
209 
210  faceId_.setSize(nFaces);
211  facePatchId_.setSize(nFaces);
212  faceSign_.setSize(nFaces);
213  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
214 
215  forAll(faceId_, facei)
216  {
217  faceId_[facei] = facei;
218  facePatchId_[facei] = patchid;
219  faceSign_[facei] = 1;
220  }
221 }
222 
223 
224 void Foam::functionObjects::fieldValues::surfaceFieldValue::sampledSurfaceFaces
225 (
226  const dictionary& dict
227 )
228 {
229  surfacePtr_ = sampledSurface::New
230  (
231  name(),
232  mesh_,
233  dict.subDict("sampledSurfaceDict")
234  );
235  surfacePtr_().update();
236  nFaces_ = returnReduce(surfacePtr_().faces().size(), sumOp<label>());
237 }
238 
239 
240 void Foam::functionObjects::fieldValues::surfaceFieldValue::combineMeshGeometry
241 (
242  faceList& faces,
243  pointField& points
244 ) const
245 {
246  List<faceList> allFaces(Pstream::nProcs());
248 
249  labelList globalFacesIs(faceId_);
250  forAll(globalFacesIs, i)
251  {
252  if (facePatchId_[i] != -1)
253  {
254  label patchi = facePatchId_[i];
255  globalFacesIs[i] += mesh_.boundaryMesh()[patchi].start();
256  }
257  }
258 
259  // Add local faces and points to the all* lists
261  (
262  IndirectList<face>(mesh_.faces(), globalFacesIs),
263  mesh_.points()
264  );
265  allFaces[Pstream::myProcNo()] = pp.localFaces();
266  allPoints[Pstream::myProcNo()] = pp.localPoints();
267 
268  Pstream::gatherList(allFaces);
269  Pstream::gatherList(allPoints);
270 
271  // Renumber and flatten
272  label nFaces = 0;
273  label nPoints = 0;
274  forAll(allFaces, proci)
275  {
276  nFaces += allFaces[proci].size();
277  nPoints += allPoints[proci].size();
278  }
279 
280  faces.setSize(nFaces);
281  points.setSize(nPoints);
282 
283  nFaces = 0;
284  nPoints = 0;
285 
286  // My own data first
287  {
288  const faceList& fcs = allFaces[Pstream::myProcNo()];
289  forAll(fcs, i)
290  {
291  const face& f = fcs[i];
292  face& newF = faces[nFaces++];
293  newF.setSize(f.size());
294  forAll(f, fp)
295  {
296  newF[fp] = f[fp] + nPoints;
297  }
298  }
299 
300  const pointField& pts = allPoints[Pstream::myProcNo()];
301  forAll(pts, i)
302  {
303  points[nPoints++] = pts[i];
304  }
305  }
306 
307  // Other proc data follows
308  forAll(allFaces, proci)
309  {
310  if (proci != Pstream::myProcNo())
311  {
312  const faceList& fcs = allFaces[proci];
313  forAll(fcs, i)
314  {
315  const face& f = fcs[i];
316  face& newF = faces[nFaces++];
317  newF.setSize(f.size());
318  forAll(f, fp)
319  {
320  newF[fp] = f[fp] + nPoints;
321  }
322  }
323 
324  const pointField& pts = allPoints[proci];
325  forAll(pts, i)
326  {
327  points[nPoints++] = pts[i];
328  }
329  }
330  }
331 
332  // Merge
333  labelList oldToNew;
334  pointField newPoints;
335  bool hasMerged = mergePoints
336  (
337  points,
338  small,
339  false,
340  oldToNew,
341  newPoints
342  );
343 
344  if (hasMerged)
345  {
346  if (debug)
347  {
348  Pout<< "Merged from " << points.size()
349  << " down to " << newPoints.size() << " points" << endl;
350  }
351 
352  points.transfer(newPoints);
353  forAll(faces, i)
354  {
355  inplaceRenumber(oldToNew, faces[i]);
356  }
357  }
358 }
359 
360 
361 void Foam::functionObjects::fieldValues::surfaceFieldValue::
362 combineSurfaceGeometry
363 (
364  faceList& faces,
365  pointField& points
366 ) const
367 {
368  if (regionType_ == regionTypes::sampledSurface)
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::surfaceFieldValue::totalArea() const
403 {
404  if (regionType_ == regionTypes::sampledSurface)
405  {
406  return gSum(surfacePtr_().magSf());
407  }
408  else
409  {
410  return gSum(filterField(mesh_.magSf()));
411  }
412 }
413 
414 
415 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
416 
418 (
419  const dictionary& dict
420 )
421 {
422  switch (regionType_)
423  {
424  case regionTypes::faceZone:
425  {
426  dict.lookup("name") >> regionName_;
427  setFaceZoneFaces();
428  break;
429  }
430  case regionTypes::patch:
431  {
432  dict.lookup("name") >> regionName_;
433  setPatchFaces();
434  break;
435  }
436  case regionTypes::sampledSurface:
437  {
438  sampledSurfaceFaces(dict);
439  regionName_ = surfacePtr_().name();
440  break;
441  }
442  default:
443  {
445  << type() << " " << name() << ": "
446  << regionTypeNames_[regionType_] << "(" << regionName_ << "):"
447  << nl << " Unknown region type. Valid region types are:"
448  << regionTypeNames_.sortedToc() << nl << exit(FatalError);
449  }
450  }
451 
452  if (nFaces_ == 0)
453  {
455  << type() << " " << name() << ": "
456  << regionTypeNames_[regionType_] << "(" << regionName_ << "):" << nl
457  << " Region has no faces" << exit(FatalError);
458  }
459 
460  if (regionType_ == regionTypes::sampledSurface)
461  {
462  surfacePtr_().update();
463  }
464 
465  totalArea_ = totalArea();
466 
467  Info<< type() << " " << name() << ":" << nl
468  << " total faces = " << nFaces_
469  << nl
470  << " total area = " << totalArea_
471  << nl;
472 
473  if (dict.readIfPresent("weightFields", weightFieldNames_))
474  {
475  Info<< name() << " " << operationTypeNames_[operation_]
476  << " weight fields " << weightFieldNames_;
477  }
478  else if (dict.found("weightField"))
479  {
480  weightFieldNames_.setSize(1);
481  dict.lookup("weightField") >> weightFieldNames_[0];
482 
483  Info<< name() << " " << operationTypeNames_[operation_]
484  << " weight field " << weightFieldNames_[0];
485  }
486 
487  if (dict.readIfPresent("scaleFactor", scaleFactor_))
488  {
489  Info<< " scale factor = " << scaleFactor_ << nl;
490  }
491 
492  Info<< nl << endl;
493 
494  if (writeFields_)
495  {
496  const word surfaceFormat(dict.lookup("surfaceFormat"));
497 
498  surfaceWriterPtr_.reset
499  (
500  surfaceWriter::New(surfaceFormat, dict).ptr()
501  );
502  }
503 }
504 
505 
507 (
508  const label i
509 )
510 {
511  if (operation_ != operationType::none)
512  {
513  writeCommented(file(), "Region type : ");
514  file() << regionTypeNames_[regionType_] << " " << regionName_ << endl;
515  writeCommented(file(), "Faces : ");
516  file() << nFaces_ << endl;
517  writeCommented(file(), "Area : ");
518  file() << totalArea_ << endl;
519 
520  writeCommented(file(), "Time");
521  if (writeArea_)
522  {
523  file() << tab << "Area";
524  }
525 
526  forAll(fields_, fieldi)
527  {
528  file()
529  << tab << operationTypeNames_[operation_]
530  << "(" << fields_[fieldi] << ")";
531  }
532 
533  file() << endl;
534  }
535 }
536 
537 
539 (
540  const Field<scalar>& values,
541  const scalarField& signs,
542  const scalarField& weights,
543  const vectorField& Sf,
544  scalar& result
545 ) const
546 {
547  switch (operation_)
548  {
549  case operationType::sumDirection:
550  {
551  const vector n(dict_.lookup("direction"));
552  result = gSum(weights*pos0(values*(Sf & n))*mag(values));
553  return true;
554  }
555  case operationType::sumDirectionBalance:
556  {
557  const vector n(dict_.lookup("direction"));
558  const scalarField nv(values*(Sf & n));
559  result = gSum(weights*(pos0(nv)*mag(values) - neg(nv)*mag(values)));
560  return true;
561  }
562  default:
563  {
564  // Fall through to same-type operations
565  return processValuesTypeType
566  (
567  values,
568  signs,
569  weights,
570  Sf,
571  result
572  );
573  }
574  }
575 }
576 
577 
579 (
580  const Field<vector>& values,
581  const scalarField& signs,
582  const scalarField& weights,
583  const vectorField& Sf,
584  scalar& result
585 ) const
586 {
587  switch (operation_)
588  {
589  case operationType::areaNormalAverage:
590  {
591  result = gSum(weights*values & Sf)/gSum(mag(weights*Sf));
592  return true;
593  }
594  case operationType::areaNormalIntegrate:
595  {
596  result = gSum(weights*values & Sf);
597  return true;
598  }
599  default:
600  {
601  return false;
602  }
603  }
604 }
605 
606 
608 (
609  const Field<vector>& values,
610  const scalarField& signs,
611  const scalarField& weights,
612  const vectorField& Sf,
613  vector& result
614 ) const
615 {
616  switch (operation_)
617  {
618  case operationType::sumDirection:
619  {
620  const vector n = normalised(dict_.lookup<vector>("direction"));
621  const scalarField nv(n & values);
622  result = gSum(weights*pos0(nv)*n*(nv));
623  return true;
624  }
625  case operationType::sumDirectionBalance:
626  {
627  const vector n = normalised(dict_.lookup<vector>("direction"));
628  const scalarField nv(n & values);
629  result = gSum(weights*pos0(nv)*n*(nv));
630  return true;
631  }
632  default:
633  {
634  // Fall through to same-type operations
635  return processValuesTypeType
636  (
637  values,
638  signs,
639  weights,
640  Sf,
641  result
642  );
643  }
644  }
645 }
646 
647 
648 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
649 
651 (
652  const word& name,
653  const Time& runTime,
654  const dictionary& dict
655 )
656 :
657  fieldValue(name, runTime, dict, typeName),
658  surfaceWriterPtr_(nullptr),
659  regionType_(regionTypeNames_.read(dict.lookup("regionType"))),
660  operation_(operationTypeNames_.read(dict.lookup("operation"))),
661  weightFieldNames_(),
662  scaleFactor_(1),
663  writeArea_(dict.lookupOrDefault("writeArea", false)),
664  nFaces_(0),
665  faceId_(),
666  facePatchId_(),
667  faceSign_()
668 {
669  read(dict);
670 }
671 
673 (
674  const word& name,
675  const objectRegistry& obr,
676  const dictionary& dict
677 )
678 :
679  fieldValue(name, obr, dict, typeName),
680  surfaceWriterPtr_(nullptr),
681  regionType_(regionTypeNames_.read(dict.lookup("regionType"))),
682  operation_(operationTypeNames_.read(dict.lookup("operation"))),
683  weightFieldNames_(),
684  scaleFactor_(1),
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 (regionType_ == regionTypes::sampledSurface)
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 (writeFields_)
744  {
745  faceList faces;
747 
748  if (regionType_ == regionTypes::sampledSurface)
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 the sign and weight fields and the surface normals
770  const scalarField signs
771  (
772  regionType_ == regionTypes::sampledSurface
773  ? scalarField(surfacePtr_().Sf().size(), 1)
774  : List<scalar>(faceSign_)
775  );
776  scalarField weights(signs.size(), 1);
777  forAll(weightFieldNames_, i)
778  {
779  weights *= getFieldValues<scalar>(weightFieldNames_[i]);
780  }
781  const vectorField Sf
782  (
783  regionType_ == regionTypes::sampledSurface
784  ? surfacePtr_().Sf()
785  : (signs*filterField(mesh_.Sf()))()
786  );
787 
788  // Process the fields
789  forAll(fields_, i)
790  {
791  const word& fieldName = fields_[i];
792  bool ok = false;
793 
794  #define writeValuesFieldType(fieldType, none) \
795  ok = \
796  ok \
797  || writeValues<fieldType> \
798  ( \
799  fieldName, \
800  signs, \
801  weights, \
802  Sf \
803  );
805  #undef writeValuesFieldType
806 
807  if (!ok)
808  {
810  << "Requested field " << fieldName
811  << " not found in database and not processed"
812  << endl;
813  }
814  }
815 
816  if (operation_ != operationType::none && Pstream::master())
817  {
818  file() << endl;
819  }
820 
821  Log << endl;
822 
823  return true;
824 }
825 
826 
827 // ************************************************************************* //
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:643
#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:156
addToRunTimeSelectionTable(functionObject, fieldValueDelta, dictionary)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
An abstract class for surfaces with sampling.
static const NamedEnum< operationType, 16 > operationTypeNames_
Operation type names.
static autoPtr< surfaceWriter > New(const word &writeType, const IOstream::streamFormat writeFormat)
Return a reference to the selected surfaceWriter.
Definition: surfaceWriter.C:45
const polyBoundaryMesh & boundaryMesh() const
Return boundaryMesh reference.
Definition: polyPatch.C:291
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
#define FOR_ALL_FIELD_TYPES(Macro,...)
Definition: fieldTypes.H:50
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.
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:982
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
#define writeValuesFieldType(fieldType, none)
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)
Form normalised(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:413
const pointField & points
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
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.
bool processValues(const Field< Type > &values, const scalarField &signs, const scalarField &weights, const vectorField &Sf, ResultType &result) const
Apply the operation to the values, and return true if successful.
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.
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:309
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:844
virtual bool write()
Write.
Definition: fieldValue.C:114
label whichFace(const label l) const
Return label of face in patch from global face label.
Definition: polyPatch.H:389