foamVtuSizing.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2016-2021 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "foamVtuSizing.H"
29 #include "foamVtkCore.H"
30 #include "polyMesh.H"
31 #include "cellShape.H"
32 
33 // Only used in this file
34 #include "foamVtuSizingImpl.C"
35 
36 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
37 
38 void Foam::vtk::vtuSizing::presizeMaps(foamVtkMeshMaps& maps) const
39 {
40  maps.cellMap().resize(this->nFieldCells());
41  maps.additionalIds().resize(this->nAddPoints());
42 }
43 
44 
45 void Foam::vtk::vtuSizing::checkSizes
46 (
47  const vtk::vtuSizing& sizing,
48 
49  const label cellTypes_size,
50  const label vertLabels_size,
51  const label vertOffset_size,
52  const label faceLabels_size,
53  const label faceOffset_size,
54 
55  const enum contentType output,
56  const label cellMap_size,
57  const label addPointsIds_size
58 )
59 {
60  label nErrors = 0;
61 
62  #undef CHECK_SIZING
63  #define CHECK_SIZING(what, sizeInput, sizeExpected) \
64  if (sizeInput != sizeExpected) \
65  { \
66  if (!nErrors++) \
67  { \
68  FatalErrorInFunction << "VTK sizing error" << nl; \
69  } \
70  FatalError \
71  << " " << what << " size=" << sizeInput \
72  << " expected " << sizeExpected << nl; \
73  }
74 
75 
76  CHECK_SIZING("cellTypes", cellTypes_size, sizing.nFieldCells());
77  CHECK_SIZING("cellMap", cellMap_size, sizing.nFieldCells());
78  CHECK_SIZING("addPointsIds", addPointsIds_size, sizing.nAddPoints());
79 
80  switch (output)
81  {
82  case contentType::LEGACY:
83  {
84  CHECK_SIZING("legacy", vertLabels_size, sizing.sizeLegacy());
85  break;
86  }
87 
88  case contentType::XML:
89  {
90  // XML uses connectivity/offset pair.
92  (
93  "connectivity",
94  vertLabels_size,
95  sizing.sizeXml(slotType::CELLS)
96  );
98  (
99  "offsets",
100  vertOffset_size,
101  sizing.sizeXml(slotType::CELLS_OFFSETS)
102  );
103  if (sizing.nFaceLabels())
104  {
106  (
107  "faces",
108  faceLabels_size,
109  sizing.sizeXml(slotType::FACES)
110  );
111 
113  (
114  "faceOffsets",
115  faceOffset_size,
116  sizing.sizeXml(slotType::FACES_OFFSETS)
117  );
118  }
119  break;
120  }
121 
122  case contentType::INTERNAL1:
123  {
124  // VTK-internal1 connectivity/offset pair.
126  (
127  "connectivity",
128  vertLabels_size,
129  sizing.sizeInternal1(slotType::CELLS)
130  );
132  (
133  "offsets",
134  vertOffset_size,
135  sizing.sizeInternal1(slotType::CELLS_OFFSETS)
136  );
137  if (sizing.nFaceLabels())
138  {
140  (
141  "faces",
142  faceLabels_size,
143  sizing.sizeInternal1(slotType::FACES)
144  );
146  (
147  "faceOffsets",
148  faceOffset_size,
149  sizing.sizeInternal1(slotType::FACES_OFFSETS)
150  );
151  }
152  break;
153  }
154 
155  case contentType::INTERNAL2:
156  {
157  // VTK-internal2 connectivity/offset pair.
159  (
160  "connectivity",
161  vertLabels_size,
162  sizing.sizeInternal2(slotType::CELLS)
163  );
165  (
166  "offsets",
167  vertOffset_size,
168  sizing.sizeInternal2(slotType::CELLS_OFFSETS)
169  );
170  if (sizing.nFaceLabels())
171  {
173  (
174  "faces",
175  faceLabels_size,
176  sizing.sizeInternal2(slotType::FACES)
177  );
179  (
180  "faceOffsets",
181  faceOffset_size,
182  sizing.sizeInternal2(slotType::FACES_OFFSETS)
183  );
184  }
185  break;
186  }
187  }
188 
189  if (nErrors)
190  {
191  FatalError
192  << nl
193  << "Total of " << nErrors << " sizing errors encountered!"
194  << exit(FatalError);
195  }
196 
197  #undef CHECK_SIZING
198 }
199 
200 
201 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
202 
204 {
205  clear();
206 }
207 
208 
210 (
211  const polyMesh& mesh,
212  const bool decompose
213 )
214 {
215  clear();
216  reset(mesh, decompose);
217 }
218 
219 
220 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
221 
223 {
224  decompose_ = false;
225  selectionMode_ = FULL_MESH;
226  nCells_ = 0;
227  nPoints_ = 0;
228  nVertLabels_ = 0;
229 
230  nFaceLabels_ = 0;
231  nCellsPoly_ = 0;
232  nVertPoly_ = 0;
233 
234  nAddCells_ = 0;
235  nAddPoints_ = 0;
236  nAddVerts_ = 0;
237 }
238 
239 
241 (
242  const polyMesh& mesh,
243  const bool decompose
244 )
245 {
246  reset(mesh, labelUList::null(), decompose);
247 }
248 
249 
251 (
252  const polyMesh& mesh,
253  const labelUList& subsetCellsIds,
254  const bool decompose
255 )
256 {
257  // References to cell shape models
260  const cellModel& prism = cellModel::ref(cellModel::PRISM);
262  const cellModel& wedge = cellModel::ref(cellModel::WEDGE);
263  const cellModel& tetWedge = cellModel::ref(cellModel::TETWEDGE);
264 
265  const cellShapeList& shapes = mesh.cellShapes();
266 
267  // Unique vertex labels per polyhedral
268  labelHashSet hashUniqId(2*256);
269 
270 
271  // Special treatment for mesh subsets.
272  const bool isSubsetMesh
273  (
274  notNull(subsetCellsIds)
275  );
276 
277  if (isSubsetMesh)
278  {
279  decompose_ = false; // Disallow decomposition for subset mode
280  selectionMode_ = selectionModeType::SUBSET_MESH;
281  }
282  else
283  {
284  decompose_ = decompose; // Disallow decomposition
285  selectionMode_ = selectionModeType::FULL_MESH;
286  }
287 
288  const label nInputCells =
289  (
290  isSubsetMesh
291  ? subsetCellsIds.size()
292  : shapes.size()
293  );
294 
295  nCells_ = nInputCells;
296  nPoints_ = mesh.nPoints();
297  nAddCells_ = 0;
298  nAddVerts_ = 0;
299 
300  nCellsPoly_ = nCells_;
301  nVertLabels_ = 0;
302  nFaceLabels_ = 0;
303  nVertPoly_ = 0;
304 
305  for (label inputi = 0; inputi < nInputCells; ++inputi)
306  {
307  const label celli(isSubsetMesh ? subsetCellsIds[inputi] : inputi);
308 
309  const cellShape& shape = shapes[celli];
310  const cellModel& model = shape.model();
311 
312  if
313  (
314  model == tet
315  || model == pyr
316  || model == prism
317  || model == hex
318  )
319  {
320  // Normal primitive - not a poly
321  --nCellsPoly_;
322  nVertLabels_ += shape.size();
323  }
324  else if (model == tetWedge && decompose_)
325  {
326  nVertLabels_ += 6; // Treat as squeezed prism (VTK_WEDGE)
327  }
328  else if (model == wedge && decompose_)
329  {
330  nVertLabels_ += 8; // Treat as squeezed hex
331  }
332  else if (decompose_)
333  {
334  // Polyhedral: Decompose into tets + pyramids.
335  ++nAddPoints_;
336 
337  // Count vertices into first decomposed cell
338  bool first = true;
339 
340  const cell& cFaces = mesh.cells()[celli];
341  for (const label facei : cFaces)
342  {
343  const face& f = mesh.faces()[facei];
344 
345  // Face decomposed into triangles and quads
346  // Tri -> Tet, Quad -> Pyr
347  label nTria = 0, nQuad = 0;
348  f.nTrianglesQuads(mesh.points(), nTria, nQuad);
349 
350  nAddCells_ += nTria + nQuad;
351  nAddVerts_ += (nTria * 4) + (nQuad * 5);
352 
353  if (first)
354  {
355  first = false;
356  --nAddCells_;
357 
358  const label nvrt = (nQuad ? 5 : 4);
359  nAddVerts_ -= nvrt;
360  nVertLabels_ += nvrt;
361  }
362  }
363  }
364  else
365  {
366  // Polyhedral: Not decomposed
367 
368  const labelList& cFaces = mesh.cells()[celli];
369 
370  // Unique node ids used (XML/INTERNAL, not needed for LEGACY)
371  hashUniqId.clear();
372 
373  // Face stream sizing:
374  // number of faces, size of each face, vertices per face
375  // [nFaces, nFace0Pts, id1, id2, ..., nFace1Pts, id1, id2, ...]
376 
377  for (const label facei : cFaces)
378  {
379  const face& f = mesh.faces()[facei];
380  nFaceLabels_ += f.size();
381 
382  hashUniqId.insert(f);
383  }
384 
385  // Legacy format only uses the face-stream.
386  // - track what *NOT* to use for legacy
387  nVertLabels_ += hashUniqId.size();
388  nVertPoly_ += hashUniqId.size();
389 
390  nFaceLabels_ += 1 + cFaces.size();
391  }
392  }
393 
394  // Requested and actually required
395  decompose_ = (decompose_ && nCellsPoly_);
396 }
397 
398 
399 // Synchronize changes here with the following:
400 // - vtuSizing::resetShapes
401 // - vtuSizing::populateArrays
402 //
404 (
405  const UList<cellShape>& shapes
406 )
407 {
410  const cellModel& prism = cellModel::ref(cellModel::PRISM);
412 
413  decompose_ = false; // Disallow decomposition
414  selectionMode_ = SHAPE_MESH;
415 
416  const label nInputCells = shapes.size();
417 
418  nCells_ = nInputCells;
419  nPoints_ = 0;
420  nAddCells_ = 0;
421  nAddVerts_ = 0;
422 
423  nCellsPoly_ = 0;
424  nVertLabels_ = 0;
425  nFaceLabels_ = 0;
426  nVertPoly_ = 0;
427 
428  label nIgnored = 0;
429 
430  for (label inputi = 0; inputi < nInputCells; ++inputi)
431  {
432  const cellShape& shape = shapes[inputi];
433  const cellModel& model = shape.model();
434 
435  if
436  (
437  model == tet
438  || model == pyr
439  || model == prism
440  || model == hex
441  )
442  {
443  nVertLabels_ += shape.size();
444 
445  // Guess for number of addressed points
446  nPoints_ = max(nPoints_, max(shape));
447  }
448  else
449  {
450  --nCells_;
451  ++nIgnored;
452  }
453  }
454 
455  if (nIgnored)
456  {
458  << "Encountered " << nIgnored << " unsupported cell shapes"
459  << " ... this is likely not good" << nl
460  << exit(FatalError);
461  }
462 
463  if (nCells_)
464  {
465  ++nPoints_;
466  }
467 }
468 
469 
471 (
472  const enum contentType output,
473  const enum slotType slot
474 ) const
475 {
476  switch (output)
477  {
478  case contentType::LEGACY:
479  {
480  switch (slot)
481  {
482  case slotType::CELLS:
483  // legacy uses connectivity for primitives, but directly
484  // stores face streams into connectivity as well.
485  // size-prefix per cell
486  return
487  (
488  nVertLabels() + nAddVerts() - nVertPoly() // primitives
489  + nFaceLabels() // face-stream (poly)
490  + nFieldCells() // nFieldCells (size prefix)
491  );
492  break;
493 
494  default:
495  break;
496  }
497  break;
498  }
499 
500  case contentType::XML:
501  {
502  switch (slot)
503  {
504  case slotType::CELLS:
505  return (nVertLabels() + nAddVerts());
506  break;
507 
508  case slotType::CELLS_OFFSETS:
509  return nFieldCells();
510  break;
511 
512  case slotType::FACES:
513  return nFaceLabels();
514  break;
515 
516  case slotType::FACES_OFFSETS:
517  return nFaceLabels() ? nFieldCells() : 0;
518  break;
519  }
520  break;
521  }
522 
523  case contentType::INTERNAL1:
524  {
525  switch (slot)
526  {
527  case slotType::CELLS:
528  // size-prefix per cell
529  return (nVertLabels() + nAddVerts() + nFieldCells());
530  break;
531 
532  case slotType::CELLS_OFFSETS:
533  return nFieldCells();
534  break;
535 
536  case slotType::FACES:
537  return nFaceLabels();
538  break;
539 
540  case slotType::FACES_OFFSETS:
541  return nFaceLabels() ? nFieldCells() : 0;
542  break;
543  }
544  break;
545  }
546 
547  case contentType::INTERNAL2:
548  {
549  switch (slot)
550  {
551  case slotType::CELLS:
552  return (nVertLabels() + nAddVerts());
553  break;
554 
555  case slotType::CELLS_OFFSETS:
556  return (nFieldCells() + 1);
557  break;
558 
559  case slotType::FACES:
560  return nFaceLabels();
561  break;
562 
563  case slotType::FACES_OFFSETS:
564  return nFaceLabels() ? nFieldCells() : 0;
565  break;
566  }
567  break;
568  }
569  }
570 
571  return 0;
572 }
573 
574 
575 // * * * * * * * * * * * * * * Populate Lists * * * * * * * * * * * * * * * //
576 
578 (
579  const polyMesh& mesh,
581  labelUList& vertLabels,
582  foamVtkMeshMaps& maps
583 ) const
584 {
585  // Leave as zero-sized so that populateArrays doesn't fill it.
586  List<label> unused;
587 
588  presizeMaps(maps);
589 
590  populateArrays
591  (
592  mesh,
593  *this,
594  cellTypes,
595  vertLabels,
596  unused, // offsets
597  unused, // faces
598  unused, // facesOffsets
599  contentType::LEGACY,
600  maps.cellMap(),
601  maps.additionalIds()
602  );
603 }
604 
605 
607 (
608  const UList<cellShape>& shapes,
610  labelUList& vertLabels,
611  foamVtkMeshMaps& maps
612 ) const
613 {
614  // Leave as zero-sized so that populateArrays doesn't fill it.
615  List<label> unused;
616 
617  presizeMaps(maps);
618 
619  populateArrays
620  (
621  shapes,
622  *this,
623  cellTypes,
624  vertLabels,
625  unused, // offsets
626  unused, // faces
627  unused, // facesOffsets
628  contentType::LEGACY,
629  maps.cellMap(),
630  maps.additionalIds()
631  );
632 }
633 
634 
636 (
637  const polyMesh& mesh,
639  labelUList& connectivity,
640  labelUList& offsets,
641  labelUList& faces,
642  labelUList& facesOffsets,
643  foamVtkMeshMaps& maps
644 ) const
645 {
646  presizeMaps(maps);
647 
648  populateArrays
649  (
650  mesh,
651  *this,
652  cellTypes,
653  connectivity,
654  offsets,
655  faces,
656  facesOffsets,
657  contentType::XML,
658  maps.cellMap(),
659  maps.additionalIds()
660  );
661 }
662 
663 
665 (
666  const UList<cellShape>& shapes,
668  labelUList& connectivity,
669  labelUList& offsets,
670  labelUList& faces,
671  labelUList& facesOffsets,
672  foamVtkMeshMaps& maps
673 ) const
674 {
675  // Leave as zero-sized so that populateArrays doesn't fill it.
676  List<label> unused;
677 
678  presizeMaps(maps);
679 
680  populateArrays
681  (
682  shapes,
683  *this,
684  cellTypes,
685  connectivity,
686  offsets,
687  unused, // faces
688  unused, // facesOffsets
689  contentType::XML,
690  maps.cellMap(),
691  maps.additionalIds()
692  );
693 }
694 
695 
696 #undef definePopulateInternalMethod
697 #define definePopulateInternalMethod(Type) \
698  \
699  void Foam::vtk::vtuSizing::populateInternal \
700  ( \
701  const polyMesh& mesh, \
702  UList<uint8_t>& cellTypes, \
703  UList<Type>& connectivity, \
704  UList<Type>& offsets, \
705  UList<Type>& faces, \
706  UList<Type>& facesOffsets, \
707  foamVtkMeshMaps& maps, \
708  const enum contentType output \
709  ) const \
710  { \
711  presizeMaps(maps); \
712  \
713  populateArrays \
714  ( \
715  mesh, \
716  *this, \
717  cellTypes, \
718  connectivity, \
719  offsets, \
720  faces, \
721  facesOffsets, \
722  output, \
723  maps.cellMap(), \
724  maps.additionalIds() \
725  ); \
726  } \
727  \
728  void Foam::vtk::vtuSizing::populateInternal \
729  ( \
730  const polyMesh& mesh, \
731  UList<uint8_t>& cellTypes, \
732  UList<Type>& connectivity, \
733  UList<Type>& offsets, \
734  UList<Type>& faces, \
735  UList<Type>& facesOffsets, \
736  labelUList& cellMap, \
737  labelUList& addPointsIds, \
738  const enum contentType output \
739  ) const \
740  { \
741  populateArrays \
742  ( \
743  mesh, \
744  *this, \
745  cellTypes, \
746  connectivity, \
747  offsets, \
748  faces, \
749  facesOffsets, \
750  output, \
751  cellMap, \
752  addPointsIds \
753  ); \
754  }
755 
756 
760 
761 
762 #undef definePopulateInternalMethod
763 
764 
765 // * * * * * * * * * * * * * * Renumber vertices * * * * * * * * * * * * * * //
766 
768 (
769  const labelUList& vertLabels,
770  const label globalPointOffset
771 )
772 {
773  if (!globalPointOffset)
774  {
775  return vertLabels;
776  }
777 
778  labelList output(vertLabels);
779  renumberVertLabelsLegacy(output, globalPointOffset);
780 
781  return output;
782 }
783 
784 
786 (
787  labelUList& vertLabels,
788  const label globalPointOffset
789 )
790 {
791  if (!globalPointOffset)
792  {
793  return;
794  }
795 
796  // LEGACY vertLabels = "cells" contains
797  // - connectivity
798  // [nLabels, vertex labels...]
799  // - face-stream
800  // [nLabels nFaces, nFace0Pts, id1,id2,..., nFace1Pts, id1,id2,...]
801 
802  // Note the simplest volume cell is a tet (4 points, 4 faces)
803  // As a poly-face stream this would have
804  // 2 for nLabels, nFaces
805  // 4 labels (size + ids) per face * 4 == 16 labels
806  //
807  // Therefore anything with 18 labels or more must be a poly
808 
809  auto iter = vertLabels.begin();
810  const auto last = vertLabels.end();
811 
812  while (iter < last)
813  {
814  label nLabels = *iter; // nLabels (for this cell)
815  ++iter;
816 
817  if (nLabels < 18)
818  {
819  // Normal primitive type
820 
821  while (nLabels--)
822  {
823  *iter += globalPointOffset;
824  ++iter;
825  }
826  }
827  else
828  {
829  // Polyhedral face-stream (explained above)
830 
831  label nFaces = *iter;
832  ++iter;
833 
834  while (nFaces--)
835  {
836  nLabels = *iter; // nLabels (for this face)
837  ++iter;
838 
839  while (nLabels--)
840  {
841  *iter += globalPointOffset;
842  ++iter;
843  }
844  }
845  }
846  }
847 }
848 
849 
851 (
852  const labelUList& vertLabels,
853  const label globalPointOffset
854 )
855 {
856  if (!globalPointOffset)
857  {
858  return vertLabels;
859  }
860 
861  labelList output(vertLabels);
862  renumberVertLabelsXml(output, globalPointOffset);
863 
864  return output;
865 }
866 
867 
869 (
870  labelUList& vertLabels,
871  const label globalPointOffset
872 )
873 {
874  if (!globalPointOffset)
875  {
876  return;
877  }
878 
879  // XML vertLabels = "connectivity" contains
880  // [cell1-verts, cell2-verts, ...]
881 
882  for (label& vertId : vertLabels)
883  {
884  vertId += globalPointOffset;
885  }
886 }
887 
888 
890 (
891  const labelUList& faceLabels,
892  const label globalPointOffset
893 )
894 {
895  if (!globalPointOffset)
896  {
897  return faceLabels;
898  }
899 
900  labelList output(faceLabels);
901  renumberFaceLabelsXml(output, globalPointOffset);
902 
903  return output;
904 }
905 
906 
908 (
909  labelUList& faceLabels,
910  const label globalPointOffset
911 )
912 {
913  if (!globalPointOffset)
914  {
915  return;
916  }
917 
918  // XML face-stream
919  // [nFaces, nFace0Pts, id1,id2,..., nFace1Pts, id1,id2,...]
920 
921  auto iter = faceLabels.begin();
922  const auto last = faceLabels.end();
923 
924  while (iter < last)
925  {
926  label nFaces = *iter;
927  ++iter;
928 
929  while (nFaces--)
930  {
931  label nLabels = *iter;
932  ++iter;
933 
934  while (nLabels--)
935  {
936  *iter += globalPointOffset;
937  ++iter;
938  }
939  }
940  }
941 }
942 
943 
945 (
946  const labelUList& faceOffsets,
947  const label prevOffset
948 )
949 {
950  if (!prevOffset)
951  {
952  return faceOffsets;
953  }
954 
955  labelList output(faceOffsets);
956  renumberFaceOffsetsXml(output, prevOffset);
957 
958  return output;
959 }
960 
961 
963 (
964  labelUList& faceOffsets,
965  const label prevOffset
966 )
967 {
968  if (!prevOffset)
969  {
970  return;
971  }
972 
973  // offsets
974  // [-1, off1, off2, ... -1, ..]
975 
976  for (label& val : faceOffsets)
977  {
978  if (val != -1)
979  {
980  val += prevOffset;
981  }
982  }
983 }
984 
985 
986 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
987 
989 {
990  os << "nFieldCells:" << nFieldCells();
991  if (nAddCells_)
992  {
993  os << " (" << nCells_ << "+" << nAddCells_ << ")";
994  }
995  else
996  {
997  os << " (poly:" << nCellsPoly_ << ")";
998  }
999 
1000  os << " nFieldPoints:" << nFieldPoints();
1001  if (nAddPoints_)
1002  {
1003  os << " (" << nPoints_ << "+" << nAddPoints_ << ")";
1004  }
1005 
1006  os << " nVertLabels:" << (nVertLabels_ + nAddVerts_);
1007  if (nAddVerts_)
1008  {
1009  os << " (" << nVertLabels_ << "+" << nAddVerts_ << ")";
1010  }
1011  else if (nVertPoly_)
1012  {
1013  os << " (poly:" << nVertPoly_ << ")";
1014  }
1015 
1016  os << " nFaceLabels:" << nFaceLabels_;
1017  os << " legacy-count:" << sizeLegacy();
1018 }
1019 
1020 
1021 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1022 
1024 {
1025  return
1026  (
1027  decompose() == rhs.decompose()
1028  // required? && pointOffset() == rhs.pointOffset()
1029  && nCells() == rhs.nCells()
1030  && nPoints() == rhs.nPoints()
1031  && nVertLabels() == rhs.nVertLabels()
1032  && nFaceLabels() == rhs.nFaceLabels()
1033  && nCellsPoly() == rhs.nCellsPoly()
1034  && nVertPoly() == rhs.nVertPoly()
1035  && nAddCells() == rhs.nAddCells()
1036  && nAddPoints() == rhs.nAddPoints()
1037  && nAddVerts() == rhs.nAddVerts()
1038  );
1039 }
1040 
1041 
1043 {
1044  return !operator==(rhs);
1045 }
1046 
1047 
1048 // ************************************************************************* //
Foam::vtk::vtuSizing::nAddCells
label nAddCells() const noexcept
Number of additional (decomposed) cells for the mesh.
Definition: foamVtuSizingI.H:81
Foam::vtk::vtuSizing::renumberFaceLabelsXml
static void renumberFaceLabelsXml(labelUList &faceLabels, const label globalPointOffset)
Renumber faces stream labels by global point offset - XML format.
Definition: foamVtuSizing.C:908
Foam::vtk::vtuSizing::nCells
label nCells() const noexcept
Number of cells for the mesh.
Definition: foamVtuSizingI.H:45
Foam::vtk::vtuSizing::nVertLabels
label nVertLabels() const noexcept
Number of vertex labels for the mesh.
Definition: foamVtuSizingI.H:57
Foam::output
static Ostream & output(Ostream &os, const IntRange< T > &range)
Definition: IntRanges.C:66
Foam::vtk::vtuSizing::resetShapes
void resetShapes(const UList< cellShape > &shapes)
Reset sizing using primitive shapes only (ADVANCED USAGE)
Definition: foamVtuSizing.C:404
Foam::UList::end
iterator end() noexcept
Return an iterator to end traversing the UList.
Definition: UListI.H:350
Foam::vtk::vtuSizing::populateLegacy
void populateLegacy(const polyMesh &mesh, UList< uint8_t > &cellTypes, labelUList &connectivity, foamVtkMeshMaps &maps) const
Populate lists for Legacy output.
Definition: foamVtuSizing.C:578
Foam::vtk::vtuSizing::nFaceLabels
label nFaceLabels() const noexcept
Number of polyhedral face labels for the mesh.
Definition: foamVtuSizingI.H:63
Foam::vtk::vtuSizing::contentType
contentType
Types of content that the storage may represent.
Definition: foamVtuSizing.H:208
Foam::cellModel::HEX
hex
Definition: cellModel.H:81
Foam::foamVtkMeshMaps::additionalIds
const labelList & additionalIds() const noexcept
Any additional (user) labels.
Definition: foamVtkMeshMapsI.H:69
Foam::vtk::vtuSizing::info
void info(Ostream &os) const
Report some information.
Definition: foamVtuSizing.C:988
Foam::vtk::vtuSizing::nAddVerts
label nAddVerts() const noexcept
Number of additional (decomposed) vertices for the mesh.
Definition: foamVtuSizingI.H:93
Foam::vtk::vtuSizing::nAddPoints
label nAddPoints() const noexcept
Number of additional (decomposed) points for the mesh.
Definition: foamVtuSizingI.H:87
Foam::vtk::dataArrayAttr::FACES
"faces"
polyMesh.H
Foam::HashSet< label, Hash< label > >
foamVtuSizingImpl.C
Foam::cellModel::TET
tet
Definition: cellModel.H:85
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
foamVtuSizing.H
Foam::vtk::vtuSizing::copyVertLabelsLegacy
static labelList copyVertLabelsLegacy(const labelUList &connectivity, const label globalPointOffset)
Copy vertex labels with a global point offset - legacy format.
Definition: foamVtuSizing.C:768
Foam::vtk::vtuSizing::operator!=
bool operator!=(const vtuSizing &rhs) const
Test inequality.
Definition: foamVtuSizing.C:1042
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::vtk::fileTag::CELLS
"Cells"
CHECK_SIZING
#define CHECK_SIZING(what, sizeInput, sizeExpected)
Foam::UList::begin
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:329
Foam::vtk::vtuSizing::populateXml
void populateXml(const polyMesh &mesh, UList< uint8_t > &cellTypes, labelUList &connectivity, labelUList &offsets, labelUList &faces, labelUList &facesOffsets, foamVtkMeshMaps &maps) const
Populate lists for XML output.
Definition: foamVtuSizing.C:636
Foam::operator==
tmp< faMatrix< Type > > operator==(const faMatrix< Type > &, const faMatrix< Type > &)
Foam::vtk::vtuSizing::nCellsPoly
label nCellsPoly() const noexcept
Number of polyhedral cells for the mesh.
Definition: foamVtuSizingI.H:69
Foam::vtk::vtuSizing::nFieldCells
label nFieldCells() const noexcept
Number of field cells = nCells + nAddCells.
Definition: foamVtuSizingI.H:99
Foam::cellModel::PRISM
prism
Definition: cellModel.H:83
Foam::cellModel::ref
static const cellModel & ref(const modelType model)
Look up reference to cellModel by enumeration. Fatal on failure.
Definition: cellModels.C:157
Foam::vtk::vtuSizing::populateShapesLegacy
void populateShapesLegacy(const UList< cellShape > &shapes, UList< uint8_t > &cellTypes, labelUList &connectivity, foamVtkMeshMaps &maps) const
Reset list for primitive shapes only (ADVANCED USAGE)
Definition: foamVtuSizing.C:607
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::vtk::vtuSizing::renumberVertLabelsLegacy
static void renumberVertLabelsLegacy(labelUList &connectivity, const label globalPointOffset)
Renumber vertex labels by global point offset - legacy format.
Definition: foamVtuSizing.C:786
Foam::notNull
bool notNull(const T *ptr)
True if ptr is not a pointer (of type T) to the nullObject.
Definition: nullObject.H:207
Foam::vtk::vtuSizing::slotType
slotType
The possible storage 'slots' that can be used.
Definition: foamVtuSizing.H:217
Foam::FatalError
error FatalError
Foam::vtk::vtuSizing::decompose
bool decompose() const noexcept
Query the decompose flag (normally off)
Definition: foamVtuSizingI.H:32
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::cellShape
An analytical geometric cellShape.
Definition: cellShape.H:69
Foam::vtk::vtuSizing::nVertPoly
label nVertPoly() const noexcept
Number of vertex labels for polyhedral cells of the mesh.
Definition: foamVtuSizingI.H:75
Foam::vtk::vtuSizing
Sizing descriptions and routines for transcribing an OpenFOAM volume mesh into a VTK unstructured gri...
Definition: foamVtuSizing.H:201
Foam::vtk::vtuSizing::renumberFaceOffsetsXml
static void renumberFaceOffsetsXml(labelUList &faceOffsets, const label prevOffset)
Renumber face offsets with an offset from previous - XML format.
Definition: foamVtuSizing.C:963
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::hex
IOstream & hex(IOstream &io)
Definition: IOstream.H:446
reset
meshPtr reset(new Foam::fvMesh(Foam::IOobject(regionName, runTime.timeName(), runTime, Foam::IOobject::MUST_READ), false))
Foam::cellModel::WEDGE
wedge
Definition: cellModel.H:82
Foam::vtk::vtuSizing::clear
void clear() noexcept
Reset all sizes to zero.
Definition: foamVtuSizing.C:222
Foam::vtk::vtuSizing::copyVertLabelsXml
static labelList copyVertLabelsXml(const labelUList &connectivity, const label globalPointOffset)
Copy vertex labels with a global point offset - XML format.
Definition: foamVtuSizing.C:851
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::cellModel::PYR
pyr
Definition: cellModel.H:84
cellTypes
const labelList & cellTypes
Definition: setCellMask.H:33
clear
patchWriters clear()
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::vtk::vtuSizing::copyFaceOffsetsXml
static labelList copyFaceOffsetsXml(const labelUList &faceOffsets, const label prevOffset)
Copy face offsets with an offset from previous - XML format.
Definition: foamVtuSizing.C:945
f
labelList f(nPoints)
Foam::List< cellShape >
Foam::vtk::vtuSizing::renumberVertLabelsXml
static void renumberVertLabelsXml(labelUList &connectivity, const label globalPointOffset)
Renumber vertex labels by global point offset - XML format.
Definition: foamVtuSizing.C:869
Foam::vtk::vtuSizing::copyFaceLabelsXml
static labelList copyFaceLabelsXml(const labelUList &faceLabels, const label globalPointOffset)
Copy faces stream labels with a global point offset - XML format.
Definition: foamVtuSizing.C:890
Foam::UList< label >
foamVtkCore.H
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:191
Foam::vtk::vtuSizing::populateShapesXml
void populateShapesXml(const UList< cellShape > &shapes, UList< uint8_t > &cellTypes, labelUList &connectivity, labelUList &offsets, labelUList &faces, labelUList &facesOffsets, foamVtkMeshMaps &maps) const
Reset list for primitive shapes only (ADVANCED USAGE)
Definition: foamVtuSizing.C:665
definePopulateInternalMethod
#define definePopulateInternalMethod(Type)
Definition: foamVtuSizing.C:697
cellShape.H
Foam::vtk::vtuSizing::vtuSizing
vtuSizing() noexcept
Default construct.
Definition: foamVtuSizing.C:203
Foam::cellModel::TETWEDGE
tetWedge
Definition: cellModel.H:87
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::cellModel
Maps a geometry to a set of cell primitives.
Definition: cellModel.H:72
Foam::UList::size
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
Foam::vtk::vtuSizing::nPoints
label nPoints() const noexcept
Number of points for the mesh.
Definition: foamVtuSizingI.H:51
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::UList::null
static const UList< T > & null()
Return a UList reference to a nullObject.
Definition: UListI.H:53
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::vtk::vtuSizing::operator==
bool operator==(const vtuSizing &rhs) const
Test equality.
Definition: foamVtuSizing.C:1023
Foam::vtk::vtuSizing::reset
void reset(const polyMesh &mesh, const bool decompose=false)
Reset sizing by analyzing the mesh.
Definition: foamVtuSizing.C:241
Foam::foamVtkMeshMaps
Bookkeeping for mesh subsetting and/or polyhedral cell decomposition. Although the main use case is f...
Definition: foamVtkMeshMaps.H:58
Foam::foamVtkMeshMaps::cellMap
const labelList & cellMap() const noexcept
Original cell ids for all cells (regular and decomposed).
Definition: foamVtkMeshMapsI.H:41
Foam::vtk::vtuSizing::sizeOf
label sizeOf(const enum contentType output, const enum slotType slot) const
Return the required size for the storage slot.
Definition: foamVtuSizing.C:471
Foam::cellShape::model
const cellModel & model() const
Model reference.
Definition: cellShapeI.H:126