xref: /petsc/src/dm/impls/plex/plexvtk.c (revision 07c2e4feb6773e78bda63e3a89d5b841667f9670)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I   "petscdmplex.h"   I*/
2552f7358SJed Brown #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h>
3552f7358SJed Brown 
DMPlexVTKGetCellType_Internal(DM dm,PetscInt dim,PetscInt corners,PetscInt * cellType)4d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexVTKGetCellType_Internal(DM dm, PetscInt dim, PetscInt corners, PetscInt *cellType)
5d71ae5a4SJacob Faibussowitsch {
6552f7358SJed Brown   PetscFunctionBegin;
7552f7358SJed Brown   *cellType = -1;
8552f7358SJed Brown   switch (dim) {
9552f7358SJed Brown   case 0:
10552f7358SJed Brown     switch (corners) {
11552f7358SJed Brown     case 1:
12552f7358SJed Brown       *cellType = 1; /* VTK_VERTEX */
13552f7358SJed Brown       break;
14d71ae5a4SJacob Faibussowitsch     default:
15d71ae5a4SJacob Faibussowitsch       break;
16552f7358SJed Brown     }
17552f7358SJed Brown     break;
18552f7358SJed Brown   case 1:
19552f7358SJed Brown     switch (corners) {
20552f7358SJed Brown     case 2:
21552f7358SJed Brown       *cellType = 3; /* VTK_LINE */
22552f7358SJed Brown       break;
23552f7358SJed Brown     case 3:
24552f7358SJed Brown       *cellType = 21; /* VTK_QUADRATIC_EDGE */
25552f7358SJed Brown       break;
26d71ae5a4SJacob Faibussowitsch     default:
27d71ae5a4SJacob Faibussowitsch       break;
28552f7358SJed Brown     }
29552f7358SJed Brown     break;
30552f7358SJed Brown   case 2:
31552f7358SJed Brown     switch (corners) {
32552f7358SJed Brown     case 3:
33552f7358SJed Brown       *cellType = 5; /* VTK_TRIANGLE */
34552f7358SJed Brown       break;
35552f7358SJed Brown     case 4:
36552f7358SJed Brown       *cellType = 9; /* VTK_QUAD */
37552f7358SJed Brown       break;
38552f7358SJed Brown     case 6:
39552f7358SJed Brown       *cellType = 22; /* VTK_QUADRATIC_TRIANGLE */
40552f7358SJed Brown       break;
41552f7358SJed Brown     case 9:
42552f7358SJed Brown       *cellType = 23; /* VTK_QUADRATIC_QUAD */
43552f7358SJed Brown       break;
44d71ae5a4SJacob Faibussowitsch     default:
45d71ae5a4SJacob Faibussowitsch       break;
46552f7358SJed Brown     }
47552f7358SJed Brown     break;
48552f7358SJed Brown   case 3:
49552f7358SJed Brown     switch (corners) {
50552f7358SJed Brown     case 4:
51552f7358SJed Brown       *cellType = 10; /* VTK_TETRA */
52552f7358SJed Brown       break;
53cf4edabeSMatthew G. Knepley     case 5:
54cf4edabeSMatthew G. Knepley       *cellType = 14; /* VTK_PYRAMID */
55cf4edabeSMatthew G. Knepley       break;
562f029394SStefano Zampini     case 6:
572f029394SStefano Zampini       *cellType = 13; /* VTK_WEDGE */
582f029394SStefano Zampini       break;
59552f7358SJed Brown     case 8:
60552f7358SJed Brown       *cellType = 12; /* VTK_HEXAHEDRON */
61552f7358SJed Brown       break;
62552f7358SJed Brown     case 10:
63552f7358SJed Brown       *cellType = 24; /* VTK_QUADRATIC_TETRA */
64552f7358SJed Brown       break;
65552f7358SJed Brown     case 27:
66552f7358SJed Brown       *cellType = 29; /* VTK_QUADRATIC_HEXAHEDRON */
67552f7358SJed Brown       break;
68d71ae5a4SJacob Faibussowitsch     default:
69d71ae5a4SJacob Faibussowitsch       break;
70552f7358SJed Brown     }
71552f7358SJed Brown   }
723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
73552f7358SJed Brown }
74552f7358SJed Brown 
75ffeef943SBarry Smith /*@
76552f7358SJed Brown   DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
77552f7358SJed Brown 
78552f7358SJed Brown   Collective
79552f7358SJed Brown 
804165533cSJose E. Roman   Input Parameters:
81a1cb98faSBarry Smith + odm    - The `DMPLEX` specifying the mesh, passed as a `PetscObject`
82a1cb98faSBarry Smith - viewer - viewer of type `PETSCVIEWERVTK`
83552f7358SJed Brown 
84552f7358SJed Brown   Level: developer
85552f7358SJed Brown 
86552f7358SJed Brown   Note:
87552f7358SJed Brown   This function is a callback used by the VTK viewer to actually write the file.
88552f7358SJed 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.
89552f7358SJed Brown   Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
90552f7358SJed Brown 
911cc06b55SBarry Smith .seealso: [](ch_unstructured), `DM`, `PETSCVIEWEREXODUSII`, `DMPLEX`, `PETSCVIEWERVTK`
92552f7358SJed Brown @*/
DMPlexVTKWriteAll(PetscObject odm,PetscViewer viewer)93d71ae5a4SJacob Faibussowitsch PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer)
94d71ae5a4SJacob Faibussowitsch {
95552f7358SJed Brown   DM dm = (DM)odm;
96552f7358SJed Brown 
97552f7358SJed Brown   PetscFunctionBegin;
98552f7358SJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
99a1cb98faSBarry Smith   PetscValidHeaderSpecificType(viewer, PETSC_VIEWER_CLASSID, 2, PETSCVIEWERVTK);
100552f7358SJed Brown   switch (viewer->format) {
101*ce78bad3SBarry Smith   case PETSC_VIEWER_DEFAULT:
102d71ae5a4SJacob Faibussowitsch   case PETSC_VIEWER_VTK_VTU:
103d71ae5a4SJacob Faibussowitsch     PetscCall(DMPlexVTKWriteAll_VTU(dm, viewer));
104d71ae5a4SJacob Faibussowitsch     break;
105d71ae5a4SJacob Faibussowitsch   default:
106d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]);
107552f7358SJed Brown   }
1083ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
109552f7358SJed Brown }
110