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-2022 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 "processorFvPatch.H"
28 #include "processorCyclicFvPatch.H"
29 #include "sampledSurface.H"
30 #include "mergePoints.H"
31 #include "indirectPrimitivePatch.H"
32 #include "PatchTools.H"
33 #include "fvMeshStitcher.H"
34 #include "polyTopoChangeMap.H"
35 #include "polyMeshMap.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace functionObjects
43 {
44 namespace fieldValues
45 {
46  defineTypeNameAndDebug(surfaceFieldValue, 0);
47  addToRunTimeSelectionTable(fieldValue, surfaceFieldValue, dictionary);
48  addToRunTimeSelectionTable(functionObject, surfaceFieldValue, dictionary);
49 }
50 }
51 }
52 
53 template<>
54 const char* Foam::NamedEnum
55 <
57  3
58 >::names[] =
59 {
60  "faceZone",
61  "patch",
62  "sampledSurface"
63 };
64 
65 template<>
66 const char* Foam::NamedEnum
67 <
69  16
70 >::names[] =
71 {
72  "none",
73  "sum",
74  "sumMag",
75  "sumDirection",
76  "sumDirectionBalance",
77  "orientedSum",
78  "average",
79  "areaAverage",
80  "areaIntegrate",
81  "min",
82  "max",
83  "minMag",
84  "maxMag",
85  "CoV",
86  "areaNormalAverage",
87  "areaNormalIntegrate"
88 };
89 
90 const Foam::NamedEnum
91 <
93  3
95 
96 const Foam::NamedEnum
97 <
99  16
101 
102 
103 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
104 
105 void Foam::functionObjects::fieldValues::surfaceFieldValue::setFaceZoneFaces()
106 {
107  label zoneId = mesh_.faceZones().findZoneID(regionName_);
108 
109  if (zoneId < 0)
110  {
112  << type() << " " << name() << ": "
113  << regionTypeNames_[regionType_] << "(" << regionName_ << "):" << nl
114  << " Unknown face zone name: " << regionName_
115  << ". Valid face zones are: " << mesh_.faceZones().names()
116  << nl << exit(FatalError);
117  }
118 
119  // Ensure addressing is built on all processes
120  mesh_.polyBFacePatches();
121  mesh_.polyBFacePatchFaces();
122 
123  const faceZone& zone = mesh_.faceZones()[zoneId];
124 
125  DynamicList<label> faceIds(zone.size());
126  DynamicList<label> facePatchIds(zone.size());
127  DynamicList<label> faceSigns(zone.size());
128 
129  forAll(zone, zoneFacei)
130  {
131  const label facei = zone[zoneFacei];
132  const label faceSign = zone.flipMap()[zoneFacei] ? -1 : 1;
133 
134  if (mesh_.isInternalFace(facei))
135  {
136  faceIds.append(facei);
137  facePatchIds.append(-1);
138  faceSigns.append(faceSign);
139  }
140  else
141  {
142  const label bFacei = facei - mesh_.nInternalFaces();
143 
144  const labelUList patches = mesh_.polyBFacePatches()[bFacei];
145  const labelUList patchFaces = mesh_.polyBFacePatchFaces()[bFacei];
146 
147  forAll(patches, i)
148  {
149  // Don't include processor patch faces twice
150  const fvPatch& fvp = mesh_.boundary()[patches[i]];
151  if
152  (
153  isType<processorFvPatch>(fvp)
154  && refCast<const processorFvPatch>(fvp).neighbour()
155  )
156  {
157  continue;
158  }
159 
160  faceIds.append(patchFaces[i]);
161  facePatchIds.append(patches[i]);
162  faceSigns.append(faceSign);
163  }
164  }
165  }
166 
167  faceId_.transfer(faceIds);
168  facePatchId_.transfer(facePatchIds);
169  faceSign_.transfer(faceSigns);
170 
171  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
172 }
173 
174 
175 void Foam::functionObjects::fieldValues::surfaceFieldValue::setPatchFaces()
176 {
177  const label patchId = mesh_.boundaryMesh().findPatchID(regionName_);
178 
179  if (patchId < 0)
180  {
182  << type() << " " << name() << ": "
183  << regionTypeNames_[regionType_] << "(" << regionName_ << "):" << nl
184  << " Unknown patch name: " << regionName_
185  << ". Valid patch names are: "
186  << mesh_.boundaryMesh().names() << nl
187  << exit(FatalError);
188  }
189 
190  const fvPatch& fvp = mesh_.boundary()[patchId];
191 
192  faceId_ = identity(fvp.size());
193  facePatchId_ = labelList(fvp.size(), patchId);
194  faceSign_ = labelList(fvp.size(), 1);
195 
196  // If we have selected a cyclic, also include any associated processor
197  // cyclic faces
198  forAll(mesh_.boundary(), patchi)
199  {
200  const fvPatch& fvp = mesh_.boundary()[patchi];
201 
202  if
203  (
204  isA<processorCyclicFvPatch>(fvp)
205  && refCast<const processorCyclicFvPatch>(fvp).referPatchID() == patchId
206  )
207  {
208  faceId_.append(identity(fvp.size()));
209  facePatchId_.append(labelList(fvp.size(), patchi));
210  faceSign_.append(labelList(fvp.size(), 1));
211  }
212  }
213 
214  nFaces_ = returnReduce(faceId_.size(), sumOp<label>());
215 }
216 
217 
218 void Foam::functionObjects::fieldValues::surfaceFieldValue::sampledSurfaceFaces
219 (
220  const dictionary& dict
221 )
222 {
223  surfacePtr_ = sampledSurface::New
224  (
225  name(),
226  mesh_,
227  dict.subDict("sampledSurfaceDict")
228  );
229  surfacePtr_().update();
230  nFaces_ = returnReduce(surfacePtr_().faces().size(), sumOp<label>());
231 }
232 
233 
234 void Foam::functionObjects::fieldValues::surfaceFieldValue::combineMeshGeometry
235 (
236  faceList& faces,
237  pointField& points
238 ) const
239 {
240  List<faceList> allFaces(Pstream::nProcs());
242 
243  labelList globalFacesIs(faceId_);
244  forAll(globalFacesIs, i)
245  {
246  if (facePatchId_[i] != -1)
247  {
248  label patchi = facePatchId_[i];
249  globalFacesIs[i] += mesh_.boundaryMesh()[patchi].start();
250  }
251  }
252 
253  // Add local faces and points to the all* lists
255  (
256  IndirectList<face>(mesh_.faces(), globalFacesIs),
257  mesh_.points()
258  );
259  allFaces[Pstream::myProcNo()] = pp.localFaces();
260  allPoints[Pstream::myProcNo()] = pp.localPoints();
261 
262  Pstream::gatherList(allFaces);
263  Pstream::gatherList(allPoints);
264 
265  // Renumber and flatten
266  label nFaces = 0;
267  label nPoints = 0;
268  forAll(allFaces, proci)
269  {
270  nFaces += allFaces[proci].size();
271  nPoints += allPoints[proci].size();
272  }
273 
274  faces.setSize(nFaces);
275  points.setSize(nPoints);
276 
277  nFaces = 0;
278  nPoints = 0;
279 
280  // My own data first
281  {
282  const faceList& fcs = allFaces[Pstream::myProcNo()];
283  forAll(fcs, i)
284  {
285  const face& f = fcs[i];
286  face& newF = faces[nFaces++];
287  newF.setSize(f.size());
288  forAll(f, fp)
289  {
290  newF[fp] = f[fp] + nPoints;
291  }
292  }
293 
294  const pointField& pts = allPoints[Pstream::myProcNo()];
295  forAll(pts, i)
296  {
297  points[nPoints++] = pts[i];
298  }
299  }
300 
301  // Other proc data follows
302  forAll(allFaces, proci)
303  {
304  if (proci != Pstream::myProcNo())
305  {
306  const faceList& fcs = allFaces[proci];
307  forAll(fcs, i)
308  {
309  const face& f = fcs[i];
310  face& newF = faces[nFaces++];
311  newF.setSize(f.size());
312  forAll(f, fp)
313  {
314  newF[fp] = f[fp] + nPoints;
315  }
316  }
317 
318  const pointField& pts = allPoints[proci];
319  forAll(pts, i)
320  {
321  points[nPoints++] = pts[i];
322  }
323  }
324  }
325 
326  // Merge
327  labelList oldToNew;
328  pointField newPoints;
329  bool hasMerged = mergePoints
330  (
331  points,
332  small,
333  false,
334  oldToNew,
335  newPoints
336  );
337 
338  if (hasMerged)
339  {
340  if (debug)
341  {
342  Pout<< "Merged from " << points.size()
343  << " down to " << newPoints.size() << " points" << endl;
344  }
345 
346  points.transfer(newPoints);
347  forAll(faces, i)
348  {
349  inplaceRenumber(oldToNew, faces[i]);
350  }
351  }
352 }
353 
354 
355 void Foam::functionObjects::fieldValues::surfaceFieldValue::
356 combineSurfaceGeometry
357 (
358  faceList& faces,
359  pointField& points
360 ) const
361 {
362  if (regionType_ == regionTypes::sampledSurface)
363  {
364  const sampledSurface& s = surfacePtr_();
365 
366  if (Pstream::parRun())
367  {
368  // Dimension as fraction of mesh bounding box
369  scalar mergeDim = 1e-10*mesh_.bounds().mag();
370 
371  labelList pointsMap;
372 
374  (
375  mergeDim,
377  (
378  SubList<face>(s.faces(), s.faces().size()),
379  s.points()
380  ),
381  points,
382  faces,
383  pointsMap
384  );
385  }
386  else
387  {
388  faces = s.faces();
389  points = s.points();
390  }
391  }
392 }
393 
394 
395 Foam::scalar
396 Foam::functionObjects::fieldValues::surfaceFieldValue::totalArea() const
397 {
398  if (regionType_ == regionTypes::sampledSurface)
399  {
400  return gSum(surfacePtr_().magSf());
401  }
402  else
403  {
404  return gSum(filterField(mesh_.magSf()));
405  }
406 }
407 
408 
409 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
410 
412 (
413  const dictionary& dict
414 )
415 {
416  switch (regionType_)
417  {
418  case regionTypes::faceZone:
419  {
420  dict.lookup("name") >> regionName_;
421  setFaceZoneFaces();
422  break;
423  }
424  case regionTypes::patch:
425  {
426  dict.lookup("name") >> regionName_;
427  setPatchFaces();
428  break;
429  }
430  case regionTypes::sampledSurface:
431  {
432  sampledSurfaceFaces(dict);
433  regionName_ = surfacePtr_().name();
434  break;
435  }
436  default:
437  {
439  << type() << " " << name() << ": "
440  << regionTypeNames_[regionType_] << "(" << regionName_ << "):"
441  << nl << " Unknown region type. Valid region types are:"
442  << regionTypeNames_.sortedToc() << nl << exit(FatalError);
443  }
444  }
445 
446  if (nFaces_ == 0 && (!mesh_.stitcher().stitches() || !mesh_.conformal()))
447  {
449  << type() << " " << name() << ": "
450  << regionTypeNames_[regionType_] << "(" << regionName_ << "):" << nl
451  << " Region has no faces" << exit(FatalError);
452  }
453 
454  if (regionType_ == regionTypes::sampledSurface)
455  {
456  surfacePtr_().update();
457  }
458 
459  totalArea_ = totalArea();
460 
461  Info<< type() << " " << name() << ":" << nl
462  << " total faces = " << nFaces_
463  << nl
464  << " total area = " << totalArea_
465  << nl;
466 
467  if (dict.readIfPresent("weightFields", weightFieldNames_))
468  {
469  Info<< name() << " " << operationTypeNames_[operation_]
470  << " weight fields " << weightFieldNames_;
471  }
472  else if (dict.found("weightField"))
473  {
474  weightFieldNames_.setSize(1);
475  dict.lookup("weightField") >> weightFieldNames_[0];
476 
477  Info<< name() << " " << operationTypeNames_[operation_]
478  << " weight field " << weightFieldNames_[0];
479  }
480 
481  if (dict.readIfPresent("scaleFactor", scaleFactor_))
482  {
483  Info<< " scale factor = " << scaleFactor_ << nl;
484  }
485 
486  Info<< nl << endl;
487 
488  if (writeFields_)
489  {
490  const word surfaceFormat(dict.lookup("surfaceFormat"));
491 
492  surfaceWriterPtr_.reset
493  (
494  surfaceWriter::New(surfaceFormat, dict).ptr()
495  );
496  }
497 }
498 
499 
501 (
502  const label i
503 )
504 {
505  if (operation_ != operationType::none)
506  {
507  writeCommented(file(), "Region type : ");
508  file() << regionTypeNames_[regionType_] << " " << regionName_ << endl;
509  writeCommented(file(), "Faces : ");
510  file() << nFaces_ << endl;
511  writeCommented(file(), "Area : ");
512  file() << totalArea_ << endl;
513 
514  writeCommented(file(), "Time");
515  if (writeArea_)
516  {
517  file() << tab << "Area";
518  }
519 
520  forAll(fields_, fieldi)
521  {
522  file()
523  << tab << operationTypeNames_[operation_]
524  << "(" << fields_[fieldi] << ")";
525  }
526 
527  file() << endl;
528  }
529 }
530 
531 
533 (
534  const Field<scalar>& values,
535  const scalarField& signs,
536  const scalarField& weights,
537  const vectorField& Sf,
538  scalar& result
539 ) const
540 {
541  switch (operation_)
542  {
543  case operationType::sumDirection:
544  {
545  const vector n(dict_.lookup("direction"));
546  result = gSum(weights*pos0(values*(Sf & n))*mag(values));
547  return true;
548  }
549  case operationType::sumDirectionBalance:
550  {
551  const vector n(dict_.lookup("direction"));
552  const scalarField nv(values*(Sf & n));
553  result = gSum(weights*(pos0(nv)*mag(values) - neg(nv)*mag(values)));
554  return true;
555  }
556  default:
557  {
558  // Fall through to same-type operations
559  return processValuesTypeType
560  (
561  values,
562  signs,
563  weights,
564  Sf,
565  result
566  );
567  }
568  }
569 }
570 
571 
573 (
574  const Field<vector>& values,
575  const scalarField& signs,
576  const scalarField& weights,
577  const vectorField& Sf,
578  scalar& result
579 ) const
580 {
581  switch (operation_)
582  {
583  case operationType::areaNormalAverage:
584  {
585  result = gSum(weights*values & Sf)/gSum(mag(weights*Sf));
586  return true;
587  }
588  case operationType::areaNormalIntegrate:
589  {
590  result = gSum(weights*values & Sf);
591  return true;
592  }
593  default:
594  {
595  return false;
596  }
597  }
598 }
599 
600 
602 (
603  const Field<vector>& values,
604  const scalarField& signs,
605  const scalarField& weights,
606  const vectorField& Sf,
607  vector& result
608 ) const
609 {
610  switch (operation_)
611  {
612  case operationType::sumDirection:
613  {
614  const vector n = normalised(dict_.lookup<vector>("direction"));
615  const scalarField nv(n & values);
616  result = gSum(weights*pos0(nv)*n*(nv));
617  return true;
618  }
619  case operationType::sumDirectionBalance:
620  {
621  const vector n = normalised(dict_.lookup<vector>("direction"));
622  const scalarField nv(n & values);
623  result = gSum(weights*pos0(nv)*n*(nv));
624  return true;
625  }
626  default:
627  {
628  // Fall through to same-type operations
629  return processValuesTypeType
630  (
631  values,
632  signs,
633  weights,
634  Sf,
635  result
636  );
637  }
638  }
639 }
640 
641 
642 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
643 
645 (
646  const word& name,
647  const Time& runTime,
648  const dictionary& dict
649 )
650 :
651  fieldValue(name, runTime, dict, typeName),
652  dict_(dict),
653  surfaceWriterPtr_(nullptr),
654  regionType_(regionTypeNames_.read(dict.lookup("regionType"))),
655  operation_(operationTypeNames_.read(dict.lookup("operation"))),
656  weightFieldNames_(),
657  scaleFactor_(1),
658  writeArea_(dict.lookupOrDefault("writeArea", false)),
659  nFaces_(0),
660  faceId_(),
661  facePatchId_(),
662  faceSign_()
663 {
664  read(dict_);
665 }
666 
668 (
669  const word& name,
670  const objectRegistry& obr,
671  const dictionary& dict
672 )
673 :
674  fieldValue(name, obr, dict, typeName),
675  dict_(dict),
676  surfaceWriterPtr_(nullptr),
677  regionType_(regionTypeNames_.read(dict.lookup("regionType"))),
678  operation_(operationTypeNames_.read(dict.lookup("operation"))),
679  weightFieldNames_(),
680  scaleFactor_(1),
681  writeArea_(dict.lookupOrDefault("writeArea", false)),
682  nFaces_(0),
683  faceId_(),
684  facePatchId_(),
685  faceSign_()
686 {
687  read(dict_);
688 }
689 
690 
691 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
692 
694 {}
695 
696 
697 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
698 
700 (
701  const dictionary& dict
702 )
703 {
704  fieldValue::read(dict);
705  initialise(dict);
706 
707  return true;
708 }
709 
710 
712 {
713  if (operation_ != operationType::none)
714  {
716  }
717 
718  if (regionType_ == regionTypes::sampledSurface)
719  {
720  surfacePtr_().update();
721  }
722 
723  if (operation_ != operationType::none && Pstream::master())
724  {
725  writeTime(file());
726  }
727 
728  if (writeArea_)
729  {
730  totalArea_ = totalArea();
731  if (operation_ != operationType::none && Pstream::master())
732  {
733  file() << tab << totalArea_;
734  }
735  Log << " total area = " << totalArea_ << endl;
736  }
737 
738  // Write the surface geometry
739  if (writeFields_)
740  {
741  faceList faces;
743 
744  if (regionType_ == regionTypes::sampledSurface)
745  {
746  combineSurfaceGeometry(faces, points);
747  }
748  else
749  {
750  combineMeshGeometry(faces, points);
751  }
752 
753  if (Pstream::master())
754  {
755  surfaceWriterPtr_->write
756  (
757  outputDir(),
758  regionTypeNames_[regionType_] + ("_" + regionName_),
759  points,
760  faces
761  );
762  }
763  }
764 
765  // Construct the sign and weight fields and the surface normals
766  const scalarField signs
767  (
768  regionType_ == regionTypes::sampledSurface
769  ? scalarField(surfacePtr_().Sf().size(), 1)
770  : List<scalar>(faceSign_)
771  );
772  scalarField weights(signs.size(), 1);
773  forAll(weightFieldNames_, i)
774  {
775  weights *= getFieldValues<scalar>(weightFieldNames_[i]);
776  }
777  const vectorField Sf
778  (
779  regionType_ == regionTypes::sampledSurface
780  ? surfacePtr_().Sf()
781  : (signs*filterField(mesh_.Sf()))()
782  );
783 
784  // Process the fields
785  forAll(fields_, i)
786  {
787  const word& fieldName = fields_[i];
788  bool ok = false;
789 
790  #define writeValuesFieldType(fieldType, none) \
791  ok = \
792  ok \
793  || writeValues<fieldType> \
794  ( \
795  fieldName, \
796  signs, \
797  weights, \
798  Sf \
799  );
801  #undef writeValuesFieldType
802 
803  if (!ok)
804  {
806  << "Requested field " << fieldName
807  << " not found in database and not processed"
808  << endl;
809  }
810  }
811 
812  if (operation_ != operationType::none && Pstream::master())
813  {
814  file() << endl;
815  }
816 
817  Log << endl;
818 
819  return true;
820 }
821 
822 
824 (
825  const polyMesh& mesh
826 )
827 {
828  if (&mesh == &mesh_)
829  {
830  // It may be necessary to reset if the mesh moves. The total area might
831  // change, as might non-conformal faces.
832  initialise(dict_);
833  }
834 }
835 
836 
838 (
839  const polyTopoChangeMap& map
840 )
841 {
842  if (&map.mesh() == &mesh_)
843  {
844  initialise(dict_);
845  }
846 }
847 
848 
850 (
851  const polyMeshMap& map
852 )
853 {
854  if (&map.mesh() == &mesh_)
855  {
856  initialise(dict_);
857  }
858 }
859 
860 
861 // ************************************************************************* //
const fvPatchList & patches
label patchId(-1)
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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:306
An abstract class for surfaces with sampling.
static const NamedEnum< operationType, 16 > operationTypeNames_
Operation type names.
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:249
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
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:63
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
labelList identity(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
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:69
virtual void movePoints(const polyMesh &)
Update for mesh point-motion.
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:1002
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 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects...
Definition: DynamicList.H:56
Base class for zones.
Definition: zone.H:57
A class for handling words, derived from string.
Definition: word.H:59
label nPoints
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void append(const T &)
Append an element at the end of the list.
Definition: ListI.H:178
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.
List< label > labelList
A List of labels.
Definition: labelList.H:56
static autoPtr< surfaceWriter > New(const word &writeType, const IOstream::streamFormat writeFormat, const IOstream::compressionType writeCompression)
Select given write options.
Definition: surfaceWriter.C:75
static autoPtr< sampledSurface > New(const word &name, const polyMesh &, const dictionary &)
Return a reference to the selected surface.
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
const polyMesh & mesh() const
Return polyMesh.
virtual label size() const
Return size.
Definition: fvPatch.H:157
dimensionedScalar pos0(const dimensionedScalar &ds)
virtual void mapMesh(const polyMeshMap &)
Update from another mesh using the given map.
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
virtual void topoChange(const polyTopoChangeMap &)
Update topology using the given map.
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
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
FOR_ALL_FIELD_TYPES(DefineFvWallInfoType)
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
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:65
static const NamedEnum< regionTypes, 3 > regionTypeNames_
region type names
Registry of regIOobjects.
Class containing mesh-to-mesh mapping information.
Definition: polyMeshMap.H:50
const polyMesh & mesh() const
Return polyMesh.
Definition: polyMeshMap.H:72
static void gatherList(const List< commsStruct > &comms, List< T > &Values, const int tag, const label comm)
Gather data but keep individual values separate.
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:864
virtual bool write()
Write.
Definition: fieldValue.C:114