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-2023 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 "meshObjects.H"
37 #include "pointMesh.H"
38 
39 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
40 
41 namespace Foam
42 {
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 
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  (
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(),
367  IOobject::AUTO_WRITE
368  ),
369  move(points)
370  ),
371  faces_
372  (
373  IOobject
374  (
375  "faces",
376  instance(),
377  meshSubDir,
378  *this,
379  io.readOpt(),
380  IOobject::AUTO_WRITE
381  ),
382  move(faces)
383  ),
384  owner_
385  (
386  IOobject
387  (
388  "owner",
389  instance(),
390  meshSubDir,
391  *this,
392  io.readOpt(),
393  IOobject::AUTO_WRITE
394  ),
395  move(owner)
396  ),
397  neighbour_
398  (
399  IOobject
400  (
401  "neighbour",
402  instance(),
403  meshSubDir,
404  *this,
405  io.readOpt(),
406  IOobject::AUTO_WRITE
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(),
420  IOobject::AUTO_WRITE
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(),
440  IOobject::NO_WRITE
441  ),
442  *this,
443  PtrList<pointZone>()
444  ),
445  faceZones_
446  (
447  IOobject
448  (
449  "faceZones",
450  instance(),
451  meshSubDir,
452  *this,
453  io.readOpt(),
454  IOobject::NO_WRITE
455  ),
456  *this,
457  PtrList<faceZone>()
458  ),
459  cellZones_
460  (
461  IOobject
462  (
463  "cellZones",
464  instance(),
465  meshSubDir,
466  *this,
467  io.readOpt(),
468  IOobject::NO_WRITE
469  ),
470  *this,
471  PtrList<cellZone>()
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,
520  IOobject::NO_READ,
521  IOobject::AUTO_WRITE
522  ),
523  move(points)
524  ),
525  faces_
526  (
527  IOobject
528  (
529  "faces",
530  instance(),
531  meshSubDir,
532  *this,
533  IOobject::NO_READ,
534  IOobject::AUTO_WRITE
535  ),
536  move(faces)
537  ),
538  owner_
539  (
540  IOobject
541  (
542  "owner",
543  instance(),
544  meshSubDir,
545  *this,
546  IOobject::NO_READ,
547  IOobject::AUTO_WRITE
548  ),
549  label(0)
550  ),
551  neighbour_
552  (
553  IOobject
554  (
555  "neighbour",
556  instance(),
557  meshSubDir,
558  *this,
559  IOobject::NO_READ,
560  IOobject::AUTO_WRITE
561  ),
562  label(0)
563  ),
564  clearedPrimitives_(false),
565  boundary_
566  (
567  IOobject
568  (
569  "boundary",
570  instance(),
571  meshSubDir,
572  *this,
573  IOobject::NO_READ,
574  IOobject::AUTO_WRITE
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,
593  IOobject::NO_READ,
594  IOobject::NO_WRITE
595  ),
596  *this,
597  0
598  ),
599  faceZones_
600  (
601  IOobject
602  (
603  "faceZones",
604  instance(),
605  meshSubDir,
606  *this,
607  IOobject::NO_READ,
608  IOobject::NO_WRITE
609  ),
610  *this,
611  0
612  ),
613  cellZones_
614  (
615  IOobject
616  (
617  "cellZones",
618  instance(),
619  meshSubDir,
620  *this,
621  IOobject::NO_READ,
622  IOobject::NO_WRITE
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  objectRegistry(move(mesh)),
674  primitiveMesh(move(mesh)),
675  points_(move(mesh.points_)),
676  faces_(move(mesh.faces_)),
677  owner_(move(mesh.owner_)),
678  neighbour_(move(mesh.neighbour_)),
679  clearedPrimitives_(mesh.clearedPrimitives_),
680  boundary_(move(mesh.boundary_)),
681  bounds_(move(mesh.bounds_)),
682  comm_(mesh.comm_),
683  geometricD_(mesh.geometricD_),
684  solutionD_(mesh.solutionD_),
685  tetBasePtIsPtr_(move(mesh.tetBasePtIsPtr_)),
686  cellTreePtr_(move(mesh.cellTreePtr_)),
687  pointZones_(move(mesh.pointZones_)),
688  faceZones_(move(mesh.faceZones_)),
689  cellZones_(move(mesh.cellZones_)),
690  globalMeshDataPtr_(move(mesh.globalMeshDataPtr_)),
691  curMotionTimeIndex_(mesh.curMotionTimeIndex_),
692  oldPointsPtr_(move(mesh.oldPointsPtr_)),
693  oldCellCentresPtr_(move(mesh.oldCellCentresPtr_)),
694  storeOldCellCentres_(mesh.storeOldCellCentres_),
695  moving_(mesh.moving_),
696  topoChanged_(mesh.topoChanged_)
697 {}
698 
699 
701 (
702  pointField&& points,
703  faceList&& faces,
704  labelList&& owner,
705  labelList&& neighbour,
706  const labelList& patchSizes,
707  const labelList& patchStarts,
708  const bool validBoundary
709 )
710 {
711  // Clear addressing. Keep geometric props and updateable props for mapping.
712  clearAddressing(true);
713 
714  // Take over new primitive data.
715  // Optimised to avoid overwriting data at all
716  if (notNull(points))
717  {
718  points_ = move(points);
719  bounds_ = boundBox(points_, validBoundary);
720  }
721 
722  if (notNull(faces))
723  {
724  faces_ = move(faces);
725  }
726 
727  if (notNull(owner))
728  {
729  owner_ = move(owner);
730  }
731 
732  if (notNull(neighbour))
733  {
734  neighbour_ = move(neighbour);
735  }
736 
737 
738  // Reset patch sizes and starts
739  forAll(boundary_, patchi)
740  {
741  boundary_[patchi] = polyPatch
742  (
743  boundary_[patchi],
744  boundary_,
745  patchi,
746  patchSizes[patchi],
747  patchStarts[patchi]
748  );
749  }
750 
751 
752  // Flags the mesh files as being changed
753  setInstance(time().name());
754 
755  // Check if the faces and cells are valid
756  forAll(faces_, facei)
757  {
758  const face& curFace = faces_[facei];
759 
760  if (min(curFace) < 0 || max(curFace) > points_.size())
761  {
763  << "Face " << facei << " contains vertex labels out of range: "
764  << curFace << " Max point index = " << points_.size()
765  << abort(FatalError);
766  }
767  }
768 
769 
770  // Set the primitive mesh from the owner_, neighbour_.
771  // Works out from patch end where the active faces stop.
772  initMesh();
773 
774 
775  if (validBoundary)
776  {
777  // Note that we assume that all the patches stay the same and are
778  // correct etc. so we can already use the patches to do
779  // processor-processor comms.
780 
781  // Calculate topology for the patches (processor-processor comms etc.)
782  boundary_.topoChange();
783 
784  // Calculate the geometry for the patches (transformation tensors etc.)
785  boundary_.calcGeometry();
786 
787  // Warn if global empty mesh
788  if
789  (
790  (returnReduce(nPoints(), sumOp<label>()) == 0)
791  || (returnReduce(nCells(), sumOp<label>()) == 0)
792  )
793  {
795  << "no points or no cells in mesh"
796  << exit(FatalError);
797  }
798  }
799 }
800 
801 
803 {
804  // Clear addressing. Keep geometric and updatable properties for mapping.
805  clearAddressing(true);
806  otherMesh.clearAddressing(true);
807 
808  // Swap the primitives
809  points_.swap(otherMesh.points_);
810  bounds_ = boundBox(points_, true);
811  faces_.swap(otherMesh.faces_);
812  owner_.swap(otherMesh.owner_);
813  neighbour_.swap(otherMesh.neighbour_);
814 
815  // Clear the boundary data
816  boundary_.clearGeom();
817  boundary_.clearAddressing();
818  otherMesh.boundary_.clearGeom();
819  otherMesh.boundary_.clearAddressing();
820 
821  // Swap the boundaries
822  auto updatePatches = []
823  (
824  const polyPatchList& otherPatches,
825  polyBoundaryMesh& boundaryMesh
826  )
827  {
828  boundaryMesh.resize(otherPatches.size());
829 
830  forAll(otherPatches, otherPatchi)
831  {
832  // Clone processor patches, as the decomposition may be different
833  // in the other mesh. Just change the size and start of other
834  // patches.
835 
836  if (isA<processorPolyPatch>(otherPatches[otherPatchi]))
837  {
838  boundaryMesh.set
839  (
840  otherPatchi,
841  otherPatches[otherPatchi].clone(boundaryMesh)
842  );
843  }
844  else
845  {
846  boundaryMesh[otherPatchi] = polyPatch
847  (
848  boundaryMesh[otherPatchi],
849  boundaryMesh,
850  otherPatchi,
851  otherPatches[otherPatchi].size(),
852  otherPatches[otherPatchi].start()
853  );
854  }
855  }
856  };
857 
858  {
859  const polyPatchList patches
860  (
861  boundary_,
862  boundary_
863  );
864  const polyPatchList otherPatches
865  (
866  otherMesh.boundary_,
867  otherMesh.boundary_
868  );
869 
870  updatePatches(otherPatches, boundary_);
871  updatePatches(patches, otherMesh.boundary_);
872  }
873 
874  // Parallel data depends on the patch ordering so force recalculation
875  globalMeshDataPtr_.clear();
876  otherMesh.globalMeshDataPtr_.clear();
877 
878  // Flags the mesh files as being changed
879  setInstance(time().name());
880  otherMesh.setInstance(time().name());
881 
882  // Check if the faces and cells are valid
883  auto checkFaces = [](const polyMesh& mesh)
884  {
885  forAll(mesh.faces_, facei)
886  {
887  const face& curFace = mesh.faces_[facei];
888 
889  if (min(curFace) < 0 || max(curFace) > mesh.points_.size())
890  {
892  << "Face " << facei << " contains vertex labels out of "
893  << "range: " << curFace << " Max point index = "
894  << mesh.points_.size() << abort(FatalError);
895  }
896  }
897  };
898 
899  checkFaces(*this);
900  checkFaces(otherMesh);
901 
902  // Set the primitive mesh from the owner_, neighbour_.
903  // Works out from patch end where the active faces stop.
904  initMesh();
905  otherMesh.initMesh();
906 
907  // Calculate topology for the patches (processor-processor comms etc.)
908  boundary_.topoChange();
909  otherMesh.boundary_.topoChange();
910 
911  // Calculate the geometry for the patches (transformation tensors etc.)
912  boundary_.calcGeometry();
913  otherMesh.boundary_.calcGeometry();
914 
915  // Update the optional pointMesh with respect to the updated polyMesh
916  if (foundObject<pointMesh>(pointMesh::typeName))
917  {
918  pointMesh::New(*this).reset();
919  }
920  if (otherMesh.foundObject<pointMesh>(pointMesh::typeName))
921  {
922  pointMesh::New(*this).reset();
923  }
924 
925  // Update point zones
926  if (pointZones_.size() == otherMesh.pointZones_.size())
927  {
928  pointZones_.clearAddressing();
929  otherMesh.pointZones_.clearAddressing();
930 
931  forAll(pointZones_, i)
932  {
933  pointZones_[i].swap(otherMesh.pointZones_[i]);
934  }
935  }
936  else
937  {
939  << "Number of pointZones in other mesh = "
940  << otherMesh.pointZones_.size()
941  << " is not the same as in the mesh = "
942  << pointZones_.size()
943  << exit(FatalError);
944  }
945 
946  // Update face zones
947  if (faceZones_.size() == otherMesh.faceZones_.size())
948  {
949  faceZones_.clearAddressing();
950  otherMesh.faceZones_.clearAddressing();
951 
952  forAll(faceZones_, i)
953  {
954  faceZones_[i].swap(otherMesh.faceZones_[i]);
955  }
956  }
957  else
958  {
960  << "Number of faceZones in other mesh = "
961  << otherMesh.faceZones_.size()
962  << " is not the same as in the mesh = "
963  << faceZones_.size()
964  << exit(FatalError);
965  }
966 
967  // Update cell zones
968  if (cellZones_.size() == otherMesh.cellZones_.size())
969  {
970  cellZones_.clearAddressing();
971  otherMesh.cellZones_.clearAddressing();
972 
973  forAll(cellZones_, i)
974  {
975  cellZones_[i].swap(otherMesh.cellZones_[i]);
976  }
977  }
978  else
979  {
981  << "Number of cellZones in other mesh = "
982  << otherMesh.cellZones_.size()
983  << " is not the same as in the mesh = "
984  << cellZones_.size()
985  << exit(FatalError);
986  }
987 }
988 
989 
990 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
991 
993 {
994  clearOut();
995  resetMotion();
996 }
997 
998 
999 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
1000 
1002 {
1003  if (objectRegistry::dbDir() == defaultRegion)
1004  {
1005  return parent().dbDir();
1006  }
1007  else
1008  {
1009  return objectRegistry::dbDir();
1010  }
1011 }
1012 
1013 
1015 {
1016  return dbDir()/meshSubDir;
1017 }
1018 
1019 
1021 {
1022  return points_.instance();
1023 }
1024 
1025 
1027 {
1028  return faces_.instance();
1029 }
1030 
1031 
1033 {
1034  return points_.writeOpt();
1035 }
1036 
1037 
1039 {
1040  return faces_.writeOpt();
1041 }
1042 
1043 
1045 {
1046  if (geometricD_.x() == 0)
1047  {
1048  calcDirections();
1049  }
1050 
1051  return geometricD_;
1052 }
1053 
1054 
1056 {
1057  return cmptSum(geometricD() + Vector<label>::one)/2;
1058 }
1059 
1060 
1062 {
1063  if (solutionD_.x() == 0)
1064  {
1065  calcDirections();
1066  }
1067 
1068  return solutionD_;
1069 }
1070 
1071 
1073 {
1074  return cmptSum(solutionD() + Vector<label>::one)/2;
1075 }
1076 
1077 
1079 {
1080  if (tetBasePtIsPtr_.empty())
1081  {
1082  if (debug)
1083  {
1085  << "Forcing storage of base points."
1086  << endl;
1087  }
1088 
1089  tetBasePtIsPtr_.reset
1090  (
1091  new labelIOList
1092  (
1093  IOobject
1094  (
1095  "tetBasePtIs",
1096  instance(),
1097  meshSubDir,
1098  *this,
1101  ),
1103  )
1104  );
1105  }
1106 
1107  return tetBasePtIsPtr_();
1108 }
1109 
1110 
1113 {
1114  if (cellTreePtr_.empty())
1115  {
1116  cellTreePtr_.reset
1117  (
1119  (
1120  treeDataCell
1121  (
1122  false, // not cache bb
1123  *this,
1124  CELL_TETS // use tet-decomposition for any inside test
1125  ),
1126  treeBoundBox(points()).extend(1e-4),
1127  8, // maxLevel
1128  10, // leafsize
1129  5.0 // duplicity
1130  )
1131  );
1132  }
1133 
1134  return cellTreePtr_();
1135 }
1136 
1137 
1139 (
1140  const List<polyPatch*>& p,
1141  const bool validBoundary
1142 )
1143 {
1144  if (boundaryMesh().size())
1145  {
1147  << "boundary already exists"
1148  << abort(FatalError);
1149  }
1150 
1151  // Reset valid directions
1152  geometricD_ = Zero;
1153  solutionD_ = Zero;
1154 
1155  boundary_.setSize(p.size());
1156 
1157  // Copy the patch pointers
1158  forAll(p, pI)
1159  {
1160  boundary_.set(pI, p[pI]);
1161  }
1162 
1163  // parallelData depends on the processorPatch ordering so force
1164  // recalculation. Problem: should really be done in removeBoundary but
1165  // there is some info in parallelData which might be interesting in between
1166  // removeBoundary and addPatches.
1167  globalMeshDataPtr_.clear();
1168 
1169  if (validBoundary)
1170  {
1171  // Calculate topology for the patches (processor-processor comms etc.)
1172  boundary_.topoChange();
1173 
1174  // Calculate the geometry for the patches (transformation tensors etc.)
1175  boundary_.calcGeometry();
1176 
1177  boundary_.checkDefinition();
1178  }
1179 }
1180 
1181 
1183 (
1184  const List<pointZone*>& pz,
1185  const List<faceZone*>& fz,
1186  const List<cellZone*>& cz
1187 )
1188 {
1189  if (pointZones().size() || faceZones().size() || cellZones().size())
1190  {
1192  << "point, face or cell zone already exists"
1193  << abort(FatalError);
1194  }
1195 
1196  // Point zones
1197  if (pz.size())
1198  {
1199  pointZones_.setSize(pz.size());
1200 
1201  // Copy the zone pointers
1202  forAll(pz, pI)
1203  {
1204  pointZones_.set(pI, pz[pI]);
1205  }
1206 
1207  pointZones_.writeOpt() = IOobject::AUTO_WRITE;
1208  }
1209 
1210  // Face zones
1211  if (fz.size())
1212  {
1213  faceZones_.setSize(fz.size());
1214 
1215  // Copy the zone pointers
1216  forAll(fz, fI)
1217  {
1218  faceZones_.set(fI, fz[fI]);
1219  }
1220 
1221  faceZones_.writeOpt() = IOobject::AUTO_WRITE;
1222  }
1223 
1224  // Cell zones
1225  if (cz.size())
1226  {
1227  cellZones_.setSize(cz.size());
1228 
1229  // Copy the zone pointers
1230  forAll(cz, cI)
1231  {
1232  cellZones_.set(cI, cz[cI]);
1233  }
1234 
1235  cellZones_.writeOpt() = IOobject::AUTO_WRITE;
1236  }
1237 }
1238 
1239 
1241 (
1242  const labelUList& newToOld,
1243  const bool validBoundary
1244 )
1245 {
1246  // Clear local fields and e.g. polyMesh parallelInfo. Do not clearGeom
1247  // so we keep PatchMeshObjects intact.
1248  boundary_.clearGeom();
1249  clearAddressing(true);
1250  // Clear all but PatchMeshObjects
1252  <
1253  polyMesh,
1256  >
1257  (
1258  *this
1259  );
1261  <
1262  pointMesh,
1265  >
1266  (
1267  *this
1268  );
1269 
1270  boundary_.reorderPatches(newToOld, validBoundary);
1271 
1272  // Warn mesh objects
1273  meshObjects::reorderPatches<polyMesh>(*this, newToOld, validBoundary);
1274  meshObjects::reorderPatches<pointMesh>(*this, newToOld, validBoundary);
1275 }
1276 
1277 
1279 (
1280  const label insertPatchi,
1281  const polyPatch& patch,
1282  const dictionary& patchFieldDict,
1283  const word& defaultPatchFieldType,
1284  const bool validBoundary
1285 )
1286 {
1287  const label sz = boundary_.size();
1288 
1289  label startFacei = nFaces();
1290  if (insertPatchi < sz)
1291  {
1292  startFacei = boundary_[insertPatchi].start();
1293  }
1294 
1295  // Create reordering list
1296  // patches before insert position stay as is
1297  // patches after insert position move one up
1298  labelList newToOld(boundary_.size()+1);
1299  for (label i = 0; i < insertPatchi; i++)
1300  {
1301  newToOld[i] = i;
1302  }
1303  for (label i = insertPatchi; i < sz; i++)
1304  {
1305  newToOld[i+1] = i;
1306  }
1307  newToOld[insertPatchi] = -1;
1308 
1309  reorderPatches(newToOld, false);
1310 
1311  // Clear local fields and e.g. polyMesh parallelInfo.
1312  //clearGeom(); // would clear out pointMesh as well
1313  boundary_.clearGeom();
1314  clearAddressing(true);
1315 
1316  // Clear all but PatchMeshObjects
1318  <
1319  polyMesh,
1322  >
1323  (
1324  *this
1325  );
1327  <
1328  pointMesh,
1331  >
1332  (
1333  *this
1334  );
1335 
1336 
1337  // Insert polyPatch
1338  boundary_.set
1339  (
1340  insertPatchi,
1341  patch.clone
1342  (
1343  boundary_,
1344  insertPatchi, // index
1345  0, // size
1346  startFacei // start
1347  )
1348  );
1349 
1350  if (validBoundary)
1351  {
1352  boundary_.topoChange();
1353  }
1354 
1355  // Warn mesh objects
1356  meshObjects::addPatch<polyMesh>(*this, insertPatchi);
1357  meshObjects::addPatch<pointMesh>(*this, insertPatchi);
1358 }
1359 
1360 
1362 {
1363  if (clearedPrimitives_)
1364  {
1366  << "points deallocated"
1367  << abort(FatalError);
1368  }
1369 
1370  return points_;
1371 }
1372 
1373 
1375 {
1376  if (clearedPrimitives_)
1377  {
1379  << "faces deallocated"
1380  << abort(FatalError);
1381  }
1382 
1383  return faces_;
1384 }
1385 
1386 
1388 {
1389  return owner_;
1390 }
1391 
1392 
1394 {
1395  return neighbour_;
1396 }
1397 
1398 
1400 {
1401  if (!moving_)
1402  {
1403  return points_;
1404  }
1405 
1406  if (oldPointsPtr_.empty())
1407  {
1409  << "Old points have not been stored"
1410  << exit(FatalError);
1411  }
1412 
1413  return oldPointsPtr_();
1414 }
1415 
1416 
1418 {
1419  storeOldCellCentres_ = true;
1420 
1421  if (!moving_)
1422  {
1423  return cellCentres();
1424  }
1425 
1426  if (oldCellCentresPtr_.empty())
1427  {
1429  << "Old cell centres have not been stored"
1430  << exit(FatalError);
1431  }
1432 
1433  return oldCellCentresPtr_();
1434 }
1435 
1436 
1438 {
1439  return io.upToDate(points_);
1440 }
1441 
1442 
1444 {
1445  io.eventNo() = points_.eventNo()+1;
1446 }
1447 
1448 
1450 {
1451  if (debug)
1452  {
1454  << "Set points for time " << time().value()
1455  << " index " << time().timeIndex() << endl;
1456  }
1457 
1459 
1460  points_ = newPoints;
1461 
1462  setPointsInstance(time().name());
1463 
1464  // Adjust parallel shared points
1465  if (globalMeshDataPtr_.valid())
1466  {
1467  globalMeshDataPtr_().movePoints(points_);
1468  }
1469 
1470  // Force recalculation of all geometric data with new points
1471 
1472  bounds_ = boundBox(points_);
1473  boundary_.movePoints(points_);
1474 
1475  pointZones_.movePoints(points_);
1476  faceZones_.movePoints(points_);
1477  cellZones_.movePoints(points_);
1478 
1479  // Cell tree might become invalid
1480  cellTreePtr_.clear();
1481 
1482  // Reset valid directions (could change with rotation)
1483  geometricD_ = Zero;
1484  solutionD_ = Zero;
1485 
1486  meshObjects::movePoints<polyMesh>(*this);
1487  meshObjects::movePoints<pointMesh>(*this);
1488 }
1489 
1490 
1492 (
1493  const pointField& newPoints
1494 )
1495 {
1496  if (debug)
1497  {
1499  << "Moving points for time " << time().value()
1500  << " index " << time().timeIndex() << endl;
1501  }
1502 
1503  // Pick up old points and cell centres
1504  if (curMotionTimeIndex_ != time().timeIndex())
1505  {
1506  oldPointsPtr_.clear();
1507  oldPointsPtr_.reset(new pointField(points_));
1508  if (storeOldCellCentres_)
1509  {
1510  oldCellCentresPtr_.clear();
1511  oldCellCentresPtr_.reset(new pointField(cellCentres()));
1512  }
1513  curMotionTimeIndex_ = time().timeIndex();
1514  }
1515 
1516  points_ = newPoints;
1517 
1518  setPointsInstance(time().name());
1519 
1521  (
1522  points_,
1523  oldPoints()
1524  );
1525 
1526  // Adjust parallel shared points
1527  if (globalMeshDataPtr_.valid())
1528  {
1529  globalMeshDataPtr_().movePoints(points_);
1530  }
1531 
1532  // Force recalculation of all geometric data with new points
1533 
1534  bounds_ = boundBox(points_);
1535  boundary_.movePoints(points_);
1536 
1537  pointZones_.movePoints(points_);
1538  faceZones_.movePoints(points_);
1539  cellZones_.movePoints(points_);
1540 
1541  // Cell tree might become invalid
1542  cellTreePtr_.clear();
1543 
1544  // Reset valid directions (could change with rotation)
1545  geometricD_ = Zero;
1546  solutionD_ = Zero;
1547 
1548  meshObjects::movePoints<polyMesh>(*this);
1549  meshObjects::movePoints<pointMesh>(*this);
1550 
1551  return sweptVols;
1552 }
1553 
1554 
1556 {
1557  curMotionTimeIndex_ = -1;
1558  oldPointsPtr_.clear();
1559  oldCellCentresPtr_.clear();
1560 }
1561 
1562 
1564 {
1565  if (globalMeshDataPtr_.empty())
1566  {
1567  if (debug)
1568  {
1569  Pout<< "polyMesh::globalData() const : "
1570  << "Constructing parallelData from processor topology"
1571  << endl;
1572  }
1573  // Construct globalMeshData using processorPatch information only.
1574  globalMeshDataPtr_.reset(new globalMeshData(*this));
1575  }
1576 
1577  return globalMeshDataPtr_();
1578 }
1579 
1580 
1582 {
1583  return comm_;
1584 }
1585 
1586 
1588 {
1589  return comm_;
1590 }
1591 
1592 
1593 void Foam::polyMesh::removeFiles(const fileName& instanceDir) const
1594 {
1595  fileName meshFilesPath = thisDb().time().path()/instanceDir/meshDir();
1596 
1597  rm(meshFilesPath/"points");
1598  rm(meshFilesPath/"faces");
1599  rm(meshFilesPath/"owner");
1600  rm(meshFilesPath/"neighbour");
1601  rm(meshFilesPath/"cells");
1602  rm(meshFilesPath/"boundary");
1603  rm(meshFilesPath/"pointZones");
1604  rm(meshFilesPath/"faceZones");
1605  rm(meshFilesPath/"cellZones");
1606  rm(meshFilesPath/"meshModifiers");
1607  rm(meshFilesPath/"parallelData");
1608 
1609  // remove subdirectories
1610  if (isDir(meshFilesPath/"sets"))
1611  {
1612  rmDir(meshFilesPath/"sets");
1613  }
1614 }
1615 
1616 
1618 {
1619  removeFiles(instance());
1620 }
1621 
1622 
1624 (
1625  const point& p,
1626  label& celli,
1627  label& tetFacei,
1628  label& tetPti
1629 ) const
1630 {
1631  celli = -1;
1632  tetFacei = -1;
1633  tetPti = -1;
1634 
1635  const indexedOctree<treeDataCell>& tree = cellTree();
1636 
1637  // Find point inside cell
1638  celli = tree.findInside(p);
1639 
1640  if (celli != -1)
1641  {
1642  // Check the nearest cell to see if the point is inside.
1643  findTetFacePt(celli, p, tetFacei, tetPti);
1644  }
1645 }
1646 
1647 
1649 (
1650  const label celli,
1651  const point& p,
1652  label& tetFacei,
1653  label& tetPti
1654 ) const
1655 {
1656  const polyMesh& mesh = *this;
1657 
1658  tetIndices tet(polyMeshTetDecomposition::findTet(mesh, celli, p));
1659  tetFacei = tet.face();
1660  tetPti = tet.tetPt();
1661 }
1662 
1663 
1665 (
1666  const point& p,
1667  label celli,
1668  const cellDecomposition decompMode
1669 ) const
1670 {
1671  switch (decompMode)
1672  {
1673  case FACE_PLANES:
1674  {
1675  return primitiveMesh::pointInCell(p, celli);
1676  }
1677  break;
1678 
1679  case FACE_CENTRE_TRIS:
1680  {
1681  // only test that point is on inside of plane defined by cell face
1682  // triangles
1683  const cell& cFaces = cells()[celli];
1684 
1685  forAll(cFaces, cFacei)
1686  {
1687  label facei = cFaces[cFacei];
1688  const face& f = faces_[facei];
1689  const point& fc = faceCentres()[facei];
1690  bool isOwn = (owner_[facei] == celli);
1691 
1692  forAll(f, fp)
1693  {
1694  label pointi;
1695  label nextPointi;
1696 
1697  if (isOwn)
1698  {
1699  pointi = f[fp];
1700  nextPointi = f.nextLabel(fp);
1701  }
1702  else
1703  {
1704  pointi = f.nextLabel(fp);
1705  nextPointi = f[fp];
1706  }
1707 
1708  triPointRef faceTri
1709  (
1710  points()[pointi],
1711  points()[nextPointi],
1712  fc
1713  );
1714 
1715  vector proj = p - faceTri.centre();
1716 
1717  if ((faceTri.area() & proj) > 0)
1718  {
1719  return false;
1720  }
1721  }
1722  }
1723  return true;
1724  }
1725  break;
1726 
1727  case FACE_DIAG_TRIS:
1728  {
1729  // only test that point is on inside of plane defined by cell face
1730  // triangles
1731  const cell& cFaces = cells()[celli];
1732 
1733  forAll(cFaces, cFacei)
1734  {
1735  label facei = cFaces[cFacei];
1736  const face& f = faces_[facei];
1737 
1738  for (label tetPti = 1; tetPti < f.size() - 1; tetPti++)
1739  {
1740  // Get tetIndices of face triangle
1741  tetIndices faceTetIs(celli, facei, tetPti);
1742 
1743  triPointRef faceTri = faceTetIs.faceTri(*this);
1744 
1745  vector proj = p - faceTri.centre();
1746 
1747  if ((faceTri.area() & proj) > 0)
1748  {
1749  return false;
1750  }
1751  }
1752  }
1753 
1754  return true;
1755  }
1756  break;
1757 
1758  case CELL_TETS:
1759  {
1760  label tetFacei;
1761  label tetPti;
1762 
1763  findTetFacePt(celli, p, tetFacei, tetPti);
1764 
1765  return tetFacei != -1;
1766  }
1767  break;
1768  }
1769 
1770  return false;
1771 }
1772 
1773 
1775 (
1776  const point& p,
1777  const cellDecomposition decompMode
1778 ) const
1779 {
1780  if
1781  (
1782  Pstream::parRun()
1783  && (decompMode == FACE_DIAG_TRIS || decompMode == CELL_TETS)
1784  )
1785  {
1786  // Force construction of face-diagonal decomposition before testing
1787  // for zero cells.
1788  //
1789  // If parallel running a local domain might have zero cells so never
1790  // construct the face-diagonal decomposition which uses parallel
1791  // transfers.
1792  (void)tetBasePtIs();
1793  }
1794 
1795  if (nCells() == 0)
1796  {
1797  return -1;
1798  }
1799 
1800  if (decompMode == CELL_TETS)
1801  {
1802  // Advanced search method utilising an octree
1803  // and tet-decomposition of the cells
1804 
1805  label celli;
1806  label tetFacei;
1807  label tetPti;
1808 
1809  findCellFacePt(p, celli, tetFacei, tetPti);
1810 
1811  return celli;
1812  }
1813  else
1814  {
1815  // Approximate search avoiding the construction of an octree
1816  // and cell decomposition
1817 
1818  // Find the nearest cell centre to this location
1819  label celli = findNearestCell(p);
1820 
1821  // If point is in the nearest cell return
1822  if (pointInCell(p, celli, decompMode))
1823  {
1824  return celli;
1825  }
1826  else
1827  {
1828  // Point is not in the nearest cell so search all cells
1829 
1830  for (label celli = 0; celli < nCells(); celli++)
1831  {
1832  if (pointInCell(p, celli, decompMode))
1833  {
1834  return celli;
1835  }
1836  }
1837 
1838  return -1;
1839  }
1840  }
1841 }
1842 
1843 
1844 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
A List of objects of type <Type> with automated input and output using a compact storage....
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const word & headerClassName() const
Return name of the class name read from header.
Definition: IOobject.H:316
writeOption
Enumeration defining the write options.
Definition: IOobject.H:126
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void clearAddressing()
Clear addressing.
Definition: MeshZones.C:387
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)=0
Reordered/removed trailing patches. If validBoundary call is parallel.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
MeshObject types:
Definition: MeshObjects.H:60
Inter-processor communications stream.
Definition: UPstream.H:59
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
virtual void topoChange(const polyTopoChangeMap &map)=0
Update topology using the given map.
static const direction nComponents
Number of components in this vector space.
Definition: VectorSpace.H:99
Templated 3D Vector derived from VectorSpace adding construction from 3 components,...
Definition: Vector.H:60
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
A subset of mesh cells.
Definition: cellZone.H:64
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:68
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
A class for handling file names.
Definition: fileName.H:82
fileName path() const
Return directory path name (part before last /)
Definition: fileName.C:265
Various mesh related information for a parallel run. Upon construction, constructs all info using par...
label findInside(const point &) const
Find shape containing point. Only implemented for certain.
static void clearUpto(objectRegistry &)
Clear all meshObjects derived from FromType up to (but not including)
Registry of regIOobjects.
const Time & time() const
Return time.
virtual const fileName & dbDir() const
Local directory path of this objectRegistry relative to the time.
bool foundObject(const word &name) const
Is the named Type in registry.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:52
void reset()
Reset pointMesh with respect to the updated polyMesh.
Definition: pointMesh.C:188
A subset of mesh points. The labels of points in the zone can be obtained from the addressing() list.
Definition: pointZone.H:65
Foam::polyBoundaryMesh.
void topoChange()
Correct polyBoundaryMesh after topology update.
void clearGeom()
Clear geometry at this level and at patches.
void clearAddressing()
Clear addressing at this level and at patches.
static labelList findFaceBasePts(const polyMesh &mesh, scalar tol=minTetQuality, bool report=false)
Find a suitable base point for each face for decomposition.
static tetIndices findTet(const polyMesh &mesh, label cI, const point &pt)
Find the tet decomposition of the cell containing the given point.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual ~polyMesh()
Destructor.
Definition: polyMesh.C:992
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition: polyMesh.C:1624
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:1026
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:1055
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:269
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:1014
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1492
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:1279
cellDecomposition
Enumeration defining the decomposition of the cell for.
Definition: polyMesh.H:101
IOobject::writeOption facesWriteOpt() const
Return the points write option.
Definition: polyMesh.C:1038
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1374
virtual const fileName & dbDir() const
Override the objectRegistry dbDir for a single-region case.
Definition: polyMesh.C:1001
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1775
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:701
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1387
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:1078
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:1649
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1563
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:266
bool pointInCell(const point &p, label celli, const cellDecomposition=CELL_TETS) const
Test if point p is in the celli.
Definition: polyMesh.C:1665
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:1072
virtual const pointField & oldPoints() const
Return old points for mesh motion.
Definition: polyMesh.C:1399
void swap(polyMesh &)
Swap mesh.
Definition: polyMesh.C:802
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:1020
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1581
polyMesh(const IOobject &io)
Construct from IOobject.
Definition: polyMesh.C:163
void clearAddressing(const bool isMeshUpdate=false)
Clear addressing.
Definition: polyMeshClear.C:77
void resetMotion() const
Reset motion.
Definition: polyMesh.C:1555
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:403
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1393
void addPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:1139
bool upToDatePoints(const regIOobject &io) const
Return true if io is up-to-date with points.
Definition: polyMesh.C:1437
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Add mesh zones.
Definition: polyMesh.C:1183
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1361
void setUpToDatePoints(regIOobject &io) const
Set io to be up-to-date with points.
Definition: polyMesh.C:1443
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition: polyMesh.C:1112
IOobject::writeOption pointsWriteOpt() const
Return the points write option.
Definition: polyMesh.C:1032
void removeFiles() const
Remove all files from mesh instance()
Definition: polyMesh.C:1617
virtual const pointField & oldCellCentres() const
Return old cell centres for mesh motion.
Definition: polyMesh.C:1417
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorder and trim existing patches. If validBoundary the new.
Definition: polyMesh.C:1241
void setInstance(const fileName &)
Set the instance for mesh files.
Definition: polyMeshIO.C:78
virtual void setPoints(const pointField &)
Reset the points.
Definition: polyMesh.C:1449
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:1061
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:1044
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
virtual autoPtr< polyPatch > clone(const polyBoundaryMesh &bm) const
Construct and return a clone, resetting the boundary mesh.
Definition: polyPatch.H:215
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:75
label nCells() const
label nPoints() const
bool pointInCell(const point &p, label celli) const
Return true if the point is in the cell.
void clearGeom()
Clear geometry.
tmp< scalarField > movePoints(const pointField &p, const pointField &oldP)
Move points, returns volumes swept by faces in motion.
const vectorField & faceAreas() const
const cellList & cells() const
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:55
bool upToDate(const regIOobject &) const
Return true if up-to-date with respect to given object.
Definition: regIOobject.C:329
virtual bool write(const bool write=true) const
Write using setting from DB.
label eventNo() const
Event number at last update.
Definition: regIOobjectI.H:89
Storage and named access for the indices of a tet which is part of the decomposition of a cell.
Definition: tetIndices.H:82
label face() const
Return the face.
Definition: tetIndicesI.H:43
label tetPt() const
Return the characterising tetPtI.
Definition: tetIndicesI.H:55
triPointRef faceTri(const polyMesh &mesh) const
Return the geometry corresponding to the tri on the face for.
Definition: tetIndicesI.H:123
A class for managing temporary objects.
Definition: tmp.H:55
Standard boundBox + extra functionality for use in octree.
Definition: treeBoundBox.H:90
Encapsulation of data needed to search in/for cells. Used to find the cell containing a point (e....
Definition: treeDataCell.H:55
A triangle primitive used to calculate face areas and swept volumes.
Definition: triangle.H:80
vector area() const
Return vector area.
Definition: triangleI.H:96
Point centre() const
Return centre (centroid)
Definition: triangleI.H:89
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
label patchi
const pointField & points
label nPoints
const cellShapeList & cells
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
#define InfoInFunction
Report an information message using Foam::Info.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
const doubleScalar e
Definition: doubleScalar.H:105
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 rm(const fileName &)
Remove a file, returning true if successful otherwise false.
Definition: POSIX.C:1017
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
bool notNull(const T &t)
Return true if t is not a reference to the nullObject of type T.
Definition: nullObjectI.H:58
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
dimensioned< scalar > mag(const dimensioned< Type > &)
bool rmDir(const fileName &)
Remove a directory and its contents.
Definition: POSIX.C:1047
defineTypeNameAndDebug(combustionModel, 0)
T clone(const T &t)
Definition: List.H:55
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:42
dimensionSet cmptMag(const dimensionSet &)
Definition: dimensionSet.C:275
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
error FatalError
Cmpt cmptSum(const VectorSpace< Form, Cmpt, Ncmpts > &vs)
Definition: VectorSpaceI.H:500
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
uint8_t direction
Definition: direction.H:45
label timeIndex
Definition: getTimeIndex.H:4
labelList f(nPoints)
volScalarField & p