xref: /petsc/src/dm/impls/plex/plexvtk.c (revision ce78bad369055609e946c9d2c25ea67a45873e27)
1552f7358SJed Brown #define PETSCDM_DLL
2af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I   "petscdmplex.h"   I*/
3552f7358SJed Brown #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h>
4552f7358SJed Brown 
5d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexVTKGetCellType_Internal(DM dm, PetscInt dim, PetscInt corners, PetscInt *cellType)
6d71ae5a4SJacob Faibussowitsch {
7552f7358SJed Brown   PetscFunctionBegin;
8552f7358SJed Brown   *cellType = -1;
9552f7358SJed Brown   switch (dim) {
10552f7358SJed Brown   case 0:
11552f7358SJed Brown     switch (corners) {
12552f7358SJed Brown     case 1:
13552f7358SJed Brown       *cellType = 1; /* VTK_VERTEX */
14552f7358SJed Brown       break;
15d71ae5a4SJacob Faibussowitsch     default:
16d71ae5a4SJacob Faibussowitsch       break;
17552f7358SJed Brown     }
18552f7358SJed Brown     break;
19552f7358SJed Brown   case 1:
20552f7358SJed Brown     switch (corners) {
21552f7358SJed Brown     case 2:
22552f7358SJed Brown       *cellType = 3; /* VTK_LINE */
23552f7358SJed Brown       break;
24552f7358SJed Brown     case 3:
25552f7358SJed Brown       *cellType = 21; /* VTK_QUADRATIC_EDGE */
26552f7358SJed Brown       break;
27d71ae5a4SJacob Faibussowitsch     default:
28d71ae5a4SJacob Faibussowitsch       break;
29552f7358SJed Brown     }
30552f7358SJed Brown     break;
31552f7358SJed Brown   case 2:
32552f7358SJed Brown     switch (corners) {
33552f7358SJed Brown     case 3:
34552f7358SJed Brown       *cellType = 5; /* VTK_TRIANGLE */
35552f7358SJed Brown       break;
36552f7358SJed Brown     case 4:
37552f7358SJed Brown       *cellType = 9; /* VTK_QUAD */
38552f7358SJed Brown       break;
39552f7358SJed Brown     case 6:
40552f7358SJed Brown       *cellType = 22; /* VTK_QUADRATIC_TRIANGLE */
41552f7358SJed Brown       break;
42552f7358SJed Brown     case 9:
43552f7358SJed Brown       *cellType = 23; /* VTK_QUADRATIC_QUAD */
44552f7358SJed Brown       break;
45d71ae5a4SJacob Faibussowitsch     default:
46d71ae5a4SJacob Faibussowitsch       break;
47552f7358SJed Brown     }
48552f7358SJed Brown     break;
49552f7358SJed Brown   case 3:
50552f7358SJed Brown     switch (corners) {
51552f7358SJed Brown     case 4:
52552f7358SJed Brown       *cellType = 10; /* VTK_TETRA */
53552f7358SJed Brown       break;
54cf4edabeSMatthew G. Knepley     case 5:
55cf4edabeSMatthew G. Knepley       *cellType = 14; /* VTK_PYRAMID */
56cf4edabeSMatthew G. Knepley       break;
572f029394SStefano Zampini     case 6:
582f029394SStefano Zampini       *cellType = 13; /* VTK_WEDGE */
592f029394SStefano Zampini       break;
60552f7358SJed Brown     case 8:
61552f7358SJed Brown       *cellType = 12; /* VTK_HEXAHEDRON */
62552f7358SJed Brown       break;
63552f7358SJed Brown     case 10:
64552f7358SJed Brown       *cellType = 24; /* VTK_QUADRATIC_TETRA */
65552f7358SJed Brown       break;
66552f7358SJed Brown     case 27:
67552f7358SJed Brown       *cellType = 29; /* VTK_QUADRATIC_HEXAHEDRON */
68552f7358SJed Brown       break;
69d71ae5a4SJacob Faibussowitsch     default:
70d71ae5a4SJacob Faibussowitsch       break;
71552f7358SJed Brown     }
72552f7358SJed Brown   }
733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
74552f7358SJed Brown }
75552f7358SJed Brown 
76ffeef943SBarry Smith /*@
77552f7358SJed Brown   DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
78552f7358SJed Brown 
79552f7358SJed Brown   Collective
80552f7358SJed Brown 
814165533cSJose E. Roman   Input Parameters:
82a1cb98faSBarry Smith + odm    - The `DMPLEX` specifying the mesh, passed as a `PetscObject`
83a1cb98faSBarry Smith - viewer - viewer of type `PETSCVIEWERVTK`
84552f7358SJed Brown 
85552f7358SJed Brown   Level: developer
86552f7358SJed Brown 
87552f7358SJed Brown   Note:
88552f7358SJed Brown   This function is a callback used by the VTK viewer to actually write the file.
89552f7358SJed Brown   The reason for this odd model is that the VTK file format does not provide any way to write one field at a time.
90552f7358SJed Brown   Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
91552f7358SJed Brown 
921cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `PETSCVIEWEREXODUSII`, `DMPLEX`, `PETSCVIEWERVTK`
93552f7358SJed Brown @*/
94d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer)
95d71ae5a4SJacob Faibussowitsch {
96552f7358SJed Brown   DM dm = (DM)odm;
97552f7358SJed Brown 
98552f7358SJed Brown   PetscFunctionBegin;
99552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
100a1cb98faSBarry Smith   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 2, PETSCVIEWERVTK);
101552f7358SJed Brown   switch (viewer->format) {
102*ce78bad3SBarry Smith   case PETSC_VIEWER_DEFAULT:
103d71ae5a4SJacob Faibussowitsch   case PETSC_VIEWER_VTK_VTU:
104d71ae5a4SJacob Faibussowitsch     PetscCall(DMPlexVTKWriteAll_VTU(dm, viewer));
105d71ae5a4SJacob Faibussowitsch     break;
106d71ae5a4SJacob Faibussowitsch   default:
107d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]);
108552f7358SJed Brown   }
1093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
110552f7358SJed Brown }
111