xref: /petsc/src/dm/impls/da/grvtk.c (revision d8e47b638cf8f604a99e9678e1df24f82d959cd7)
1ffeef943SBarry Smith #include <petsc/private/dmdaimpl.h> /*I   "petscdmda.h"    I*/
294fbd55eSBarry Smith /*
394fbd55eSBarry Smith    Note that the API for using PETSCVIEWERVTK is totally wrong since its use requires
494fbd55eSBarry Smith    including the private vtkvimpl.h file. The code should be refactored.
594fbd55eSBarry Smith */
65c6c1daeSBarry Smith #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h>
74061b8bfSJed Brown 
8fb5bd1c2SPatrick Sanan /* Helper function which determines if any DMDA fields are named.  This is used
9fb5bd1c2SPatrick Sanan    as a proxy for the user's intention to use DMDA fields as distinct
10fb5bd1c2SPatrick Sanan    scalar-valued fields as opposed to a single vector-valued field */
DMDAGetFieldsNamed(DM da,PetscBool * fieldsnamed)11d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAGetFieldsNamed(DM da, PetscBool *fieldsnamed)
12d71ae5a4SJacob Faibussowitsch {
13fb5bd1c2SPatrick Sanan   PetscInt f, bs;
14fb5bd1c2SPatrick Sanan 
15fb5bd1c2SPatrick Sanan   PetscFunctionBegin;
16fb5bd1c2SPatrick Sanan   *fieldsnamed = PETSC_FALSE;
179566063dSJacob Faibussowitsch   PetscCall(DMDAGetDof(da, &bs));
18fb5bd1c2SPatrick Sanan   for (f = 0; f < bs; ++f) {
19fb5bd1c2SPatrick Sanan     const char *fieldname;
209566063dSJacob Faibussowitsch     PetscCall(DMDAGetFieldName(da, f, &fieldname));
21fb5bd1c2SPatrick Sanan     if (fieldname) {
22fb5bd1c2SPatrick Sanan       *fieldsnamed = PETSC_TRUE;
23fb5bd1c2SPatrick Sanan       break;
24fb5bd1c2SPatrick Sanan     }
25fb5bd1c2SPatrick Sanan   }
263ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27fb5bd1c2SPatrick Sanan }
28fb5bd1c2SPatrick Sanan 
DMDAVTKWriteAll_VTS(DM da,PetscViewer viewer)29d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAVTKWriteAll_VTS(DM da, PetscViewer viewer)
30d71ae5a4SJacob Faibussowitsch {
3130815ce0SLisandro Dalcin   const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian";
324061b8bfSJed Brown #if defined(PETSC_USE_REAL_SINGLE)
334061b8bfSJed Brown   const char precision[] = "Float32";
344061b8bfSJed Brown #elif defined(PETSC_USE_REAL_DOUBLE)
354061b8bfSJed Brown   const char precision[] = "Float64";
364061b8bfSJed Brown #else
374061b8bfSJed Brown   const char precision[] = "UnknownPrecision";
384061b8bfSJed Brown #endif
39ce94432eSBarry Smith   MPI_Comm                 comm;
402eaa9ef4SLisandro Dalcin   Vec                      Coords;
414061b8bfSJed Brown   PetscViewer_VTK         *vtk = (PetscViewer_VTK *)viewer->data;
424061b8bfSJed Brown   PetscViewerVTKObjectLink link;
434061b8bfSJed Brown   FILE                    *fp;
444061b8bfSJed Brown   PetscMPIInt              rank, size, tag;
454061b8bfSJed Brown   DMDALocalInfo            info;
46*6497c311SBarry Smith   PetscInt                 dim, mx, my, mz, cdim, bs, maxnnodes, maxbs, i, j, k;
470f2609c8SJed Brown   PetscInt64               boffset;
480298fd71SBarry Smith   PetscInt                 rloc[6], (*grloc)[6] = NULL;
494061b8bfSJed Brown   PetscScalar             *array, *array2;
504061b8bfSJed Brown 
514061b8bfSJed Brown   PetscFunctionBegin;
529566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
531dca8a05SBarry Smith   PetscCheck(!PetscDefined(USE_COMPLEX), PETSC_COMM_SELF, PETSC_ERR_SUP, "Complex values not supported");
549566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
559566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
569566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, &mx, &my, &mz, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL));
579566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
589566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(da, &Coords));
592eaa9ef4SLisandro Dalcin   if (Coords) {
602eaa9ef4SLisandro Dalcin     PetscInt csize;
619566063dSJacob Faibussowitsch     PetscCall(VecGetSize(Coords, &csize));
6208401ef6SPierre Jolivet     PetscCheck(csize % (mx * my * mz) == 0, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Coordinate vector size mismatch");
632eaa9ef4SLisandro Dalcin     cdim = csize / (mx * my * mz);
641dca8a05SBarry Smith     PetscCheck(cdim >= dim && cdim <= 3, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Coordinate vector size mismatch");
652eaa9ef4SLisandro Dalcin   } else {
662eaa9ef4SLisandro Dalcin     cdim = dim;
672eaa9ef4SLisandro Dalcin   }
684061b8bfSJed Brown 
699566063dSJacob Faibussowitsch   PetscCall(PetscFOpen(comm, vtk->filename, "wb", &fp));
709566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "<?xml version=\"1.0\"?>\n"));
714a8c153fSJed Brown   PetscCall(PetscFPrintf(comm, fp, "<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"%s\" header_type=\"UInt64\">\n", byte_order));
7263a3b9bcSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "  <StructuredGrid WholeExtent=\"%d %" PetscInt_FMT " %d %" PetscInt_FMT " %d %" PetscInt_FMT "\">\n", 0, mx - 1, 0, my - 1, 0, mz - 1));
734061b8bfSJed Brown 
749566063dSJacob Faibussowitsch   if (rank == 0) PetscCall(PetscMalloc1(size * 6, &grloc));
754061b8bfSJed Brown   rloc[0] = info.xs;
764061b8bfSJed Brown   rloc[1] = info.xm;
774061b8bfSJed Brown   rloc[2] = info.ys;
784061b8bfSJed Brown   rloc[3] = info.ym;
794061b8bfSJed Brown   rloc[4] = info.zs;
804061b8bfSJed Brown   rloc[5] = info.zm;
81810441c8SPierre Jolivet   PetscCallMPI(MPI_Gather(rloc, 6, MPIU_INT, rank == 0 ? grloc[0] : NULL, 6, MPIU_INT, 0, comm));
824061b8bfSJed Brown 
834061b8bfSJed Brown   /* Write XML header */
844061b8bfSJed Brown   maxnnodes = 0; /* Used for the temporary array size on rank 0 */
85ea2d7708SPatrick Sanan   maxbs     = 0; /* Used for the temporary array size on rank 0 */
864061b8bfSJed Brown   boffset   = 0; /* Offset into binary file */
87*6497c311SBarry Smith   for (PetscMPIInt r = 0; r < size; r++) {
884061b8bfSJed Brown     PetscInt xs = -1, xm = -1, ys = -1, ym = -1, zs = -1, zm = -1, nnodes = 0;
89*6497c311SBarry Smith 
90dd400576SPatrick Sanan     if (rank == 0) {
914061b8bfSJed Brown       xs     = grloc[r][0];
924061b8bfSJed Brown       xm     = grloc[r][1];
934061b8bfSJed Brown       ys     = grloc[r][2];
944061b8bfSJed Brown       ym     = grloc[r][3];
954061b8bfSJed Brown       zs     = grloc[r][4];
964061b8bfSJed Brown       zm     = grloc[r][5];
974061b8bfSJed Brown       nnodes = xm * ym * zm;
984061b8bfSJed Brown     }
994061b8bfSJed Brown     maxnnodes = PetscMax(maxnnodes, nnodes);
10063a3b9bcSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "    <Piece Extent=\"%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\">\n", xs, xs + xm - 1, ys, ys + ym - 1, zs, zs + zm - 1));
1019566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "      <Points>\n"));
1020f2609c8SJed Brown     PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, boffset));
1034a8c153fSJed Brown     boffset += 3 * nnodes * sizeof(PetscScalar) + sizeof(PetscInt64);
1049566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "      </Points>\n"));
1054061b8bfSJed Brown 
1069566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "      <PointData Scalars=\"ScalarPointData\">\n"));
1074061b8bfSJed Brown     for (link = vtk->link; link; link = link->next) {
1084061b8bfSJed Brown       Vec         X = (Vec)link->vec;
109fb5bd1c2SPatrick Sanan       PetscInt    bs, f;
110ea2d7708SPatrick Sanan       DM          daCurr;
111fb5bd1c2SPatrick Sanan       PetscBool   fieldsnamed;
112fb5bd1c2SPatrick Sanan       const char *vecname = "Unnamed";
113ea2d7708SPatrick Sanan 
1149566063dSJacob Faibussowitsch       PetscCall(VecGetDM(X, &daCurr));
1159566063dSJacob Faibussowitsch       PetscCall(DMDAGetInfo(daCurr, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL));
116ea2d7708SPatrick Sanan       maxbs = PetscMax(maxbs, bs);
117ea2d7708SPatrick Sanan 
11848a46eb9SPierre Jolivet       if (((PetscObject)X)->name || link != vtk->link) PetscCall(PetscObjectGetName((PetscObject)X, &vecname));
119fb5bd1c2SPatrick Sanan 
120fb5bd1c2SPatrick Sanan       /* If any fields are named, add scalar fields. Otherwise, add a vector field */
1219566063dSJacob Faibussowitsch       PetscCall(DMDAGetFieldsNamed(daCurr, &fieldsnamed));
122fb5bd1c2SPatrick Sanan       if (fieldsnamed) {
123fb5bd1c2SPatrick Sanan         for (f = 0; f < bs; f++) {
124fb5bd1c2SPatrick Sanan           char        buf[256];
125fb5bd1c2SPatrick Sanan           const char *fieldname;
1269566063dSJacob Faibussowitsch           PetscCall(DMDAGetFieldName(daCurr, f, &fieldname));
127fb5bd1c2SPatrick Sanan           if (!fieldname) {
12863a3b9bcSJacob Faibussowitsch             PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, f));
129fb5bd1c2SPatrick Sanan             fieldname = buf;
130fb5bd1c2SPatrick Sanan           }
1310f2609c8SJed Brown           PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset));
1324a8c153fSJed Brown           boffset += nnodes * sizeof(PetscScalar) + sizeof(PetscInt64);
133fb5bd1c2SPatrick Sanan         }
134fb5bd1c2SPatrick Sanan       } else {
1350f2609c8SJed Brown         PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%" PetscInt_FMT "\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, bs, boffset));
1364a8c153fSJed Brown         boffset += bs * nnodes * sizeof(PetscScalar) + sizeof(PetscInt64);
1374061b8bfSJed Brown       }
138fb5bd1c2SPatrick Sanan     }
1399566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "      </PointData>\n"));
1409566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "    </Piece>\n"));
1414061b8bfSJed Brown   }
1429566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "  </StructuredGrid>\n"));
1439566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "  <AppendedData encoding=\"raw\">\n"));
1449566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "_"));
1454061b8bfSJed Brown 
1464061b8bfSJed Brown   /* Now write the arrays. */
1474061b8bfSJed Brown   tag = ((PetscObject)viewer)->tag;
1489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(maxnnodes * PetscMax(3, maxbs), &array, maxnnodes * PetscMax(3, maxbs), &array2));
149*6497c311SBarry Smith   for (PetscMPIInt r = 0; r < size; r++) {
1504061b8bfSJed Brown     MPI_Status status;
1514061b8bfSJed Brown     PetscInt   xs = -1, xm = -1, ys = -1, ym = -1, zs = -1, zm = -1, nnodes = 0;
152*6497c311SBarry Smith 
153dd400576SPatrick Sanan     if (rank == 0) {
1544061b8bfSJed Brown       xs     = grloc[r][0];
1554061b8bfSJed Brown       xm     = grloc[r][1];
1564061b8bfSJed Brown       ys     = grloc[r][2];
1574061b8bfSJed Brown       ym     = grloc[r][3];
1584061b8bfSJed Brown       zs     = grloc[r][4];
1594061b8bfSJed Brown       zm     = grloc[r][5];
1604061b8bfSJed Brown       nnodes = xm * ym * zm;
1614061b8bfSJed Brown     } else if (r == rank) {
1624061b8bfSJed Brown       nnodes = info.xm * info.ym * info.zm;
1634061b8bfSJed Brown     }
1644061b8bfSJed Brown 
1652eaa9ef4SLisandro Dalcin     /* Write the coordinates */
1664061b8bfSJed Brown     if (Coords) {
1674061b8bfSJed Brown       const PetscScalar *coords;
168*6497c311SBarry Smith 
1699566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(Coords, &coords));
170dd400576SPatrick Sanan       if (rank == 0) {
1714061b8bfSJed Brown         if (r) {
1726622f924SJed Brown           PetscMPIInt nn;
173*6497c311SBarry Smith 
174*6497c311SBarry Smith           PetscCallMPI(MPIU_Recv(array, nnodes * cdim, MPIU_SCALAR, r, tag, comm, &status));
1759566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn));
17608401ef6SPierre Jolivet           PetscCheck(nn == nnodes * cdim, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch");
1771baa6e33SBarry Smith         } else PetscCall(PetscArraycpy(array, coords, nnodes * cdim));
1784061b8bfSJed Brown         /* Transpose coordinates to VTK (C-style) ordering */
1794061b8bfSJed Brown         for (k = 0; k < zm; k++) {
1804061b8bfSJed Brown           for (j = 0; j < ym; j++) {
1814061b8bfSJed Brown             for (i = 0; i < xm; i++) {
1824061b8bfSJed Brown               PetscInt Iloc        = i + xm * (j + ym * k);
1832eaa9ef4SLisandro Dalcin               array2[Iloc * 3 + 0] = array[Iloc * cdim + 0];
1842eaa9ef4SLisandro Dalcin               array2[Iloc * 3 + 1] = cdim > 1 ? array[Iloc * cdim + 1] : 0.0;
1852eaa9ef4SLisandro Dalcin               array2[Iloc * 3 + 2] = cdim > 2 ? array[Iloc * cdim + 2] : 0.0;
1864061b8bfSJed Brown             }
1874061b8bfSJed Brown           }
1884061b8bfSJed Brown         }
1894061b8bfSJed Brown       } else if (r == rank) {
190*6497c311SBarry Smith         PetscCallMPI(MPIU_Send((void *)coords, nnodes * cdim, MPIU_SCALAR, 0, tag, comm));
1914061b8bfSJed Brown       }
1929566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(Coords, &coords));
1934061b8bfSJed Brown     } else { /* Fabricate some coordinates using grid index */
1944061b8bfSJed Brown       for (k = 0; k < zm; k++) {
1954061b8bfSJed Brown         for (j = 0; j < ym; j++) {
1964061b8bfSJed Brown           for (i = 0; i < xm; i++) {
1974061b8bfSJed Brown             PetscInt Iloc        = i + xm * (j + ym * k);
1984061b8bfSJed Brown             array2[Iloc * 3 + 0] = xs + i;
1994061b8bfSJed Brown             array2[Iloc * 3 + 1] = ys + j;
2004061b8bfSJed Brown             array2[Iloc * 3 + 2] = zs + k;
2014061b8bfSJed Brown           }
2024061b8bfSJed Brown         }
2034061b8bfSJed Brown       }
2044061b8bfSJed Brown     }
2059566063dSJacob Faibussowitsch     PetscCall(PetscViewerVTKFWrite(viewer, fp, array2, nnodes * 3, MPIU_SCALAR));
2064061b8bfSJed Brown 
2074061b8bfSJed Brown     /* Write each of the objects queued up for this file */
2084061b8bfSJed Brown     for (link = vtk->link; link; link = link->next) {
2094061b8bfSJed Brown       Vec                X = (Vec)link->vec;
2104061b8bfSJed Brown       const PetscScalar *x;
211fb5bd1c2SPatrick Sanan       PetscInt           bs, f;
212ea2d7708SPatrick Sanan       DM                 daCurr;
213fb5bd1c2SPatrick Sanan       PetscBool          fieldsnamed;
2149566063dSJacob Faibussowitsch       PetscCall(VecGetDM(X, &daCurr));
2159566063dSJacob Faibussowitsch       PetscCall(DMDAGetInfo(daCurr, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL));
2169566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(X, &x));
217dd400576SPatrick Sanan       if (rank == 0) {
2184061b8bfSJed Brown         if (r) {
2196622f924SJed Brown           PetscMPIInt nn;
220*6497c311SBarry Smith 
221*6497c311SBarry Smith           PetscCallMPI(MPIU_Recv(array, nnodes * bs, MPIU_SCALAR, r, tag, comm, &status));
2229566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn));
223*6497c311SBarry Smith           PetscCheck(nn == nnodes * bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch receiving from rank %d", r);
2241baa6e33SBarry Smith         } else PetscCall(PetscArraycpy(array, x, nnodes * bs));
225fb5bd1c2SPatrick Sanan 
226fb5bd1c2SPatrick Sanan         /* If any fields are named, add scalar fields. Otherwise, add a vector field */
2279566063dSJacob Faibussowitsch         PetscCall(DMDAGetFieldsNamed(daCurr, &fieldsnamed));
228fb5bd1c2SPatrick Sanan         if (fieldsnamed) {
229fb5bd1c2SPatrick Sanan           for (f = 0; f < bs; f++) {
230fb5bd1c2SPatrick Sanan             /* Extract and transpose the f'th field */
231fb5bd1c2SPatrick Sanan             for (k = 0; k < zm; k++) {
232fb5bd1c2SPatrick Sanan               for (j = 0; j < ym; j++) {
233fb5bd1c2SPatrick Sanan                 for (i = 0; i < xm; i++) {
234fb5bd1c2SPatrick Sanan                   PetscInt Iloc = i + xm * (j + ym * k);
235fb5bd1c2SPatrick Sanan                   array2[Iloc]  = array[Iloc * bs + f];
236fb5bd1c2SPatrick Sanan                 }
237fb5bd1c2SPatrick Sanan               }
238fb5bd1c2SPatrick Sanan             }
2399566063dSJacob Faibussowitsch             PetscCall(PetscViewerVTKFWrite(viewer, fp, array2, nnodes, MPIU_SCALAR));
240fb5bd1c2SPatrick Sanan           }
2411baa6e33SBarry Smith         } else PetscCall(PetscViewerVTKFWrite(viewer, fp, array, bs * nnodes, MPIU_SCALAR));
242*6497c311SBarry Smith       } else if (r == rank) PetscCallMPI(MPIU_Send((void *)x, nnodes * bs, MPIU_SCALAR, 0, tag, comm));
2439566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(X, &x));
2444061b8bfSJed Brown     }
2454061b8bfSJed Brown   }
2469566063dSJacob Faibussowitsch   PetscCall(PetscFree2(array, array2));
2479566063dSJacob Faibussowitsch   PetscCall(PetscFree(grloc));
2484061b8bfSJed Brown 
2499566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "\n </AppendedData>\n"));
2509566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "</VTKFile>\n"));
2519566063dSJacob Faibussowitsch   PetscCall(PetscFClose(comm, fp));
2523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2534061b8bfSJed Brown }
2544061b8bfSJed Brown 
DMDAVTKWriteAll_VTR(DM da,PetscViewer viewer)255d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDAVTKWriteAll_VTR(DM da, PetscViewer viewer)
256d71ae5a4SJacob Faibussowitsch {
25730815ce0SLisandro Dalcin   const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian";
258a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_REAL_SINGLE)
259a13bc4e3SShao-Ching Huang   const char precision[] = "Float32";
260a13bc4e3SShao-Ching Huang #elif defined(PETSC_USE_REAL_DOUBLE)
261a13bc4e3SShao-Ching Huang   const char precision[] = "Float64";
262a13bc4e3SShao-Ching Huang #else
263a13bc4e3SShao-Ching Huang   const char precision[] = "UnknownPrecision";
264a13bc4e3SShao-Ching Huang #endif
265a13bc4e3SShao-Ching Huang   MPI_Comm                 comm;
266a13bc4e3SShao-Ching Huang   PetscViewer_VTK         *vtk = (PetscViewer_VTK *)viewer->data;
267a13bc4e3SShao-Ching Huang   PetscViewerVTKObjectLink link;
268a13bc4e3SShao-Ching Huang   FILE                    *fp;
269a13bc4e3SShao-Ching Huang   PetscMPIInt              rank, size, tag;
270a13bc4e3SShao-Ching Huang   DMDALocalInfo            info;
271*6497c311SBarry Smith   PetscInt                 dim, mx, my, mz, maxnnodes, maxbs, i, j, k;
2720f2609c8SJed Brown   PetscInt64               boffset;
273a13bc4e3SShao-Ching Huang   PetscInt                 rloc[6], (*grloc)[6] = NULL;
274a13bc4e3SShao-Ching Huang   PetscScalar             *array, *array2;
275a13bc4e3SShao-Ching Huang 
276a13bc4e3SShao-Ching Huang   PetscFunctionBegin;
2779566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)da, &comm));
2781dca8a05SBarry Smith   PetscCheck(!PetscDefined(USE_COMPLEX), PETSC_COMM_SELF, PETSC_ERR_SUP, "Complex values not supported");
2799566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
2809566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
2819566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, &mx, &my, &mz, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
2829566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
2839566063dSJacob Faibussowitsch   PetscCall(PetscFOpen(comm, vtk->filename, "wb", &fp));
2849566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "<?xml version=\"1.0\"?>\n"));
2854a8c153fSJed Brown   PetscCall(PetscFPrintf(comm, fp, "<VTKFile type=\"RectilinearGrid\" version=\"0.1\" byte_order=\"%s\" header_type=\"UInt64\">\n", byte_order));
28663a3b9bcSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "  <RectilinearGrid WholeExtent=\"%d %" PetscInt_FMT " %d %" PetscInt_FMT " %d %" PetscInt_FMT "\">\n", 0, mx - 1, 0, my - 1, 0, mz - 1));
287a13bc4e3SShao-Ching Huang 
2889566063dSJacob Faibussowitsch   if (rank == 0) PetscCall(PetscMalloc1(size * 6, &grloc));
289a13bc4e3SShao-Ching Huang   rloc[0] = info.xs;
290a13bc4e3SShao-Ching Huang   rloc[1] = info.xm;
291a13bc4e3SShao-Ching Huang   rloc[2] = info.ys;
292a13bc4e3SShao-Ching Huang   rloc[3] = info.ym;
293a13bc4e3SShao-Ching Huang   rloc[4] = info.zs;
294a13bc4e3SShao-Ching Huang   rloc[5] = info.zm;
295810441c8SPierre Jolivet   PetscCallMPI(MPI_Gather(rloc, 6, MPIU_INT, rank == 0 ? grloc[0] : NULL, 6, MPIU_INT, 0, comm));
296a13bc4e3SShao-Ching Huang 
297a13bc4e3SShao-Ching Huang   /* Write XML header */
298a13bc4e3SShao-Ching Huang   maxnnodes = 0; /* Used for the temporary array size on rank 0 */
299ea2d7708SPatrick Sanan   maxbs     = 0; /* Used for the temporary array size on rank 0 */
300a13bc4e3SShao-Ching Huang   boffset   = 0; /* Offset into binary file */
301*6497c311SBarry Smith   for (PetscMPIInt r = 0; r < size; r++) {
302a13bc4e3SShao-Ching Huang     PetscInt xs = -1, xm = -1, ys = -1, ym = -1, zs = -1, zm = -1, nnodes = 0;
303*6497c311SBarry Smith 
304dd400576SPatrick Sanan     if (rank == 0) {
305a13bc4e3SShao-Ching Huang       xs     = grloc[r][0];
306a13bc4e3SShao-Ching Huang       xm     = grloc[r][1];
307a13bc4e3SShao-Ching Huang       ys     = grloc[r][2];
308a13bc4e3SShao-Ching Huang       ym     = grloc[r][3];
309a13bc4e3SShao-Ching Huang       zs     = grloc[r][4];
310a13bc4e3SShao-Ching Huang       zm     = grloc[r][5];
311a13bc4e3SShao-Ching Huang       nnodes = xm * ym * zm;
312a13bc4e3SShao-Ching Huang     }
313a13bc4e3SShao-Ching Huang     maxnnodes = PetscMax(maxnnodes, nnodes);
31463a3b9bcSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "    <Piece Extent=\"%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\">\n", xs, xs + xm - 1, ys, ys + ym - 1, zs, zs + zm - 1));
3159566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "      <Coordinates>\n"));
3160f2609c8SJed Brown     PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"Xcoord\"  format=\"appended\"  offset=\"%" PetscInt64_FMT "\" />\n", precision, boffset));
3174a8c153fSJed Brown     boffset += xm * sizeof(PetscScalar) + sizeof(PetscInt64);
3180f2609c8SJed Brown     PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"Ycoord\"  format=\"appended\"  offset=\"%" PetscInt64_FMT "\" />\n", precision, boffset));
3194a8c153fSJed Brown     boffset += ym * sizeof(PetscScalar) + sizeof(PetscInt64);
3200f2609c8SJed Brown     PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"Zcoord\"  format=\"appended\"  offset=\"%" PetscInt64_FMT "\" />\n", precision, boffset));
3214a8c153fSJed Brown     boffset += zm * sizeof(PetscScalar) + sizeof(PetscInt64);
3229566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "      </Coordinates>\n"));
3239566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "      <PointData Scalars=\"ScalarPointData\">\n"));
324a13bc4e3SShao-Ching Huang     for (link = vtk->link; link; link = link->next) {
325a13bc4e3SShao-Ching Huang       Vec         X = (Vec)link->vec;
326fb5bd1c2SPatrick Sanan       PetscInt    bs, f;
327ea2d7708SPatrick Sanan       DM          daCurr;
328fb5bd1c2SPatrick Sanan       PetscBool   fieldsnamed;
329fb5bd1c2SPatrick Sanan       const char *vecname = "Unnamed";
330ea2d7708SPatrick Sanan 
3319566063dSJacob Faibussowitsch       PetscCall(VecGetDM(X, &daCurr));
3329566063dSJacob Faibussowitsch       PetscCall(DMDAGetInfo(daCurr, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL));
333ea2d7708SPatrick Sanan       maxbs = PetscMax(maxbs, bs);
33448a46eb9SPierre Jolivet       if (((PetscObject)X)->name || link != vtk->link) PetscCall(PetscObjectGetName((PetscObject)X, &vecname));
335ea2d7708SPatrick Sanan 
336fb5bd1c2SPatrick Sanan       /* If any fields are named, add scalar fields. Otherwise, add a vector field */
3379566063dSJacob Faibussowitsch       PetscCall(DMDAGetFieldsNamed(daCurr, &fieldsnamed));
338fb5bd1c2SPatrick Sanan       if (fieldsnamed) {
339fb5bd1c2SPatrick Sanan         for (f = 0; f < bs; f++) {
340fb5bd1c2SPatrick Sanan           char        buf[256];
341fb5bd1c2SPatrick Sanan           const char *fieldname;
3429566063dSJacob Faibussowitsch           PetscCall(DMDAGetFieldName(daCurr, f, &fieldname));
343fb5bd1c2SPatrick Sanan           if (!fieldname) {
34463a3b9bcSJacob Faibussowitsch             PetscCall(PetscSNPrintf(buf, sizeof(buf), "%" PetscInt_FMT, f));
345fb5bd1c2SPatrick Sanan             fieldname = buf;
346fb5bd1c2SPatrick Sanan           }
3470f2609c8SJed Brown           PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, fieldname, boffset));
3484a8c153fSJed Brown           boffset += nnodes * sizeof(PetscScalar) + sizeof(PetscInt64);
349fb5bd1c2SPatrick Sanan         }
350fb5bd1c2SPatrick Sanan       } else {
3510f2609c8SJed Brown         PetscCall(PetscFPrintf(comm, fp, "        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%" PetscInt_FMT "\" format=\"appended\" offset=\"%" PetscInt64_FMT "\" />\n", precision, vecname, bs, boffset));
3524a8c153fSJed Brown         boffset += bs * nnodes * sizeof(PetscScalar) + sizeof(PetscInt64);
353a13bc4e3SShao-Ching Huang       }
354fb5bd1c2SPatrick Sanan     }
3559566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "      </PointData>\n"));
3569566063dSJacob Faibussowitsch     PetscCall(PetscFPrintf(comm, fp, "    </Piece>\n"));
357a13bc4e3SShao-Ching Huang   }
3589566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "  </RectilinearGrid>\n"));
3599566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "  <AppendedData encoding=\"raw\">\n"));
3609566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "_"));
361a13bc4e3SShao-Ching Huang 
362a13bc4e3SShao-Ching Huang   /* Now write the arrays. */
363a13bc4e3SShao-Ching Huang   tag = ((PetscObject)viewer)->tag;
3649566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(PetscMax(maxnnodes * maxbs, info.xm + info.ym + info.zm), &array, PetscMax(maxnnodes * maxbs, info.xm + info.ym + info.zm), &array2));
365*6497c311SBarry Smith   for (PetscMPIInt r = 0; r < size; r++) {
366a13bc4e3SShao-Ching Huang     MPI_Status status;
367a13bc4e3SShao-Ching Huang     PetscInt   xs = -1, xm = -1, ys = -1, ym = -1, zs = -1, zm = -1, nnodes = 0;
368*6497c311SBarry Smith 
369dd400576SPatrick Sanan     if (rank == 0) {
370a13bc4e3SShao-Ching Huang       xs     = grloc[r][0];
371a13bc4e3SShao-Ching Huang       xm     = grloc[r][1];
372a13bc4e3SShao-Ching Huang       ys     = grloc[r][2];
373a13bc4e3SShao-Ching Huang       ym     = grloc[r][3];
374a13bc4e3SShao-Ching Huang       zs     = grloc[r][4];
375a13bc4e3SShao-Ching Huang       zm     = grloc[r][5];
376a13bc4e3SShao-Ching Huang       nnodes = xm * ym * zm;
377a13bc4e3SShao-Ching Huang     } else if (r == rank) {
378a13bc4e3SShao-Ching Huang       nnodes = info.xm * info.ym * info.zm;
379a13bc4e3SShao-Ching Huang     }
380a13bc4e3SShao-Ching Huang     { /* Write the coordinates */
381a13bc4e3SShao-Ching Huang       Vec Coords;
382ea2d7708SPatrick Sanan 
3839566063dSJacob Faibussowitsch       PetscCall(DMGetCoordinates(da, &Coords));
384a13bc4e3SShao-Ching Huang       if (Coords) {
385a13bc4e3SShao-Ching Huang         const PetscScalar *coords;
3869566063dSJacob Faibussowitsch         PetscCall(VecGetArrayRead(Coords, &coords));
387dd400576SPatrick Sanan         if (rank == 0) {
388a13bc4e3SShao-Ching Huang           if (r) {
389a13bc4e3SShao-Ching Huang             PetscMPIInt nn;
390*6497c311SBarry Smith 
391*6497c311SBarry Smith             PetscCallMPI(MPIU_Recv(array, xm + ym + zm, MPIU_SCALAR, r, tag, comm, &status));
3929566063dSJacob Faibussowitsch             PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn));
39308401ef6SPierre Jolivet             PetscCheck(nn == xm + ym + zm, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch");
394a13bc4e3SShao-Ching Huang           } else {
395a13bc4e3SShao-Ching Huang             /* x coordinates */
396a13bc4e3SShao-Ching Huang             for (j = 0, k = 0, i = 0; i < xm; i++) {
397a13bc4e3SShao-Ching Huang               PetscInt Iloc = i + xm * (j + ym * k);
398*6497c311SBarry Smith 
399a13bc4e3SShao-Ching Huang               array[i] = coords[Iloc * dim + 0];
400a13bc4e3SShao-Ching Huang             }
401a13bc4e3SShao-Ching Huang             /* y coordinates */
402a13bc4e3SShao-Ching Huang             for (i = 0, k = 0, j = 0; j < ym; j++) {
403a13bc4e3SShao-Ching Huang               PetscInt Iloc = i + xm * (j + ym * k);
404*6497c311SBarry Smith 
405a13bc4e3SShao-Ching Huang               array[j + xm] = dim > 1 ? coords[Iloc * dim + 1] : 0;
406a13bc4e3SShao-Ching Huang             }
407a13bc4e3SShao-Ching Huang             /* z coordinates */
408a13bc4e3SShao-Ching Huang             for (i = 0, j = 0, k = 0; k < zm; k++) {
409a13bc4e3SShao-Ching Huang               PetscInt Iloc = i + xm * (j + ym * k);
410*6497c311SBarry Smith 
411a13bc4e3SShao-Ching Huang               array[k + xm + ym] = dim > 2 ? coords[Iloc * dim + 2] : 0;
412a13bc4e3SShao-Ching Huang             }
413a13bc4e3SShao-Ching Huang           }
414a13bc4e3SShao-Ching Huang         } else if (r == rank) {
415a13bc4e3SShao-Ching Huang           xm = info.xm;
416a13bc4e3SShao-Ching Huang           ym = info.ym;
417a13bc4e3SShao-Ching Huang           zm = info.zm;
418a13bc4e3SShao-Ching Huang           for (j = 0, k = 0, i = 0; i < xm; i++) {
419a13bc4e3SShao-Ching Huang             PetscInt Iloc = i + xm * (j + ym * k);
420*6497c311SBarry Smith 
421a13bc4e3SShao-Ching Huang             array2[i] = coords[Iloc * dim + 0];
422a13bc4e3SShao-Ching Huang           }
423a13bc4e3SShao-Ching Huang           for (i = 0, k = 0, j = 0; j < ym; j++) {
424a13bc4e3SShao-Ching Huang             PetscInt Iloc = i + xm * (j + ym * k);
425*6497c311SBarry Smith 
426a13bc4e3SShao-Ching Huang             array2[j + xm] = dim > 1 ? coords[Iloc * dim + 1] : 0;
427a13bc4e3SShao-Ching Huang           }
428a13bc4e3SShao-Ching Huang           for (i = 0, j = 0, k = 0; k < zm; k++) {
429a13bc4e3SShao-Ching Huang             PetscInt Iloc = i + xm * (j + ym * k);
430*6497c311SBarry Smith 
431a13bc4e3SShao-Ching Huang             array2[k + xm + ym] = dim > 2 ? coords[Iloc * dim + 2] : 0;
432a13bc4e3SShao-Ching Huang           }
433*6497c311SBarry Smith           PetscCallMPI(MPIU_Send((void *)array2, xm + ym + zm, MPIU_SCALAR, 0, tag, comm));
434a13bc4e3SShao-Ching Huang         }
4359566063dSJacob Faibussowitsch         PetscCall(VecRestoreArrayRead(Coords, &coords));
436a13bc4e3SShao-Ching Huang       } else { /* Fabricate some coordinates using grid index */
437ad540459SPierre Jolivet         for (i = 0; i < xm; i++) array[i] = xs + i;
438ad540459SPierre Jolivet         for (j = 0; j < ym; j++) array[j + xm] = ys + j;
439ad540459SPierre Jolivet         for (k = 0; k < zm; k++) array[k + xm + ym] = zs + k;
440a13bc4e3SShao-Ching Huang       }
441dd400576SPatrick Sanan       if (rank == 0) {
442f4f49eeaSPierre Jolivet         PetscCall(PetscViewerVTKFWrite(viewer, fp, &array[0], xm, MPIU_SCALAR));
443f4f49eeaSPierre Jolivet         PetscCall(PetscViewerVTKFWrite(viewer, fp, &array[xm], ym, MPIU_SCALAR));
444f4f49eeaSPierre Jolivet         PetscCall(PetscViewerVTKFWrite(viewer, fp, &array[xm + ym], zm, MPIU_SCALAR));
445a13bc4e3SShao-Ching Huang       }
446a13bc4e3SShao-Ching Huang     }
447a13bc4e3SShao-Ching Huang 
448a13bc4e3SShao-Ching Huang     /* Write each of the objects queued up for this file */
449a13bc4e3SShao-Ching Huang     for (link = vtk->link; link; link = link->next) {
450a13bc4e3SShao-Ching Huang       Vec                X = (Vec)link->vec;
451a13bc4e3SShao-Ching Huang       const PetscScalar *x;
452fb5bd1c2SPatrick Sanan       PetscInt           bs, f;
453ea2d7708SPatrick Sanan       DM                 daCurr;
454fb5bd1c2SPatrick Sanan       PetscBool          fieldsnamed;
455*6497c311SBarry Smith 
4569566063dSJacob Faibussowitsch       PetscCall(VecGetDM(X, &daCurr));
4579566063dSJacob Faibussowitsch       PetscCall(DMDAGetInfo(daCurr, NULL, NULL, NULL, NULL, NULL, NULL, NULL, &bs, NULL, NULL, NULL, NULL, NULL));
458a13bc4e3SShao-Ching Huang 
4599566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(X, &x));
460dd400576SPatrick Sanan       if (rank == 0) {
461a13bc4e3SShao-Ching Huang         if (r) {
462a13bc4e3SShao-Ching Huang           PetscMPIInt nn;
463*6497c311SBarry Smith 
464*6497c311SBarry Smith           PetscCallMPI(MPIU_Recv(array, nnodes * bs, MPIU_SCALAR, r, tag, comm, &status));
4659566063dSJacob Faibussowitsch           PetscCallMPI(MPI_Get_count(&status, MPIU_SCALAR, &nn));
466*6497c311SBarry Smith           PetscCheck(nn == nnodes * bs, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Array size mismatch receiving from rank %d", r);
4671baa6e33SBarry Smith         } else PetscCall(PetscArraycpy(array, x, nnodes * bs));
468fb5bd1c2SPatrick Sanan         /* If any fields are named, add scalar fields. Otherwise, add a vector field */
4699566063dSJacob Faibussowitsch         PetscCall(DMDAGetFieldsNamed(daCurr, &fieldsnamed));
470fb5bd1c2SPatrick Sanan         if (fieldsnamed) {
471fb5bd1c2SPatrick Sanan           for (f = 0; f < bs; f++) {
472fb5bd1c2SPatrick Sanan             /* Extract and transpose the f'th field */
473fb5bd1c2SPatrick Sanan             for (k = 0; k < zm; k++) {
474fb5bd1c2SPatrick Sanan               for (j = 0; j < ym; j++) {
475fb5bd1c2SPatrick Sanan                 for (i = 0; i < xm; i++) {
476fb5bd1c2SPatrick Sanan                   PetscInt Iloc = i + xm * (j + ym * k);
477fb5bd1c2SPatrick Sanan                   array2[Iloc]  = array[Iloc * bs + f];
478fb5bd1c2SPatrick Sanan                 }
479fb5bd1c2SPatrick Sanan               }
480fb5bd1c2SPatrick Sanan             }
4819566063dSJacob Faibussowitsch             PetscCall(PetscViewerVTKFWrite(viewer, fp, array2, nnodes, MPIU_SCALAR));
482fb5bd1c2SPatrick Sanan           }
483fb5bd1c2SPatrick Sanan         }
4849566063dSJacob Faibussowitsch         PetscCall(PetscViewerVTKFWrite(viewer, fp, array, nnodes * bs, MPIU_SCALAR));
485ea2d7708SPatrick Sanan 
486a13bc4e3SShao-Ching Huang       } else if (r == rank) {
487*6497c311SBarry Smith         PetscCallMPI(MPIU_Send((void *)x, nnodes * bs, MPIU_SCALAR, 0, tag, comm));
488a13bc4e3SShao-Ching Huang       }
4899566063dSJacob Faibussowitsch       PetscCall(VecRestoreArrayRead(X, &x));
490a13bc4e3SShao-Ching Huang     }
491a13bc4e3SShao-Ching Huang   }
4929566063dSJacob Faibussowitsch   PetscCall(PetscFree2(array, array2));
4939566063dSJacob Faibussowitsch   PetscCall(PetscFree(grloc));
494a13bc4e3SShao-Ching Huang 
4959566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "\n </AppendedData>\n"));
4969566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(comm, fp, "</VTKFile>\n"));
4979566063dSJacob Faibussowitsch   PetscCall(PetscFClose(comm, fp));
4983ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
499a13bc4e3SShao-Ching Huang }
500a13bc4e3SShao-Ching Huang 
501ffeef943SBarry Smith /*@
5024061b8bfSJed Brown   DMDAVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
5034061b8bfSJed Brown 
5044061b8bfSJed Brown   Collective
5054061b8bfSJed Brown 
5064165533cSJose E. Roman   Input Parameters:
507dce8aebaSBarry Smith + odm    - `DMDA` specifying the grid layout, passed as a `PetscObject`
508dce8aebaSBarry Smith - viewer - viewer of type `PETSCVIEWERVTK`
5094061b8bfSJed Brown 
5104061b8bfSJed Brown   Level: developer
5114061b8bfSJed Brown 
512fb5bd1c2SPatrick Sanan   Notes:
51312b4a537SBarry Smith   This function is a callback used by the `PETSCVIEWERVTK` viewer to actually write the file.
5144061b8bfSJed 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.
5154061b8bfSJed Brown   Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
5164061b8bfSJed Brown 
51712b4a537SBarry Smith   If any fields have been named (see e.g. `DMDASetFieldName()`), then individual scalar
518fb5bd1c2SPatrick Sanan   fields are written. Otherwise, a single multi-dof (vector) field is written.
519fb5bd1c2SPatrick Sanan 
52012b4a537SBarry Smith .seealso: [](sec_struct), `DMDA`, `DM`, `PETSCVIEWERVTK`, `DMDASetFieldName()`
5214061b8bfSJed Brown @*/
DMDAVTKWriteAll(PetscObject odm,PetscViewer viewer)522d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVTKWriteAll(PetscObject odm, PetscViewer viewer)
523d71ae5a4SJacob Faibussowitsch {
5244061b8bfSJed Brown   DM        dm = (DM)odm;
5254061b8bfSJed Brown   PetscBool isvtk;
5264061b8bfSJed Brown 
5274061b8bfSJed Brown   PetscFunctionBegin;
5284061b8bfSJed Brown   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5294061b8bfSJed Brown   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
5309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERVTK, &isvtk));
5317a8be351SBarry Smith   PetscCheck(isvtk, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use viewer type %s", ((PetscObject)viewer)->type_name);
5324061b8bfSJed Brown   switch (viewer->format) {
533d71ae5a4SJacob Faibussowitsch   case PETSC_VIEWER_VTK_VTS:
534d71ae5a4SJacob Faibussowitsch     PetscCall(DMDAVTKWriteAll_VTS(dm, viewer));
535d71ae5a4SJacob Faibussowitsch     break;
536d71ae5a4SJacob Faibussowitsch   case PETSC_VIEWER_VTK_VTR:
537d71ae5a4SJacob Faibussowitsch     PetscCall(DMDAVTKWriteAll_VTR(dm, viewer));
538d71ae5a4SJacob Faibussowitsch     break;
539d71ae5a4SJacob Faibussowitsch   default:
540d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "No support for format '%s'", PetscViewerFormats[viewer->format]);
5414061b8bfSJed Brown   }
5423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5434061b8bfSJed Brown }
544