xref: /petsc/src/dm/impls/da/grvtk.c (revision ea2d7708fa6bebc714ed9a4ec0380dd837201800)
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 
84061b8bfSJed Brown static PetscErrorCode DMDAVTKWriteAll_VTS(DM da,PetscViewer viewer)
94061b8bfSJed Brown {
104061b8bfSJed Brown #if defined(PETSC_USE_REAL_SINGLE)
114061b8bfSJed Brown   const char precision[] = "Float32";
124061b8bfSJed Brown #elif defined(PETSC_USE_REAL_DOUBLE)
134061b8bfSJed Brown   const char precision[] = "Float64";
144061b8bfSJed Brown #else
154061b8bfSJed Brown   const char precision[] = "UnknownPrecision";
164061b8bfSJed Brown #endif
17ce94432eSBarry Smith   MPI_Comm                 comm;
182eaa9ef4SLisandro Dalcin   Vec                      Coords;
194061b8bfSJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
204061b8bfSJed Brown   PetscViewerVTKObjectLink link;
214061b8bfSJed Brown   FILE                     *fp;
224061b8bfSJed Brown   PetscMPIInt              rank,size,tag;
234061b8bfSJed Brown   DMDALocalInfo            info;
24*ea2d7708SPatrick Sanan   PetscInt                 dim,mx,my,mz,cdim,bs,boffset,maxnnodes,maxbs,i,j,k,r;
250298fd71SBarry Smith   PetscInt                 rloc[6],(*grloc)[6] = NULL;
264061b8bfSJed Brown   PetscScalar              *array,*array2;
274061b8bfSJed Brown   PetscReal                gmin[3],gmax[3];
284061b8bfSJed Brown   PetscErrorCode           ierr;
294061b8bfSJed Brown 
304061b8bfSJed Brown   PetscFunctionBegin;
31094921d9SJed Brown   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
324061b8bfSJed Brown #if defined(PETSC_USE_COMPLEX)
333bf036e2SBarry Smith   SETERRQ(comm,PETSC_ERR_SUP,"Complex values not supported");
344061b8bfSJed Brown #endif
354061b8bfSJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
364061b8bfSJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
374061b8bfSJed Brown   ierr = DMDAGetInfo(da,&dim,&mx,&my,&mz,0,0,0,&bs,0,0,0,0,0);CHKERRQ(ierr);
384061b8bfSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
394061b8bfSJed Brown   ierr = DMDAGetBoundingBox(da,gmin,gmax);CHKERRQ(ierr);
402eaa9ef4SLisandro Dalcin   ierr = DMGetCoordinates(da,&Coords);CHKERRQ(ierr);
412eaa9ef4SLisandro Dalcin   if (Coords) {
422eaa9ef4SLisandro Dalcin     PetscInt csize;
432eaa9ef4SLisandro Dalcin     ierr = VecGetSize(Coords,&csize);CHKERRQ(ierr);
442eaa9ef4SLisandro Dalcin     if (csize % (mx*my*mz)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Coordinate vector size mismatch");
452eaa9ef4SLisandro Dalcin     cdim = csize/(mx*my*mz);
462eaa9ef4SLisandro Dalcin     if (cdim < dim || cdim > 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Coordinate vector size mismatch");
472eaa9ef4SLisandro Dalcin   } else {
482eaa9ef4SLisandro Dalcin     cdim = dim;
492eaa9ef4SLisandro Dalcin   }
504061b8bfSJed Brown 
514061b8bfSJed Brown   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
524061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
53519f805aSKarl Rupp #if defined(PETSC_WORDS_BIGENDIAN)
544061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr);
554061b8bfSJed Brown #else
564061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr);
574061b8bfSJed Brown #endif
584061b8bfSJed 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);
594061b8bfSJed Brown 
60785e854fSJed Brown   if (!rank) {ierr = PetscMalloc1(size*6,&grloc);CHKERRQ(ierr);}
614061b8bfSJed Brown   rloc[0] = info.xs;
624061b8bfSJed Brown   rloc[1] = info.xm;
634061b8bfSJed Brown   rloc[2] = info.ys;
644061b8bfSJed Brown   rloc[3] = info.ym;
654061b8bfSJed Brown   rloc[4] = info.zs;
664061b8bfSJed Brown   rloc[5] = info.zm;
674061b8bfSJed Brown   ierr    = MPI_Gather(rloc,6,MPIU_INT,&grloc[0][0],6,MPIU_INT,0,comm);CHKERRQ(ierr);
684061b8bfSJed Brown 
694061b8bfSJed Brown   /* Write XML header */
704061b8bfSJed Brown   maxnnodes = 0;                /* Used for the temporary array size on rank 0 */
71*ea2d7708SPatrick Sanan   maxbs     = 0;                /* Used for the temporary array size on rank 0 */
724061b8bfSJed Brown   boffset   = 0;                /* Offset into binary file */
734061b8bfSJed Brown   for (r=0; r<size; r++) {
744061b8bfSJed Brown     PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
754061b8bfSJed Brown     if (!rank) {
764061b8bfSJed Brown       xs     = grloc[r][0];
774061b8bfSJed Brown       xm     = grloc[r][1];
784061b8bfSJed Brown       ys     = grloc[r][2];
794061b8bfSJed Brown       ym     = grloc[r][3];
804061b8bfSJed Brown       zs     = grloc[r][4];
814061b8bfSJed Brown       zm     = grloc[r][5];
824061b8bfSJed Brown       nnodes = xm*ym*zm;
834061b8bfSJed Brown     }
844061b8bfSJed Brown     maxnnodes = PetscMax(maxnnodes,nnodes);
854061b8bfSJed 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);
864061b8bfSJed Brown     ierr     = PetscFPrintf(comm,fp,"      <Points>\n");CHKERRQ(ierr);
874061b8bfSJed Brown     ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
884061b8bfSJed Brown     boffset += 3*nnodes*sizeof(PetscScalar) + sizeof(int);
894061b8bfSJed Brown     ierr     = PetscFPrintf(comm,fp,"      </Points>\n");CHKERRQ(ierr);
904061b8bfSJed Brown 
914061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"      <PointData Scalars=\"ScalarPointData\">\n");CHKERRQ(ierr);
924061b8bfSJed Brown     for (link=vtk->link; link; link=link->next) {
934061b8bfSJed Brown       Vec        X = (Vec)link->vec;
94*ea2d7708SPatrick Sanan       PetscInt   bs;
95*ea2d7708SPatrick Sanan       DM         daCurr;
96*ea2d7708SPatrick Sanan       const char *vecname = "Unnamed Vec data";
97*ea2d7708SPatrick Sanan 
98*ea2d7708SPatrick Sanan       ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr);
99*ea2d7708SPatrick Sanan       ierr = DMDAGetInfo(daCurr,0,0,0,0,0,0,0,&bs,0,0,0,0,0);CHKERRQ(ierr);
100*ea2d7708SPatrick Sanan       maxbs = PetscMax(maxbs,bs);
101*ea2d7708SPatrick Sanan 
102*ea2d7708SPatrick Sanan       if (((PetscObject)X)->name || link != vtk->link) {
1034061b8bfSJed Brown         ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
1044061b8bfSJed Brown       }
105*ea2d7708SPatrick Sanan       ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%D\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,bs,boffset);CHKERRQ(ierr);
106*ea2d7708SPatrick Sanan       boffset += bs*nnodes*sizeof(PetscScalar) + sizeof(int);
1074061b8bfSJed Brown     }
1084061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"      </PointData>\n");CHKERRQ(ierr);
1094061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"    </Piece>\n");CHKERRQ(ierr);
1104061b8bfSJed Brown   }
1114061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"  </StructuredGrid>\n");CHKERRQ(ierr);
1124061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
1134061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
1144061b8bfSJed Brown 
1154061b8bfSJed Brown   /* Now write the arrays. */
1164061b8bfSJed Brown   tag  = ((PetscObject)viewer)->tag;
117*ea2d7708SPatrick Sanan   ierr = PetscMalloc2(maxnnodes*PetscMax(3,maxbs),&array,maxnnodes*3,&array2);CHKERRQ(ierr);
1184061b8bfSJed Brown   for (r=0; r<size; r++) {
1194061b8bfSJed Brown     MPI_Status status;
1204061b8bfSJed Brown     PetscInt   xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
1214061b8bfSJed Brown     if (!rank) {
1224061b8bfSJed Brown       xs     = grloc[r][0];
1234061b8bfSJed Brown       xm     = grloc[r][1];
1244061b8bfSJed Brown       ys     = grloc[r][2];
1254061b8bfSJed Brown       ym     = grloc[r][3];
1264061b8bfSJed Brown       zs     = grloc[r][4];
1274061b8bfSJed Brown       zm     = grloc[r][5];
1284061b8bfSJed Brown       nnodes = xm*ym*zm;
1294061b8bfSJed Brown     } else if (r == rank) {
1304061b8bfSJed Brown       nnodes = info.xm*info.ym*info.zm;
1314061b8bfSJed Brown     }
1324061b8bfSJed Brown 
1332eaa9ef4SLisandro Dalcin     /* Write the coordinates */
1344061b8bfSJed Brown     if (Coords) {
1354061b8bfSJed Brown       const PetscScalar *coords;
1364061b8bfSJed Brown       ierr = VecGetArrayRead(Coords,&coords);CHKERRQ(ierr);
1374061b8bfSJed Brown       if (!rank) {
1384061b8bfSJed Brown         if (r) {
1396622f924SJed Brown           PetscMPIInt nn;
1402eaa9ef4SLisandro Dalcin           ierr = MPI_Recv(array,nnodes*cdim,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr);
1414061b8bfSJed Brown           ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr);
1422eaa9ef4SLisandro Dalcin           if (nn != nnodes*cdim) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
1434061b8bfSJed Brown         } else {
1442eaa9ef4SLisandro Dalcin           ierr = PetscMemcpy(array,coords,nnodes*cdim*sizeof(PetscScalar));CHKERRQ(ierr);
1454061b8bfSJed Brown         }
1464061b8bfSJed Brown         /* Transpose coordinates to VTK (C-style) ordering */
1474061b8bfSJed Brown         for (k=0; k<zm; k++) {
1484061b8bfSJed Brown           for (j=0; j<ym; j++) {
1494061b8bfSJed Brown             for (i=0; i<xm; i++) {
1504061b8bfSJed Brown               PetscInt Iloc = i+xm*(j+ym*k);
1512eaa9ef4SLisandro Dalcin               array2[Iloc*3+0] = array[Iloc*cdim + 0];
1522eaa9ef4SLisandro Dalcin               array2[Iloc*3+1] = cdim > 1 ? array[Iloc*cdim + 1] : 0.0;
1532eaa9ef4SLisandro Dalcin               array2[Iloc*3+2] = cdim > 2 ? array[Iloc*cdim + 2] : 0.0;
1544061b8bfSJed Brown             }
1554061b8bfSJed Brown           }
1564061b8bfSJed Brown         }
1574061b8bfSJed Brown       } else if (r == rank) {
1582eaa9ef4SLisandro Dalcin         ierr = MPI_Send((void*)coords,nnodes*cdim,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr);
1594061b8bfSJed Brown       }
1604061b8bfSJed Brown       ierr = VecRestoreArrayRead(Coords,&coords);CHKERRQ(ierr);
1614061b8bfSJed Brown     } else {       /* Fabricate some coordinates using grid index */
1624061b8bfSJed Brown       for (k=0; k<zm; k++) {
1634061b8bfSJed Brown         for (j=0; j<ym; j++) {
1644061b8bfSJed Brown           for (i=0; i<xm; i++) {
1654061b8bfSJed Brown             PetscInt Iloc = i+xm*(j+ym*k);
1664061b8bfSJed Brown             array2[Iloc*3+0] = xs+i;
1674061b8bfSJed Brown             array2[Iloc*3+1] = ys+j;
1684061b8bfSJed Brown             array2[Iloc*3+2] = zs+k;
1694061b8bfSJed Brown           }
1704061b8bfSJed Brown         }
1714061b8bfSJed Brown       }
1724061b8bfSJed Brown     }
17394fbd55eSBarry Smith     ierr = PetscViewerVTKFWrite(viewer,fp,array2,nnodes*3,MPIU_SCALAR);CHKERRQ(ierr);
1744061b8bfSJed Brown 
1754061b8bfSJed Brown     /* Write each of the objects queued up for this file */
1764061b8bfSJed Brown     for (link=vtk->link; link; link=link->next) {
1774061b8bfSJed Brown       Vec               X = (Vec)link->vec;
1784061b8bfSJed Brown       const PetscScalar *x;
179*ea2d7708SPatrick Sanan       PetscInt          bs;
180*ea2d7708SPatrick Sanan       DM                daCurr;
181*ea2d7708SPatrick Sanan       ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr);
182*ea2d7708SPatrick Sanan       ierr = DMDAGetInfo(daCurr,0,0,0,0, 0,0,0,&bs,0,0,0,0,0);CHKERRQ(ierr);
1834061b8bfSJed Brown       ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
1844061b8bfSJed Brown       if (!rank) {
1854061b8bfSJed Brown         if (r) {
1866622f924SJed Brown           PetscMPIInt nn;
1874061b8bfSJed Brown           ierr = MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr);
1884061b8bfSJed Brown           ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr);
1894061b8bfSJed Brown           if (nn != nnodes*bs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %D",r);
1904061b8bfSJed Brown         } else {
1914061b8bfSJed Brown           ierr = PetscMemcpy(array,x,nnodes*bs*sizeof(PetscScalar));CHKERRQ(ierr);
1924061b8bfSJed Brown         }
193*ea2d7708SPatrick Sanan         ierr = PetscViewerVTKFWrite(viewer,fp,array,bs*nnodes,MPIU_SCALAR);CHKERRQ(ierr);
1944061b8bfSJed Brown       } else if (r == rank) {
1954061b8bfSJed Brown         ierr = MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr);
1964061b8bfSJed Brown       }
1974061b8bfSJed Brown       ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
1984061b8bfSJed Brown     }
1994061b8bfSJed Brown   }
2004061b8bfSJed Brown   ierr = PetscFree2(array,array2);CHKERRQ(ierr);
2014061b8bfSJed Brown   ierr = PetscFree(grloc);CHKERRQ(ierr);
2024061b8bfSJed Brown 
2034061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr);
2044061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
2054061b8bfSJed Brown   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
2064061b8bfSJed Brown   PetscFunctionReturn(0);
2074061b8bfSJed Brown }
2084061b8bfSJed Brown 
2094061b8bfSJed Brown 
210a13bc4e3SShao-Ching Huang static PetscErrorCode DMDAVTKWriteAll_VTR(DM da,PetscViewer viewer)
211a13bc4e3SShao-Ching Huang {
212a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_REAL_SINGLE)
213a13bc4e3SShao-Ching Huang   const char precision[] = "Float32";
214a13bc4e3SShao-Ching Huang #elif defined(PETSC_USE_REAL_DOUBLE)
215a13bc4e3SShao-Ching Huang   const char precision[] = "Float64";
216a13bc4e3SShao-Ching Huang #else
217a13bc4e3SShao-Ching Huang   const char precision[] = "UnknownPrecision";
218a13bc4e3SShao-Ching Huang #endif
219a13bc4e3SShao-Ching Huang   MPI_Comm                 comm;
220a13bc4e3SShao-Ching Huang   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
221a13bc4e3SShao-Ching Huang   PetscViewerVTKObjectLink link;
222a13bc4e3SShao-Ching Huang   FILE                     *fp;
223a13bc4e3SShao-Ching Huang   PetscMPIInt              rank,size,tag;
224a13bc4e3SShao-Ching Huang   DMDALocalInfo            info;
225*ea2d7708SPatrick Sanan   PetscInt                 dim,mx,my,mz,boffset,maxnnodes,maxbs,i,j,k,r;
226a13bc4e3SShao-Ching Huang   PetscInt                 rloc[6],(*grloc)[6] = NULL;
227a13bc4e3SShao-Ching Huang   PetscScalar              *array,*array2;
228a13bc4e3SShao-Ching Huang   PetscReal                gmin[3],gmax[3];
229a13bc4e3SShao-Ching Huang   PetscErrorCode           ierr;
230a13bc4e3SShao-Ching Huang 
231a13bc4e3SShao-Ching Huang   PetscFunctionBegin;
232a13bc4e3SShao-Ching Huang   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
233a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_COMPLEX)
234a13bc4e3SShao-Ching Huang   SETERRQ(comm,PETSC_ERR_SUP,"Complex values not supported");
235a13bc4e3SShao-Ching Huang #endif
236a13bc4e3SShao-Ching Huang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
237a13bc4e3SShao-Ching Huang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
238*ea2d7708SPatrick Sanan   ierr = DMDAGetInfo(da,&dim,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr);
239a13bc4e3SShao-Ching Huang   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
240a13bc4e3SShao-Ching Huang   ierr = DMDAGetBoundingBox(da,gmin,gmax);CHKERRQ(ierr);
241a13bc4e3SShao-Ching Huang   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
242a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
243a13bc4e3SShao-Ching Huang #if defined(PETSC_WORDS_BIGENDIAN)
244a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"RectilinearGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr);
245a13bc4e3SShao-Ching Huang #else
246a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"RectilinearGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr);
247a13bc4e3SShao-Ching Huang #endif
248a13bc4e3SShao-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);
249a13bc4e3SShao-Ching Huang 
250785e854fSJed Brown   if (!rank) {ierr = PetscMalloc1(size*6,&grloc);CHKERRQ(ierr);}
251a13bc4e3SShao-Ching Huang   rloc[0] = info.xs;
252a13bc4e3SShao-Ching Huang   rloc[1] = info.xm;
253a13bc4e3SShao-Ching Huang   rloc[2] = info.ys;
254a13bc4e3SShao-Ching Huang   rloc[3] = info.ym;
255a13bc4e3SShao-Ching Huang   rloc[4] = info.zs;
256a13bc4e3SShao-Ching Huang   rloc[5] = info.zm;
257a13bc4e3SShao-Ching Huang   ierr    = MPI_Gather(rloc,6,MPIU_INT,&grloc[0][0],6,MPIU_INT,0,comm);CHKERRQ(ierr);
258a13bc4e3SShao-Ching Huang 
259a13bc4e3SShao-Ching Huang   /* Write XML header */
260a13bc4e3SShao-Ching Huang   maxnnodes = 0;                /* Used for the temporary array size on rank 0 */
261*ea2d7708SPatrick Sanan   maxbs     = 0;                /* Used for the temporary array size on rank 0 */
262a13bc4e3SShao-Ching Huang   boffset   = 0;                /* Offset into binary file */
263a13bc4e3SShao-Ching Huang   for (r=0; r<size; r++) {
264a13bc4e3SShao-Ching Huang     PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
265a13bc4e3SShao-Ching Huang     if (!rank) {
266a13bc4e3SShao-Ching Huang       xs     = grloc[r][0];
267a13bc4e3SShao-Ching Huang       xm     = grloc[r][1];
268a13bc4e3SShao-Ching Huang       ys     = grloc[r][2];
269a13bc4e3SShao-Ching Huang       ym     = grloc[r][3];
270a13bc4e3SShao-Ching Huang       zs     = grloc[r][4];
271a13bc4e3SShao-Ching Huang       zm     = grloc[r][5];
272a13bc4e3SShao-Ching Huang       nnodes = xm*ym*zm;
273a13bc4e3SShao-Ching Huang     }
274a13bc4e3SShao-Ching Huang     maxnnodes = PetscMax(maxnnodes,nnodes);
275a13bc4e3SShao-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);
276a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"      <Coordinates>\n");CHKERRQ(ierr);
277a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Xcoord\"  format=\"appended\"  offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
278a13bc4e3SShao-Ching Huang     boffset += xm*sizeof(PetscScalar) + sizeof(int);
279a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Ycoord\"  format=\"appended\"  offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
280a13bc4e3SShao-Ching Huang     boffset += ym*sizeof(PetscScalar) + sizeof(int);
281a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Zcoord\"  format=\"appended\"  offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
282a13bc4e3SShao-Ching Huang     boffset += zm*sizeof(PetscScalar) + sizeof(int);
283a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"      </Coordinates>\n");CHKERRQ(ierr);
284a13bc4e3SShao-Ching Huang     ierr = PetscFPrintf(comm,fp,"      <PointData Scalars=\"ScalarPointData\">\n");CHKERRQ(ierr);
285a13bc4e3SShao-Ching Huang     for (link=vtk->link; link; link=link->next) {
286a13bc4e3SShao-Ching Huang       Vec        X = (Vec)link->vec;
287*ea2d7708SPatrick Sanan       PetscInt   bs;
288*ea2d7708SPatrick Sanan       DM         daCurr;
289*ea2d7708SPatrick Sanan       const char *vecname = "Unnamed Vec data";
290*ea2d7708SPatrick Sanan 
291*ea2d7708SPatrick Sanan       ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr);
292*ea2d7708SPatrick Sanan       ierr = DMDAGetInfo(daCurr,0,0,0,0,0,0,0,&bs,0,0,0,0,0);CHKERRQ(ierr);
293*ea2d7708SPatrick Sanan       maxbs = PetscMax(maxbs,bs);
294*ea2d7708SPatrick Sanan       if (((PetscObject)X)->name || link != vtk->link) {
295a13bc4e3SShao-Ching Huang         ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
296a13bc4e3SShao-Ching Huang       }
297*ea2d7708SPatrick Sanan 
298*ea2d7708SPatrick Sanan      ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%D\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,bs,boffset);CHKERRQ(ierr);
299*ea2d7708SPatrick Sanan      boffset += bs*nnodes*sizeof(PetscScalar) + sizeof(int);
300a13bc4e3SShao-Ching Huang     }
301a13bc4e3SShao-Ching Huang     ierr = PetscFPrintf(comm,fp,"      </PointData>\n");CHKERRQ(ierr);
302a13bc4e3SShao-Ching Huang     ierr = PetscFPrintf(comm,fp,"    </Piece>\n");CHKERRQ(ierr);
303a13bc4e3SShao-Ching Huang   }
304a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"  </RectilinearGrid>\n");CHKERRQ(ierr);
305a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
306a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
307a13bc4e3SShao-Ching Huang 
308a13bc4e3SShao-Ching Huang   /* Now write the arrays. */
309a13bc4e3SShao-Ching Huang   tag  = ((PetscObject)viewer)->tag;
310*ea2d7708SPatrick Sanan   ierr = PetscMalloc2(PetscMax(maxnnodes*maxbs,info.xm+info.ym+info.zm),&array,info.xm+info.ym+info.zm,&array2);CHKERRQ(ierr);
311a13bc4e3SShao-Ching Huang   for (r=0; r<size; r++) {
312a13bc4e3SShao-Ching Huang     MPI_Status status;
313a13bc4e3SShao-Ching Huang     PetscInt   xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
314a13bc4e3SShao-Ching Huang     if (!rank) {
315a13bc4e3SShao-Ching Huang       xs     = grloc[r][0];
316a13bc4e3SShao-Ching Huang       xm     = grloc[r][1];
317a13bc4e3SShao-Ching Huang       ys     = grloc[r][2];
318a13bc4e3SShao-Ching Huang       ym     = grloc[r][3];
319a13bc4e3SShao-Ching Huang       zs     = grloc[r][4];
320a13bc4e3SShao-Ching Huang       zm     = grloc[r][5];
321a13bc4e3SShao-Ching Huang       nnodes = xm*ym*zm;
322a13bc4e3SShao-Ching Huang     } else if (r == rank) {
323a13bc4e3SShao-Ching Huang       nnodes = info.xm*info.ym*info.zm;
324a13bc4e3SShao-Ching Huang     }
325a13bc4e3SShao-Ching Huang     {                           /* Write the coordinates */
326a13bc4e3SShao-Ching Huang       Vec Coords;
327*ea2d7708SPatrick Sanan 
328a13bc4e3SShao-Ching Huang       ierr = DMGetCoordinates(da,&Coords);CHKERRQ(ierr);
329a13bc4e3SShao-Ching Huang       if (Coords) {
330a13bc4e3SShao-Ching Huang         const PetscScalar *coords;
331a13bc4e3SShao-Ching Huang         ierr = VecGetArrayRead(Coords,&coords);CHKERRQ(ierr);
332a13bc4e3SShao-Ching Huang         if (!rank) {
333a13bc4e3SShao-Ching Huang           if (r) {
334a13bc4e3SShao-Ching Huang             PetscMPIInt nn;
335a13bc4e3SShao-Ching Huang             ierr = MPI_Recv(array,xm+ym+zm,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr);
336a13bc4e3SShao-Ching Huang             ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr);
337a13bc4e3SShao-Ching Huang             if (nn != xm+ym+zm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
338a13bc4e3SShao-Ching Huang           } else {
339a13bc4e3SShao-Ching Huang             /* x coordinates */
340a13bc4e3SShao-Ching Huang             for (j=0, k=0, i=0; i<xm; i++) {
341a13bc4e3SShao-Ching Huang               PetscInt Iloc = i+xm*(j+ym*k);
342a13bc4e3SShao-Ching Huang               array[i] = coords[Iloc*dim + 0];
343a13bc4e3SShao-Ching Huang             }
344a13bc4e3SShao-Ching Huang             /* y coordinates */
345a13bc4e3SShao-Ching Huang             for (i=0, k=0, j=0; j<ym; j++) {
346a13bc4e3SShao-Ching Huang               PetscInt Iloc = i+xm*(j+ym*k);
347a13bc4e3SShao-Ching Huang               array[j+xm] = dim > 1 ? coords[Iloc*dim + 1] : 0;
348a13bc4e3SShao-Ching Huang             }
349a13bc4e3SShao-Ching Huang             /* z coordinates */
350a13bc4e3SShao-Ching Huang             for (i=0, j=0, k=0; k<zm; k++) {
351a13bc4e3SShao-Ching Huang               PetscInt Iloc = i+xm*(j+ym*k);
352a13bc4e3SShao-Ching Huang               array[k+xm+ym] = dim > 2 ? coords[Iloc*dim + 2] : 0;
353a13bc4e3SShao-Ching Huang             }
354a13bc4e3SShao-Ching Huang           }
355a13bc4e3SShao-Ching Huang         } else if (r == rank) {
356a13bc4e3SShao-Ching Huang           xm = info.xm;
357a13bc4e3SShao-Ching Huang           ym = info.ym;
358a13bc4e3SShao-Ching Huang           zm = info.zm;
359a13bc4e3SShao-Ching Huang           for (j=0, k=0, i=0; i<xm; i++) {
360a13bc4e3SShao-Ching Huang             PetscInt Iloc = i+xm*(j+ym*k);
361a13bc4e3SShao-Ching Huang             array2[i] = coords[Iloc*dim + 0];
362a13bc4e3SShao-Ching Huang           }
363a13bc4e3SShao-Ching Huang           for (i=0, k=0, j=0; j<ym; j++) {
364a13bc4e3SShao-Ching Huang             PetscInt Iloc = i+xm*(j+ym*k);
365a13bc4e3SShao-Ching Huang             array2[j+xm] = dim > 1 ? coords[Iloc*dim + 1] : 0;
366a13bc4e3SShao-Ching Huang           }
367a13bc4e3SShao-Ching Huang           for (i=0, j=0, k=0; k<zm; k++) {
368a13bc4e3SShao-Ching Huang             PetscInt Iloc = i+xm*(j+ym*k);
369a13bc4e3SShao-Ching Huang             array2[k+xm+ym] = dim > 2 ? coords[Iloc*dim + 2] : 0;
370a13bc4e3SShao-Ching Huang           }
371a13bc4e3SShao-Ching Huang           ierr = MPI_Send((void*)array2,xm+ym+zm,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr);
372a13bc4e3SShao-Ching Huang         }
373a13bc4e3SShao-Ching Huang         ierr = VecRestoreArrayRead(Coords,&coords);CHKERRQ(ierr);
374a13bc4e3SShao-Ching Huang       } else {       /* Fabricate some coordinates using grid index */
3753d3eaba7SBarry Smith         for (i=0; i<xm; i++) {
376a13bc4e3SShao-Ching Huang           array[i] = xs+i;
377a13bc4e3SShao-Ching Huang         }
3783d3eaba7SBarry Smith         for (j=0; j<ym; j++) {
379a13bc4e3SShao-Ching Huang           array[j+xm] = ys+j;
380a13bc4e3SShao-Ching Huang         }
3813d3eaba7SBarry Smith         for (k=0; k<zm; k++) {
382a13bc4e3SShao-Ching Huang           array[k+xm+ym] = zs+k;
383a13bc4e3SShao-Ching Huang         }
384a13bc4e3SShao-Ching Huang       }
385a13bc4e3SShao-Ching Huang       if (!rank) {
38694fbd55eSBarry Smith         ierr = PetscViewerVTKFWrite(viewer,fp,&(array[0]),xm,MPIU_SCALAR);CHKERRQ(ierr);
38794fbd55eSBarry Smith         ierr = PetscViewerVTKFWrite(viewer,fp,&(array[xm]),ym,MPIU_SCALAR);CHKERRQ(ierr);
38894fbd55eSBarry Smith         ierr = PetscViewerVTKFWrite(viewer,fp,&(array[xm+ym]),zm,MPIU_SCALAR);CHKERRQ(ierr);
389a13bc4e3SShao-Ching Huang       }
390a13bc4e3SShao-Ching Huang     }
391a13bc4e3SShao-Ching Huang 
392a13bc4e3SShao-Ching Huang     /* Write each of the objects queued up for this file */
393a13bc4e3SShao-Ching Huang     for (link=vtk->link; link; link=link->next) {
394a13bc4e3SShao-Ching Huang       Vec               X = (Vec)link->vec;
395a13bc4e3SShao-Ching Huang       const PetscScalar *x;
396*ea2d7708SPatrick Sanan       PetscInt          bs;
397*ea2d7708SPatrick Sanan       DM                daCurr;
398*ea2d7708SPatrick Sanan       ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr);
399*ea2d7708SPatrick Sanan       ierr = DMDAGetInfo(daCurr,0,0,0,0,0,0,0,&bs,0,0,0,0,0);CHKERRQ(ierr);
400a13bc4e3SShao-Ching Huang 
401a13bc4e3SShao-Ching Huang       ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
402a13bc4e3SShao-Ching Huang       if (!rank) {
403a13bc4e3SShao-Ching Huang         if (r) {
404a13bc4e3SShao-Ching Huang           PetscMPIInt nn;
405a13bc4e3SShao-Ching Huang           ierr = MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr);
406a13bc4e3SShao-Ching Huang           ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr);
407a13bc4e3SShao-Ching Huang           if (nn != nnodes*bs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %D",r);
408a13bc4e3SShao-Ching Huang         } else {
409a13bc4e3SShao-Ching Huang           ierr = PetscMemcpy(array,x,nnodes*bs*sizeof(PetscScalar));CHKERRQ(ierr);
410a13bc4e3SShao-Ching Huang         }
411*ea2d7708SPatrick Sanan         ierr = PetscViewerVTKFWrite(viewer,fp,array,nnodes*bs,MPIU_SCALAR);CHKERRQ(ierr);
412*ea2d7708SPatrick Sanan 
413a13bc4e3SShao-Ching Huang       } else if (r == rank) {
414a13bc4e3SShao-Ching Huang         ierr = MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr);
415a13bc4e3SShao-Ching Huang       }
416a13bc4e3SShao-Ching Huang       ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
417a13bc4e3SShao-Ching Huang     }
418a13bc4e3SShao-Ching Huang   }
419a13bc4e3SShao-Ching Huang   ierr = PetscFree2(array,array2);CHKERRQ(ierr);
420a13bc4e3SShao-Ching Huang   ierr = PetscFree(grloc);CHKERRQ(ierr);
421a13bc4e3SShao-Ching Huang 
422a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr);
423a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
424a13bc4e3SShao-Ching Huang   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
425a13bc4e3SShao-Ching Huang   PetscFunctionReturn(0);
426a13bc4e3SShao-Ching Huang }
427a13bc4e3SShao-Ching Huang 
4284061b8bfSJed Brown /*@C
4294061b8bfSJed Brown    DMDAVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
4304061b8bfSJed Brown 
4314061b8bfSJed Brown    Collective
4324061b8bfSJed Brown 
4334061b8bfSJed Brown    Input Arguments:
4344061b8bfSJed Brown    odm - DM specifying the grid layout, passed as a PetscObject
4354061b8bfSJed Brown    viewer - viewer of type VTK
4364061b8bfSJed Brown 
4374061b8bfSJed Brown    Level: developer
4384061b8bfSJed Brown 
4394061b8bfSJed Brown    Note:
4404061b8bfSJed Brown    This function is a callback used by the VTK viewer to actually write the file.
4414061b8bfSJed 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.
4424061b8bfSJed Brown    Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
4434061b8bfSJed Brown 
4444061b8bfSJed Brown .seealso: PETSCVIEWERVTK
4454061b8bfSJed Brown @*/
4464061b8bfSJed Brown PetscErrorCode DMDAVTKWriteAll(PetscObject odm,PetscViewer viewer)
4474061b8bfSJed Brown {
4484061b8bfSJed Brown   DM             dm = (DM)odm;
4494061b8bfSJed Brown   PetscBool      isvtk;
4504061b8bfSJed Brown   PetscErrorCode ierr;
4514061b8bfSJed Brown 
4524061b8bfSJed Brown   PetscFunctionBegin;
4534061b8bfSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4544061b8bfSJed Brown   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
455251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr);
456ce94432eSBarry Smith   if (!isvtk) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_INCOMP,"Cannot use viewer type %s",((PetscObject)viewer)->type_name);
4574061b8bfSJed Brown   switch (viewer->format) {
4584061b8bfSJed Brown   case PETSC_VIEWER_VTK_VTS:
4594061b8bfSJed Brown     ierr = DMDAVTKWriteAll_VTS(dm,viewer);CHKERRQ(ierr);
4604061b8bfSJed Brown     break;
461a13bc4e3SShao-Ching Huang   case PETSC_VIEWER_VTK_VTR:
462a13bc4e3SShao-Ching Huang     ierr = DMDAVTKWriteAll_VTR(dm,viewer);CHKERRQ(ierr);
463a13bc4e3SShao-Ching Huang     break;
464ce94432eSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No support for format '%s'",PetscViewerFormats[viewer->format]);
4654061b8bfSJed Brown   }
4664061b8bfSJed Brown   PetscFunctionReturn(0);
4674061b8bfSJed Brown }
468