polyMesh.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017, 2020 OpenFOAM Foundation
9  Copyright (C) 2016-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "polyMesh.H"
30 #include "Time.H"
31 #include "cellIOList.H"
32 #include "wedgePolyPatch.H"
33 #include "emptyPolyPatch.H"
34 #include "globalMeshData.H"
35 #include "processorPolyPatch.H"
37 #include "indexedOctree.H"
38 #include "treeDataCell.H"
39 #include "MeshObject.H"
40 #include "pointMesh.H"
41 
42 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
43 
44 namespace Foam
45 {
46  defineTypeNameAndDebug(polyMesh, 0);
47 }
48 
50 
52 
53 
54 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
55 
56 void Foam::polyMesh::calcDirections() const
57 {
58  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
59  {
60  solutionD_[cmpt] = 1;
61  }
62 
63  // Knock out empty and wedge directions. Note:they will be present on all
64  // domains.
65 
66  bool hasEmptyPatches = false;
67  bool hasWedgePatches = false;
68 
69  vector emptyDirVec = Zero;
70  vector wedgeDirVec = Zero;
71 
72  forAll(boundaryMesh(), patchi)
73  {
74  const polyPatch& pp = boundaryMesh()[patchi];
75  if (isA<emptyPolyPatch>(pp))
76  {
77  // Force calculation of geometric properties, independent of
78  // size. This avoids parallel synchronisation problems.
79  const vectorField::subField fa(pp.faceAreas());
80 
81  if (pp.size())
82  {
83  hasEmptyPatches = true;
84  emptyDirVec += sum(cmptMag(fa));
85  }
86  }
87  else if (isA<wedgePolyPatch>(pp))
88  {
89  const wedgePolyPatch& wpp = refCast<const wedgePolyPatch>(pp);
90 
91  // Force calculation of geometric properties, independent of
92  // size. This avoids parallel synchronisation problems.
93  (void)wpp.faceNormals();
94 
95  if (pp.size())
96  {
97  hasWedgePatches = true;
98  wedgeDirVec += cmptMag(wpp.centreNormal());
99  }
100  }
101  }
102 
103 
104  if (returnReduceOr(hasEmptyPatches))
105  {
106  reduce(emptyDirVec, sumOp<vector>());
107 
108  emptyDirVec.normalise();
109 
110  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
111  {
112  if (emptyDirVec[cmpt] > 1e-6)
113  {
114  solutionD_[cmpt] = -1;
115  }
116  else
117  {
118  solutionD_[cmpt] = 1;
119  }
120  }
121  }
122 
123 
124  // Knock out wedge directions
125 
126  geometricD_ = solutionD_;
127 
128  if (returnReduceOr(hasWedgePatches))
129  {
130  reduce(wedgeDirVec, sumOp<vector>());
131 
132  wedgeDirVec.normalise();
133 
134  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
135  {
136  if (wedgeDirVec[cmpt] > 1e-6)
137  {
138  geometricD_[cmpt] = -1;
139  }
140  else
141  {
142  geometricD_[cmpt] = 1;
143  }
144  }
145  }
146 }
147 
148 
149 Foam::autoPtr<Foam::labelIOList> Foam::polyMesh::readTetBasePtIs() const
150 {
151  IOobject io
152  (
153  "tetBasePtIs",
154  instance(),
155  meshSubDir,
156  *this,
160  );
161 
162  if (io.typeHeaderOk<labelIOList>(true))
163  {
165  }
167  return nullptr;
168 }
169 
170 
171 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
172 
173 Foam::polyMesh::polyMesh(const IOobject& io, const bool doInit)
174 :
176  primitiveMesh(),
177  data_(static_cast<const objectRegistry&>(*this)),
178  points_
179  (
180  IOobject
181  (
182  "points",
183  time().findInstance(meshDir(), "points"),
184  meshSubDir,
185  *this,
186  IOobject::MUST_READ,
187  IOobject::NO_WRITE
188  )
189  ),
190  faces_
191  (
192  IOobject
193  (
194  "faces",
195  time().findInstance(meshDir(), "faces"),
196  meshSubDir,
197  *this,
198  IOobject::MUST_READ,
199  IOobject::NO_WRITE
200  )
201  ),
202  owner_
203  (
204  IOobject
205  (
206  "owner",
207  faces_.instance(),
208  meshSubDir,
209  *this,
210  IOobject::READ_IF_PRESENT,
211  IOobject::NO_WRITE
212  )
213  ),
214  neighbour_
215  (
216  IOobject
217  (
218  "neighbour",
219  faces_.instance(),
220  meshSubDir,
221  *this,
222  IOobject::READ_IF_PRESENT,
223  IOobject::NO_WRITE
224  )
225  ),
226  clearedPrimitives_(false),
227  boundary_
228  (
229  IOobject
230  (
231  "boundary",
232  time().findInstance // allow 'newer' boundary file
233  (
234  meshDir(),
235  "boundary",
236  IOobject::MUST_READ,
237  faces_.instance()
238  ),
239  meshSubDir,
240  *this,
241  IOobject::MUST_READ,
242  IOobject::NO_WRITE
243  ),
244  *this
245  ),
246  bounds_(points_),
247  comm_(UPstream::worldComm),
248  geometricD_(Zero),
249  solutionD_(Zero),
250  tetBasePtIsPtr_(nullptr),
251  cellTreePtr_(nullptr),
252  pointZones_
253  (
254  IOobject
255  (
256  "pointZones",
257  faces_.instance(),
258  meshSubDir,
259  *this,
260  IOobject::READ_IF_PRESENT,
261  IOobject::NO_WRITE
262  ),
263  *this,
264  PtrList<entry>()
265  ),
266  faceZones_
267  (
268  IOobject
269  (
270  "faceZones",
271  faces_.instance(),
272  meshSubDir,
273  *this,
274  IOobject::READ_IF_PRESENT,
275  IOobject::NO_WRITE
276  ),
277  *this,
278  PtrList<entry>()
279  ),
280  cellZones_
281  (
282  IOobject
283  (
284  "cellZones",
285  faces_.instance(),
286  meshSubDir,
287  *this,
288  IOobject::READ_IF_PRESENT,
289  IOobject::NO_WRITE
290  ),
291  *this,
292  PtrList<entry>()
293  ),
294  globalMeshDataPtr_(nullptr),
295  moving_(false),
296  topoChanging_(false),
297  storeOldCellCentres_(false),
298  curMotionTimeIndex_(time().timeIndex()),
299  oldPointsPtr_(nullptr),
300  oldCellCentresPtr_(nullptr)
301 {
302  if (owner_.hasHeaderClass())
303  {
304  initMesh();
305  }
306  else
307  {
308  cellCompactIOList cLst
309  (
310  IOobject
311  (
312  "cells",
313  time().findInstance(meshDir(), "cells"),
314  meshSubDir,
315  *this,
318  )
319  );
320 
321  // Set the primitive mesh
322  initMesh(cLst);
323 
324  owner_.write();
325  neighbour_.write();
326  }
327 
328  if (returnReduceOr(boundary_.empty()))
329  {
331  << "Missing mesh boundary on one or more domains" << endl;
332 
333  // Warn if global empty mesh
334  if (returnReduceAnd(!nPoints()))
335  {
337  << "No points in mesh" << endl;
338  }
339  if (returnReduceAnd(!nCells()))
340  {
342  << "No cells in mesh" << endl;
343  }
344  }
345 
346  if (doInit)
347  {
348  polyMesh::init(false); // do not init lower levels
349  }
350 }
351 
352 
353 bool Foam::polyMesh::init(const bool doInit)
354 {
355  if (doInit)
356  {
357  primitiveMesh::init(doInit);
358  }
359 
360  // Calculate topology for the patches (processor-processor comms etc.)
361  boundary_.updateMesh();
362 
363  // Calculate the geometry for the patches (transformation tensors etc.)
364  boundary_.calcGeometry();
365 
366  // Initialise demand-driven data
367  calcDirections();
368 
369  return false;
370 }
371 
372 
373 Foam::polyMesh::polyMesh
374 (
375  const IOobject& io,
376  pointField&& points,
377  faceList&& faces,
378  labelList&& owner,
379  labelList&& neighbour,
380  const bool syncPar
381 )
382 :
384  primitiveMesh(),
385  data_(static_cast<const objectRegistry&>(*this)),
386  points_
387  (
388  IOobject
389  (
390  "points",
391  instance(),
392  meshSubDir,
393  *this,
394  IOobject::NO_READ, //io.readOpt(),
395  io.writeOpt()
396  ),
397  std::move(points)
398  ),
399  faces_
400  (
401  IOobject
402  (
403  "faces",
404  instance(),
405  meshSubDir,
406  *this,
407  IOobject::NO_READ, //io.readOpt(),
408  io.writeOpt()
409  ),
410  std::move(faces)
411  ),
412  owner_
413  (
414  IOobject
415  (
416  "owner",
417  instance(),
418  meshSubDir,
419  *this,
420  IOobject::NO_READ, //io.readOpt(),
421  io.writeOpt()
422  ),
423  std::move(owner)
424  ),
425  neighbour_
426  (
427  IOobject
428  (
429  "neighbour",
430  instance(),
431  meshSubDir,
432  *this,
433  IOobject::NO_READ, //io.readOpt(),
434  io.writeOpt()
435  ),
436  std::move(neighbour)
437  ),
438  clearedPrimitives_(false),
439  boundary_
440  (
441  IOobject
442  (
443  "boundary",
444  instance(),
445  meshSubDir,
446  *this,
447  IOobject::NO_READ, // ignore since no alternative can be supplied
448  io.writeOpt()
449  ),
450  *this,
451  polyPatchList()
452  ),
453  bounds_(points_, syncPar),
454  comm_(UPstream::worldComm),
455  geometricD_(Zero),
456  solutionD_(Zero),
457  tetBasePtIsPtr_(nullptr),
458  cellTreePtr_(nullptr),
459  pointZones_
460  (
461  IOobject
462  (
463  "pointZones",
464  instance(),
465  meshSubDir,
466  *this,
467  IOobject::NO_READ, // ignore since no alternative can be supplied
468  IOobject::NO_WRITE
469  ),
470  *this,
471  Foam::zero{}
472  ),
473  faceZones_
474  (
475  IOobject
476  (
477  "faceZones",
478  instance(),
479  meshSubDir,
480  *this,
481  IOobject::NO_READ, // ignore since no alternative can be supplied
483  ),
484  *this,
485  Foam::zero{}
486  ),
487  cellZones_
488  (
489  IOobject
490  (
491  "cellZones",
492  instance(),
493  meshSubDir,
494  *this,
495  IOobject::NO_READ, // ignore since no alternative can be supplied
497  ),
498  *this,
499  Foam::zero{}
500  ),
501  globalMeshDataPtr_(nullptr),
502  moving_(false),
503  topoChanging_(false),
504  storeOldCellCentres_(false),
505  curMotionTimeIndex_(time().timeIndex()),
506  oldPointsPtr_(nullptr),
507  oldCellCentresPtr_(nullptr)
508 {
509  // Check if the faces and cells are valid
510  forAll(faces_, facei)
511  {
512  const face& curFace = faces_[facei];
513 
514  if (min(curFace) < 0 || max(curFace) > points_.size())
515  {
517  << "Face " << facei << "contains vertex labels out of range: "
518  << curFace << " Max point index = " << points_.size()
519  << abort(FatalError);
520  }
521  }
523  // Set the primitive mesh
524  initMesh();
525 }
526 
527 
528 Foam::polyMesh::polyMesh
529 (
530  const IOobject& io,
531  pointField&& points,
532  faceList&& faces,
533  cellList&& cells,
534  const bool syncPar
535 )
536 :
538  primitiveMesh(),
539  data_(static_cast<const objectRegistry&>(*this)),
540  points_
541  (
542  IOobject
543  (
544  "points",
545  instance(),
546  meshSubDir,
547  *this,
548  IOobject::NO_READ,
549  io.writeOpt()
550  ),
551  std::move(points)
552  ),
553  faces_
554  (
555  IOobject
556  (
557  "faces",
558  instance(),
559  meshSubDir,
560  *this,
561  IOobject::NO_READ,
562  io.writeOpt()
563  ),
564  std::move(faces)
565  ),
566  owner_
567  (
568  IOobject
569  (
570  "owner",
571  instance(),
572  meshSubDir,
573  *this,
574  IOobject::NO_READ,
575  io.writeOpt()
576  ),
577  Foam::zero{}
578  ),
579  neighbour_
580  (
581  IOobject
582  (
583  "neighbour",
584  instance(),
585  meshSubDir,
586  *this,
588  io.writeOpt()
589  ),
590  Foam::zero{}
591  ),
592  clearedPrimitives_(false),
593  boundary_
594  (
595  IOobject
596  (
597  "boundary",
598  instance(),
599  meshSubDir,
600  *this,
602  io.writeOpt()
603  ),
604  *this,
605  Foam::zero{}
606  ),
607  bounds_(points_, syncPar),
608  comm_(UPstream::worldComm),
609  geometricD_(Zero),
610  solutionD_(Zero),
611  tetBasePtIsPtr_(nullptr),
612  cellTreePtr_(nullptr),
613  pointZones_
614  (
615  IOobject
616  (
617  "pointZones",
618  instance(),
619  meshSubDir,
620  *this,
623  ),
624  *this,
625  Foam::zero{}
626  ),
627  faceZones_
628  (
629  IOobject
630  (
631  "faceZones",
632  instance(),
633  meshSubDir,
634  *this,
637  ),
638  *this,
639  Foam::zero{}
640  ),
641  cellZones_
642  (
643  IOobject
644  (
645  "cellZones",
646  instance(),
647  meshSubDir,
648  *this,
651  ),
652  *this,
653  Foam::zero{}
654  ),
655  globalMeshDataPtr_(nullptr),
656  moving_(false),
657  topoChanging_(false),
658  storeOldCellCentres_(false),
659  curMotionTimeIndex_(time().timeIndex()),
660  oldPointsPtr_(nullptr),
661  oldCellCentresPtr_(nullptr)
662 {
663  // Check if faces are valid
664  forAll(faces_, facei)
665  {
666  const face& curFace = faces_[facei];
667 
668  if (min(curFace) < 0 || max(curFace) > points_.size())
669  {
671  << "Face " << facei << "contains vertex labels out of range: "
672  << curFace << " Max point index = " << points_.size()
673  << abort(FatalError);
674  }
675  }
676 
677  // Transfer in cell list
678  cellList cLst(std::move(cells));
679 
680  // Check if cells are valid
681  forAll(cLst, celli)
682  {
683  const cell& curCell = cLst[celli];
684 
685  if (min(curCell) < 0 || max(curCell) > faces_.size())
686  {
688  << "Cell " << celli << "contains face labels out of range: "
689  << curCell << " Max face index = " << faces_.size()
690  << abort(FatalError);
691  }
692  }
694  // Set the primitive mesh
695  initMesh(cLst);
696 }
697 
698 
699 Foam::polyMesh::polyMesh
700 (
701  const IOobject& io,
702  const Foam::zero,
703  const bool syncPar
704 )
705 :
706  polyMesh(io, pointField(), faceList(), labelList(), labelList(), syncPar)
707 {}
708 
709 
711 (
713  autoPtr<faceList>&& faces,
714  autoPtr<labelList>&& owner,
715  autoPtr<labelList>&& neighbour,
716  const labelUList& patchSizes,
717  const labelUList& patchStarts,
718  const bool validBoundary
719 )
720 {
721  // Clear addressing. Keep geometric props and updateable props for mapping.
722  clearAddressing(true);
723 
724  // Take over new primitive data.
725  // Optimized to avoid overwriting data at all
726  if (points)
727  {
728  points_.transfer(*points);
729  bounds_ = boundBox(points_, validBoundary);
730  }
731 
732  if (faces)
733  {
734  faces_.transfer(*faces);
735  }
736 
737  if (owner)
738  {
739  owner_.transfer(*owner);
740  }
741 
742  if (neighbour)
743  {
744  neighbour_.transfer(*neighbour);
745  }
746 
747 
748  // Reset patch sizes and starts
749  forAll(boundary_, patchi)
750  {
751  boundary_[patchi] = polyPatch
752  (
753  boundary_[patchi],
754  boundary_,
755  patchi,
756  patchSizes[patchi],
757  patchStarts[patchi]
758  );
759  }
760 
761 
762  // Flags the mesh files as being changed
763  setInstance(time().timeName());
764 
765  // Check if the faces and cells are valid
766  forAll(faces_, facei)
767  {
768  const face& curFace = faces_[facei];
769 
770  if (min(curFace) < 0 || max(curFace) > points_.size())
771  {
773  << "Face " << facei << " contains vertex labels out of range: "
774  << curFace << " Max point index = " << points_.size()
775  << abort(FatalError);
776  }
777  }
778 
779 
780  // Set the primitive mesh from the owner_, neighbour_.
781  // Works out from patch end where the active faces stop.
782  initMesh();
783 
784 
785  if (validBoundary)
786  {
787  // Note that we assume that all the patches stay the same and are
788  // correct etc. so we can already use the patches to do
789  // processor-processor comms.
790 
791  // Calculate topology for the patches (processor-processor comms etc.)
792  boundary_.updateMesh();
793 
794  // Calculate the geometry for the patches (transformation tensors etc.)
795  boundary_.calcGeometry();
796 
797  // Warn if global empty mesh
798  if (returnReduceAnd(!nPoints()) || returnReduceAnd(!nCells()))
799  {
801  << "No points or no cells in mesh" << endl;
802  }
803  }
804 }
805 
806 
807 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
808 
810 {
811  clearOut();
812  resetMotion();
813 }
814 
815 
816 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
818 const Foam::word& Foam::polyMesh::regionName(const word& region)
819 {
820  return (region == polyMesh::defaultRegion ? word::null : region);
821 }
822 
823 
825 {
826  if (region.empty() || region == polyMesh::defaultRegion)
827  {
828  return polyMesh::meshSubDir;
829  }
831  return (region / polyMesh::meshSubDir);
832 }
833 
834 
835 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
836 
838 {
840  {
841  return parent().dbDir();
842  }
843 
844  return objectRegistry::dbDir();
845 }
846 
849 {
850  return dbDir()/meshSubDir;
851 }
852 
855 {
857 }
858 
861 {
862  return points_.instance();
863 }
864 
867 {
868  return faces_.instance();
869 }
870 
871 
873 {
874  if (geometricD_.x() == 0)
875  {
876  calcDirections();
877  }
878 
879  return geometricD_;
880 }
881 
883 Foam::label Foam::polyMesh::nGeometricD() const
884 {
885  return cmptSum(geometricD() + Vector<label>::one)/2;
886 }
887 
888 
890 {
891  if (solutionD_.x() == 0)
892  {
893  calcDirections();
894  }
895 
896  return solutionD_;
897 }
898 
900 Foam::label Foam::polyMesh::nSolutionD() const
901 {
902  return cmptSum(solutionD() + Vector<label>::one)/2;
903 }
904 
905 
907 {
908  if (!tetBasePtIsPtr_)
909  {
910  if (debug)
911  {
913  << "Forcing storage of base points."
914  << endl;
915  }
916 
917  labelList basePts
918  (
920  );
921 
922  tetBasePtIsPtr_.reset
923  (
924  new labelIOList
925  (
926  IOobject
927  (
928  "tetBasePtIs",
929  instance(),
930  meshSubDir,
931  *this,
935  ),
936  std::move(basePts)
937  )
938  );
939  }
940 
941  return *tetBasePtIsPtr_;
942 }
943 
944 
947 {
948  if (!cellTreePtr_)
949  {
950  Random rndGen(261782);
951 
952  treeBoundBox overallBb(points());
953  overallBb.inflate(rndGen, 1e-4, ROOTVSMALL);
954 
955  cellTreePtr_.reset
956  (
957  new indexedOctree<treeDataCell>
958  (
959  treeDataCell
960  (
961  false, // not cache bb
962  *this,
963  CELL_TETS // use tet-decomposition for any inside test
964  ),
965  overallBb,
966  8, // maxLevel
967  10, // leafsize
968  5.0 // duplicity
969  )
970  );
971  }
972 
973  return *cellTreePtr_;
974 }
975 
976 
978 (
979  polyPatchList& plist,
980  const bool validBoundary
981 )
982 {
983  if (boundaryMesh().size())
984  {
986  << "boundary already exists"
987  << abort(FatalError);
988  }
989 
990  // Reset valid directions
991  geometricD_ = Zero;
992  solutionD_ = Zero;
993 
994  boundary_.transfer(plist);
995 
996  // parallelData depends on the processorPatch ordering so force
997  // recalculation. Problem: should really be done in removeBoundary but
998  // there is some info in parallelData which might be interesting inbetween
999  // removeBoundary and addPatches.
1000  globalMeshDataPtr_.reset(nullptr);
1001 
1002  if (validBoundary)
1003  {
1004  // Calculate topology for the patches (processor-processor comms etc.)
1005  boundary_.updateMesh();
1006 
1007  // Calculate the geometry for the patches (transformation tensors etc.)
1008  boundary_.calcGeometry();
1010  boundary_.checkDefinition();
1011  }
1012 }
1013 
1014 
1016 (
1017  PtrList<pointZone>&& pz,
1018  PtrList<faceZone>&& fz,
1019  PtrList<cellZone>&& cz
1020 )
1021 {
1022  if (pointZones_.size() || faceZones_.size() || cellZones_.size())
1023  {
1025  << "point, face or cell zone already exists"
1026  << abort(FatalError);
1027  }
1028 
1029  // Point zones - take ownership of the pointers
1030  if (pz.size())
1031  {
1032  pointZones_.clear();
1033  pointZones_.transfer(pz);
1034  pointZones_.writeOpt(IOobject::AUTO_WRITE);
1035  }
1036 
1037  // Face zones - take ownership of the pointers
1038  if (fz.size())
1039  {
1040  faceZones_.clear();
1041  faceZones_.transfer(fz);
1042  faceZones_.writeOpt(IOobject::AUTO_WRITE);
1043  }
1044 
1045  // Cell zones - take ownership of the pointers
1046  if (cz.size())
1047  {
1048  cellZones_.clear();
1049  cellZones_.transfer(cz);
1050  cellZones_.writeOpt(IOobject::AUTO_WRITE);
1051  }
1052 }
1053 
1054 
1056 (
1057  const List<polyPatch*>& p,
1058  const bool validBoundary
1059 )
1060 {
1061  // Acquire ownership of the pointers
1062  polyPatchList plist(const_cast<List<polyPatch*>&>(p));
1063 
1064  addPatches(plist, validBoundary);
1065 }
1066 
1067 
1069 (
1070  const List<pointZone*>& pz,
1071  const List<faceZone*>& fz,
1072  const List<cellZone*>& cz
1073 )
1074 {
1075  // Acquire ownership of the pointers
1076  addZones
1077  (
1079  PtrList<faceZone>(const_cast<List<faceZone*>&>(fz)),
1080  PtrList<cellZone>(const_cast<List<cellZone*>&>(cz))
1081  );
1082 }
1083 
1084 
1086 {
1087  if (clearedPrimitives_)
1088  {
1090  << "points deallocated"
1092  }
1093 
1094  return points_;
1095 }
1096 
1099 {
1100  return io.upToDate(points_);
1101 }
1102 
1105 {
1106  io.eventNo() = points_.eventNo()+1;
1107 }
1108 
1109 
1111 {
1112  if (clearedPrimitives_)
1113  {
1115  << "faces deallocated"
1117  }
1118 
1119  return faces_;
1120 }
1121 
1124 {
1125  return owner_;
1126 }
1127 
1130 {
1131  return neighbour_;
1132 }
1133 
1134 
1136 {
1137  if (!moving_)
1138  {
1139  return points_;
1140  }
1141 
1142  if (!oldPointsPtr_)
1143  {
1144  if (debug)
1145  {
1147  }
1148 
1149  oldPointsPtr_.reset(new pointField(points_));
1150  curMotionTimeIndex_ = time().timeIndex();
1151  }
1152 
1153  return *oldPointsPtr_;
1154 }
1155 
1156 
1158 {
1159  storeOldCellCentres_ = true;
1160 
1161  if (!moving_)
1162  {
1163  return cellCentres();
1164  }
1165 
1166  if (!oldCellCentresPtr_)
1167  {
1168  oldCellCentresPtr_.reset(new pointField(cellCentres()));
1169  }
1170 
1171  return *oldCellCentresPtr_;
1172 }
1173 
1174 
1175 void Foam::polyMesh::movePoints(const pointField& newPoints)
1176 {
1178  << "Moving points for time " << time().value()
1179  << " index " << time().timeIndex() << endl;
1180 
1181  if (newPoints.size() != points_.size())
1182  {
1184  << "Size of newPoints " << newPoints.size()
1185  << " does not correspond to current mesh points size "
1186  << points_.size()
1187  << exit(FatalError);
1188  }
1189 
1190 
1191  moving(true);
1192 
1193  // Pick up old points
1194  if (curMotionTimeIndex_ != time().timeIndex())
1195  {
1196  if (debug)
1197  {
1198  Info<< "tmp<scalarField> polyMesh::movePoints(const pointField&) : "
1199  << " Storing current points for time " << time().value()
1200  << " index " << time().timeIndex() << endl;
1201  }
1202 
1203  if (storeOldCellCentres_)
1204  {
1205  oldCellCentresPtr_.clear();
1206  oldCellCentresPtr_.reset(new pointField(cellCentres()));
1207  }
1208 
1209  // Mesh motion in the new time step
1210  oldPointsPtr_.clear();
1211  oldPointsPtr_.reset(new pointField(points_));
1212  curMotionTimeIndex_ = time().timeIndex();
1213  }
1214 
1215  points_ = newPoints;
1216 
1217  bool moveError = false;
1218  if (debug)
1219  {
1220  // Check mesh motion
1221  if (checkMeshMotion(points_, true))
1222  {
1223  moveError = true;
1224 
1226  << "Moving the mesh with given points will "
1227  << "invalidate the mesh." << nl
1228  << "Mesh motion should not be executed." << endl;
1229  }
1230  }
1231 
1232  points_.writeOpt(IOobject::AUTO_WRITE);
1233  points_.instance() = time().timeName();
1234  points_.eventNo() = getEvent();
1235 
1236  if (tetBasePtIsPtr_)
1237  {
1238  tetBasePtIsPtr_->writeOpt(IOobject::AUTO_WRITE);
1239  tetBasePtIsPtr_->instance() = time().timeName();
1240  tetBasePtIsPtr_->eventNo() = getEvent();
1241  }
1242 
1243  // Currently a no-op; earlier versions set meshPhi and call
1244  // primitiveMesh::clearGeom
1245  (void)primitiveMesh::movePoints(points_, oldPoints());
1246 
1247  // Update the mesh geometry (via fvGeometryScheme)
1248  // - updateGeom is virtual -> calls fvMesh::updateGeom (or higher)
1249  // - fvMesh::updateGeom defers to surfaceInterpolation::updateGeom(),
1250  // which defers to fvGeometryScheme::movePoints()
1251  // - set the mesh flux
1252  // - clear out/recalculate stale geometry
1253  updateGeom();
1254 
1255  // Adjust parallel shared points
1256  if (globalMeshDataPtr_)
1257  {
1258  globalMeshDataPtr_->movePoints(points_);
1259  }
1260 
1261  // Force recalculation of all geometric data with new points
1262 
1263  bounds_ = boundBox(points_);
1264  boundary_.movePoints(points_);
1265 
1266  pointZones_.movePoints(points_);
1267  faceZones_.movePoints(points_);
1268  cellZones_.movePoints(points_);
1269 
1270  // Reset cell tree - it gets built from mesh geometry so might have
1271  // wrong boxes. It is correct as long as none of the cells leaves
1272  // the boxes it is in which most likely is almost never the case except
1273  // for tiny displacements. An alternative is to check the displacements
1274  // to see if they are tiny - imagine a big windtunnel with a small rotating
1275  // object. In this case the processors without the rotating object wouldn't
1276  // have to clear any geometry. However your critical path still stays the
1277  // same so no time would be gained (unless the decomposition gets weighted).
1278  // Small benefit for lots of scope for problems so not done.
1279  cellTreePtr_.clear();
1280 
1281  // Reset valid directions (could change with rotation)
1282  geometricD_ = Zero;
1283  solutionD_ = Zero;
1284 
1285  // Note: tet-base decomposition does not get cleared. Ideally your face
1286  // decomposition should not change during mesh motion ...
1287 
1288 
1289  meshObject::movePoints<polyMesh>(*this);
1290  meshObject::movePoints<pointMesh>(*this);
1291 
1292  const_cast<Time&>(time()).functionObjects().movePoints(*this);
1293 
1294 
1295  if (debug && moveError)
1296  {
1297  // Write mesh to ease debugging. Note we want to avoid calling
1298  // e.g. fvMesh::write since meshPhi not yet complete.
1299  polyMesh::write();
1300  }
1301 }
1302 
1303 
1304 void Foam::polyMesh::resetMotion() const
1306  curMotionTimeIndex_ = 0;
1307  oldPointsPtr_.clear();
1308  oldCellCentresPtr_.clear();
1309 }
1310 
1313 {
1314  return bool(globalMeshDataPtr_);
1315 }
1316 
1317 
1319 {
1320  if (!globalMeshDataPtr_)
1321  {
1322  if (debug)
1323  {
1324  Pout<< "polyMesh::globalData() const : "
1325  << "Constructing parallelData from processor topology"
1326  << endl;
1327  }
1328  // Construct globalMeshData using processorPatch information only.
1329  globalMeshDataPtr_.reset(new globalMeshData(*this));
1330  }
1331 
1332  return *globalMeshDataPtr_;
1333 }
1334 
1335 
1336 void Foam::polyMesh::removeFiles(const fileName& instanceDir) const
1337 {
1338  fileName meshFilesPath = thisDb().time().path()/instanceDir/meshDir();
1339 
1340  rm(meshFilesPath/"points");
1341  rm(meshFilesPath/"faces");
1342  rm(meshFilesPath/"owner");
1343  rm(meshFilesPath/"neighbour");
1344  rm(meshFilesPath/"cells");
1345  rm(meshFilesPath/"boundary");
1346  rm(meshFilesPath/"pointZones");
1347  rm(meshFilesPath/"faceZones");
1348  rm(meshFilesPath/"cellZones");
1349  rm(meshFilesPath/"meshModifiers");
1350  rm(meshFilesPath/"parallelData");
1351 
1352  // remove subdirectories
1353  if (isDir(meshFilesPath/"sets"))
1354  {
1355  rmDir(meshFilesPath/"sets");
1356  }
1357 }
1358 
1359 
1361 {
1362  removeFiles(instance());
1363 }
1364 
1365 
1367 (
1368  const point& p,
1369  label& celli,
1370  label& tetFacei,
1371  label& tetPti
1372 ) const
1373 {
1374  celli = -1;
1375  tetFacei = -1;
1376  tetPti = -1;
1377 
1378  const indexedOctree<treeDataCell>& tree = cellTree();
1379 
1380  // Find point inside cell
1381  celli = tree.findInside(p);
1382 
1383  if (celli != -1)
1384  {
1385  // Check the nearest cell to see if the point is inside.
1386  findTetFacePt(celli, p, tetFacei, tetPti);
1387  }
1388 }
1389 
1390 
1392 (
1393  const label celli,
1394  const point& p,
1395  label& tetFacei,
1396  label& tetPti
1397 ) const
1398 {
1399  const polyMesh& mesh = *this;
1400 
1402  tetFacei = tet.face();
1403  tetPti = tet.tetPt();
1404 }
1405 
1406 
1408 (
1409  const point& p,
1410  label celli,
1411  const cellDecomposition decompMode
1412 ) const
1413 {
1414  switch (decompMode)
1415  {
1416  case FACE_PLANES:
1417  {
1418  return primitiveMesh::pointInCell(p, celli);
1419  }
1420  break;
1421 
1422  case FACE_CENTRE_TRIS:
1423  {
1424  // only test that point is on inside of plane defined by cell face
1425  // triangles
1426  const cell& cFaces = cells()[celli];
1427 
1428  forAll(cFaces, cFacei)
1429  {
1430  label facei = cFaces[cFacei];
1431  const face& f = faces_[facei];
1432  const point& fc = faceCentres()[facei];
1433  bool isOwn = (owner_[facei] == celli);
1434 
1435  forAll(f, fp)
1436  {
1437  label pointi;
1438  label nextPointi;
1439 
1440  if (isOwn)
1441  {
1442  pointi = f[fp];
1443  nextPointi = f.nextLabel(fp);
1444  }
1445  else
1446  {
1447  pointi = f.nextLabel(fp);
1448  nextPointi = f[fp];
1449  }
1450 
1451  triPointRef faceTri
1452  (
1453  points()[pointi],
1454  points()[nextPointi],
1455  fc
1456  );
1457 
1458  vector proj = p - faceTri.centre();
1459 
1460  if ((faceTri.areaNormal() & proj) > 0)
1461  {
1462  return false;
1463  }
1464  }
1465  }
1466  return true;
1467  }
1468  break;
1469 
1470  case FACE_DIAG_TRIS:
1471  {
1472  // only test that point is on inside of plane defined by cell face
1473  // triangles
1474  const cell& cFaces = cells()[celli];
1475 
1476  forAll(cFaces, cFacei)
1477  {
1478  label facei = cFaces[cFacei];
1479  const face& f = faces_[facei];
1480 
1481  for (label tetPti = 1; tetPti < f.size() - 1; tetPti++)
1482  {
1483  // Get tetIndices of face triangle
1484  tetIndices faceTetIs(celli, facei, tetPti);
1485 
1486  triPointRef faceTri = faceTetIs.faceTri(*this);
1487 
1488  vector proj = p - faceTri.centre();
1489 
1490  if ((faceTri.areaNormal() & proj) > 0)
1491  {
1492  return false;
1493  }
1494  }
1495  }
1496 
1497  return true;
1498  }
1499  break;
1500 
1501  case CELL_TETS:
1502  {
1503  label tetFacei;
1504  label tetPti;
1505 
1506  findTetFacePt(celli, p, tetFacei, tetPti);
1507 
1508  return tetFacei != -1;
1509  }
1510  break;
1511  }
1512 
1513  return false;
1514 }
1515 
1516 
1517 Foam::label Foam::polyMesh::findCell
1518 (
1519  const point& p,
1520  const cellDecomposition decompMode
1521 ) const
1522 {
1523  if
1524  (
1525  Pstream::parRun()
1526  && (decompMode == FACE_DIAG_TRIS || decompMode == CELL_TETS)
1527  )
1528  {
1529  // Force construction of face-diagonal decomposition before testing
1530  // for zero cells.
1531  //
1532  // If parallel running a local domain might have zero cells so never
1533  // construct the face-diagonal decomposition which uses parallel
1534  // transfers.
1535  (void)tetBasePtIs();
1536  }
1537 
1538  if (nCells() == 0)
1539  {
1540  return -1;
1541  }
1542 
1543  if (decompMode == CELL_TETS)
1544  {
1545  // Advanced search method utilizing an octree
1546  // and tet-decomposition of the cells
1547 
1548  label celli;
1549  label tetFacei;
1550  label tetPti;
1551 
1552  findCellFacePt(p, celli, tetFacei, tetPti);
1553 
1554  return celli;
1555  }
1556  else
1557  {
1558  // Approximate search avoiding the construction of an octree
1559  // and cell decomposition
1560 
1561  if (Pstream::parRun() && decompMode == FACE_DIAG_TRIS)
1562  {
1563  // Force construction of face-diagonal decomposition before testing
1564  // for zero cells. If parallel running a local domain might have
1565  // zero cells so never construct the face-diagonal decomposition
1566  // (which uses parallel transfers)
1567  (void)tetBasePtIs();
1568  }
1569 
1570  // Find the nearest cell centre to this location
1571  label celli = findNearestCell(p);
1572 
1573  // If point is in the nearest cell return
1574  if (pointInCell(p, celli, decompMode))
1575  {
1576  return celli;
1577  }
1578  else
1579  {
1580  // Point is not in the nearest cell so search all cells
1581 
1582  for (label celli = 0; celli < nCells(); celli++)
1583  {
1584  if (pointInCell(p, celli, decompMode))
1585  {
1586  return celli;
1587  }
1588  }
1589 
1590  return -1;
1591  }
1592  }
1593 }
1594 
1595 
1596 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
1597 
1599 (
1600  IOstreamOption streamOpt,
1601  const bool writeOnProc
1602 ) const
1603 {
1604  // Currently no special treatment. Just write the objects
1605 
1606  return objectRegistry::writeObject(streamOpt, writeOnProc);
1607 }
1608 
1609 
1610 // ************************************************************************* //
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
Definition: polyMesh.C:346
uint8_t direction
Definition: direction.H:46
A class for handling file names.
Definition: fileName.H:72
static labelList findFaceBasePts(const polyMesh &mesh, scalar tol=minTetQuality, bool report=false)
Find a suitable base point for each face for decomposition into tets.
List< cell > cellList
List of cell.
Definition: cellListFwd.H:39
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
virtual void movePoints(const pointField &)
Move points.
Definition: polyMesh.C:1168
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:859
label nPoints() const noexcept
Number of mesh points.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
Cmpt cmptSum(const SphericalTensor< Cmpt > &st)
Return the sum of components of a SphericalTensor.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:899
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:847
virtual void setUpToDatePoints(regIOobject &io) const
Set io to be up-to-date with points.
Definition: polyMesh.C:1097
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
void resetMotion() const
Reset motion.
Definition: polyMesh.C:1297
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:882
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Random rndGen
Definition: createFields.H:23
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:410
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:1049
void resetPrimitives(autoPtr< pointField > &&points, autoPtr< faceList > &&faces, autoPtr< labelList > &&owner, autoPtr< labelList > &&neighbour, const labelUList &patchSizes, const labelUList &patchStarts, const bool validBoundary=true)
Reset mesh primitive data. Assumes all patch info correct.
Definition: polyMesh.C:704
A bounding box defined in terms of min/max extrema points.
Definition: boundBox.H:63
bool pointInCell(const point &p, label celli) const
Return true if the point is in the cell.
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
virtual const fileName & dbDir() const
Override the objectRegistry dbDir for a single-region case.
Definition: polyMesh.C:830
Ignore writing from objectRegistry::writeObject()
static label worldComm
Communicator for all ranks. May differ from commGlobal() if local worlds are in use.
Definition: UPstream.H:409
bool pointInCell(const point &p, label celli, const cellDecomposition=CELL_TETS) const
Test if point p is in the celli.
Definition: polyMesh.C:1401
void movePoints(const pointField &p, const pointField &oldP)
Move points.
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:876
bool hasHeaderClass() const noexcept
True if headerClassName() is non-empty (after reading)
Definition: IOobjectI.H:251
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1078
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc) const
Write the objects using stream options.
static tetIndices findTet(const polyMesh &mesh, label cI, const point &pt)
Find the tet decomposition of the cell containing the given point.
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:865
bool isDir(const fileName &name, const bool followLink=true)
Does the name exist as a DIRECTORY in the file system?
Definition: POSIX.C:860
PtrList< polyPatch > polyPatchList
Store lists of polyPatch as a PtrList.
Definition: polyPatch.H:56
word timeName
Definition: getTimeIndex.H:3
writeOption writeOpt() const noexcept
Get the write option.
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:853
IOList< label > labelIOList
IO for a List of label.
Definition: labelIOList.H:32
dynamicFvMesh & mesh
void removeFiles() const
Remove all files from mesh instance()
Definition: polyMesh.C:1353
const cellShapeList & cells
const pointField & points
bool hasGlobalData() const noexcept
Is demand-driven parallel info available?
Definition: polyMesh.C:1305
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
A class for handling words, derived from Foam::string.
Definition: word.H:63
label nPoints
const Time & time() const noexcept
Return time registry.
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition: polyMesh.C:1360
#define DebugInFunction
Report an information message using Foam::Info.
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:405
Tree tree(triangles.begin(), triangles.end())
void addZones(PtrList< pointZone > &&pz, PtrList< faceZone > &&fz, PtrList< cellZone > &&cz)
Add mesh zones.
Definition: polyMesh.C:1009
virtual const pointField & oldPoints() const
Return old points (mesh motion)
Definition: polyMesh.C:1128
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static const word null
An empty word.
Definition: word.H:84
const globalMeshData & globalData() const
Return parallel info (demand-driven)
Definition: polyMesh.C:1311
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:78
Vector< scalar > vector
Definition: vector.H:57
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1103
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
errorManip< error > abort(error &err)
Definition: errorManip.H:139
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info.
bool rmDir(const fileName &directory, const bool silent=false, const bool emptyOnly=false)
Remove a directory and its contents recursively,.
Definition: POSIX.C:1433
Templated 3D Vector derived from VectorSpace adding construction from 3 components, element access using x(), y() and z() member functions and the inner-product (dot-product) and cross-product operators.
Definition: Vector.H:58
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
const direction noexcept
Definition: Scalar.H:258
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
CompactIOList< cell, label > cellCompactIOList
Compact IO for a List of cell.
Definition: cellIOList.H:35
labelList f(nPoints)
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:893
virtual bool upToDatePoints(const regIOobject &io) const
Return true if io is up-to-date with points.
Definition: polyMesh.C:1091
bool empty() const noexcept
True if the list is empty (ie, size() is zero)
Definition: UPtrListI.H:99
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition: polyMesh.C:939
virtual const pointField & oldCellCentres() const
Return old cellCentres (mesh motion)
Definition: polyMesh.C:1150
vector point
Point is a vector.
Definition: point.H:37
#define WarningInFunction
Report a warning using Foam::Warning.
label nCells() const noexcept
Number of mesh cells.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:841
Nothing to be read.
Automatically write from objectRegistry::writeObject()
virtual bool init(const bool doInit)
Initialise all non-demand-driven data.
void findTetFacePt(const label celli, const point &p, label &tetFacei, label &tetPti) const
Find the tetFacei and tetPti for point p in celli.
Definition: polyMesh.C:1385
triangle< point, const point & > triPointRef
A triangle using referred points.
Definition: triangleFwd.H:39
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
Reduce inplace (cf. MPI Allreduce) using specified communication schedule.
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1511
virtual ~polyMesh()
Destructor.
Definition: polyMesh.C:802
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:66
messageStream Info
Information stream (stdout output on master, null elsewhere)
A class representing the concept of 0 (zero) that can be used to avoid manipulating objects known to ...
Definition: zero.H:57
SubField< vector > subField
Declare type of subField.
Definition: Field.H:128
virtual const fileName & dbDir() const
Local directory path of this objectRegistry relative to the time.
Reading is optional [identical to READ_IF_PRESENT].
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
Registry of regIOobjects.
static autoPtr< T > New(Args &&... args)
Construct autoPtr with forwarding arguments.
Definition: autoPtr.H:178
void addPatches(polyPatchList &plist, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:971
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Inter-processor communications stream.
Definition: UPstream.H:60
Do not request registration (bool: false)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
label timeIndex
Definition: getTimeIndex.H:24
A keyword and a list of tokens is an &#39;entry&#39;.
Definition: entry.H:63
bool rm(const fileName &file)
Remove a file (or its gz equivalent), returning true if successful.
Definition: POSIX.C:1404
#define InfoInFunction
Report an information message using Foam::Info.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
virtual bool writeObject(IOstreamOption streamOpt, const bool writeOnProc=true) const
Write items held in the objectRegistry. Normally includes mesh components (points, faces, etc) and any registered fields.
Definition: polyMesh.C:1592