polyMesh.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 "polyMesh.H"
27 #include "Time.H"
28 #include "cellIOList.H"
29 #include "wedgePolyPatch.H"
30 #include "emptyPolyPatch.H"
31 #include "globalMeshData.H"
32 #include "processorPolyPatch.H"
34 #include "indexedOctree.H"
35 #include "treeDataCell.H"
36 #include "MeshObject.H"
37 #include "pointMesh.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
43  defineTypeNameAndDebug(polyMesh, 0);
44 
45  word polyMesh::defaultRegion = "region0";
46  word polyMesh::meshSubDir = "polyMesh";
47 }
48 
49 
50 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51 
52 void Foam::polyMesh::calcDirections() const
53 {
54  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
55  {
56  solutionD_[cmpt] = 1;
57  }
58 
59  // Knock out empty and wedge directions. Note:they will be present on all
60  // domains.
61 
62  label nEmptyPatches = 0;
63  label nWedgePatches = 0;
64 
65  vector emptyDirVec = Zero;
66  vector wedgeDirVec = Zero;
67 
68  forAll(boundaryMesh(), patchi)
69  {
70  if (boundaryMesh()[patchi].size())
71  {
72  if (isA<emptyPolyPatch>(boundaryMesh()[patchi]))
73  {
74  nEmptyPatches++;
75  emptyDirVec += sum(cmptMag(boundaryMesh()[patchi].faceAreas()));
76  }
77  else if (isA<wedgePolyPatch>(boundaryMesh()[patchi]))
78  {
79  const wedgePolyPatch& wpp = refCast<const wedgePolyPatch>
80  (
81  boundaryMesh()[patchi]
82  );
83 
84  nWedgePatches++;
85  wedgeDirVec += cmptMag(wpp.centreNormal());
86  }
87  }
88  }
89 
90  reduce(nEmptyPatches, maxOp<label>());
91  reduce(nWedgePatches, maxOp<label>());
92 
93  if (nEmptyPatches)
94  {
95  reduce(emptyDirVec, sumOp<vector>());
96 
97  emptyDirVec /= mag(emptyDirVec);
98 
99  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
100  {
101  if (emptyDirVec[cmpt] > 1e-6)
102  {
103  solutionD_[cmpt] = -1;
104  }
105  else
106  {
107  solutionD_[cmpt] = 1;
108  }
109  }
110  }
111 
112 
113  // Knock out wedge directions
114 
115  geometricD_ = solutionD_;
116 
117  if (nWedgePatches)
118  {
119  reduce(wedgeDirVec, sumOp<vector>());
120 
121  wedgeDirVec /= mag(wedgeDirVec);
122 
123  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
124  {
125  if (wedgeDirVec[cmpt] > 1e-6)
126  {
127  geometricD_[cmpt] = -1;
128  }
129  else
130  {
131  geometricD_[cmpt] = 1;
132  }
133  }
134  }
135 }
136 
137 
138 Foam::autoPtr<Foam::labelIOList> Foam::polyMesh::readTetBasePtIs() const
139 {
140  typeIOobject<labelIOList> io
141  (
142  "tetBasePtIs",
143  instance(),
144  meshSubDir,
145  *this,
148  );
149 
150  if (io.headerOk())
151  {
152  return autoPtr<labelIOList>(new labelIOList(io));
153  }
154  else
155  {
156  return autoPtr<labelIOList>(nullptr);
157  }
158 }
159 
160 
161 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
162 
164 :
165  objectRegistry(io),
166  primitiveMesh(),
167  points_
168  (
169  IOobject
170  (
171  "points",
172  time().findInstance(meshDir(), "points"),
173  meshSubDir,
174  *this,
175  IOobject::MUST_READ,
176  IOobject::NO_WRITE
177  )
178  ),
179  faces_
180  (
181  IOobject
182  (
183  "faces",
184  time().findInstance(meshDir(), "faces"),
185  meshSubDir,
186  *this,
187  IOobject::MUST_READ,
188  IOobject::NO_WRITE
189  )
190  ),
191  owner_
192  (
193  IOobject
194  (
195  "owner",
196  faces_.instance(),
197  meshSubDir,
198  *this,
199  IOobject::READ_IF_PRESENT,
200  IOobject::NO_WRITE
201  )
202  ),
203  neighbour_
204  (
205  IOobject
206  (
207  "neighbour",
208  faces_.instance(),
209  meshSubDir,
210  *this,
211  IOobject::READ_IF_PRESENT,
212  IOobject::NO_WRITE
213  )
214  ),
215  clearedPrimitives_(false),
216  boundary_
217  (
218  IOobject
219  (
220  "boundary",
221  time().findInstance(meshDir(), "boundary"),
222  meshSubDir,
223  *this,
224  IOobject::MUST_READ,
225  IOobject::NO_WRITE
226  ),
227  *this
228  ),
229  bounds_(points_),
230  comm_(UPstream::worldComm),
231  geometricD_(Zero),
232  solutionD_(Zero),
233  tetBasePtIsPtr_(readTetBasePtIs()),
234  cellTreePtr_(nullptr),
235  pointZones_
236  (
237  IOobject
238  (
239  "pointZones",
240  time().findInstance
241  (
242  meshDir(),
243  "pointZones",
244  IOobject::READ_IF_PRESENT
245  ),
246  meshSubDir,
247  *this,
248  IOobject::READ_IF_PRESENT,
249  IOobject::NO_WRITE
250  ),
251  *this
252  ),
253  faceZones_
254  (
255  IOobject
256  (
257  "faceZones",
258  time().findInstance
259  (
260  meshDir(),
261  "faceZones",
262  IOobject::READ_IF_PRESENT
263  ),
264  meshSubDir,
265  *this,
266  IOobject::READ_IF_PRESENT,
267  IOobject::NO_WRITE
268  ),
269  *this
270  ),
271  cellZones_
272  (
273  IOobject
274  (
275  "cellZones",
276  time().findInstance
277  (
278  meshDir(),
279  "cellZones",
280  IOobject::READ_IF_PRESENT
281  ),
282  meshSubDir,
283  *this,
284  IOobject::READ_IF_PRESENT,
285  IOobject::NO_WRITE
286  ),
287  *this
288  ),
289  globalMeshDataPtr_(nullptr),
290  curMotionTimeIndex_(-1),
291  oldPointsPtr_(nullptr),
292  oldCellCentresPtr_(nullptr),
293  storeOldCellCentres_(false),
294  moving_(false),
295  topoChanged_(false)
296 {
297  if (!owner_.headerClassName().empty())
298  {
299  initMesh();
300  }
301  else
302  {
303  cellCompactIOList cLst
304  (
305  IOobject
306  (
307  "cells",
308  time().findInstance(meshDir(), "cells"),
309  meshSubDir,
310  *this,
313  )
314  );
315 
316  // Set the primitive mesh
317  initMesh(cLst);
318 
319  owner_.write();
320  neighbour_.write();
321  }
322 
323  // Calculate topology for the patches (processor-processor comms etc.)
324  boundary_.topoChange();
325 
326  // Calculate the geometry for the patches (transformation tensors etc.)
327  boundary_.calcGeometry();
328 
329  // Warn if global empty mesh
330  if (returnReduce(nPoints(), sumOp<label>()) == 0)
331  {
333  << "no points in mesh" << endl;
334  }
335  if (returnReduce(nCells(), sumOp<label>()) == 0)
336  {
338  << "no cells in mesh" << endl;
339  }
340 
341  // Initialise demand-driven data
342  calcDirections();
343 }
344 
345 
347 (
348  const IOobject& io,
349  pointField&& points,
350  faceList&& faces,
351  labelList&& owner,
352  labelList&& neighbour,
353  const bool syncPar
354 )
355 :
356  objectRegistry(io),
357  primitiveMesh(),
358  points_
359  (
360  IOobject
361  (
362  "points",
363  instance(),
364  meshSubDir,
365  *this,
366  io.readOpt(),
368  ),
369  move(points)
370  ),
371  faces_
372  (
373  IOobject
374  (
375  "faces",
376  instance(),
377  meshSubDir,
378  *this,
379  io.readOpt(),
381  ),
382  move(faces)
383  ),
384  owner_
385  (
386  IOobject
387  (
388  "owner",
389  instance(),
390  meshSubDir,
391  *this,
392  io.readOpt(),
394  ),
395  move(owner)
396  ),
397  neighbour_
398  (
399  IOobject
400  (
401  "neighbour",
402  instance(),
403  meshSubDir,
404  *this,
405  io.readOpt(),
407  ),
408  move(neighbour)
409  ),
410  clearedPrimitives_(false),
411  boundary_
412  (
413  IOobject
414  (
415  "boundary",
416  instance(),
417  meshSubDir,
418  *this,
419  io.readOpt(),
421  ),
422  *this,
423  polyPatchList()
424  ),
425  bounds_(points_, syncPar),
426  comm_(UPstream::worldComm),
427  geometricD_(Zero),
428  solutionD_(Zero),
429  tetBasePtIsPtr_(readTetBasePtIs()),
430  cellTreePtr_(nullptr),
431  pointZones_
432  (
433  IOobject
434  (
435  "pointZones",
436  instance(),
437  meshSubDir,
438  *this,
439  io.readOpt(),
441  ),
442  *this,
444  ),
445  faceZones_
446  (
447  IOobject
448  (
449  "faceZones",
450  instance(),
451  meshSubDir,
452  *this,
453  io.readOpt(),
455  ),
456  *this,
458  ),
459  cellZones_
460  (
461  IOobject
462  (
463  "cellZones",
464  instance(),
465  meshSubDir,
466  *this,
467  io.readOpt(),
469  ),
470  *this,
472  ),
473  globalMeshDataPtr_(nullptr),
474  curMotionTimeIndex_(-1),
475  oldPointsPtr_(nullptr),
476  oldCellCentresPtr_(nullptr),
477  storeOldCellCentres_(false),
478  moving_(false),
479  topoChanged_(false)
480 {
481  // Check if the faces and cells are valid
482  forAll(faces_, facei)
483  {
484  const face& curFace = faces_[facei];
485 
486  if (min(curFace) < 0 || max(curFace) > points_.size())
487  {
489  << "Face " << facei << "contains vertex labels out of range: "
490  << curFace << " Max point index = " << points_.size()
491  << abort(FatalError);
492  }
493  }
494 
495  // Set the primitive mesh
496  initMesh();
497 }
498 
499 
500 
502 (
503  const IOobject& io,
504  pointField&& points,
505  faceList&& faces,
506  cellList&& cells,
507  const bool syncPar
508 )
509 :
510  objectRegistry(io),
511  primitiveMesh(),
512  points_
513  (
514  IOobject
515  (
516  "points",
517  instance(),
518  meshSubDir,
519  *this,
522  ),
523  move(points)
524  ),
525  faces_
526  (
527  IOobject
528  (
529  "faces",
530  instance(),
531  meshSubDir,
532  *this,
535  ),
536  move(faces)
537  ),
538  owner_
539  (
540  IOobject
541  (
542  "owner",
543  instance(),
544  meshSubDir,
545  *this,
548  ),
549  0
550  ),
551  neighbour_
552  (
553  IOobject
554  (
555  "neighbour",
556  instance(),
557  meshSubDir,
558  *this,
561  ),
562  0
563  ),
564  clearedPrimitives_(false),
565  boundary_
566  (
567  IOobject
568  (
569  "boundary",
570  instance(),
571  meshSubDir,
572  *this,
575  ),
576  *this,
577  0
578  ),
579  bounds_(points_, syncPar),
580  comm_(UPstream::worldComm),
581  geometricD_(Zero),
582  solutionD_(Zero),
583  tetBasePtIsPtr_(readTetBasePtIs()),
584  cellTreePtr_(nullptr),
585  pointZones_
586  (
587  IOobject
588  (
589  "pointZones",
590  instance(),
591  meshSubDir,
592  *this,
595  ),
596  *this,
597  0
598  ),
599  faceZones_
600  (
601  IOobject
602  (
603  "faceZones",
604  instance(),
605  meshSubDir,
606  *this,
609  ),
610  *this,
611  0
612  ),
613  cellZones_
614  (
615  IOobject
616  (
617  "cellZones",
618  instance(),
619  meshSubDir,
620  *this,
623  ),
624  *this,
625  0
626  ),
627  globalMeshDataPtr_(nullptr),
628  curMotionTimeIndex_(-1),
629  oldPointsPtr_(nullptr),
630  oldCellCentresPtr_(nullptr),
631  storeOldCellCentres_(false),
632  moving_(false),
633  topoChanged_(false)
634 {
635  // Check if faces are valid
636  forAll(faces_, facei)
637  {
638  const face& curFace = faces_[facei];
639 
640  if (min(curFace) < 0 || max(curFace) > points_.size())
641  {
643  << "Face " << facei << "contains vertex labels out of range: "
644  << curFace << " Max point index = " << points_.size()
645  << abort(FatalError);
646  }
647  }
648 
649  // transfer in cell list
650  cellList cLst(move(cells));
651 
652  // Check if cells are valid
653  forAll(cLst, celli)
654  {
655  const cell& curCell = cLst[celli];
656 
657  if (min(curCell) < 0 || max(curCell) > faces_.size())
658  {
660  << "Cell " << celli << "contains face labels out of range: "
661  << curCell << " Max face index = " << faces_.size()
662  << abort(FatalError);
663  }
664  }
665 
666  // Set the primitive mesh
667  initMesh(cLst);
668 }
669 
670 
672 (
673  pointField&& points,
674  faceList&& faces,
675  labelList&& owner,
676  labelList&& neighbour,
677  const labelList& patchSizes,
678  const labelList& patchStarts,
679  const bool validBoundary
680 )
681 {
682  // Clear addressing. Keep geometric props and updateable props for mapping.
683  clearAddressing(true);
684 
685  // Take over new primitive data.
686  // Optimised to avoid overwriting data at all
687  if (notNull(points))
688  {
689  points_ = move(points);
690  bounds_ = boundBox(points_, validBoundary);
691  }
692 
693  if (notNull(faces))
694  {
695  faces_ = move(faces);
696  }
697 
698  if (notNull(owner))
699  {
700  owner_ = move(owner);
701  }
702 
703  if (notNull(neighbour))
704  {
705  neighbour_ = move(neighbour);
706  }
707 
708 
709  // Reset patch sizes and starts
710  forAll(boundary_, patchi)
711  {
712  boundary_[patchi] = polyPatch
713  (
714  boundary_[patchi],
715  boundary_,
716  patchi,
717  patchSizes[patchi],
718  patchStarts[patchi]
719  );
720  }
721 
722 
723  // Flags the mesh files as being changed
725 
726  // Check if the faces and cells are valid
727  forAll(faces_, facei)
728  {
729  const face& curFace = faces_[facei];
730 
731  if (min(curFace) < 0 || max(curFace) > points_.size())
732  {
734  << "Face " << facei << " contains vertex labels out of range: "
735  << curFace << " Max point index = " << points_.size()
736  << abort(FatalError);
737  }
738  }
739 
740 
741  // Set the primitive mesh from the owner_, neighbour_.
742  // Works out from patch end where the active faces stop.
743  initMesh();
744 
745 
746  if (validBoundary)
747  {
748  // Note that we assume that all the patches stay the same and are
749  // correct etc. so we can already use the patches to do
750  // processor-processor comms.
751 
752  // Calculate topology for the patches (processor-processor comms etc.)
753  boundary_.topoChange();
754 
755  // Calculate the geometry for the patches (transformation tensors etc.)
756  boundary_.calcGeometry();
757 
758  // Warn if global empty mesh
759  if
760  (
761  (returnReduce(nPoints(), sumOp<label>()) == 0)
762  || (returnReduce(nCells(), sumOp<label>()) == 0)
763  )
764  {
766  << "no points or no cells in mesh" << endl;
767  }
768  }
769 }
770 
771 
772 void Foam::polyMesh::reset(const polyMesh& newMesh)
773 {
774  // Clear addressing. Keep geometric props and updateable props for mapping.
775  clearAddressing(true);
776 
777  points_ = newMesh.points();
778  bounds_ = boundBox(points_, true);
779  faces_ = newMesh.faces();
780  owner_ = newMesh.faceOwner();
781  neighbour_ = newMesh.faceNeighbour();
782 
783  const polyPatchList& patches(newMesh.boundaryMesh());
784 
785  boundary_.clearGeom();
786  boundary_.clearAddressing();
787 
788  // Reset the number of patches in case the decomposition changed
789  boundary_.setSize(patches.size());
790 
791  forAll(boundary_, patchi)
792  {
793  // Construct new processor patches in case the decomposition changed
794  if (!isA<processorPolyPatch>(patches[patchi]))
795  {
796  boundary_[patchi] = polyPatch
797  (
798  boundary_[patchi],
799  boundary_,
800  patchi,
801  patches[patchi].size(),
802  patches[patchi].start()
803  );
804  }
805  else
806  {
807  boundary_.set(patchi, patches[patchi].clone(boundary_));
808  }
809  }
810 
811  // parallelData depends on the processorPatch ordering so force
812  // recalculation. Problem: should really be done in removeBoundary but
813  // there is some info in parallelData which might be interesting in between
814  // removeBoundary and addPatches.
815  globalMeshDataPtr_.clear();
816 
817  // Flags the mesh files as being changed
819 
820  // Check if the faces and cells are valid
821  forAll(faces_, facei)
822  {
823  const face& curFace = faces_[facei];
824 
825  if (min(curFace) < 0 || max(curFace) > points_.size())
826  {
828  << "Face " << facei << " contains vertex labels out of range: "
829  << curFace << " Max point index = " << points_.size()
830  << abort(FatalError);
831  }
832  }
833 
834  // Set the primitive mesh from the owner_, neighbour_.
835  // Works out from patch end where the active faces stop.
836  initMesh();
837 
838  // Calculate topology for the patches (processor-processor comms etc.)
839  boundary_.topoChange();
840 
841  // Calculate the geometry for the patches (transformation tensors etc.)
842  boundary_.calcGeometry();
843 
844  // Update the optional pointMesh with respect to the updated polyMesh
845  if (foundObject<pointMesh>(pointMesh::typeName))
846  {
847  pointMesh::New(*this).reset();
848  }
849 }
850 
851 
852 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
853 
855 {
856  clearOut();
857  resetMotion();
858 }
859 
860 
861 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
862 
864 {
866  {
867  return parent().dbDir();
868  }
869  else
870  {
871  return objectRegistry::dbDir();
872  }
873 }
874 
875 
877 {
878  return dbDir()/meshSubDir;
879 }
880 
881 
883 {
884  return points_.instance();
885 }
886 
887 
889 {
890  return faces_.instance();
891 }
892 
893 
895 {
896  if (geometricD_.x() == 0)
897  {
898  calcDirections();
899  }
900 
901  return geometricD_;
902 }
903 
904 
906 {
907  return cmptSum(geometricD() + Vector<label>::one)/2;
908 }
909 
910 
912 {
913  if (solutionD_.x() == 0)
914  {
915  calcDirections();
916  }
917 
918  return solutionD_;
919 }
920 
921 
923 {
924  return cmptSum(solutionD() + Vector<label>::one)/2;
925 }
926 
927 
929 {
930  if (tetBasePtIsPtr_.empty())
931  {
932  if (debug)
933  {
935  << "Forcing storage of base points."
936  << endl;
937  }
938 
939  tetBasePtIsPtr_.reset
940  (
941  new labelIOList
942  (
943  IOobject
944  (
945  "tetBasePtIs",
946  instance(),
947  meshSubDir,
948  *this,
951  ),
953  )
954  );
955  }
956 
957  return tetBasePtIsPtr_();
958 }
959 
960 
963 {
964  if (cellTreePtr_.empty())
965  {
966  cellTreePtr_.reset
967  (
969  (
971  (
972  false, // not cache bb
973  *this,
974  CELL_TETS // use tet-decomposition for any inside test
975  ),
976  treeBoundBox(points()).extend(1e-4),
977  8, // maxLevel
978  10, // leafsize
979  5.0 // duplicity
980  )
981  );
982  }
983 
984  return cellTreePtr_();
985 }
986 
987 
989 (
990  const List<polyPatch*>& p,
991  const bool validBoundary
992 )
993 {
994  if (boundaryMesh().size())
995  {
997  << "boundary already exists"
998  << abort(FatalError);
999  }
1000 
1001  // Reset valid directions
1002  geometricD_ = Zero;
1003  solutionD_ = Zero;
1004 
1005  boundary_.setSize(p.size());
1006 
1007  // Copy the patch pointers
1008  forAll(p, pI)
1009  {
1010  boundary_.set(pI, p[pI]);
1011  }
1012 
1013  // parallelData depends on the processorPatch ordering so force
1014  // recalculation. Problem: should really be done in removeBoundary but
1015  // there is some info in parallelData which might be interesting in between
1016  // removeBoundary and addPatches.
1017  globalMeshDataPtr_.clear();
1018 
1019  if (validBoundary)
1020  {
1021  // Calculate topology for the patches (processor-processor comms etc.)
1022  boundary_.topoChange();
1023 
1024  // Calculate the geometry for the patches (transformation tensors etc.)
1025  boundary_.calcGeometry();
1026 
1027  boundary_.checkDefinition();
1028  }
1029 }
1030 
1031 
1034  const List<pointZone*>& pz,
1035  const List<faceZone*>& fz,
1036  const List<cellZone*>& cz
1037 )
1038 {
1039  if (pointZones().size() || faceZones().size() || cellZones().size())
1040  {
1042  << "point, face or cell zone already exists"
1043  << abort(FatalError);
1044  }
1045 
1046  // Point zones
1047  if (pz.size())
1048  {
1049  pointZones_.setSize(pz.size());
1050 
1051  // Copy the zone pointers
1052  forAll(pz, pI)
1053  {
1054  pointZones_.set(pI, pz[pI]);
1055  }
1056 
1057  pointZones_.writeOpt() = IOobject::AUTO_WRITE;
1058  }
1059 
1060  // Face zones
1061  if (fz.size())
1062  {
1063  faceZones_.setSize(fz.size());
1064 
1065  // Copy the zone pointers
1066  forAll(fz, fI)
1067  {
1068  faceZones_.set(fI, fz[fI]);
1069  }
1070 
1071  faceZones_.writeOpt() = IOobject::AUTO_WRITE;
1072  }
1073 
1074  // Cell zones
1075  if (cz.size())
1076  {
1077  cellZones_.setSize(cz.size());
1078 
1079  // Copy the zone pointers
1080  forAll(cz, cI)
1081  {
1082  cellZones_.set(cI, cz[cI]);
1083  }
1084 
1085  cellZones_.writeOpt() = IOobject::AUTO_WRITE;
1086  }
1087 }
1088 
1089 
1092  const labelUList& newToOld,
1093  const bool validBoundary
1094 )
1095 {
1096  // Clear local fields and e.g. polyMesh parallelInfo. Do not clearGeom
1097  // so we keep PatchMeshObjects intact.
1098  boundary_.clearGeom();
1099  clearAddressing(true);
1100  // Clear all but PatchMeshObjects
1102  <
1103  polyMesh,
1106  >
1107  (
1108  *this
1109  );
1111  <
1112  pointMesh,
1113  TopologicalMeshObject,
1115  >
1116  (
1117  *this
1118  );
1119 
1120  boundary_.reorderPatches(newToOld, validBoundary);
1121 
1122  // Warn mesh objects
1123  meshObject::reorderPatches<polyMesh>(*this, newToOld, validBoundary);
1124  meshObject::reorderPatches<pointMesh>(*this, newToOld, validBoundary);
1125 }
1126 
1127 
1130  const label insertPatchi,
1131  const polyPatch& patch,
1132  const dictionary& patchFieldDict,
1133  const word& defaultPatchFieldType,
1134  const bool validBoundary
1135 )
1136 {
1137  const label sz = boundary_.size();
1138 
1139  label startFacei = nFaces();
1140  if (insertPatchi < sz)
1141  {
1142  startFacei = boundary_[insertPatchi].start();
1143  }
1144 
1145  // Create reordering list
1146  // patches before insert position stay as is
1147  // patches after insert position move one up
1148  labelList newToOld(boundary_.size()+1);
1149  for (label i = 0; i < insertPatchi; i++)
1150  {
1151  newToOld[i] = i;
1152  }
1153  for (label i = insertPatchi; i < sz; i++)
1154  {
1155  newToOld[i+1] = i;
1156  }
1157  newToOld[insertPatchi] = -1;
1158 
1159  reorderPatches(newToOld, false);
1160 
1161  // Clear local fields and e.g. polyMesh parallelInfo.
1162  //clearGeom(); // would clear out pointMesh as well
1163  boundary_.clearGeom();
1164  clearAddressing(true);
1165 
1166  // Clear all but PatchMeshObjects
1168  <
1169  polyMesh,
1172  >
1173  (
1174  *this
1175  );
1177  <
1178  pointMesh,
1179  TopologicalMeshObject,
1181  >
1182  (
1183  *this
1184  );
1185 
1186 
1187  // Insert polyPatch
1188  boundary_.set
1189  (
1190  insertPatchi,
1191  patch.clone
1192  (
1193  boundary_,
1194  insertPatchi, // index
1195  0, // size
1196  startFacei // start
1197  )
1198  );
1199 
1200  if (validBoundary)
1201  {
1202  boundary_.topoChange();
1203  }
1204 
1205  // Warn mesh objects
1206  meshObject::addPatch<polyMesh>(*this, insertPatchi);
1207  meshObject::addPatch<pointMesh>(*this, insertPatchi);
1208 }
1209 
1210 
1212 {
1213  if (clearedPrimitives_)
1214  {
1216  << "points deallocated"
1217  << abort(FatalError);
1218  }
1219 
1220  return points_;
1221 }
1222 
1223 
1225 {
1226  if (clearedPrimitives_)
1227  {
1229  << "faces deallocated"
1230  << abort(FatalError);
1231  }
1232 
1233  return faces_;
1234 }
1235 
1236 
1238 {
1239  return owner_;
1240 }
1241 
1242 
1244 {
1245  return neighbour_;
1246 }
1247 
1248 
1250 {
1251  if (!moving_)
1252  {
1253  return points_;
1254  }
1255 
1256  if (oldPointsPtr_.empty())
1257  {
1259  << "Old points have not been stored"
1260  << exit(FatalError);
1261  }
1262 
1263  return oldPointsPtr_();
1264 }
1265 
1266 
1268 {
1269  storeOldCellCentres_ = true;
1270 
1271  if (!moving_)
1272  {
1273  return cellCentres();
1274  }
1275 
1276  if (oldCellCentresPtr_.empty())
1277  {
1279  << "Old cell centres have not been stored"
1280  << exit(FatalError);
1281  }
1282 
1283  return oldCellCentresPtr_();
1284 }
1285 
1286 
1288 {
1289  return io.upToDate(points_);
1290 }
1291 
1292 
1294 {
1295  io.eventNo() = points_.eventNo()+1;
1296 }
1297 
1298 
1300 {
1301  if (debug)
1302  {
1304  << "Set points for time " << time().value()
1305  << " index " << time().timeIndex() << endl;
1306  }
1307 
1309 
1310  points_ = newPoints;
1311 
1313 
1314  // Adjust parallel shared points
1315  if (globalMeshDataPtr_.valid())
1316  {
1317  globalMeshDataPtr_().movePoints(points_);
1318  }
1319 
1320  // Force recalculation of all geometric data with new points
1321 
1322  bounds_ = boundBox(points_);
1323  boundary_.movePoints(points_);
1324 
1325  pointZones_.movePoints(points_);
1326  faceZones_.movePoints(points_);
1327  cellZones_.movePoints(points_);
1328 
1329  // Cell tree might become invalid
1330  cellTreePtr_.clear();
1331 
1332  // Reset valid directions (could change with rotation)
1333  geometricD_ = Zero;
1334  solutionD_ = Zero;
1335 
1336  meshObject::movePoints<polyMesh>(*this);
1337  meshObject::movePoints<pointMesh>(*this);
1338 }
1339 
1340 
1343  const pointField& newPoints
1344 )
1345 {
1346  if (debug)
1347  {
1349  << "Moving points for time " << time().value()
1350  << " index " << time().timeIndex() << endl;
1351  }
1352 
1353  // Pick up old points and cell centres
1354  if (curMotionTimeIndex_ != time().timeIndex())
1355  {
1356  oldPointsPtr_.clear();
1357  oldPointsPtr_.reset(new pointField(points_));
1358  if (storeOldCellCentres_)
1359  {
1360  oldCellCentresPtr_.clear();
1361  oldCellCentresPtr_.reset(new pointField(cellCentres()));
1362  }
1363  curMotionTimeIndex_ = time().timeIndex();
1364  }
1365 
1366  points_ = newPoints;
1367 
1368  bool moveError = false;
1369  if (debug)
1370  {
1371  // Check mesh motion
1372  if (checkMeshMotion(points_, true))
1373  {
1374  moveError = true;
1375 
1377  << "Moving the mesh with given points will "
1378  << "invalidate the mesh." << nl
1379  << "Mesh motion should not be executed." << endl;
1380  }
1381  }
1382 
1384 
1386  (
1387  points_,
1388  oldPoints()
1389  );
1390 
1391  // Adjust parallel shared points
1392  if (globalMeshDataPtr_.valid())
1393  {
1394  globalMeshDataPtr_().movePoints(points_);
1395  }
1396 
1397  // Force recalculation of all geometric data with new points
1398 
1399  bounds_ = boundBox(points_);
1400  boundary_.movePoints(points_);
1401 
1402  pointZones_.movePoints(points_);
1403  faceZones_.movePoints(points_);
1404  cellZones_.movePoints(points_);
1405 
1406  // Cell tree might become invalid
1407  cellTreePtr_.clear();
1408 
1409  // Reset valid directions (could change with rotation)
1410  geometricD_ = Zero;
1411  solutionD_ = Zero;
1412 
1413  meshObject::movePoints<polyMesh>(*this);
1414  meshObject::movePoints<pointMesh>(*this);
1415 
1416  if (debug && moveError)
1417  {
1418  // Write mesh to ease debugging. Note we want to avoid calling
1419  // e.g. fvMesh::write since meshPhi not yet complete.
1420  polyMesh::write();
1421  }
1422 
1423  return sweptVols;
1424 }
1425 
1426 
1428 {
1429  curMotionTimeIndex_ = -1;
1430  oldPointsPtr_.clear();
1431  oldCellCentresPtr_.clear();
1432 }
1433 
1434 
1436 {
1437  if (globalMeshDataPtr_.empty())
1438  {
1439  if (debug)
1440  {
1441  Pout<< "polyMesh::globalData() const : "
1442  << "Constructing parallelData from processor topology"
1443  << endl;
1444  }
1445  // Construct globalMeshData using processorPatch information only.
1446  globalMeshDataPtr_.reset(new globalMeshData(*this));
1447  }
1448 
1449  return globalMeshDataPtr_();
1450 }
1451 
1452 
1454 {
1455  return comm_;
1456 }
1457 
1458 
1460 {
1461  return comm_;
1462 }
1463 
1464 
1465 void Foam::polyMesh::removeFiles(const fileName& instanceDir) const
1466 {
1467  fileName meshFilesPath = thisDb().time().path()/instanceDir/meshDir();
1468 
1469  rm(meshFilesPath/"points");
1470  rm(meshFilesPath/"faces");
1471  rm(meshFilesPath/"owner");
1472  rm(meshFilesPath/"neighbour");
1473  rm(meshFilesPath/"cells");
1474  rm(meshFilesPath/"boundary");
1475  rm(meshFilesPath/"pointZones");
1476  rm(meshFilesPath/"faceZones");
1477  rm(meshFilesPath/"cellZones");
1478  rm(meshFilesPath/"meshModifiers");
1479  rm(meshFilesPath/"parallelData");
1480 
1481  // remove subdirectories
1482  if (isDir(meshFilesPath/"sets"))
1483  {
1484  rmDir(meshFilesPath/"sets");
1485  }
1486 }
1487 
1488 
1490 {
1491  removeFiles(instance());
1492 }
1493 
1494 
1497  const point& p,
1498  label& celli,
1499  label& tetFacei,
1500  label& tetPti
1501 ) const
1502 {
1503  celli = -1;
1504  tetFacei = -1;
1505  tetPti = -1;
1506 
1507  const indexedOctree<treeDataCell>& tree = cellTree();
1508 
1509  // Find point inside cell
1510  celli = tree.findInside(p);
1511 
1512  if (celli != -1)
1513  {
1514  // Check the nearest cell to see if the point is inside.
1515  findTetFacePt(celli, p, tetFacei, tetPti);
1516  }
1517 }
1518 
1519 
1522  const label celli,
1523  const point& p,
1524  label& tetFacei,
1525  label& tetPti
1526 ) const
1527 {
1528  const polyMesh& mesh = *this;
1529 
1530  tetIndices tet(polyMeshTetDecomposition::findTet(mesh, celli, p));
1531  tetFacei = tet.face();
1532  tetPti = tet.tetPt();
1533 }
1534 
1535 
1538  const point& p,
1539  label celli,
1540  const cellDecomposition decompMode
1541 ) const
1542 {
1543  switch (decompMode)
1544  {
1545  case FACE_PLANES:
1546  {
1547  return primitiveMesh::pointInCell(p, celli);
1548  }
1549  break;
1550 
1551  case FACE_CENTRE_TRIS:
1552  {
1553  // only test that point is on inside of plane defined by cell face
1554  // triangles
1555  const cell& cFaces = cells()[celli];
1556 
1557  forAll(cFaces, cFacei)
1558  {
1559  label facei = cFaces[cFacei];
1560  const face& f = faces_[facei];
1561  const point& fc = faceCentres()[facei];
1562  bool isOwn = (owner_[facei] == celli);
1563 
1564  forAll(f, fp)
1565  {
1566  label pointi;
1567  label nextPointi;
1568 
1569  if (isOwn)
1570  {
1571  pointi = f[fp];
1572  nextPointi = f.nextLabel(fp);
1573  }
1574  else
1575  {
1576  pointi = f.nextLabel(fp);
1577  nextPointi = f[fp];
1578  }
1579 
1580  triPointRef faceTri
1581  (
1582  points()[pointi],
1583  points()[nextPointi],
1584  fc
1585  );
1586 
1587  vector proj = p - faceTri.centre();
1588 
1589  if ((faceTri.area() & proj) > 0)
1590  {
1591  return false;
1592  }
1593  }
1594  }
1595  return true;
1596  }
1597  break;
1598 
1599  case FACE_DIAG_TRIS:
1600  {
1601  // only test that point is on inside of plane defined by cell face
1602  // triangles
1603  const cell& cFaces = cells()[celli];
1604 
1605  forAll(cFaces, cFacei)
1606  {
1607  label facei = cFaces[cFacei];
1608  const face& f = faces_[facei];
1609 
1610  for (label tetPti = 1; tetPti < f.size() - 1; tetPti++)
1611  {
1612  // Get tetIndices of face triangle
1613  tetIndices faceTetIs(celli, facei, tetPti);
1614 
1615  triPointRef faceTri = faceTetIs.faceTri(*this);
1616 
1617  vector proj = p - faceTri.centre();
1618 
1619  if ((faceTri.area() & proj) > 0)
1620  {
1621  return false;
1622  }
1623  }
1624  }
1625 
1626  return true;
1627  }
1628  break;
1629 
1630  case CELL_TETS:
1631  {
1632  label tetFacei;
1633  label tetPti;
1634 
1635  findTetFacePt(celli, p, tetFacei, tetPti);
1636 
1637  return tetFacei != -1;
1638  }
1639  break;
1640  }
1641 
1642  return false;
1643 }
1644 
1645 
1648  const point& p,
1649  const cellDecomposition decompMode
1650 ) const
1651 {
1652  if
1653  (
1654  Pstream::parRun()
1655  && (decompMode == FACE_DIAG_TRIS || decompMode == CELL_TETS)
1656  )
1657  {
1658  // Force construction of face-diagonal decomposition before testing
1659  // for zero cells.
1660  //
1661  // If parallel running a local domain might have zero cells so never
1662  // construct the face-diagonal decomposition which uses parallel
1663  // transfers.
1664  (void)tetBasePtIs();
1665  }
1666 
1667  if (nCells() == 0)
1668  {
1669  return -1;
1670  }
1671 
1672  if (decompMode == CELL_TETS)
1673  {
1674  // Advanced search method utilising an octree
1675  // and tet-decomposition of the cells
1676 
1677  label celli;
1678  label tetFacei;
1679  label tetPti;
1680 
1681  findCellFacePt(p, celli, tetFacei, tetPti);
1682 
1683  return celli;
1684  }
1685  else
1686  {
1687  // Approximate search avoiding the construction of an octree
1688  // and cell decomposition
1689 
1690  // Find the nearest cell centre to this location
1691  label celli = findNearestCell(p);
1692 
1693  // If point is in the nearest cell return
1694  if (pointInCell(p, celli, decompMode))
1695  {
1696  return celli;
1697  }
1698  else
1699  {
1700  // Point is not in the nearest cell so search all cells
1701 
1702  for (label celli = 0; celli < nCells(); celli++)
1703  {
1704  if (pointInCell(p, celli, decompMode))
1705  {
1706  return celli;
1707  }
1708  }
1709 
1710  return -1;
1711  }
1712  }
1713 }
1714 
1715 
1716 // ************************************************************************* //
const fvPatchList & patches
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:453
PtrList< polyPatch > polyPatchList
container classes for polyPatch
Definition: polyPatchList.H:45
bool moving_
Member data pending transfer to fvMesh.
Definition: polyMesh.H:307
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorder and trim existing patches. If validBoundary the new.
Definition: polyMesh.C:1091
bool upToDate(const regIOobject &) const
Return true if up-to-date with respect to given object.
Definition: regIOobject.C:329
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:288
label tetPt() const
Return the characterising tetPtI.
Definition: tetIndicesI.H:55
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1342
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void clearAddressing()
Clear topological data.
A triangle primitive used to calculate face areas and swept volumes.
Definition: triangle.H:57
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
A class for handling file names.
Definition: fileName.H:79
static labelList findFaceBasePts(const polyMesh &mesh, scalar tol=minTetQuality, bool report=false)
Find a suitable base point for each face for decomposition.
bool set(const label) const
Is element set.
Definition: PtrListI.H:65
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:888
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
const meshCellZones & cellZones() const
Return cell zones.
Definition: polyMesh.H:501
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1243
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:928
void setUpToDatePoints(regIOobject &io) const
Set io to be up-to-date with points.
Definition: polyMesh.C:1293
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:74
void resetMotion() const
Reset motion.
Definition: polyMesh.C:1427
label nFaces() const
triPointRef faceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the face for.
Definition: tetIndicesI.H:123
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
void clearOut()
Clear all geometry and addressing unnecessary for CFD.
uint8_t direction
Definition: direction.H:45
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:328
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:911
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
label nCells() const
Cmpt cmptSum(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:500
void reset(const polyMesh &)
Reset mesh.
Definition: polyMesh.C:772
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:325
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:100
void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorders the patches. Ordering does not have to be done in.
primitiveMesh()
Construct null.
Definition: primitiveMesh.C:39
void clearGeom()
Clear geometry at this level and at patches.
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:58
void setPointsInstance(const fileName &)
Set the instance for the points files.
Definition: polyMeshIO.C:57
const cellList & cells() const
bool pointInCell(const point &p, label celli) const
Return true if the point is in the cell.
const objectRegistry & parent() const
Return the parent objectRegistry.
virtual const fileName & dbDir() const
Override the objectRegistry dbDir for a single-region case.
Definition: polyMesh.C:863
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
static label worldComm
Default communicator (all processors)
Definition: UPstream.H:278
bool pointInCell(const point &p, label celli, const cellDecomposition=CELL_TETS) const
Test if point p is in the celli.
Definition: polyMesh.C:1537
fvMesh & mesh
void movePoints(const pointField &)
Correct zones after moving points.
Definition: MeshZones.C:505
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:100
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:905
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
const meshPointZones & pointZones() const
Return point zones.
Definition: polyMesh.H:489
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1211
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
vector area() const
Return vector area.
Definition: triangleI.H:96
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:894
objectRegistry(const Time &db, const label nIoObjects=128)
Construct the time objectRegistry given an initial estimate.
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: polyPatch.H:225
void reset()
Reset pointMesh with respect to the updated polyMesh.
Definition: pointMesh.C:188
bool topoChanged_
Is the mesh topology changing.
Definition: polyMesh.H:310
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:99
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:882
void clearGeom()
Clear geometry.
tmp< scalarField > movePoints(const pointField &p, const pointField &oldP)
Move points, returns volumes swept by faces in motion.
void topoChange()
Correct polyBoundaryMesh after topology update.
void removeFiles() const
Remove all files from mesh instance()
Definition: polyMesh.C:1489
void movePoints(const pointField &)
Correct polyBoundaryMesh after moving points.
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
A class for handling words, derived from string.
Definition: word.H:59
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition: polyMesh.C:1496
bool checkDefinition(const bool report=false) const
Check boundary definition. Return true if in error.
static void clearUpto(objectRegistry &)
Clear all meshObject derived from FromType up to (but not including)
Definition: MeshObject.C:540
virtual void addPatch(const label insertPatchi, const polyPatch &patch, const dictionary &patchFieldDict, const word &defaultPatchFieldType, const bool validBoundary)
Add/insert single patch. If validBoundary the new situation.
Definition: polyMesh.C:1129
virtual const pointField & oldPoints() const
Return old points for mesh motion.
Definition: polyMesh.C:1249
dimensionSet cmptMag(const dimensionSet &)
Definition: dimensionSet.C:275
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1237
const Type & value() const
Return const reference to value.
Point centre() const
Return centre (centroid)
Definition: triangleI.H:89
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
label eventNo() const
Event number at last update.
Definition: regIOobjectI.H:89
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1435
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
Storage and named access for the indices of a tet which is part of the decomposition of a cell...
Definition: tetIndices.H:81
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1224
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1453
void addZones(const List< pointZone *> &pz, const List< faceZone *> &fz, const List< cellZone *> &cz)
Add mesh zones.
Definition: polyMesh.C:1033
const vectorField & cellCentres() const
errorManip< error > abort(error &err)
Definition: errorManip.H:131
void addPatches(const List< polyPatch *> &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:989
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
label face() const
Return the face.
Definition: tetIndicesI.H:43
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
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:57
bool rmDir(const fileName &)
Remove a directory and its contents.
Definition: POSIX.C:1047
static const char nl
Definition: Ostream.H:260
label findInside(const point &) const
Find shape containing point. Only implemented for certain.
void resetPrimitives(pointField &&points, faceList &&faces, labelList &&owner, labelList &&neighbour, const labelList &patchSizes, const labelList &patchStarts, const bool validBoundary=true)
Reset mesh primitive data. Assumes all patch info correct.
Definition: polyMesh.C:672
const Time & time() const
Return time.
defineTypeNameAndDebug(combustionModel, 0)
labelList f(nPoints)
A List of objects of type <T> with automated input and output using a compact storage. Behaves like IOList except when binary output in case it writes a CompactListList.
Definition: CompactIOList.H:54
Encapsulation of data needed to search in/for cells. Used to find the cell containing a point (e...
Definition: treeDataCell.H:54
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:77
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:52
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:922
bool upToDatePoints(const regIOobject &io) const
Return true if io is up-to-date with points.
Definition: polyMesh.C:1287
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
const vectorField & faceCentres() const
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition: polyMesh.C:962
writeOption writeOpt() const
Definition: IOobject.H:375
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:28
virtual const pointField & oldCellCentres() const
Return old cell centres for mesh motion.
Definition: polyMesh.C:1267
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
void clearAddressing()
Clear addressing at this level and at patches.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:57
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
const meshFaceZones & faceZones() const
Return face zones.
Definition: polyMesh.H:495
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:876
label timeIndex
Definition: getTimeIndex.H:4
label findNearestCell(const point &location) const
Find the cell with the nearest cell centre to location.
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:1521
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1647
virtual ~polyMesh()
Destructor.
Definition: polyMesh.C:854
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:52
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
dimensioned< scalar > mag(const dimensioned< Type > &)
bool headerOk()
Read header of local object without type-checking.
virtual const fileName & dbDir() const
Local directory path of this objectRegistry relative to the time.
label nPoints() const
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
const objectRegistry & thisDb() const
Return the object registry.
Definition: polyMesh.H:516
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
fileName & instance() const
Return the instance directory, constant, system, <time> etc.
Definition: IOobject.C:355
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
virtual bool write(const bool write=true) const
Write using setting from DB.
A class for managing temporary objects.
Definition: PtrList.H:53
bool rm(const fileName &)
Remove a file, returning true if successful otherwise false.
Definition: POSIX.C:1017
Registry of regIOobjects.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
virtual void setPoints(const pointField &)
Reset the points.
Definition: polyMesh.C:1299
readOption readOpt() const
Definition: IOobject.H:365
fileName path() const
Return path.
Definition: TimePaths.H:139
const word & headerClassName() const
Return name of the class name read from header.
Definition: IOobject.H:321
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
Inter-processor communications stream.
Definition: UPstream.H:58
virtual bool checkMeshMotion(const pointField &newPoints, const bool report=false, const bool detailedReport=false) const
Check mesh motion for correctness given motion points.
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:42
Namespace for OpenFOAM.
polyMesh(const IOobject &io)
Construct from IOobject.
Definition: polyMesh.C:163
#define InfoInFunction
Report an information message using Foam::Info.