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-2019 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  IOobject io
141  (
142  "tetBasePtIs",
143  instance(),
144  meshSubDir,
145  *this,
148  );
149 
150  if (io.typeHeaderOk<labelIOList>())
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  moving_(false),
291  topoChanging_(false),
292  curMotionTimeIndex_(-1),
293  oldPointsPtr_(nullptr),
294  oldCellCentresPtr_(nullptr),
295  storeOldCellCentres_(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_.updateMesh();
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  moving_(false),
475  topoChanging_(false),
476  curMotionTimeIndex_(-1),
477  oldPointsPtr_(nullptr),
478  oldCellCentresPtr_(nullptr),
479  storeOldCellCentres_(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  moving_(false),
629  topoChanging_(false),
630  curMotionTimeIndex_(-1),
631  oldPointsPtr_(nullptr),
632  oldCellCentresPtr_(nullptr),
633  storeOldCellCentres_(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  // Optimized 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_.updateMesh();
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 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
773 
775 {
776  clearOut();
777  resetMotion();
778 }
779 
780 
781 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
782 
784 {
786  {
787  return parent().dbDir();
788  }
789  else
790  {
791  return objectRegistry::dbDir();
792  }
793 }
794 
795 
797 {
798  return dbDir()/meshSubDir;
799 }
800 
801 
803 {
804  return points_.instance();
805 }
806 
807 
809 {
810  return faces_.instance();
811 }
812 
813 
815 {
816  if (geometricD_.x() == 0)
817  {
818  calcDirections();
819  }
820 
821  return geometricD_;
822 }
823 
824 
826 {
827  return cmptSum(geometricD() + Vector<label>::one)/2;
828 }
829 
830 
832 {
833  if (solutionD_.x() == 0)
834  {
835  calcDirections();
836  }
837 
838  return solutionD_;
839 }
840 
841 
843 {
844  return cmptSum(solutionD() + Vector<label>::one)/2;
845 }
846 
847 
849 {
850  if (tetBasePtIsPtr_.empty())
851  {
852  if (debug)
853  {
855  << "Forcing storage of base points."
856  << endl;
857  }
858 
859  tetBasePtIsPtr_.reset
860  (
861  new labelIOList
862  (
863  IOobject
864  (
865  "tetBasePtIs",
866  instance(),
867  meshSubDir,
868  *this,
871  ),
873  )
874  );
875  }
876 
877  return tetBasePtIsPtr_();
878 }
879 
880 
883 {
884  if (cellTreePtr_.empty())
885  {
886  cellTreePtr_.reset
887  (
889  (
891  (
892  false, // not cache bb
893  *this,
894  CELL_TETS // use tet-decomposition for any inside test
895  ),
896  treeBoundBox(points()).extend(1e-4),
897  8, // maxLevel
898  10, // leafsize
899  5.0 // duplicity
900  )
901  );
902  }
903 
904  return cellTreePtr_();
905 }
906 
907 
909 (
910  const List<polyPatch*>& p,
911  const bool validBoundary
912 )
913 {
914  if (boundaryMesh().size())
915  {
917  << "boundary already exists"
918  << abort(FatalError);
919  }
920 
921  // Reset valid directions
922  geometricD_ = Zero;
923  solutionD_ = Zero;
924 
925  boundary_.setSize(p.size());
926 
927  // Copy the patch pointers
928  forAll(p, pI)
929  {
930  boundary_.set(pI, p[pI]);
931  }
932 
933  // parallelData depends on the processorPatch ordering so force
934  // recalculation. Problem: should really be done in removeBoundary but
935  // there is some info in parallelData which might be interesting in between
936  // removeBoundary and addPatches.
937  globalMeshDataPtr_.clear();
938 
939  if (validBoundary)
940  {
941  // Calculate topology for the patches (processor-processor comms etc.)
942  boundary_.updateMesh();
943 
944  // Calculate the geometry for the patches (transformation tensors etc.)
945  boundary_.calcGeometry();
946 
947  boundary_.checkDefinition();
948  }
949 }
950 
951 
953 (
954  const List<pointZone*>& pz,
955  const List<faceZone*>& fz,
956  const List<cellZone*>& cz
957 )
958 {
959  if (pointZones().size() || faceZones().size() || cellZones().size())
960  {
962  << "point, face or cell zone already exists"
963  << abort(FatalError);
964  }
965 
966  // Point zones
967  if (pz.size())
968  {
969  pointZones_.setSize(pz.size());
970 
971  // Copy the zone pointers
972  forAll(pz, pI)
973  {
974  pointZones_.set(pI, pz[pI]);
975  }
976 
977  pointZones_.writeOpt() = IOobject::AUTO_WRITE;
978  }
979 
980  // Face zones
981  if (fz.size())
982  {
983  faceZones_.setSize(fz.size());
984 
985  // Copy the zone pointers
986  forAll(fz, fI)
987  {
988  faceZones_.set(fI, fz[fI]);
989  }
990 
991  faceZones_.writeOpt() = IOobject::AUTO_WRITE;
992  }
993 
994  // Cell zones
995  if (cz.size())
996  {
997  cellZones_.setSize(cz.size());
998 
999  // Copy the zone pointers
1000  forAll(cz, cI)
1001  {
1002  cellZones_.set(cI, cz[cI]);
1003  }
1004 
1005  cellZones_.writeOpt() = IOobject::AUTO_WRITE;
1006  }
1007 }
1008 
1009 
1012  const labelUList& newToOld,
1013  const bool validBoundary
1014 )
1015 {
1016  // Clear local fields and e.g. polyMesh parallelInfo. Do not clearGeom
1017  // so we keep PatchMeshObjects intact.
1018  boundary_.clearGeom();
1019  clearAddressing(true);
1020  // Clear all but PatchMeshObjects
1022  <
1023  polyMesh,
1026  >
1027  (
1028  *this
1029  );
1031  <
1032  pointMesh,
1033  TopologicalMeshObject,
1035  >
1036  (
1037  *this
1038  );
1039 
1040  boundary_.shuffle(newToOld, validBoundary);
1041 
1042  // Warn mesh objects
1043  meshObject::reorderPatches<polyMesh>(*this, newToOld, validBoundary);
1044  meshObject::reorderPatches<pointMesh>(*this, newToOld, validBoundary);
1045 }
1046 
1047 
1050  const label insertPatchi,
1051  const polyPatch& patch,
1052  const dictionary& patchFieldDict,
1053  const word& defaultPatchFieldType,
1054  const bool validBoundary
1055 )
1056 {
1057  const label sz = boundary_.size();
1058 
1059  label startFacei = nFaces();
1060  if (insertPatchi < sz)
1061  {
1062  startFacei = boundary_[insertPatchi].start();
1063  }
1064 
1065  // Create reordering list
1066  // patches before insert position stay as is
1067  // patches after insert position move one up
1068  labelList newToOld(boundary_.size()+1);
1069  for (label i = 0; i < insertPatchi; i++)
1070  {
1071  newToOld[i] = i;
1072  }
1073  for (label i = insertPatchi; i < sz; i++)
1074  {
1075  newToOld[i+1] = i;
1076  }
1077  newToOld[insertPatchi] = -1;
1078 
1079  reorderPatches(newToOld, false);
1080 
1081  // Clear local fields and e.g. polyMesh parallelInfo.
1082  //clearGeom(); // would clear out pointMesh as well
1083  boundary_.clearGeom();
1084  clearAddressing(true);
1085 
1086  // Clear all but PatchMeshObjects
1088  <
1089  polyMesh,
1092  >
1093  (
1094  *this
1095  );
1097  <
1098  pointMesh,
1099  TopologicalMeshObject,
1101  >
1102  (
1103  *this
1104  );
1105 
1106 
1107  // Insert polyPatch
1108  boundary_.set
1109  (
1110  insertPatchi,
1111  patch.clone
1112  (
1113  boundary_,
1114  insertPatchi, // index
1115  0, // size
1116  startFacei // start
1117  )
1118  );
1119 
1120  if (validBoundary)
1121  {
1122  boundary_.updateMesh();
1123  }
1124 
1125  // Warn mesh objects
1126  meshObject::addPatch<polyMesh>(*this, insertPatchi);
1127  meshObject::addPatch<pointMesh>(*this, insertPatchi);
1128 }
1129 
1130 
1132 {
1133  if (clearedPrimitives_)
1134  {
1136  << "points deallocated"
1137  << abort(FatalError);
1138  }
1139 
1140  return points_;
1141 }
1142 
1143 
1145 {
1146  return io.upToDate(points_);
1147 }
1148 
1149 
1151 {
1152  io.eventNo() = points_.eventNo()+1;
1153 }
1154 
1155 
1157 {
1158  if (clearedPrimitives_)
1159  {
1161  << "faces deallocated"
1162  << abort(FatalError);
1163  }
1164 
1165  return faces_;
1166 }
1167 
1168 
1170 {
1171  return owner_;
1172 }
1173 
1174 
1176 {
1177  return neighbour_;
1178 }
1179 
1180 
1182 {
1183  if (!moving_)
1184  {
1185  return points_;
1186  }
1187 
1188  if (oldPointsPtr_.empty())
1189  {
1191  << "Old points have not been stored"
1192  << exit(FatalError);
1193  }
1194 
1195  return oldPointsPtr_();
1196 }
1197 
1198 
1200 {
1201  storeOldCellCentres_ = true;
1202 
1203  if (!moving_)
1204  {
1205  return cellCentres();
1206  }
1207 
1208  if (oldCellCentresPtr_.empty())
1209  {
1211  << "Old cell centres have not been stored"
1212  << exit(FatalError);
1213  }
1214 
1215  return oldCellCentresPtr_();
1216 }
1217 
1218 
1221  const polyMesh& mesh
1222 )
1223 {
1224  const word instance
1225  (
1226  mesh.time().findInstance
1227  (
1228  mesh.meshDir(),
1229  "points0",
1231  )
1232  );
1233 
1234  if (instance != mesh.time().constant())
1235  {
1236  // Points0 written to a time folder
1237 
1238  return IOobject
1239  (
1240  "points0",
1241  instance,
1243  mesh,
1246  false
1247  );
1248  }
1249  else
1250  {
1251  // Check that points0 are actually in constant directory
1252 
1253  IOobject io
1254  (
1255  "points0",
1256  instance,
1258  mesh,
1261  false
1262  );
1263 
1264  if (io.typeHeaderOk<pointIOField>())
1265  {
1266  return io;
1267  }
1268  else
1269  {
1270  // Copy of original mesh points
1271  return IOobject
1272  (
1273  "points",
1274  instance,
1276  mesh,
1279  false
1280  );
1281  }
1282  }
1283 }
1284 
1285 
1288  const pointField& newPoints
1289 )
1290 {
1291  if (debug)
1292  {
1294  << "Moving points for time " << time().value()
1295  << " index " << time().timeIndex() << endl;
1296  }
1297 
1298  moving(true);
1299 
1300  // Pick up old points and cell centres
1301  if (curMotionTimeIndex_ != time().timeIndex())
1302  {
1303  oldPointsPtr_.clear();
1304  oldPointsPtr_.reset(new pointField(points_));
1305  if (storeOldCellCentres_)
1306  {
1307  oldCellCentresPtr_.clear();
1308  oldCellCentresPtr_.reset(new pointField(cellCentres()));
1309  }
1310  curMotionTimeIndex_ = time().timeIndex();
1311  }
1312 
1313  points_ = newPoints;
1314 
1315  bool moveError = false;
1316  if (debug)
1317  {
1318  // Check mesh motion
1319  if (checkMeshMotion(points_, true))
1320  {
1321  moveError = true;
1322 
1324  << "Moving the mesh with given points will "
1325  << "invalidate the mesh." << nl
1326  << "Mesh motion should not be executed." << endl;
1327  }
1328  }
1329 
1330  points_.writeOpt() = IOobject::AUTO_WRITE;
1331  points_.instance() = time().timeName();
1332  points_.eventNo() = getEvent();
1333 
1334  if (tetBasePtIsPtr_.valid())
1335  {
1336  tetBasePtIsPtr_().writeOpt() = IOobject::AUTO_WRITE;
1337  tetBasePtIsPtr_().instance() = time().timeName();
1338  tetBasePtIsPtr_().eventNo() = getEvent();
1339  }
1340 
1342  (
1343  points_,
1344  oldPoints()
1345  );
1346 
1347  // Adjust parallel shared points
1348  if (globalMeshDataPtr_.valid())
1349  {
1350  globalMeshDataPtr_().movePoints(points_);
1351  }
1352 
1353  // Force recalculation of all geometric data with new points
1354 
1355  bounds_ = boundBox(points_);
1356  boundary_.movePoints(points_);
1357 
1358  pointZones_.movePoints(points_);
1359  faceZones_.movePoints(points_);
1360  cellZones_.movePoints(points_);
1361 
1362  // Cell tree might become invalid
1363  cellTreePtr_.clear();
1364 
1365  // Reset valid directions (could change with rotation)
1366  geometricD_ = Zero;
1367  solutionD_ = Zero;
1368 
1369  meshObject::movePoints<polyMesh>(*this);
1370  meshObject::movePoints<pointMesh>(*this);
1371 
1372  const_cast<Time&>(time()).functionObjects().movePoints(*this);
1373 
1374 
1375  if (debug && moveError)
1376  {
1377  // Write mesh to ease debugging. Note we want to avoid calling
1378  // e.g. fvMesh::write since meshPhi not yet complete.
1379  polyMesh::write();
1380  }
1381 
1382  return sweptVols;
1383 }
1384 
1385 
1387 {
1388  curMotionTimeIndex_ = -1;
1389  oldPointsPtr_.clear();
1390  oldCellCentresPtr_.clear();
1391 }
1392 
1393 
1395 {
1396  if (globalMeshDataPtr_.empty())
1397  {
1398  if (debug)
1399  {
1400  Pout<< "polyMesh::globalData() const : "
1401  << "Constructing parallelData from processor topology"
1402  << endl;
1403  }
1404  // Construct globalMeshData using processorPatch information only.
1405  globalMeshDataPtr_.reset(new globalMeshData(*this));
1406  }
1407 
1408  return globalMeshDataPtr_();
1409 }
1410 
1411 
1413 {
1414  return comm_;
1415 }
1416 
1417 
1419 {
1420  return comm_;
1421 }
1422 
1423 
1424 void Foam::polyMesh::removeFiles(const fileName& instanceDir) const
1425 {
1426  fileName meshFilesPath = thisDb().time().path()/instanceDir/meshDir();
1427 
1428  rm(meshFilesPath/"points");
1429  rm(meshFilesPath/"faces");
1430  rm(meshFilesPath/"owner");
1431  rm(meshFilesPath/"neighbour");
1432  rm(meshFilesPath/"cells");
1433  rm(meshFilesPath/"boundary");
1434  rm(meshFilesPath/"pointZones");
1435  rm(meshFilesPath/"faceZones");
1436  rm(meshFilesPath/"cellZones");
1437  rm(meshFilesPath/"meshModifiers");
1438  rm(meshFilesPath/"parallelData");
1439 
1440  // remove subdirectories
1441  if (isDir(meshFilesPath/"sets"))
1442  {
1443  rmDir(meshFilesPath/"sets");
1444  }
1445 }
1446 
1447 
1449 {
1450  removeFiles(instance());
1451 }
1452 
1453 
1456  const point& p,
1457  label& celli,
1458  label& tetFacei,
1459  label& tetPti
1460 ) const
1461 {
1462  celli = -1;
1463  tetFacei = -1;
1464  tetPti = -1;
1465 
1466  const indexedOctree<treeDataCell>& tree = cellTree();
1467 
1468  // Find point inside cell
1469  celli = tree.findInside(p);
1470 
1471  if (celli != -1)
1472  {
1473  // Check the nearest cell to see if the point is inside.
1474  findTetFacePt(celli, p, tetFacei, tetPti);
1475  }
1476 }
1477 
1478 
1481  const label celli,
1482  const point& p,
1483  label& tetFacei,
1484  label& tetPti
1485 ) const
1486 {
1487  const polyMesh& mesh = *this;
1488 
1489  tetIndices tet(polyMeshTetDecomposition::findTet(mesh, celli, p));
1490  tetFacei = tet.face();
1491  tetPti = tet.tetPt();
1492 }
1493 
1494 
1497  const point& p,
1498  label celli,
1499  const cellDecomposition decompMode
1500 ) const
1501 {
1502  switch (decompMode)
1503  {
1504  case FACE_PLANES:
1505  {
1506  return primitiveMesh::pointInCell(p, celli);
1507  }
1508  break;
1509 
1510  case FACE_CENTRE_TRIS:
1511  {
1512  // only test that point is on inside of plane defined by cell face
1513  // triangles
1514  const cell& cFaces = cells()[celli];
1515 
1516  forAll(cFaces, cFacei)
1517  {
1518  label facei = cFaces[cFacei];
1519  const face& f = faces_[facei];
1520  const point& fc = faceCentres()[facei];
1521  bool isOwn = (owner_[facei] == celli);
1522 
1523  forAll(f, fp)
1524  {
1525  label pointi;
1526  label nextPointi;
1527 
1528  if (isOwn)
1529  {
1530  pointi = f[fp];
1531  nextPointi = f.nextLabel(fp);
1532  }
1533  else
1534  {
1535  pointi = f.nextLabel(fp);
1536  nextPointi = f[fp];
1537  }
1538 
1539  triPointRef faceTri
1540  (
1541  points()[pointi],
1542  points()[nextPointi],
1543  fc
1544  );
1545 
1546  vector proj = p - faceTri.centre();
1547 
1548  if ((faceTri.area() & proj) > 0)
1549  {
1550  return false;
1551  }
1552  }
1553  }
1554  return true;
1555  }
1556  break;
1557 
1558  case FACE_DIAG_TRIS:
1559  {
1560  // only test that point is on inside of plane defined by cell face
1561  // triangles
1562  const cell& cFaces = cells()[celli];
1563 
1564  forAll(cFaces, cFacei)
1565  {
1566  label facei = cFaces[cFacei];
1567  const face& f = faces_[facei];
1568 
1569  for (label tetPti = 1; tetPti < f.size() - 1; tetPti++)
1570  {
1571  // Get tetIndices of face triangle
1572  tetIndices faceTetIs(celli, facei, tetPti);
1573 
1574  triPointRef faceTri = faceTetIs.faceTri(*this);
1575 
1576  vector proj = p - faceTri.centre();
1577 
1578  if ((faceTri.area() & proj) > 0)
1579  {
1580  return false;
1581  }
1582  }
1583  }
1584 
1585  return true;
1586  }
1587  break;
1588 
1589  case CELL_TETS:
1590  {
1591  label tetFacei;
1592  label tetPti;
1593 
1594  findTetFacePt(celli, p, tetFacei, tetPti);
1595 
1596  return tetFacei != -1;
1597  }
1598  break;
1599  }
1600 
1601  return false;
1602 }
1603 
1604 
1607  const point& p,
1608  const cellDecomposition decompMode
1609 ) const
1610 {
1611  if
1612  (
1613  Pstream::parRun()
1614  && (decompMode == FACE_DIAG_TRIS || decompMode == CELL_TETS)
1615  )
1616  {
1617  // Force construction of face-diagonal decomposition before testing
1618  // for zero cells.
1619  //
1620  // If parallel running a local domain might have zero cells so never
1621  // construct the face-diagonal decomposition which uses parallel
1622  // transfers.
1623  (void)tetBasePtIs();
1624  }
1625 
1626  if (nCells() == 0)
1627  {
1628  return -1;
1629  }
1630 
1631  if (decompMode == CELL_TETS)
1632  {
1633  // Advanced search method utilizing an octree
1634  // and tet-decomposition of the cells
1635 
1636  label celli;
1637  label tetFacei;
1638  label tetPti;
1639 
1640  findCellFacePt(p, celli, tetFacei, tetPti);
1641 
1642  return celli;
1643  }
1644  else
1645  {
1646  // Approximate search avoiding the construction of an octree
1647  // and cell decomposition
1648 
1649  // Find the nearest cell centre to this location
1650  label celli = findNearestCell(p);
1651 
1652  // If point is in the nearest cell return
1653  if (pointInCell(p, celli, decompMode))
1654  {
1655  return celli;
1656  }
1657  else
1658  {
1659  // Point is not in the nearest cell so search all cells
1660 
1661  for (label celli = 0; celli < nCells(); celli++)
1662  {
1663  if (pointInCell(p, celli, decompMode))
1664  {
1665  return celli;
1666  }
1667  }
1668 
1669  return -1;
1670  }
1671  }
1672 }
1673 
1674 
1675 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
PtrList< polyPatch > polyPatchList
container classes for polyPatch
Definition: polyPatchList.H:45
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorder and trim existing patches. If validBoundary the new.
Definition: polyMesh.C:1011
bool upToDate(const regIOobject &) const
Return true if up-to-date with respect to given object.
Definition: regIOobject.C:322
label tetPt() const
Return the characterising tetPtI.
Definition: tetIndicesI.H:52
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1287
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void clearAddressing()
Clear topological data.
fileName path() const
Return path.
Definition: Time.H:266
A triangle primitive used to calculate face areas and swept volumes.
Definition: triangle.H:57
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
bool typeHeaderOk(const bool checkType=true)
Read header (uses typeFilePath to find file) and check header.
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
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:476
bool moving() const
Is mesh moving.
Definition: polyMesh.H:512
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:808
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:75
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1175
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:848
virtual void setUpToDatePoints(regIOobject &io) const
Set io to be up-to-date with points.
Definition: polyMesh.C:1150
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:74
void resetMotion() const
Reset motion.
Definition: polyMesh.C:1386
label nFaces() const
triPointRef faceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the face for.
Definition: tetIndicesI.H:113
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:312
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:831
void shuffle(const labelUList &newToOld, const bool validBoundary)
Reorders patches. Ordering does not have to be done in.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:209
label nCells() const
Cmpt cmptSum(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:500
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
word findInstance(const fileName &dir, const word &name=word::null, const IOobject::readOption rOpt=IOobject::MUST_READ, const word &stopInstance=word::null) const
Return the location of "dir" containing the file "name".
Definition: Time.C:649
label nextLabel(const label i) const
Next vertex on face.
Definition: faceI.H:117
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
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:783
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
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:1496
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:98
void movePoints(const pointField &)
Correct zone mesh after moving points.
Definition: ZoneMesh.C:505
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:825
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1131
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:814
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:222
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:802
tmp< scalarField > movePoints(const pointField &p, const pointField &oldP)
Move points, returns volumes swept by faces in motion.
void removeFiles() const
Remove all files from mesh instance()
Definition: polyMesh.C:1448
void movePoints(const pointField &)
Correct polyBoundaryMesh after moving points.
const pointZoneMesh & pointZones() const
Return point zone mesh.
Definition: polyMesh.H:470
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
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:482
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition: polyMesh.C:1455
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:494
label getEvent() const
Return new event number.
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:1049
const word & constant() const
Return constant name.
Definition: TimePaths.H:124
virtual const pointField & oldPoints() const
Return old points for mesh motion.
Definition: polyMesh.C:1181
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1169
const Type & value() const
Return const reference to value.
Point centre() const
Return centre (centroid)
Definition: triangleI.H:89
label eventNo() const
Event number at last update.
Definition: regIOobjectI.H:83
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1394
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:1156
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1412
void addZones(const List< pointZone *> &pz, const List< faceZone *> &fz, const List< cellZone *> &cz)
Add mesh zones.
Definition: polyMesh.C:953
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:909
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
const Cmpt & x() const
Definition: VectorI.H:75
label face() const
Return the face.
Definition: tetIndicesI.H:40
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:1051
void cmptMag(FieldField< Field, Type > &cf, const FieldField< Field, Type > &f)
static const char nl
Definition: Ostream.H:265
label findInside(const point &) const
Find shape containing point. Only implemented for certain.
void updateMesh()
Correct polyBoundaryMesh after topology update.
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
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
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:32
static IOobject points0IO(const polyMesh &mesh)
Return IO object for points0.
Definition: polyMesh.C:1220
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:842
virtual bool upToDatePoints(const regIOobject &io) const
Return true if io is up-to-date with points.
Definition: polyMesh.C:1144
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:882
writeOption writeOpt() const
Definition: IOobject.H:355
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
const fileName & instance() const
Definition: IOobject.H:378
label timeIndex() const
Return current time index.
Definition: TimeStateI.H:35
virtual const pointField & oldCellCentres() const
Return old points for mesh motion.
Definition: polyMesh.C:1199
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:56
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:796
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:1480
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1606
virtual ~polyMesh()
Destructor.
Definition: polyMesh.C:774
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:55
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:87
dimensioned< scalar > mag(const dimensioned< Type > &)
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:497
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:98
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
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:1021
Registry of regIOobjects.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
readOption readOpt() const
Definition: IOobject.H:345
const word & headerClassName() const
Return name of the class name read from header.
Definition: IOobject.H:301
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
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.