xref: /petsc/src/dm/impls/da/grvtk.c (revision fb5bd1c2e8c91420af3fc7b19c5e1a447f8f7160)
1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>
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 
8*fb5bd1c2SPatrick Sanan /* Helper function which determines if any DMDA fields are named.  This is used
9*fb5bd1c2SPatrick Sanan    as a proxy for the user's intention to use DMDA fields as distinct
10*fb5bd1c2SPatrick Sanan    scalar-valued fields as opposed to a single vector-valued field */
11*fb5bd1c2SPatrick Sanan static PetscErrorCode DMDAGetFieldsNamed(DM da,PetscBool *fieldsnamed)
12*fb5bd1c2SPatrick Sanan {
13*fb5bd1c2SPatrick Sanan   PetscErrorCode ierr;
14*fb5bd1c2SPatrick Sanan   PetscInt       f,bs;
15*fb5bd1c2SPatrick Sanan 
16*fb5bd1c2SPatrick Sanan   PetscFunctionBegin;
17*fb5bd1c2SPatrick Sanan   ierr = DMDAGetDof(da,&bs);CHKERRQ(ierr);
18*fb5bd1c2SPatrick Sanan   *fieldsnamed = PETSC_FALSE;
19*fb5bd1c2SPatrick Sanan   for (f=0; f<bs; ++f) {
20*fb5bd1c2SPatrick Sanan     const char * fieldname;
21*fb5bd1c2SPatrick Sanan     ierr = DMDAGetFieldName(da,f,&fieldname);CHKERRQ(ierr);
22*fb5bd1c2SPatrick Sanan     if (fieldname) {
23*fb5bd1c2SPatrick Sanan       *fieldsnamed = PETSC_TRUE;
24*fb5bd1c2SPatrick Sanan       break;
25*fb5bd1c2SPatrick Sanan     }
26*fb5bd1c2SPatrick Sanan   }
27*fb5bd1c2SPatrick Sanan   PetscFunctionReturn(0);
28*fb5bd1c2SPatrick Sanan }
29*fb5bd1c2SPatrick Sanan 
304061b8bfSJed Brown static PetscErrorCode DMDAVTKWriteAll_VTS(DM da,PetscViewer viewer)
314061b8bfSJed Brown {
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;
46ea2d7708SPatrick Sanan   PetscInt                 dim,mx,my,mz,cdim,bs,boffset,maxnnodes,maxbs,i,j,k,r;
470298fd71SBarry Smith   PetscInt                 rloc[6],(*grloc)[6] = NULL;
484061b8bfSJed Brown   PetscScalar              *array,*array2;
494061b8bfSJed Brown   PetscReal                gmin[3],gmax[3];
504061b8bfSJed Brown   PetscErrorCode           ierr;
514061b8bfSJed Brown 
524061b8bfSJed Brown   PetscFunctionBegin;
53094921d9SJed Brown   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
544061b8bfSJed Brown #if defined(PETSC_USE_COMPLEX)
553bf036e2SBarry Smith   SETERRQ(comm,PETSC_ERR_SUP,"Complex values not supported");
564061b8bfSJed Brown #endif
574061b8bfSJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
584061b8bfSJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
594061b8bfSJed Brown   ierr = DMDAGetInfo(da,&dim,&mx,&my,&mz,0,0,0,&bs,0,0,0,0,0);CHKERRQ(ierr);
604061b8bfSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
614061b8bfSJed Brown   ierr = DMDAGetBoundingBox(da,gmin,gmax);CHKERRQ(ierr);
622eaa9ef4SLisandro Dalcin   ierr = DMGetCoordinates(da,&Coords);CHKERRQ(ierr);
632eaa9ef4SLisandro Dalcin   if (Coords) {
642eaa9ef4SLisandro Dalcin     PetscInt csize;
652eaa9ef4SLisandro Dalcin     ierr = VecGetSize(Coords,&csize);CHKERRQ(ierr);
662eaa9ef4SLisandro Dalcin     if (csize % (mx*my*mz)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Coordinate vector size mismatch");
672eaa9ef4SLisandro Dalcin     cdim = csize/(mx*my*mz);
682eaa9ef4SLisandro Dalcin     if (cdim < dim || cdim > 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Coordinate vector size mismatch");
692eaa9ef4SLisandro Dalcin   } else {
702eaa9ef4SLisandro Dalcin     cdim = dim;
712eaa9ef4SLisandro Dalcin   }
724061b8bfSJed Brown 
734061b8bfSJed Brown   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
744061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
75519f805aSKarl Rupp #if defined(PETSC_WORDS_BIGENDIAN)
764061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr);
774061b8bfSJed Brown #else
784061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr);
794061b8bfSJed Brown #endif
804061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"  <StructuredGrid WholeExtent=\"%D %D %D %D %D %D\">\n",0,mx-1,0,my-1,0,mz-1);CHKERRQ(ierr);
814061b8bfSJed Brown 
82785e854fSJed Brown   if (!rank) {ierr = PetscMalloc1(size*6,&grloc);CHKERRQ(ierr);}
834061b8bfSJed Brown   rloc[0] = info.xs;
844061b8bfSJed Brown   rloc[1] = info.xm;
854061b8bfSJed Brown   rloc[2] = info.ys;
864061b8bfSJed Brown   rloc[3] = info.ym;
874061b8bfSJed Brown   rloc[4] = info.zs;
884061b8bfSJed Brown   rloc[5] = info.zm;
894061b8bfSJed Brown   ierr    = MPI_Gather(rloc,6,MPIU_INT,&grloc[0][0],6,MPIU_INT,0,comm);CHKERRQ(ierr);
904061b8bfSJed Brown 
914061b8bfSJed Brown   /* Write XML header */
924061b8bfSJed Brown   maxnnodes = 0;                /* Used for the temporary array size on rank 0 */
93ea2d7708SPatrick Sanan   maxbs     = 0;                /* Used for the temporary array size on rank 0 */
944061b8bfSJed Brown   boffset   = 0;                /* Offset into binary file */
954061b8bfSJed Brown   for (r=0; r<size; r++) {
964061b8bfSJed Brown     PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
974061b8bfSJed Brown     if (!rank) {
984061b8bfSJed Brown       xs     = grloc[r][0];
994061b8bfSJed Brown       xm     = grloc[r][1];
1004061b8bfSJed Brown       ys     = grloc[r][2];
1014061b8bfSJed Brown       ym     = grloc[r][3];
1024061b8bfSJed Brown       zs     = grloc[r][4];
1034061b8bfSJed Brown       zm     = grloc[r][5];
1044061b8bfSJed Brown       nnodes = xm*ym*zm;
1054061b8bfSJed Brown     }
1064061b8bfSJed Brown     maxnnodes = PetscMax(maxnnodes,nnodes);
1074061b8bfSJed Brown     ierr     = PetscFPrintf(comm,fp,"    <Piece Extent=\"%D %D %D %D %D %D\">\n",xs,xs+xm-1,ys,ys+ym-1,zs,zs+zm-1);CHKERRQ(ierr);
1084061b8bfSJed Brown     ierr     = PetscFPrintf(comm,fp,"      <Points>\n");CHKERRQ(ierr);
1094061b8bfSJed Brown     ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
1104061b8bfSJed Brown     boffset += 3*nnodes*sizeof(PetscScalar) + sizeof(int);
1114061b8bfSJed Brown     ierr     = PetscFPrintf(comm,fp,"      </Points>\n");CHKERRQ(ierr);
1124061b8bfSJed Brown 
1134061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"      <PointData Scalars=\"ScalarPointData\">\n");CHKERRQ(ierr);
1144061b8bfSJed Brown     for (link=vtk->link; link; link=link->next) {
1154061b8bfSJed Brown       Vec        X = (Vec)link->vec;
116*fb5bd1c2SPatrick Sanan       PetscInt   bs,f;
117ea2d7708SPatrick Sanan       DM         daCurr;
118*fb5bd1c2SPatrick Sanan       PetscBool  fieldsnamed;
119*fb5bd1c2SPatrick Sanan       const char *vecname = "Unnamed";
120ea2d7708SPatrick Sanan 
121ea2d7708SPatrick Sanan       ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr);
122ea2d7708SPatrick Sanan       ierr = DMDAGetInfo(daCurr,0,0,0,0,0,0,0,&bs,0,0,0,0,0);CHKERRQ(ierr);
123ea2d7708SPatrick Sanan       maxbs = PetscMax(maxbs,bs);
124ea2d7708SPatrick Sanan 
125ea2d7708SPatrick Sanan       if (((PetscObject)X)->name || link != vtk->link) {
1264061b8bfSJed Brown         ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
1274061b8bfSJed Brown       }
128*fb5bd1c2SPatrick Sanan 
129*fb5bd1c2SPatrick Sanan       /* If any fields are named, add scalar fields. Otherwise, add a vector field */
130*fb5bd1c2SPatrick Sanan       ierr = DMDAGetFieldsNamed(daCurr,&fieldsnamed);CHKERRQ(ierr);
131*fb5bd1c2SPatrick Sanan       if (fieldsnamed) {
132*fb5bd1c2SPatrick Sanan         for (f=0; f<bs; f++) {
133*fb5bd1c2SPatrick Sanan           char       buf[256];
134*fb5bd1c2SPatrick Sanan           const char *fieldname;
135*fb5bd1c2SPatrick Sanan           ierr = DMDAGetFieldName(daCurr,f,&fieldname);CHKERRQ(ierr);
136*fb5bd1c2SPatrick Sanan           if (!fieldname) {
137*fb5bd1c2SPatrick Sanan             ierr      = PetscSNPrintf(buf,sizeof(buf),"%D",f);CHKERRQ(ierr);
138*fb5bd1c2SPatrick Sanan             fieldname = buf;
139*fb5bd1c2SPatrick Sanan           }
140*fb5bd1c2SPatrick Sanan           ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
141*fb5bd1c2SPatrick Sanan           boffset += nnodes*sizeof(PetscScalar) + sizeof(int);
142*fb5bd1c2SPatrick Sanan         }
143*fb5bd1c2SPatrick Sanan       } else {
144ea2d7708SPatrick Sanan         ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%D\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,bs,boffset);CHKERRQ(ierr);
145ea2d7708SPatrick Sanan         boffset += bs*nnodes*sizeof(PetscScalar) + sizeof(int);
1464061b8bfSJed Brown       }
147*fb5bd1c2SPatrick Sanan     }
1484061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"      </PointData>\n");CHKERRQ(ierr);
1494061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"    </Piece>\n");CHKERRQ(ierr);
1504061b8bfSJed Brown   }
1514061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"  </StructuredGrid>\n");CHKERRQ(ierr);
1524061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
1534061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
1544061b8bfSJed Brown 
1554061b8bfSJed Brown   /* Now write the arrays. */
1564061b8bfSJed Brown   tag  = ((PetscObject)viewer)->tag;
157*fb5bd1c2SPatrick Sanan   ierr = PetscMalloc2(maxnnodes*PetscMax(3,maxbs),&array,maxnnodes*PetscMax(3,maxbs),&array2);CHKERRQ(ierr);
1584061b8bfSJed Brown   for (r=0; r<size; r++) {
1594061b8bfSJed Brown     MPI_Status status;
1604061b8bfSJed Brown     PetscInt   xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
1614061b8bfSJed Brown     if (!rank) {
1624061b8bfSJed Brown       xs     = grloc[r][0];
1634061b8bfSJed Brown       xm     = grloc[r][1];
1644061b8bfSJed Brown       ys     = grloc[r][2];
1654061b8bfSJed Brown       ym     = grloc[r][3];
1664061b8bfSJed Brown       zs     = grloc[r][4];
1674061b8bfSJed Brown       zm     = grloc[r][5];
1684061b8bfSJed Brown       nnodes = xm*ym*zm;
1694061b8bfSJed Brown     } else if (r == rank) {
1704061b8bfSJed Brown       nnodes = info.xm*info.ym*info.zm;
1714061b8bfSJed Brown     }
1724061b8bfSJed Brown 
1732eaa9ef4SLisandro Dalcin     /* Write the coordinates */
1744061b8bfSJed Brown     if (Coords) {
1754061b8bfSJed Brown       const PetscScalar *coords;
1764061b8bfSJed Brown       ierr = VecGetArrayRead(Coords,&coords);CHKERRQ(ierr);
1774061b8bfSJed Brown       if (!rank) {
1784061b8bfSJed Brown         if (r) {
1796622f924SJed Brown           PetscMPIInt nn;
1802eaa9ef4SLisandro Dalcin           ierr = MPI_Recv(array,nnodes*cdim,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr);
1814061b8bfSJed Brown           ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr);
1822eaa9ef4SLisandro Dalcin           if (nn != nnodes*cdim) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
1834061b8bfSJed Brown         } else {
1842eaa9ef4SLisandro Dalcin           ierr = PetscMemcpy(array,coords,nnodes*cdim*sizeof(PetscScalar));CHKERRQ(ierr);
1854061b8bfSJed Brown         }
1864061b8bfSJed Brown         /* Transpose coordinates to VTK (C-style) ordering */
1874061b8bfSJed Brown         for (k=0; k<zm; k++) {
1884061b8bfSJed Brown           for (j=0; j<ym; j++) {
1894061b8bfSJed Brown             for (i=0; i<xm; i++) {
1904061b8bfSJed Brown               PetscInt Iloc = i+xm*(j+ym*k);
1912eaa9ef4SLisandro Dalcin               array2[Iloc*3+0] = array[Iloc*cdim + 0];
1922eaa9ef4SLisandro Dalcin               array2[Iloc*3+1] = cdim > 1 ? array[Iloc*cdim + 1] : 0.0;
1932eaa9ef4SLisandro Dalcin               array2[Iloc*3+2] = cdim > 2 ? array[Iloc*cdim + 2] : 0.0;
1944061b8bfSJed Brown             }
1954061b8bfSJed Brown           }
1964061b8bfSJed Brown         }
1974061b8bfSJed Brown       } else if (r == rank) {
1982eaa9ef4SLisandro Dalcin         ierr = MPI_Send((void*)coords,nnodes*cdim,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr);
1994061b8bfSJed Brown       }
2004061b8bfSJed Brown       ierr = VecRestoreArrayRead(Coords,&coords);CHKERRQ(ierr);
2014061b8bfSJed Brown     } else {       /* Fabricate some coordinates using grid index */
2024061b8bfSJed Brown       for (k=0; k<zm; k++) {
2034061b8bfSJed Brown         for (j=0; j<ym; j++) {
2044061b8bfSJed Brown           for (i=0; i<xm; i++) {
2054061b8bfSJed Brown             PetscInt Iloc = i+xm*(j+ym*k);
2064061b8bfSJed Brown             array2[Iloc*3+0] = xs+i;
2074061b8bfSJed Brown             array2[Iloc*3+1] = ys+j;
2084061b8bfSJed Brown             array2[Iloc*3+2] = zs+k;
2094061b8bfSJed Brown           }
2104061b8bfSJed Brown         }
2114061b8bfSJed Brown       }
2124061b8bfSJed Brown     }
21394fbd55eSBarry Smith     ierr = PetscViewerVTKFWrite(viewer,fp,array2,nnodes*3,MPIU_SCALAR);CHKERRQ(ierr);
2144061b8bfSJed Brown 
2154061b8bfSJed Brown     /* Write each of the objects queued up for this file */
2164061b8bfSJed Brown     for (link=vtk->link; link; link=link->next) {
2174061b8bfSJed Brown       Vec               X = (Vec)link->vec;
2184061b8bfSJed Brown       const PetscScalar *x;
219*fb5bd1c2SPatrick Sanan       PetscInt          bs,f;
220ea2d7708SPatrick Sanan       DM                daCurr;
221*fb5bd1c2SPatrick Sanan       PetscBool         fieldsnamed;
222ea2d7708SPatrick Sanan       ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr);
223ea2d7708SPatrick Sanan       ierr = DMDAGetInfo(daCurr,0,0,0,0, 0,0,0,&bs,0,0,0,0,0);CHKERRQ(ierr);
2244061b8bfSJed Brown       ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
2254061b8bfSJed Brown       if (!rank) {
2264061b8bfSJed Brown         if (r) {
2276622f924SJed Brown           PetscMPIInt nn;
2284061b8bfSJed Brown           ierr = MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr);
2294061b8bfSJed Brown           ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr);
2304061b8bfSJed Brown           if (nn != nnodes*bs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %D",r);
2314061b8bfSJed Brown         } else {
2324061b8bfSJed Brown           ierr = PetscMemcpy(array,x,nnodes*bs*sizeof(PetscScalar));CHKERRQ(ierr);
2334061b8bfSJed Brown         }
234*fb5bd1c2SPatrick Sanan 
235*fb5bd1c2SPatrick Sanan         /* If any fields are named, add scalar fields. Otherwise, add a vector field */
236*fb5bd1c2SPatrick Sanan         ierr = DMDAGetFieldsNamed(daCurr,&fieldsnamed);CHKERRQ(ierr);
237*fb5bd1c2SPatrick Sanan         if (fieldsnamed) {
238*fb5bd1c2SPatrick Sanan           for (f=0; f<bs; f++) {
239*fb5bd1c2SPatrick Sanan             /* Extract and transpose the f'th field */
240*fb5bd1c2SPatrick Sanan             for (k=0; k<zm; k++) {
241*fb5bd1c2SPatrick Sanan               for (j=0; j<ym; j++) {
242*fb5bd1c2SPatrick Sanan                 for (i=0; i<xm; i++) {
243*fb5bd1c2SPatrick Sanan                   PetscInt Iloc = i+xm*(j+ym*k);
244*fb5bd1c2SPatrick Sanan                   array2[Iloc] = array[Iloc*bs + f];
245*fb5bd1c2SPatrick Sanan                 }
246*fb5bd1c2SPatrick Sanan               }
247*fb5bd1c2SPatrick Sanan             }
248*fb5bd1c2SPatrick Sanan             ierr = PetscViewerVTKFWrite(viewer,fp,array2,nnodes,MPIU_SCALAR);CHKERRQ(ierr);
249*fb5bd1c2SPatrick Sanan           }
250*fb5bd1c2SPatrick Sanan         } else {
251ea2d7708SPatrick Sanan           ierr = PetscViewerVTKFWrite(viewer,fp,array,bs*nnodes,MPIU_SCALAR);CHKERRQ(ierr);
252*fb5bd1c2SPatrick Sanan         }
2534061b8bfSJed Brown       } else if (r == rank) {
2544061b8bfSJed Brown         ierr = MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr);
2554061b8bfSJed Brown       }
2564061b8bfSJed Brown       ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
2574061b8bfSJed Brown     }
2584061b8bfSJed Brown   }
2594061b8bfSJed Brown   ierr = PetscFree2(array,array2);CHKERRQ(ierr);
2604061b8bfSJed Brown   ierr = PetscFree(grloc);CHKERRQ(ierr);
2614061b8bfSJed Brown 
2624061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr);
2634061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
2644061b8bfSJed Brown   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
2654061b8bfSJed Brown   PetscFunctionReturn(0);
2664061b8bfSJed Brown }
2674061b8bfSJed Brown 
2684061b8bfSJed Brown 
269a13bc4e3SShao-Ching Huang static PetscErrorCode DMDAVTKWriteAll_VTR(DM da,PetscViewer viewer)
270a13bc4e3SShao-Ching Huang {
271a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_REAL_SINGLE)
272a13bc4e3SShao-Ching Huang   const char precision[] = "Float32";
273a13bc4e3SShao-Ching Huang #elif defined(PETSC_USE_REAL_DOUBLE)
274a13bc4e3SShao-Ching Huang   const char precision[] = "Float64";
275a13bc4e3SShao-Ching Huang #else
276a13bc4e3SShao-Ching Huang   const char precision[] = "UnknownPrecision";
277a13bc4e3SShao-Ching Huang #endif
278a13bc4e3SShao-Ching Huang   MPI_Comm                 comm;
279a13bc4e3SShao-Ching Huang   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
280a13bc4e3SShao-Ching Huang   PetscViewerVTKObjectLink link;
281a13bc4e3SShao-Ching Huang   FILE                     *fp;
282a13bc4e3SShao-Ching Huang   PetscMPIInt              rank,size,tag;
283a13bc4e3SShao-Ching Huang   DMDALocalInfo            info;
284ea2d7708SPatrick Sanan   PetscInt                 dim,mx,my,mz,boffset,maxnnodes,maxbs,i,j,k,r;
285a13bc4e3SShao-Ching Huang   PetscInt                 rloc[6],(*grloc)[6] = NULL;
286a13bc4e3SShao-Ching Huang   PetscScalar              *array,*array2;
287a13bc4e3SShao-Ching Huang   PetscReal                gmin[3],gmax[3];
288a13bc4e3SShao-Ching Huang   PetscErrorCode           ierr;
289a13bc4e3SShao-Ching Huang 
290a13bc4e3SShao-Ching Huang   PetscFunctionBegin;
291a13bc4e3SShao-Ching Huang   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
292a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_COMPLEX)
293a13bc4e3SShao-Ching Huang   SETERRQ(comm,PETSC_ERR_SUP,"Complex values not supported");
294a13bc4e3SShao-Ching Huang #endif
295a13bc4e3SShao-Ching Huang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
296a13bc4e3SShao-Ching Huang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
297ea2d7708SPatrick Sanan   ierr = DMDAGetInfo(da,&dim,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
298a13bc4e3SShao-Ching Huang   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
299a13bc4e3SShao-Ching Huang   ierr = DMDAGetBoundingBox(da,gmin,gmax);CHKERRQ(ierr);
300a13bc4e3SShao-Ching Huang   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
301a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
302a13bc4e3SShao-Ching Huang #if defined(PETSC_WORDS_BIGENDIAN)
303a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"RectilinearGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr);
304a13bc4e3SShao-Ching Huang #else
305a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"RectilinearGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr);
306a13bc4e3SShao-Ching Huang #endif
307a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"  <RectilinearGrid WholeExtent=\"%D %D %D %D %D %D\">\n",0,mx-1,0,my-1,0,mz-1);CHKERRQ(ierr);
308a13bc4e3SShao-Ching Huang 
309785e854fSJed Brown   if (!rank) {ierr = PetscMalloc1(size*6,&grloc);CHKERRQ(ierr);}
310a13bc4e3SShao-Ching Huang   rloc[0] = info.xs;
311a13bc4e3SShao-Ching Huang   rloc[1] = info.xm;
312a13bc4e3SShao-Ching Huang   rloc[2] = info.ys;
313a13bc4e3SShao-Ching Huang   rloc[3] = info.ym;
314a13bc4e3SShao-Ching Huang   rloc[4] = info.zs;
315a13bc4e3SShao-Ching Huang   rloc[5] = info.zm;
316a13bc4e3SShao-Ching Huang   ierr    = MPI_Gather(rloc,6,MPIU_INT,&grloc[0][0],6,MPIU_INT,0,comm);CHKERRQ(ierr);
317a13bc4e3SShao-Ching Huang 
318a13bc4e3SShao-Ching Huang   /* Write XML header */
319a13bc4e3SShao-Ching Huang   maxnnodes = 0;                /* Used for the temporary array size on rank 0 */
320ea2d7708SPatrick Sanan   maxbs     = 0;                /* Used for the temporary array size on rank 0 */
321a13bc4e3SShao-Ching Huang   boffset   = 0;                /* Offset into binary file */
322a13bc4e3SShao-Ching Huang   for (r=0; r<size; r++) {
323a13bc4e3SShao-Ching Huang     PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
324a13bc4e3SShao-Ching Huang     if (!rank) {
325a13bc4e3SShao-Ching Huang       xs     = grloc[r][0];
326a13bc4e3SShao-Ching Huang       xm     = grloc[r][1];
327a13bc4e3SShao-Ching Huang       ys     = grloc[r][2];
328a13bc4e3SShao-Ching Huang       ym     = grloc[r][3];
329a13bc4e3SShao-Ching Huang       zs     = grloc[r][4];
330a13bc4e3SShao-Ching Huang       zm     = grloc[r][5];
331a13bc4e3SShao-Ching Huang       nnodes = xm*ym*zm;
332a13bc4e3SShao-Ching Huang     }
333a13bc4e3SShao-Ching Huang     maxnnodes = PetscMax(maxnnodes,nnodes);
334a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"    <Piece Extent=\"%D %D %D %D %D %D\">\n",xs,xs+xm-1,ys,ys+ym-1,zs,zs+zm-1);CHKERRQ(ierr);
335a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"      <Coordinates>\n");CHKERRQ(ierr);
336a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Xcoord\"  format=\"appended\"  offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
337a13bc4e3SShao-Ching Huang     boffset += xm*sizeof(PetscScalar) + sizeof(int);
338a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Ycoord\"  format=\"appended\"  offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
339a13bc4e3SShao-Ching Huang     boffset += ym*sizeof(PetscScalar) + sizeof(int);
340a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Zcoord\"  format=\"appended\"  offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
341a13bc4e3SShao-Ching Huang     boffset += zm*sizeof(PetscScalar) + sizeof(int);
342a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"      </Coordinates>\n");CHKERRQ(ierr);
343a13bc4e3SShao-Ching Huang     ierr = PetscFPrintf(comm,fp,"      <PointData Scalars=\"ScalarPointData\">\n");CHKERRQ(ierr);
344a13bc4e3SShao-Ching Huang     for (link=vtk->link; link; link=link->next) {
345a13bc4e3SShao-Ching Huang       Vec        X = (Vec)link->vec;
346*fb5bd1c2SPatrick Sanan       PetscInt   bs,f;
347ea2d7708SPatrick Sanan       DM         daCurr;
348*fb5bd1c2SPatrick Sanan       PetscBool  fieldsnamed;
349*fb5bd1c2SPatrick Sanan       const char *vecname = "Unnamed";
350ea2d7708SPatrick Sanan 
351ea2d7708SPatrick Sanan       ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr);
352ea2d7708SPatrick Sanan       ierr = DMDAGetInfo(daCurr,0,0,0,0,0,0,0,&bs,0,0,0,0,0);CHKERRQ(ierr);
353ea2d7708SPatrick Sanan       maxbs = PetscMax(maxbs,bs);
354ea2d7708SPatrick Sanan       if (((PetscObject)X)->name || link != vtk->link) {
355a13bc4e3SShao-Ching Huang         ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
356a13bc4e3SShao-Ching Huang       }
357ea2d7708SPatrick Sanan 
358*fb5bd1c2SPatrick Sanan       /* If any fields are named, add scalar fields. Otherwise, add a vector field */
359*fb5bd1c2SPatrick Sanan       ierr = DMDAGetFieldsNamed(daCurr,&fieldsnamed);CHKERRQ(ierr);
360*fb5bd1c2SPatrick Sanan       if (fieldsnamed) {
361*fb5bd1c2SPatrick Sanan         for (f=0; f<bs; f++) {
362*fb5bd1c2SPatrick Sanan           char       buf[256];
363*fb5bd1c2SPatrick Sanan           const char *fieldname;
364*fb5bd1c2SPatrick Sanan           ierr = DMDAGetFieldName(daCurr,f,&fieldname);CHKERRQ(ierr);
365*fb5bd1c2SPatrick Sanan           if (!fieldname) {
366*fb5bd1c2SPatrick Sanan             ierr      = PetscSNPrintf(buf,sizeof(buf),"%D",f);CHKERRQ(ierr);
367*fb5bd1c2SPatrick Sanan             fieldname = buf;
368*fb5bd1c2SPatrick Sanan           }
369*fb5bd1c2SPatrick Sanan           ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
370*fb5bd1c2SPatrick Sanan           boffset += nnodes*sizeof(PetscScalar) + sizeof(int);
371*fb5bd1c2SPatrick Sanan         }
372*fb5bd1c2SPatrick Sanan       } else {
373ea2d7708SPatrick Sanan         ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%D\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,bs,boffset);CHKERRQ(ierr);
374ea2d7708SPatrick Sanan         boffset += bs*nnodes*sizeof(PetscScalar) + sizeof(int);
375a13bc4e3SShao-Ching Huang       }
376*fb5bd1c2SPatrick Sanan     }
377a13bc4e3SShao-Ching Huang     ierr = PetscFPrintf(comm,fp,"      </PointData>\n");CHKERRQ(ierr);
378a13bc4e3SShao-Ching Huang     ierr = PetscFPrintf(comm,fp,"    </Piece>\n");CHKERRQ(ierr);
379a13bc4e3SShao-Ching Huang   }
380a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"  </RectilinearGrid>\n");CHKERRQ(ierr);
381a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
382a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
383a13bc4e3SShao-Ching Huang 
384a13bc4e3SShao-Ching Huang   /* Now write the arrays. */
385a13bc4e3SShao-Ching Huang   tag  = ((PetscObject)viewer)->tag;
386*fb5bd1c2SPatrick Sanan   ierr = PetscMalloc2(PetscMax(maxnnodes*maxbs,info.xm+info.ym+info.zm),&array,PetscMax(maxnnodes*maxbs,info.xm+info.ym+info.zm),&array2);CHKERRQ(ierr);
387a13bc4e3SShao-Ching Huang   for (r=0; r<size; r++) {
388a13bc4e3SShao-Ching Huang     MPI_Status status;
389a13bc4e3SShao-Ching Huang     PetscInt   xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
390a13bc4e3SShao-Ching Huang     if (!rank) {
391a13bc4e3SShao-Ching Huang       xs     = grloc[r][0];
392a13bc4e3SShao-Ching Huang       xm     = grloc[r][1];
393a13bc4e3SShao-Ching Huang       ys     = grloc[r][2];
394a13bc4e3SShao-Ching Huang       ym     = grloc[r][3];
395a13bc4e3SShao-Ching Huang       zs     = grloc[r][4];
396a13bc4e3SShao-Ching Huang       zm     = grloc[r][5];
397a13bc4e3SShao-Ching Huang       nnodes = xm*ym*zm;
398a13bc4e3SShao-Ching Huang     } else if (r == rank) {
399a13bc4e3SShao-Ching Huang       nnodes = info.xm*info.ym*info.zm;
400a13bc4e3SShao-Ching Huang     }
401a13bc4e3SShao-Ching Huang     {                           /* Write the coordinates */
402a13bc4e3SShao-Ching Huang       Vec Coords;
403ea2d7708SPatrick Sanan 
404a13bc4e3SShao-Ching Huang       ierr = DMGetCoordinates(da,&Coords);CHKERRQ(ierr);
405a13bc4e3SShao-Ching Huang       if (Coords) {
406a13bc4e3SShao-Ching Huang         const PetscScalar *coords;
407a13bc4e3SShao-Ching Huang         ierr = VecGetArrayRead(Coords,&coords);CHKERRQ(ierr);
408a13bc4e3SShao-Ching Huang         if (!rank) {
409a13bc4e3SShao-Ching Huang           if (r) {
410a13bc4e3SShao-Ching Huang             PetscMPIInt nn;
411a13bc4e3SShao-Ching Huang             ierr = MPI_Recv(array,xm+ym+zm,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr);
412a13bc4e3SShao-Ching Huang             ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr);
413a13bc4e3SShao-Ching Huang             if (nn != xm+ym+zm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
414a13bc4e3SShao-Ching Huang           } else {
415a13bc4e3SShao-Ching Huang             /* x coordinates */
416a13bc4e3SShao-Ching Huang             for (j=0, k=0, i=0; i<xm; i++) {
417a13bc4e3SShao-Ching Huang               PetscInt Iloc = i+xm*(j+ym*k);
418a13bc4e3SShao-Ching Huang               array[i] = coords[Iloc*dim + 0];
419a13bc4e3SShao-Ching Huang             }
420a13bc4e3SShao-Ching Huang             /* y coordinates */
421a13bc4e3SShao-Ching Huang             for (i=0, k=0, j=0; j<ym; j++) {
422a13bc4e3SShao-Ching Huang               PetscInt Iloc = i+xm*(j+ym*k);
423a13bc4e3SShao-Ching Huang               array[j+xm] = dim > 1 ? coords[Iloc*dim + 1] : 0;
424a13bc4e3SShao-Ching Huang             }
425a13bc4e3SShao-Ching Huang             /* z coordinates */
426a13bc4e3SShao-Ching Huang             for (i=0, j=0, k=0; k<zm; k++) {
427a13bc4e3SShao-Ching Huang               PetscInt Iloc = i+xm*(j+ym*k);
428a13bc4e3SShao-Ching Huang               array[k+xm+ym] = dim > 2 ? coords[Iloc*dim + 2] : 0;
429a13bc4e3SShao-Ching Huang             }
430a13bc4e3SShao-Ching Huang           }
431a13bc4e3SShao-Ching Huang         } else if (r == rank) {
432a13bc4e3SShao-Ching Huang           xm = info.xm;
433a13bc4e3SShao-Ching Huang           ym = info.ym;
434a13bc4e3SShao-Ching Huang           zm = info.zm;
435a13bc4e3SShao-Ching Huang           for (j=0, k=0, i=0; i<xm; i++) {
436a13bc4e3SShao-Ching Huang             PetscInt Iloc = i+xm*(j+ym*k);
437a13bc4e3SShao-Ching Huang             array2[i] = coords[Iloc*dim + 0];
438a13bc4e3SShao-Ching Huang           }
439a13bc4e3SShao-Ching Huang           for (i=0, k=0, j=0; j<ym; j++) {
440a13bc4e3SShao-Ching Huang             PetscInt Iloc = i+xm*(j+ym*k);
441a13bc4e3SShao-Ching Huang             array2[j+xm] = dim > 1 ? coords[Iloc*dim + 1] : 0;
442a13bc4e3SShao-Ching Huang           }
443a13bc4e3SShao-Ching Huang           for (i=0, j=0, k=0; k<zm; k++) {
444a13bc4e3SShao-Ching Huang             PetscInt Iloc = i+xm*(j+ym*k);
445a13bc4e3SShao-Ching Huang             array2[k+xm+ym] = dim > 2 ? coords[Iloc*dim + 2] : 0;
446a13bc4e3SShao-Ching Huang           }
447a13bc4e3SShao-Ching Huang           ierr = MPI_Send((void*)array2,xm+ym+zm,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr);
448a13bc4e3SShao-Ching Huang         }
449a13bc4e3SShao-Ching Huang         ierr = VecRestoreArrayRead(Coords,&coords);CHKERRQ(ierr);
450a13bc4e3SShao-Ching Huang       } else {       /* Fabricate some coordinates using grid index */
4513d3eaba7SBarry Smith         for (i=0; i<xm; i++) {
452a13bc4e3SShao-Ching Huang           array[i] = xs+i;
453a13bc4e3SShao-Ching Huang         }
4543d3eaba7SBarry Smith         for (j=0; j<ym; j++) {
455a13bc4e3SShao-Ching Huang           array[j+xm] = ys+j;
456a13bc4e3SShao-Ching Huang         }
4573d3eaba7SBarry Smith         for (k=0; k<zm; k++) {
458a13bc4e3SShao-Ching Huang           array[k+xm+ym] = zs+k;
459a13bc4e3SShao-Ching Huang         }
460a13bc4e3SShao-Ching Huang       }
461a13bc4e3SShao-Ching Huang       if (!rank) {
46294fbd55eSBarry Smith         ierr = PetscViewerVTKFWrite(viewer,fp,&(array[0]),xm,MPIU_SCALAR);CHKERRQ(ierr);
46394fbd55eSBarry Smith         ierr = PetscViewerVTKFWrite(viewer,fp,&(array[xm]),ym,MPIU_SCALAR);CHKERRQ(ierr);
46494fbd55eSBarry Smith         ierr = PetscViewerVTKFWrite(viewer,fp,&(array[xm+ym]),zm,MPIU_SCALAR);CHKERRQ(ierr);
465a13bc4e3SShao-Ching Huang       }
466a13bc4e3SShao-Ching Huang     }
467a13bc4e3SShao-Ching Huang 
468a13bc4e3SShao-Ching Huang     /* Write each of the objects queued up for this file */
469a13bc4e3SShao-Ching Huang     for (link=vtk->link; link; link=link->next) {
470a13bc4e3SShao-Ching Huang       Vec               X = (Vec)link->vec;
471a13bc4e3SShao-Ching Huang       const PetscScalar *x;
472*fb5bd1c2SPatrick Sanan       PetscInt          bs,f;
473ea2d7708SPatrick Sanan       DM                daCurr;
474*fb5bd1c2SPatrick Sanan       PetscBool         fieldsnamed;
475ea2d7708SPatrick Sanan       ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr);
476ea2d7708SPatrick Sanan       ierr = DMDAGetInfo(daCurr,0,0,0,0,0,0,0,&bs,0,0,0,0,0);CHKERRQ(ierr);
477a13bc4e3SShao-Ching Huang 
478a13bc4e3SShao-Ching Huang       ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
479a13bc4e3SShao-Ching Huang       if (!rank) {
480a13bc4e3SShao-Ching Huang         if (r) {
481a13bc4e3SShao-Ching Huang           PetscMPIInt nn;
482a13bc4e3SShao-Ching Huang           ierr = MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr);
483a13bc4e3SShao-Ching Huang           ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr);
484a13bc4e3SShao-Ching Huang           if (nn != nnodes*bs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %D",r);
485a13bc4e3SShao-Ching Huang         } else {
486a13bc4e3SShao-Ching Huang           ierr = PetscMemcpy(array,x,nnodes*bs*sizeof(PetscScalar));CHKERRQ(ierr);
487a13bc4e3SShao-Ching Huang         }
488*fb5bd1c2SPatrick Sanan         /* If any fields are named, add scalar fields. Otherwise, add a vector field */
489*fb5bd1c2SPatrick Sanan         ierr = DMDAGetFieldsNamed(daCurr,&fieldsnamed);CHKERRQ(ierr);
490*fb5bd1c2SPatrick Sanan         if (fieldsnamed) {
491*fb5bd1c2SPatrick Sanan           for (f=0; f<bs; f++) {
492*fb5bd1c2SPatrick Sanan             /* Extract and transpose the f'th field */
493*fb5bd1c2SPatrick Sanan             for (k=0; k<zm; k++) {
494*fb5bd1c2SPatrick Sanan               for (j=0; j<ym; j++) {
495*fb5bd1c2SPatrick Sanan                 for (i=0; i<xm; i++) {
496*fb5bd1c2SPatrick Sanan                   PetscInt Iloc = i+xm*(j+ym*k);
497*fb5bd1c2SPatrick Sanan                   array2[Iloc] = array[Iloc*bs + f];
498*fb5bd1c2SPatrick Sanan                 }
499*fb5bd1c2SPatrick Sanan               }
500*fb5bd1c2SPatrick Sanan             }
501*fb5bd1c2SPatrick Sanan             ierr = PetscViewerVTKFWrite(viewer,fp,array2,nnodes,MPIU_SCALAR);CHKERRQ(ierr);
502*fb5bd1c2SPatrick Sanan           }
503*fb5bd1c2SPatrick Sanan         }
504ea2d7708SPatrick Sanan         ierr = PetscViewerVTKFWrite(viewer,fp,array,nnodes*bs,MPIU_SCALAR);CHKERRQ(ierr);
505ea2d7708SPatrick Sanan 
506a13bc4e3SShao-Ching Huang       } else if (r == rank) {
507a13bc4e3SShao-Ching Huang         ierr = MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr);
508a13bc4e3SShao-Ching Huang       }
509a13bc4e3SShao-Ching Huang       ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
510a13bc4e3SShao-Ching Huang     }
511a13bc4e3SShao-Ching Huang   }
512a13bc4e3SShao-Ching Huang   ierr = PetscFree2(array,array2);CHKERRQ(ierr);
513a13bc4e3SShao-Ching Huang   ierr = PetscFree(grloc);CHKERRQ(ierr);
514a13bc4e3SShao-Ching Huang 
515a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr);
516a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
517a13bc4e3SShao-Ching Huang   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
518a13bc4e3SShao-Ching Huang   PetscFunctionReturn(0);
519a13bc4e3SShao-Ching Huang }
520a13bc4e3SShao-Ching Huang 
5214061b8bfSJed Brown /*@C
5224061b8bfSJed Brown    DMDAVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
5234061b8bfSJed Brown 
5244061b8bfSJed Brown    Collective
5254061b8bfSJed Brown 
5264061b8bfSJed Brown    Input Arguments:
5274061b8bfSJed Brown    odm - DM specifying the grid layout, passed as a PetscObject
5284061b8bfSJed Brown    viewer - viewer of type VTK
5294061b8bfSJed Brown 
5304061b8bfSJed Brown    Level: developer
5314061b8bfSJed Brown 
532*fb5bd1c2SPatrick Sanan    Notes:
5334061b8bfSJed Brown    This function is a callback used by the VTK viewer to actually write the file.
5344061b8bfSJed 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.
5354061b8bfSJed Brown    Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
5364061b8bfSJed Brown 
537*fb5bd1c2SPatrick Sanan    If any fields have been named (see e.g. DMDASetFieldName()), then individual scalar
538*fb5bd1c2SPatrick Sanan    fields are written. Otherwise, a single multi-dof (vector) field is written.
539*fb5bd1c2SPatrick Sanan 
5404061b8bfSJed Brown .seealso: PETSCVIEWERVTK
5414061b8bfSJed Brown @*/
5424061b8bfSJed Brown PetscErrorCode DMDAVTKWriteAll(PetscObject odm,PetscViewer viewer)
5434061b8bfSJed Brown {
5444061b8bfSJed Brown   DM             dm = (DM)odm;
5454061b8bfSJed Brown   PetscBool      isvtk;
5464061b8bfSJed Brown   PetscErrorCode ierr;
5474061b8bfSJed Brown 
5484061b8bfSJed Brown   PetscFunctionBegin;
5494061b8bfSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5504061b8bfSJed Brown   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
551251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr);
552ce94432eSBarry Smith   if (!isvtk) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_INCOMP,"Cannot use viewer type %s",((PetscObject)viewer)->type_name);
5534061b8bfSJed Brown   switch (viewer->format) {
5544061b8bfSJed Brown   case PETSC_VIEWER_VTK_VTS:
5554061b8bfSJed Brown     ierr = DMDAVTKWriteAll_VTS(dm,viewer);CHKERRQ(ierr);
5564061b8bfSJed Brown     break;
557a13bc4e3SShao-Ching Huang   case PETSC_VIEWER_VTK_VTR:
558a13bc4e3SShao-Ching Huang     ierr = DMDAVTKWriteAll_VTR(dm,viewer);CHKERRQ(ierr);
559a13bc4e3SShao-Ching Huang     break;
560ce94432eSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No support for format '%s'",PetscViewerFormats[viewer->format]);
5614061b8bfSJed Brown   }
5624061b8bfSJed Brown   PetscFunctionReturn(0);
5634061b8bfSJed Brown }
564