Actual source code: plexvtk.c

  1: #define PETSCDM_DLL
  2: #include <petsc/private/dmpleximpl.h>
  3: #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h>

  5: PetscErrorCode DMPlexVTKGetCellType_Internal(DM dm, PetscInt dim, PetscInt corners, PetscInt *cellType)
  6: {
  7:   PetscFunctionBegin;
  8:   *cellType = -1;
  9:   switch (dim) {
 10:   case 0:
 11:     switch (corners) {
 12:     case 1:
 13:       *cellType = 1; /* VTK_VERTEX */
 14:       break;
 15:     default:
 16:       break;
 17:     }
 18:     break;
 19:   case 1:
 20:     switch (corners) {
 21:     case 2:
 22:       *cellType = 3; /* VTK_LINE */
 23:       break;
 24:     case 3:
 25:       *cellType = 21; /* VTK_QUADRATIC_EDGE */
 26:       break;
 27:     default:
 28:       break;
 29:     }
 30:     break;
 31:   case 2:
 32:     switch (corners) {
 33:     case 3:
 34:       *cellType = 5; /* VTK_TRIANGLE */
 35:       break;
 36:     case 4:
 37:       *cellType = 9; /* VTK_QUAD */
 38:       break;
 39:     case 6:
 40:       *cellType = 22; /* VTK_QUADRATIC_TRIANGLE */
 41:       break;
 42:     case 9:
 43:       *cellType = 23; /* VTK_QUADRATIC_QUAD */
 44:       break;
 45:     default:
 46:       break;
 47:     }
 48:     break;
 49:   case 3:
 50:     switch (corners) {
 51:     case 4:
 52:       *cellType = 10; /* VTK_TETRA */
 53:       break;
 54:     case 5:
 55:       *cellType = 14; /* VTK_PYRAMID */
 56:       break;
 57:     case 6:
 58:       *cellType = 13; /* VTK_WEDGE */
 59:       break;
 60:     case 8:
 61:       *cellType = 12; /* VTK_HEXAHEDRON */
 62:       break;
 63:     case 10:
 64:       *cellType = 24; /* VTK_QUADRATIC_TETRA */
 65:       break;
 66:     case 27:
 67:       *cellType = 29; /* VTK_QUADRATIC_HEXAHEDRON */
 68:       break;
 69:     default:
 70:       break;
 71:     }
 72:   }
 73:   PetscFunctionReturn(PETSC_SUCCESS);
 74: }

 76: /*@
 77:   DMPlexVTKWriteAll - Write a file containing all the fields that have been provided to the viewer

 79:   Collective

 81:   Input Parameters:
 82: + odm    - The `DMPLEX` specifying the mesh, passed as a `PetscObject`
 83: - viewer - viewer of type `PETSCVIEWERVTK`

 85:   Level: developer

 87:   Note:
 88:   This function is a callback used by the VTK viewer to actually write the file.
 89:   The reason for this odd model is that the VTK file format does not provide any way to write one field at a time.
 90:   Instead, metadata for the entire file needs to be available up-front before you can start writing the file.

 92: .seealso: [](ch_unstructured), `DM`, `PETSCVIEWEREXODUSII`, `DMPLEX`, `PETSCVIEWERVTK`
 93: @*/
 94: PetscErrorCode DMPlexVTKWriteAll(PetscObject odm, PetscViewer viewer)
 95: {
 96:   DM dm = (DM)odm;

 98:   PetscFunctionBegin;
101:   switch (viewer->format) {
102:   case PETSC_VIEWER_DEFAULT:
103:   case PETSC_VIEWER_VTK_VTU:
104:     PetscCall(DMPlexVTKWriteAll_VTU(dm, viewer));
105:     break;
106:   default:
107:     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]);
108:   }
109:   PetscFunctionReturn(PETSC_SUCCESS);
110: }