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