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 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
51 void Foam::polyMesh::calcDirections() const
52 {
53  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
54  {
55  solutionD_[cmpt] = 1;
56  }
57 
58  // Knock out empty and wedge directions. Note:they will be present on all
59  // domains.
60 
61  label nEmptyPatches = 0;
62  label nWedgePatches = 0;
63 
64  vector emptyDirVec = Zero;
65  vector wedgeDirVec = Zero;
66 
68  {
69  if (boundaryMesh()[patchi].size())
70  {
71  if (isA<emptyPolyPatch>(boundaryMesh()[patchi]))
72  {
73  nEmptyPatches++;
74  emptyDirVec += sum(cmptMag(boundaryMesh()[patchi].faceAreas()));
75  }
76  else if (isA<wedgePolyPatch>(boundaryMesh()[patchi]))
77  {
78  const wedgePolyPatch& wpp = refCast<const wedgePolyPatch>
79  (
81  );
82 
83  nWedgePatches++;
84  wedgeDirVec += cmptMag(wpp.centreNormal());
85  }
86  }
87  }
88 
89  reduce(nEmptyPatches, maxOp<label>());
90  reduce(nWedgePatches, maxOp<label>());
91 
92  if (nEmptyPatches)
93  {
94  reduce(emptyDirVec, sumOp<vector>());
95 
96  emptyDirVec /= mag(emptyDirVec);
97 
98  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
99  {
100  if (emptyDirVec[cmpt] > 1e-6)
101  {
102  solutionD_[cmpt] = -1;
103  }
104  else
105  {
106  solutionD_[cmpt] = 1;
107  }
108  }
109  }
110 
111 
112  // Knock out wedge directions
113 
114  geometricD_ = solutionD_;
115 
116  if (nWedgePatches)
117  {
118  reduce(wedgeDirVec, sumOp<vector>());
119 
120  wedgeDirVec /= mag(wedgeDirVec);
121 
122  for (direction cmpt=0; cmpt<vector::nComponents; cmpt++)
123  {
124  if (wedgeDirVec[cmpt] > 1e-6)
125  {
126  geometricD_[cmpt] = -1;
127  }
128  else
129  {
130  geometricD_[cmpt] = 1;
131  }
132  }
133  }
134 }
135 
136 
137 Foam::autoPtr<Foam::labelIOList> Foam::polyMesh::readTetBasePtIs() const
138 {
139  typeIOobject<labelIOList> io
140  (
141  "tetBasePtIs",
142  instance(),
143  meshSubDir,
144  *this,
147  );
148 
149  if (io.headerOk())
150  {
151  return autoPtr<labelIOList>(new labelIOList(io));
152  }
153  else
154  {
155  return autoPtr<labelIOList>(nullptr);
156  }
157 }
158 
159 
160 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
161 
163 :
164  objectRegistry(io),
165  primitiveMesh(),
166  points_
167  (
168  IOobject
169  (
170  "points",
171  time().findInstance(meshDir(), "points"),
172  meshSubDir,
173  *this,
174  IOobject::MUST_READ,
175  IOobject::NO_WRITE
176  )
177  ),
178  faces_
179  (
180  IOobject
181  (
182  "faces",
183  time().findInstance(meshDir(), "faces"),
184  meshSubDir,
185  *this,
186  IOobject::MUST_READ,
187  IOobject::NO_WRITE
188  )
189  ),
190  owner_
191  (
192  IOobject
193  (
194  "owner",
195  faces_.instance(),
196  meshSubDir,
197  *this,
198  IOobject::READ_IF_PRESENT,
199  IOobject::NO_WRITE
200  )
201  ),
202  neighbour_
203  (
204  IOobject
205  (
206  "neighbour",
207  faces_.instance(),
208  meshSubDir,
209  *this,
210  IOobject::READ_IF_PRESENT,
211  IOobject::NO_WRITE
212  )
213  ),
214  clearedPrimitives_(false),
215  boundary_
216  (
217  IOobject
218  (
219  "boundary",
220  time().findInstance(meshDir(), "boundary"),
221  meshSubDir,
222  *this,
223  IOobject::MUST_READ,
224  IOobject::NO_WRITE
225  ),
226  *this
227  ),
228  bounds_(points_),
229  comm_(UPstream::worldComm),
230  geometricD_(Zero),
231  solutionD_(Zero),
232  tetBasePtIsPtr_(readTetBasePtIs()),
233  cellTreePtr_(nullptr),
234  pointZones_
235  (
236  IOobject
237  (
238  "pointZones",
239  time().findInstance
240  (
241  meshDir(),
242  "pointZones",
243  IOobject::READ_IF_PRESENT
244  ),
245  meshSubDir,
246  *this,
247  IOobject::READ_IF_PRESENT,
248  IOobject::NO_WRITE
249  ),
250  *this
251  ),
252  faceZones_
253  (
254  IOobject
255  (
256  "faceZones",
257  time().findInstance
258  (
259  meshDir(),
260  "faceZones",
261  IOobject::READ_IF_PRESENT
262  ),
263  meshSubDir,
264  *this,
265  IOobject::READ_IF_PRESENT,
266  IOobject::NO_WRITE
267  ),
268  *this
269  ),
270  cellZones_
271  (
272  IOobject
273  (
274  "cellZones",
275  time().findInstance
276  (
277  meshDir(),
278  "cellZones",
279  IOobject::READ_IF_PRESENT
280  ),
281  meshSubDir,
282  *this,
283  IOobject::READ_IF_PRESENT,
284  IOobject::NO_WRITE
285  ),
286  *this
287  ),
288  globalMeshDataPtr_(nullptr),
289  curMotionTimeIndex_(-1),
290  oldPointsPtr_(nullptr),
291  oldCellCentresPtr_(nullptr),
292  storeOldCellCentres_(false),
293  moving_(false),
294  topoChanged_(false)
295 {
296  if (!owner_.headerClassName().empty())
297  {
298  initMesh();
299  }
300  else
301  {
302  cellCompactIOList cLst
303  (
304  IOobject
305  (
306  "cells",
307  time().findInstance(meshDir(), "cells"),
308  meshSubDir,
309  *this,
312  )
313  );
314 
315  // Set the primitive mesh
316  initMesh(cLst);
317 
318  owner_.write();
319  neighbour_.write();
320  }
321 
322  // Calculate topology for the patches (processor-processor comms etc.)
323  boundary_.topoChange();
324 
325  // Calculate the geometry for the patches (transformation tensors etc.)
326  boundary_.calcGeometry();
327 
328  // Warn if global empty mesh
329  if (returnReduce(nPoints(), sumOp<label>()) == 0)
330  {
332  << "no points in mesh" << endl;
333  }
334  if (returnReduce(nCells(), sumOp<label>()) == 0)
335  {
337  << "no cells in mesh" << endl;
338  }
339 
340  // Initialise demand-driven data
341  calcDirections();
342 }
343 
344 
346 (
347  const IOobject& io,
348  pointField&& points,
349  faceList&& faces,
350  labelList&& owner,
351  labelList&& neighbour,
352  const bool syncPar
353 )
354 :
355  objectRegistry(io),
356  primitiveMesh(),
357  points_
358  (
359  IOobject
360  (
361  "points",
362  instance(),
363  meshSubDir,
364  *this,
365  io.readOpt(),
366  IOobject::AUTO_WRITE
367  ),
368  move(points)
369  ),
370  faces_
371  (
372  IOobject
373  (
374  "faces",
375  instance(),
376  meshSubDir,
377  *this,
378  io.readOpt(),
379  IOobject::AUTO_WRITE
380  ),
381  move(faces)
382  ),
383  owner_
384  (
385  IOobject
386  (
387  "owner",
388  instance(),
389  meshSubDir,
390  *this,
391  io.readOpt(),
392  IOobject::AUTO_WRITE
393  ),
394  move(owner)
395  ),
396  neighbour_
397  (
398  IOobject
399  (
400  "neighbour",
401  instance(),
402  meshSubDir,
403  *this,
404  io.readOpt(),
405  IOobject::AUTO_WRITE
406  ),
407  move(neighbour)
408  ),
409  clearedPrimitives_(false),
410  boundary_
411  (
412  IOobject
413  (
414  "boundary",
415  instance(),
416  meshSubDir,
417  *this,
418  io.readOpt(),
419  IOobject::AUTO_WRITE
420  ),
421  *this,
422  polyPatchList()
423  ),
424  bounds_(points_, syncPar),
425  comm_(UPstream::worldComm),
426  geometricD_(Zero),
427  solutionD_(Zero),
428  tetBasePtIsPtr_(readTetBasePtIs()),
429  cellTreePtr_(nullptr),
430  pointZones_
431  (
432  IOobject
433  (
434  "pointZones",
435  instance(),
436  meshSubDir,
437  *this,
438  io.readOpt(),
439  IOobject::NO_WRITE
440  ),
441  *this,
442  PtrList<pointZone>()
443  ),
444  faceZones_
445  (
446  IOobject
447  (
448  "faceZones",
449  instance(),
450  meshSubDir,
451  *this,
452  io.readOpt(),
453  IOobject::NO_WRITE
454  ),
455  *this,
456  PtrList<faceZone>()
457  ),
458  cellZones_
459  (
460  IOobject
461  (
462  "cellZones",
463  instance(),
464  meshSubDir,
465  *this,
466  io.readOpt(),
467  IOobject::NO_WRITE
468  ),
469  *this,
470  PtrList<cellZone>()
471  ),
472  globalMeshDataPtr_(nullptr),
473  curMotionTimeIndex_(-1),
474  oldPointsPtr_(nullptr),
475  oldCellCentresPtr_(nullptr),
476  storeOldCellCentres_(false),
477  moving_(false),
478  topoChanged_(false)
479 {
480  // Check if the faces and cells are valid
481  forAll(faces_, facei)
482  {
483  const face& curFace = faces_[facei];
484 
485  if (min(curFace) < 0 || max(curFace) > points_.size())
486  {
488  << "Face " << facei << "contains vertex labels out of range: "
489  << curFace << " Max point index = " << points_.size()
490  << abort(FatalError);
491  }
492  }
493 
494  // Set the primitive mesh
495  initMesh();
496 }
497 
498 
499 
501 (
502  const IOobject& io,
503  pointField&& points,
504  faceList&& faces,
505  cellList&& cells,
506  const bool syncPar
507 )
508 :
509  objectRegistry(io),
510  primitiveMesh(),
511  points_
512  (
513  IOobject
514  (
515  "points",
516  instance(),
517  meshSubDir,
518  *this,
519  IOobject::NO_READ,
520  IOobject::AUTO_WRITE
521  ),
522  move(points)
523  ),
524  faces_
525  (
526  IOobject
527  (
528  "faces",
529  instance(),
530  meshSubDir,
531  *this,
532  IOobject::NO_READ,
533  IOobject::AUTO_WRITE
534  ),
535  move(faces)
536  ),
537  owner_
538  (
539  IOobject
540  (
541  "owner",
542  instance(),
543  meshSubDir,
544  *this,
545  IOobject::NO_READ,
546  IOobject::AUTO_WRITE
547  ),
548  label(0)
549  ),
550  neighbour_
551  (
552  IOobject
553  (
554  "neighbour",
555  instance(),
556  meshSubDir,
557  *this,
558  IOobject::NO_READ,
559  IOobject::AUTO_WRITE
560  ),
561  label(0)
562  ),
563  clearedPrimitives_(false),
564  boundary_
565  (
566  IOobject
567  (
568  "boundary",
569  instance(),
570  meshSubDir,
571  *this,
572  IOobject::NO_READ,
573  IOobject::AUTO_WRITE
574  ),
575  *this,
576  0
577  ),
578  bounds_(points_, syncPar),
579  comm_(UPstream::worldComm),
580  geometricD_(Zero),
581  solutionD_(Zero),
582  tetBasePtIsPtr_(readTetBasePtIs()),
583  cellTreePtr_(nullptr),
584  pointZones_
585  (
586  IOobject
587  (
588  "pointZones",
589  instance(),
590  meshSubDir,
591  *this,
592  IOobject::NO_READ,
593  IOobject::NO_WRITE
594  ),
595  *this,
596  0
597  ),
598  faceZones_
599  (
600  IOobject
601  (
602  "faceZones",
603  instance(),
604  meshSubDir,
605  *this,
606  IOobject::NO_READ,
607  IOobject::NO_WRITE
608  ),
609  *this,
610  0
611  ),
612  cellZones_
613  (
614  IOobject
615  (
616  "cellZones",
617  instance(),
618  meshSubDir,
619  *this,
620  IOobject::NO_READ,
621  IOobject::NO_WRITE
622  ),
623  *this,
624  0
625  ),
626  globalMeshDataPtr_(nullptr),
627  curMotionTimeIndex_(-1),
628  oldPointsPtr_(nullptr),
629  oldCellCentresPtr_(nullptr),
630  storeOldCellCentres_(false),
631  moving_(false),
632  topoChanged_(false)
633 {
634  // Check if faces are valid
635  forAll(faces_, facei)
636  {
637  const face& curFace = faces_[facei];
638 
639  if (min(curFace) < 0 || max(curFace) > points_.size())
640  {
642  << "Face " << facei << "contains vertex labels out of range: "
643  << curFace << " Max point index = " << points_.size()
644  << abort(FatalError);
645  }
646  }
647 
648  // transfer in cell list
649  cellList cLst(move(cells));
650 
651  // Check if cells are valid
652  forAll(cLst, celli)
653  {
654  const cell& curCell = cLst[celli];
655 
656  if (min(curCell) < 0 || max(curCell) > faces_.size())
657  {
659  << "Cell " << celli << "contains face labels out of range: "
660  << curCell << " Max face index = " << faces_.size()
661  << abort(FatalError);
662  }
663  }
664 
665  // Set the primitive mesh
666  initMesh(cLst);
667 }
668 
669 
671 :
672  objectRegistry(move(mesh)),
673  primitiveMesh(move(mesh)),
674  points_(move(mesh.points_)),
675  faces_(move(mesh.faces_)),
676  owner_(move(mesh.owner_)),
677  neighbour_(move(mesh.neighbour_)),
678  clearedPrimitives_(mesh.clearedPrimitives_),
679  boundary_(move(mesh.boundary_)),
680  bounds_(move(mesh.bounds_)),
681  comm_(mesh.comm_),
682  geometricD_(mesh.geometricD_),
683  solutionD_(mesh.solutionD_),
684  tetBasePtIsPtr_(move(mesh.tetBasePtIsPtr_)),
685  cellTreePtr_(move(mesh.cellTreePtr_)),
686  pointZones_(move(mesh.pointZones_)),
687  faceZones_(move(mesh.faceZones_)),
688  cellZones_(move(mesh.cellZones_)),
689  globalMeshDataPtr_(move(mesh.globalMeshDataPtr_)),
690  curMotionTimeIndex_(mesh.curMotionTimeIndex_),
691  oldPointsPtr_(move(mesh.oldPointsPtr_)),
692  oldCellCentresPtr_(move(mesh.oldCellCentresPtr_)),
693  storeOldCellCentres_(mesh.storeOldCellCentres_),
694  moving_(mesh.moving_),
695  topoChanged_(mesh.topoChanged_)
696 {}
697 
698 
700 (
701  pointField&& points,
702  faceList&& faces,
703  labelList&& owner,
704  labelList&& neighbour,
705  const labelList& patchSizes,
706  const labelList& patchStarts,
707  const bool validBoundary
708 )
709 {
710  // Clear addressing. Keep geometric props and updateable props for mapping.
711  clearAddressing(true);
712 
713  // Take over new primitive data.
714  // Optimised to avoid overwriting data at all
715  if (notNull(points))
716  {
717  points_ = move(points);
718  bounds_ = boundBox(points_, validBoundary);
719  }
720 
721  if (notNull(faces))
722  {
723  faces_ = move(faces);
724  }
725 
726  if (notNull(owner))
727  {
728  owner_ = move(owner);
729  }
730 
731  if (notNull(neighbour))
732  {
733  neighbour_ = move(neighbour);
734  }
735 
736 
737  // Reset patch sizes and starts
738  forAll(boundary_, patchi)
739  {
740  boundary_[patchi] = polyPatch
741  (
742  boundary_[patchi],
743  boundary_,
744  patchi,
745  patchSizes[patchi],
746  patchStarts[patchi]
747  );
748  }
749 
750 
751  // Flags the mesh files as being changed
752  setInstance(time().name());
753 
754  // Check if the faces and cells are valid
755  forAll(faces_, facei)
756  {
757  const face& curFace = faces_[facei];
758 
759  if (min(curFace) < 0 || max(curFace) > points_.size())
760  {
762  << "Face " << facei << " contains vertex labels out of range: "
763  << curFace << " Max point index = " << points_.size()
764  << abort(FatalError);
765  }
766  }
767 
768 
769  // Set the primitive mesh from the owner_, neighbour_.
770  // Works out from patch end where the active faces stop.
771  initMesh();
772 
773 
774  if (validBoundary)
775  {
776  // Note that we assume that all the patches stay the same and are
777  // correct etc. so we can already use the patches to do
778  // processor-processor comms.
779 
780  // Calculate topology for the patches (processor-processor comms etc.)
781  boundary_.topoChange();
782 
783  // Calculate the geometry for the patches (transformation tensors etc.)
784  boundary_.calcGeometry();
785 
786  // Warn if global empty mesh
787  if
788  (
789  (returnReduce(nPoints(), sumOp<label>()) == 0)
790  || (returnReduce(nCells(), sumOp<label>()) == 0)
791  )
792  {
794  << "no points or no cells in mesh"
795  << exit(FatalError);
796  }
797  }
798 }
799 
800 
802 {
803  // Clear addressing. Keep geometric and updatable properties for mapping.
804  clearAddressing(true);
805  otherMesh.clearAddressing(true);
806 
807  // Swap the primitives
808  points_.swap(otherMesh.points_);
809  bounds_ = boundBox(points_, true);
810  faces_.swap(otherMesh.faces_);
811  owner_.swap(otherMesh.owner_);
812  neighbour_.swap(otherMesh.neighbour_);
813 
814  // Clear the boundary data
815  boundary_.clearGeom();
816  boundary_.clearAddressing();
817  otherMesh.boundary_.clearGeom();
818  otherMesh.boundary_.clearAddressing();
819 
820  // Swap the boundaries
821  auto updatePatches = []
822  (
823  const polyPatchList& otherPatches,
824  polyBoundaryMesh& boundaryMesh
825  )
826  {
827  boundaryMesh.resize(otherPatches.size());
828 
829  forAll(otherPatches, otherPatchi)
830  {
831  // Clone processor patches, as the decomposition may be different
832  // in the other mesh. Just change the size and start of other
833  // patches.
834 
835  if (isA<processorPolyPatch>(otherPatches[otherPatchi]))
836  {
837  boundaryMesh.set
838  (
839  otherPatchi,
840  otherPatches[otherPatchi].clone(boundaryMesh)
841  );
842  }
843  else
844  {
845  boundaryMesh[otherPatchi] = polyPatch
846  (
847  boundaryMesh[otherPatchi],
848  boundaryMesh,
849  otherPatchi,
850  otherPatches[otherPatchi].size(),
851  otherPatches[otherPatchi].start()
852  );
853  }
854  }
855  };
856 
857  {
858  const polyPatchList patches
859  (
860  boundary_,
861  boundary_
862  );
863  const polyPatchList otherPatches
864  (
865  otherMesh.boundary_,
866  otherMesh.boundary_
867  );
868 
869  updatePatches(otherPatches, boundary_);
870  updatePatches(patches, otherMesh.boundary_);
871  }
872 
873  // Parallel data depends on the patch ordering so force recalculation
874  globalMeshDataPtr_.clear();
875  otherMesh.globalMeshDataPtr_.clear();
876 
877  // Flags the mesh files as being changed
878  setInstance(time().name());
879  otherMesh.setInstance(time().name());
880 
881  // Check if the faces and cells are valid
882  auto checkFaces = [](const polyMesh& mesh)
883  {
884  forAll(mesh.faces_, facei)
885  {
886  const face& curFace = mesh.faces_[facei];
887 
888  if (min(curFace) < 0 || max(curFace) > mesh.points_.size())
889  {
891  << "Face " << facei << " contains vertex labels out of "
892  << "range: " << curFace << " Max point index = "
893  << mesh.points_.size() << abort(FatalError);
894  }
895  }
896  };
897 
898  checkFaces(*this);
899  checkFaces(otherMesh);
900 
901  // Set the primitive mesh from the owner_, neighbour_.
902  // Works out from patch end where the active faces stop.
903  initMesh();
904  otherMesh.initMesh();
905 
906  // Calculate topology for the patches (processor-processor comms etc.)
907  boundary_.topoChange();
908  otherMesh.boundary_.topoChange();
909 
910  // Calculate the geometry for the patches (transformation tensors etc.)
911  boundary_.calcGeometry();
912  otherMesh.boundary_.calcGeometry();
913 
914  // Update the optional pointMesh with respect to the updated polyMesh
915  if (foundObject<pointMesh>(pointMesh::typeName))
916  {
917  pointMesh::New(*this).reset();
918  }
919  if (otherMesh.foundObject<pointMesh>(pointMesh::typeName))
920  {
921  pointMesh::New(*this).reset();
922  }
923 
924  // Update point zones
925  if (pointZones_.size() == otherMesh.pointZones_.size())
926  {
927  pointZones_.clearAddressing();
928  otherMesh.pointZones_.clearAddressing();
929 
930  forAll(pointZones_, i)
931  {
932  pointZones_[i].swap(otherMesh.pointZones_[i]);
933  }
934  }
935  else
936  {
938  << "Number of pointZones in other mesh = "
939  << otherMesh.pointZones_.size()
940  << " is not the same as in the mesh = "
941  << pointZones_.size()
942  << exit(FatalError);
943  }
944 
945  // Update face zones
946  if (faceZones_.size() == otherMesh.faceZones_.size())
947  {
948  faceZones_.clearAddressing();
949  otherMesh.faceZones_.clearAddressing();
950 
951  forAll(faceZones_, i)
952  {
953  faceZones_[i].swap(otherMesh.faceZones_[i]);
954  }
955  }
956  else
957  {
959  << "Number of faceZones in other mesh = "
960  << otherMesh.faceZones_.size()
961  << " is not the same as in the mesh = "
962  << faceZones_.size()
963  << exit(FatalError);
964  }
965 
966  // Update cell zones
967  if (cellZones_.size() == otherMesh.cellZones_.size())
968  {
969  cellZones_.clearAddressing();
970  otherMesh.cellZones_.clearAddressing();
971 
972  forAll(cellZones_, i)
973  {
974  cellZones_[i].swap(otherMesh.cellZones_[i]);
975  }
976  }
977  else
978  {
980  << "Number of cellZones in other mesh = "
981  << otherMesh.cellZones_.size()
982  << " is not the same as in the mesh = "
983  << cellZones_.size()
984  << exit(FatalError);
985  }
986 }
987 
988 
989 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
990 
992 {
993  clearOut();
994  resetMotion();
995 }
996 
997 
998 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
999 
1001 {
1002  if (objectRegistry::dbDir() == defaultRegion)
1003  {
1004  return parent().dbDir();
1005  }
1006  else
1007  {
1008  return objectRegistry::dbDir();
1009  }
1010 }
1011 
1012 
1014 {
1015  return dbDir()/meshSubDir;
1016 }
1017 
1018 
1020 {
1021  return points_.instance();
1022 }
1023 
1024 
1026 {
1027  return faces_.instance();
1028 }
1029 
1030 
1032 {
1033  return points_.writeOpt();
1034 }
1035 
1036 
1038 {
1039  return faces_.writeOpt();
1040 }
1041 
1042 
1044 {
1045  if (geometricD_.x() == 0)
1046  {
1047  calcDirections();
1048  }
1049 
1050  return geometricD_;
1051 }
1052 
1053 
1055 {
1056  return cmptSum(geometricD() + Vector<label>::one)/2;
1057 }
1058 
1059 
1061 {
1062  if (solutionD_.x() == 0)
1063  {
1064  calcDirections();
1065  }
1066 
1067  return solutionD_;
1068 }
1069 
1070 
1072 {
1073  return cmptSum(solutionD() + Vector<label>::one)/2;
1074 }
1075 
1076 
1078 {
1079  if (tetBasePtIsPtr_.empty())
1080  {
1081  if (debug)
1082  {
1084  << "Forcing storage of base points."
1085  << endl;
1086  }
1087 
1088  tetBasePtIsPtr_.reset
1089  (
1090  new labelIOList
1091  (
1092  IOobject
1093  (
1094  "tetBasePtIs",
1095  instance(),
1096  meshSubDir,
1097  *this,
1100  ),
1102  )
1103  );
1104  }
1105 
1106  return tetBasePtIsPtr_();
1107 }
1108 
1109 
1112 {
1113  if (cellTreePtr_.empty())
1114  {
1115  cellTreePtr_.reset
1116  (
1118  (
1119  treeDataCell
1120  (
1121  false, // not cache bb
1122  *this,
1123  CELL_TETS // use tet-decomposition for any inside test
1124  ),
1125  treeBoundBox(points()).extend(1e-4),
1126  8, // maxLevel
1127  10, // leafsize
1128  5.0 // duplicity
1129  )
1130  );
1131  }
1132 
1133  return cellTreePtr_();
1134 }
1135 
1136 
1138 (
1139  const List<polyPatch*>& p,
1140  const bool validBoundary
1141 )
1142 {
1143  if (boundaryMesh().size())
1144  {
1146  << "boundary already exists"
1147  << abort(FatalError);
1148  }
1149 
1150  // Reset valid directions
1151  geometricD_ = Zero;
1152  solutionD_ = Zero;
1153 
1154  boundary_.setSize(p.size());
1155 
1156  // Copy the patch pointers
1157  forAll(p, pI)
1158  {
1159  boundary_.set(pI, p[pI]);
1160  }
1161 
1162  // parallelData depends on the processorPatch ordering so force
1163  // recalculation. Problem: should really be done in removeBoundary but
1164  // there is some info in parallelData which might be interesting in between
1165  // removeBoundary and addPatches.
1166  globalMeshDataPtr_.clear();
1167 
1168  if (validBoundary)
1169  {
1170  // Calculate topology for the patches (processor-processor comms etc.)
1171  boundary_.topoChange();
1172 
1173  // Calculate the geometry for the patches (transformation tensors etc.)
1174  boundary_.calcGeometry();
1175 
1176  boundary_.checkDefinition();
1177  }
1178 }
1179 
1180 
1182 (
1183  const List<pointZone*>& pz,
1184  const List<faceZone*>& fz,
1185  const List<cellZone*>& cz
1186 )
1187 {
1188  if (pointZones().size() || faceZones().size() || cellZones().size())
1189  {
1191  << "point, face or cell zone already exists"
1192  << abort(FatalError);
1193  }
1194 
1195  // Point zones
1196  if (pz.size())
1197  {
1198  pointZones_.setSize(pz.size());
1199 
1200  // Copy the zone pointers
1201  forAll(pz, pI)
1202  {
1203  pointZones_.set(pI, pz[pI]);
1204  }
1205 
1206  pointZones_.writeOpt() = IOobject::AUTO_WRITE;
1207  }
1208 
1209  // Face zones
1210  if (fz.size())
1211  {
1212  faceZones_.setSize(fz.size());
1213 
1214  // Copy the zone pointers
1215  forAll(fz, fI)
1216  {
1217  faceZones_.set(fI, fz[fI]);
1218  }
1219 
1220  faceZones_.writeOpt() = IOobject::AUTO_WRITE;
1221  }
1222 
1223  // Cell zones
1224  if (cz.size())
1225  {
1226  cellZones_.setSize(cz.size());
1227 
1228  // Copy the zone pointers
1229  forAll(cz, cI)
1230  {
1231  cellZones_.set(cI, cz[cI]);
1232  }
1233 
1234  cellZones_.writeOpt() = IOobject::AUTO_WRITE;
1235  }
1236 }
1237 
1238 
1240 (
1241  const labelUList& newToOld,
1242  const bool validBoundary
1243 )
1244 {
1245  // Clear local fields and e.g. polyMesh parallelInfo. Do not clearGeom
1246  // so we keep PatchMeshObjects intact.
1247  boundary_.clearGeom();
1248  clearAddressing(true);
1249  // Clear all but PatchMeshObjects
1251  <
1252  polyMesh,
1255  >
1256  (
1257  *this
1258  );
1260  <
1261  pointMesh,
1264  >
1265  (
1266  *this
1267  );
1268 
1269  boundary_.reorderPatches(newToOld, validBoundary);
1270 
1271  // Warn mesh objects
1272  meshObjects::reorderPatches<polyMesh>(*this, newToOld, validBoundary);
1273  meshObjects::reorderPatches<pointMesh>(*this, newToOld, validBoundary);
1274 }
1275 
1276 
1278 (
1279  const label insertPatchi,
1280  const polyPatch& patch,
1281  const dictionary& patchFieldDict,
1282  const word& defaultPatchFieldType,
1283  const bool validBoundary
1284 )
1285 {
1286  const label sz = boundary_.size();
1287 
1288  label startFacei = nFaces();
1289  if (insertPatchi < sz)
1290  {
1291  startFacei = boundary_[insertPatchi].start();
1292  }
1293 
1294  // Create reordering list
1295  // patches before insert position stay as is
1296  // patches after insert position move one up
1297  labelList newToOld(boundary_.size()+1);
1298  for (label i = 0; i < insertPatchi; i++)
1299  {
1300  newToOld[i] = i;
1301  }
1302  for (label i = insertPatchi; i < sz; i++)
1303  {
1304  newToOld[i+1] = i;
1305  }
1306  newToOld[insertPatchi] = -1;
1307 
1308  reorderPatches(newToOld, false);
1309 
1310  // Clear local fields and e.g. polyMesh parallelInfo.
1311  //clearGeom(); // would clear out pointMesh as well
1312  boundary_.clearGeom();
1313  clearAddressing(true);
1314 
1315  // Clear all but PatchMeshObjects
1317  <
1318  polyMesh,
1321  >
1322  (
1323  *this
1324  );
1326  <
1327  pointMesh,
1330  >
1331  (
1332  *this
1333  );
1334 
1335 
1336  // Insert polyPatch
1337  boundary_.set
1338  (
1339  insertPatchi,
1340  patch.clone
1341  (
1342  boundary_,
1343  insertPatchi, // index
1344  0, // size
1345  startFacei // start
1346  )
1347  );
1348 
1349  if (validBoundary)
1350  {
1351  boundary_.topoChange();
1352  }
1353 
1354  // Warn mesh objects
1355  meshObjects::addPatch<polyMesh>(*this, insertPatchi);
1356  meshObjects::addPatch<pointMesh>(*this, insertPatchi);
1357 }
1358 
1359 
1361 {
1362  if (clearedPrimitives_)
1363  {
1365  << "points deallocated"
1366  << abort(FatalError);
1367  }
1368 
1369  return points_;
1370 }
1371 
1372 
1374 {
1375  if (clearedPrimitives_)
1376  {
1378  << "faces deallocated"
1379  << abort(FatalError);
1380  }
1381 
1382  return faces_;
1383 }
1384 
1385 
1387 {
1388  return owner_;
1389 }
1390 
1391 
1393 {
1394  return neighbour_;
1395 }
1396 
1397 
1399 {
1400  if (!moving_)
1401  {
1402  return points_;
1403  }
1404 
1405  if (oldPointsPtr_.empty())
1406  {
1408  << "Old points have not been stored"
1409  << exit(FatalError);
1410  }
1411 
1412  return oldPointsPtr_();
1413 }
1414 
1415 
1417 {
1418  storeOldCellCentres_ = true;
1419 
1420  if (!moving_)
1421  {
1422  return cellCentres();
1423  }
1424 
1425  if (oldCellCentresPtr_.empty())
1426  {
1428  << "Old cell centres have not been stored"
1429  << exit(FatalError);
1430  }
1431 
1432  return oldCellCentresPtr_();
1433 }
1434 
1435 
1437 {
1438  return io.upToDate(points_);
1439 }
1440 
1441 
1443 {
1444  io.eventNo() = points_.eventNo()+1;
1445 }
1446 
1447 
1449 {
1450  if (debug)
1451  {
1453  << "Set points for time " << time().value()
1454  << " index " << time().timeIndex() << endl;
1455  }
1456 
1458 
1459  points_ = newPoints;
1460 
1461  setPointsInstance(time().name());
1462 
1463  // Adjust parallel shared points
1464  if (globalMeshDataPtr_.valid())
1465  {
1466  globalMeshDataPtr_().movePoints(points_);
1467  }
1468 
1469  // Force recalculation of all geometric data with new points
1470 
1471  bounds_ = boundBox(points_);
1472  boundary_.movePoints(points_);
1473 
1474  pointZones_.movePoints(points_);
1475  faceZones_.movePoints(points_);
1476  cellZones_.movePoints(points_);
1477 
1478  // Cell tree might become invalid
1479  cellTreePtr_.clear();
1480 
1481  // Reset valid directions (could change with rotation)
1482  geometricD_ = Zero;
1483  solutionD_ = Zero;
1484 
1485  meshObjects::movePoints<polyMesh>(*this);
1486  meshObjects::movePoints<pointMesh>(*this);
1487 }
1488 
1489 
1491 (
1492  const pointField& newPoints
1493 )
1494 {
1495  if (debug)
1496  {
1498  << "Moving points for time " << time().value()
1499  << " index " << time().timeIndex() << endl;
1500  }
1501 
1502  // Pick up old points and cell centres
1503  if (curMotionTimeIndex_ != time().timeIndex())
1504  {
1505  oldPointsPtr_.clear();
1506  oldPointsPtr_.reset(new pointField(points_));
1507  if (storeOldCellCentres_)
1508  {
1509  oldCellCentresPtr_.clear();
1510  oldCellCentresPtr_.reset(new pointField(cellCentres()));
1511  }
1512  curMotionTimeIndex_ = time().timeIndex();
1513  }
1514 
1515  points_ = newPoints;
1516 
1517  setPointsInstance(time().name());
1518 
1520  (
1521  points_,
1522  oldPoints()
1523  );
1524 
1525  // Adjust parallel shared points
1526  if (globalMeshDataPtr_.valid())
1527  {
1528  globalMeshDataPtr_().movePoints(points_);
1529  }
1530 
1531  // Force recalculation of all geometric data with new points
1532 
1533  bounds_ = boundBox(points_);
1534  boundary_.movePoints(points_);
1535 
1536  pointZones_.movePoints(points_);
1537  faceZones_.movePoints(points_);
1538  cellZones_.movePoints(points_);
1539 
1540  // Cell tree might become invalid
1541  cellTreePtr_.clear();
1542 
1543  // Reset valid directions (could change with rotation)
1544  geometricD_ = Zero;
1545  solutionD_ = Zero;
1546 
1547  meshObjects::movePoints<polyMesh>(*this);
1548  meshObjects::movePoints<pointMesh>(*this);
1549 
1550  return sweptVols;
1551 }
1552 
1553 
1555 {
1556  curMotionTimeIndex_ = -1;
1557  oldPointsPtr_.clear();
1558  oldCellCentresPtr_.clear();
1559 }
1560 
1561 
1563 {
1564  if (globalMeshDataPtr_.empty())
1565  {
1566  if (debug)
1567  {
1568  Pout<< "polyMesh::globalData() const : "
1569  << "Constructing parallelData from processor topology"
1570  << endl;
1571  }
1572  // Construct globalMeshData using processorPatch information only.
1573  globalMeshDataPtr_.reset(new globalMeshData(*this));
1574  }
1575 
1576  return globalMeshDataPtr_();
1577 }
1578 
1579 
1581 {
1582  return comm_;
1583 }
1584 
1585 
1587 {
1588  return comm_;
1589 }
1590 
1591 
1592 void Foam::polyMesh::removeFiles(const fileName& instanceDir) const
1593 {
1594  fileName meshFilesPath = thisDb().time().path()/instanceDir/meshDir();
1595 
1596  rm(meshFilesPath/"points");
1597  rm(meshFilesPath/"faces");
1598  rm(meshFilesPath/"owner");
1599  rm(meshFilesPath/"neighbour");
1600  rm(meshFilesPath/"cells");
1601  rm(meshFilesPath/"boundary");
1602  rm(meshFilesPath/"pointZones");
1603  rm(meshFilesPath/"faceZones");
1604  rm(meshFilesPath/"cellZones");
1605  rm(meshFilesPath/"meshModifiers");
1606  rm(meshFilesPath/"parallelData");
1607 
1608  // remove subdirectories
1609  if (isDir(meshFilesPath/"sets"))
1610  {
1611  rmDir(meshFilesPath/"sets");
1612  }
1613 }
1614 
1615 
1617 {
1618  removeFiles(instance());
1619 }
1620 
1621 
1623 (
1624  const point& p,
1625  label& celli,
1626  label& tetFacei,
1627  label& tetPti
1628 ) const
1629 {
1630  celli = -1;
1631  tetFacei = -1;
1632  tetPti = -1;
1633 
1634  const indexedOctree<treeDataCell>& tree = cellTree();
1635 
1636  // Find point inside cell
1637  celli = tree.findInside(p);
1638 
1639  if (celli != -1)
1640  {
1641  // Check the nearest cell to see if the point is inside.
1642  findTetFacePt(celli, p, tetFacei, tetPti);
1643  }
1644 }
1645 
1646 
1648 (
1649  const label celli,
1650  const point& p,
1651  label& tetFacei,
1652  label& tetPti
1653 ) const
1654 {
1655  const polyMesh& mesh = *this;
1656 
1657  tetIndices tet(polyMeshTetDecomposition::findTet(mesh, celli, p));
1658  tetFacei = tet.face();
1659  tetPti = tet.tetPt();
1660 }
1661 
1662 
1664 (
1665  const point& p,
1666  label celli,
1667  const cellDecomposition decompMode
1668 ) const
1669 {
1670  switch (decompMode)
1671  {
1672  case FACE_PLANES:
1673  {
1674  return primitiveMesh::pointInCell(p, celli);
1675  }
1676  break;
1677 
1678  case FACE_CENTRE_TRIS:
1679  {
1680  // only test that point is on inside of plane defined by cell face
1681  // triangles
1682  const cell& cFaces = cells()[celli];
1683 
1684  forAll(cFaces, cFacei)
1685  {
1686  label facei = cFaces[cFacei];
1687  const face& f = faces_[facei];
1688  const point& fc = faceCentres()[facei];
1689  bool isOwn = (owner_[facei] == celli);
1690 
1691  forAll(f, fp)
1692  {
1693  label pointi;
1694  label nextPointi;
1695 
1696  if (isOwn)
1697  {
1698  pointi = f[fp];
1699  nextPointi = f.nextLabel(fp);
1700  }
1701  else
1702  {
1703  pointi = f.nextLabel(fp);
1704  nextPointi = f[fp];
1705  }
1706 
1707  triPointRef faceTri
1708  (
1709  points()[pointi],
1710  points()[nextPointi],
1711  fc
1712  );
1713 
1714  vector proj = p - faceTri.centre();
1715 
1716  if ((faceTri.area() & proj) > 0)
1717  {
1718  return false;
1719  }
1720  }
1721  }
1722  return true;
1723  }
1724  break;
1725 
1726  case FACE_DIAG_TRIS:
1727  {
1728  // only test that point is on inside of plane defined by cell face
1729  // triangles
1730  const cell& cFaces = cells()[celli];
1731 
1732  forAll(cFaces, cFacei)
1733  {
1734  label facei = cFaces[cFacei];
1735  const face& f = faces_[facei];
1736 
1737  for (label tetPti = 1; tetPti < f.size() - 1; tetPti++)
1738  {
1739  // Get tetIndices of face triangle
1740  tetIndices faceTetIs(celli, facei, tetPti);
1741 
1742  triPointRef faceTri = faceTetIs.faceTri(*this);
1743 
1744  vector proj = p - faceTri.centre();
1745 
1746  if ((faceTri.area() & proj) > 0)
1747  {
1748  return false;
1749  }
1750  }
1751  }
1752 
1753  return true;
1754  }
1755  break;
1756 
1757  case CELL_TETS:
1758  {
1759  label tetFacei;
1760  label tetPti;
1761 
1762  findTetFacePt(celli, p, tetFacei, tetPti);
1763 
1764  return tetFacei != -1;
1765  }
1766  break;
1767  }
1768 
1769  return false;
1770 }
1771 
1772 
1774 (
1775  const point& p,
1776  const cellDecomposition decompMode
1777 ) const
1778 {
1779  if
1780  (
1781  Pstream::parRun()
1782  && (decompMode == FACE_DIAG_TRIS || decompMode == CELL_TETS)
1783  )
1784  {
1785  // Force construction of face-diagonal decomposition before testing
1786  // for zero cells.
1787  //
1788  // If parallel running a local domain might have zero cells so never
1789  // construct the face-diagonal decomposition which uses parallel
1790  // transfers.
1791  (void)tetBasePtIs();
1792  }
1793 
1794  if (nCells() == 0)
1795  {
1796  return -1;
1797  }
1798 
1799  if (decompMode == CELL_TETS)
1800  {
1801  // Advanced search method utilising an octree
1802  // and tet-decomposition of the cells
1803 
1804  label celli;
1805  label tetFacei;
1806  label tetPti;
1807 
1808  findCellFacePt(p, celli, tetFacei, tetPti);
1809 
1810  return celli;
1811  }
1812  else
1813  {
1814  // Approximate search avoiding the construction of an octree
1815  // and cell decomposition
1816 
1817  // Find the nearest cell centre to this location
1818  label celli = findNearestCell(p);
1819 
1820  // If point is in the nearest cell return
1821  if (pointInCell(p, celli, decompMode))
1822  {
1823  return celli;
1824  }
1825  else
1826  {
1827  // Point is not in the nearest cell so search all cells
1828 
1829  for (label celli = 0; celli < nCells(); celli++)
1830  {
1831  if (pointInCell(p, celli, decompMode))
1832  {
1833  return celli;
1834  }
1835  }
1836 
1837  return -1;
1838  }
1839  }
1840 }
1841 
1842 
1843 // ************************************************************************* //
#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:190
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:991
void findCellFacePt(const point &p, label &celli, label &tetFacei, label &tetPti) const
Find the cell, tetFacei and tetPti for point p.
Definition: polyMesh.C:1623
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:1025
label nGeometricD() const
Return the number of valid geometric dimensions in the mesh.
Definition: polyMesh.C:1054
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:268
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:1013
virtual tmp< scalarField > movePoints(const pointField &)
Move points, returns volumes swept by faces in motion.
Definition: polyMesh.C:1491
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:1278
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:1037
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1373
virtual const fileName & dbDir() const
Override the objectRegistry dbDir for a single-region case.
Definition: polyMesh.C:1000
label findCell(const point &p, const cellDecomposition=CELL_TETS) const
Find cell enclosing this location and return index.
Definition: polyMesh.C:1774
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:700
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1386
const labelIOList & tetBasePtIs() const
Return the tetBasePtIs.
Definition: polyMesh.C:1077
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:1648
const globalMeshData & globalData() const
Return parallel info.
Definition: polyMesh.C:1562
bool pointInCell(const point &p, label celli, const cellDecomposition=CELL_TETS) const
Test if point p is in the celli.
Definition: polyMesh.C:1664
label nSolutionD() const
Return the number of valid solved-for dimensions in the mesh.
Definition: polyMesh.C:1071
virtual const pointField & oldPoints() const
Return old points for mesh motion.
Definition: polyMesh.C:1398
void swap(polyMesh &)
Swap mesh.
Definition: polyMesh.C:801
const fileName & pointsInstance() const
Return the current instance directory for points.
Definition: polyMesh.C:1019
label comm() const
Return communicator used for parallel communication.
Definition: polyMesh.C:1580
polyMesh(const IOobject &io)
Construct from IOobject.
Definition: polyMesh.C:162
void clearAddressing(const bool isMeshUpdate=false)
Clear addressing.
Definition: polyMeshClear.C:77
void resetMotion() const
Reset motion.
Definition: polyMesh.C:1554
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:405
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1392
void addPatches(const List< polyPatch * > &, const bool validBoundary=true)
Add boundary patches.
Definition: polyMesh.C:1138
bool upToDatePoints(const regIOobject &io) const
Return true if io is up-to-date with points.
Definition: polyMesh.C:1436
void addZones(const List< pointZone * > &pz, const List< faceZone * > &fz, const List< cellZone * > &cz)
Add mesh zones.
Definition: polyMesh.C:1182
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1360
void setUpToDatePoints(regIOobject &io) const
Set io to be up-to-date with points.
Definition: polyMesh.C:1442
const indexedOctree< treeDataCell > & cellTree() const
Return the cell search tree.
Definition: polyMesh.C:1111
IOobject::writeOption pointsWriteOpt() const
Return the points write option.
Definition: polyMesh.C:1031
void removeFiles() const
Remove all files from mesh instance()
Definition: polyMesh.C:1616
virtual const pointField & oldCellCentres() const
Return old cell centres for mesh motion.
Definition: polyMesh.C:1416
virtual void reorderPatches(const labelUList &newToOld, const bool validBoundary)
Reorder and trim existing patches. If validBoundary the new.
Definition: polyMesh.C:1240
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:271
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:1448
const Vector< label > & solutionD() const
Return the vector of solved-for directions in mesh.
Definition: polyMesh.C:1060
const Vector< label > & geometricD() const
Return the vector of geometric directions in mesh.
Definition: polyMesh.C:1043
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
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
uint8_t direction
Definition: direction.H:45
label timeIndex
Definition: getTimeIndex.H:4
labelList f(nPoints)
volScalarField & p