xref: /petsc/src/dm/impls/da/grvtk.c (revision 3d3eaba7093fd3d75f596ec754341dee3aba7588)
1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>
25c6c1daeSBarry Smith #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h>
34061b8bfSJed Brown 
44061b8bfSJed Brown #undef __FUNCT__
54061b8bfSJed Brown #define __FUNCT__ "DMDAVTKWriteAll_VTS"
64061b8bfSJed Brown static PetscErrorCode DMDAVTKWriteAll_VTS(DM da,PetscViewer viewer)
74061b8bfSJed Brown {
84061b8bfSJed Brown #if defined(PETSC_USE_REAL_SINGLE)
94061b8bfSJed Brown   const char precision[] = "Float32";
104061b8bfSJed Brown #elif defined(PETSC_USE_REAL_DOUBLE)
114061b8bfSJed Brown   const char precision[] = "Float64";
124061b8bfSJed Brown #else
134061b8bfSJed Brown   const char precision[] = "UnknownPrecision";
144061b8bfSJed Brown #endif
15ce94432eSBarry Smith   MPI_Comm                 comm;
162eaa9ef4SLisandro Dalcin   Vec                      Coords;
174061b8bfSJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
184061b8bfSJed Brown   PetscViewerVTKObjectLink link;
194061b8bfSJed Brown   FILE                     *fp;
204061b8bfSJed Brown   PetscMPIInt              rank,size,tag;
214061b8bfSJed Brown   DMDALocalInfo            info;
222eaa9ef4SLisandro Dalcin   PetscInt                 dim,mx,my,mz,cdim,bs,boffset,maxnnodes,i,j,k,f,r;
230298fd71SBarry Smith   PetscInt                 rloc[6],(*grloc)[6] = NULL;
244061b8bfSJed Brown   PetscScalar              *array,*array2;
254061b8bfSJed Brown   PetscReal                gmin[3],gmax[3];
264061b8bfSJed Brown   PetscErrorCode           ierr;
274061b8bfSJed Brown 
284061b8bfSJed Brown   PetscFunctionBegin;
29094921d9SJed Brown   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
304061b8bfSJed Brown #if defined(PETSC_USE_COMPLEX)
313bf036e2SBarry Smith   SETERRQ(comm,PETSC_ERR_SUP,"Complex values not supported");
324061b8bfSJed Brown #endif
334061b8bfSJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
344061b8bfSJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
354061b8bfSJed Brown   ierr = DMDAGetInfo(da,&dim, &mx,&my,&mz, 0,0,0, &bs,0,0,0,0,0);CHKERRQ(ierr);
364061b8bfSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
374061b8bfSJed Brown   ierr = DMDAGetBoundingBox(da,gmin,gmax);CHKERRQ(ierr);
382eaa9ef4SLisandro Dalcin   ierr = DMGetCoordinates(da,&Coords);CHKERRQ(ierr);
392eaa9ef4SLisandro Dalcin   if (Coords) {
402eaa9ef4SLisandro Dalcin     PetscInt csize;
412eaa9ef4SLisandro Dalcin     ierr = VecGetSize(Coords,&csize);CHKERRQ(ierr);
422eaa9ef4SLisandro Dalcin     if (csize % (mx*my*mz)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Coordinate vector size mismatch");
432eaa9ef4SLisandro Dalcin     cdim = csize/(mx*my*mz);
442eaa9ef4SLisandro Dalcin     if (cdim < dim || cdim > 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Coordinate vector size mismatch");
452eaa9ef4SLisandro Dalcin   } else {
462eaa9ef4SLisandro Dalcin     cdim = dim;
472eaa9ef4SLisandro Dalcin   }
484061b8bfSJed Brown 
494061b8bfSJed Brown   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
504061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
51519f805aSKarl Rupp #if defined(PETSC_WORDS_BIGENDIAN)
524061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr);
534061b8bfSJed Brown #else
544061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr);
554061b8bfSJed Brown #endif
564061b8bfSJed 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);
574061b8bfSJed Brown 
58785e854fSJed Brown   if (!rank) {ierr = PetscMalloc1(size*6,&grloc);CHKERRQ(ierr);}
594061b8bfSJed Brown   rloc[0] = info.xs;
604061b8bfSJed Brown   rloc[1] = info.xm;
614061b8bfSJed Brown   rloc[2] = info.ys;
624061b8bfSJed Brown   rloc[3] = info.ym;
634061b8bfSJed Brown   rloc[4] = info.zs;
644061b8bfSJed Brown   rloc[5] = info.zm;
654061b8bfSJed Brown   ierr    = MPI_Gather(rloc,6,MPIU_INT,&grloc[0][0],6,MPIU_INT,0,comm);CHKERRQ(ierr);
664061b8bfSJed Brown 
674061b8bfSJed Brown   /* Write XML header */
684061b8bfSJed Brown   maxnnodes = 0;                /* Used for the temporary array size on rank 0 */
694061b8bfSJed Brown   boffset   = 0;                /* Offset into binary file */
704061b8bfSJed Brown   for (r=0; r<size; r++) {
714061b8bfSJed Brown     PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
724061b8bfSJed Brown     if (!rank) {
734061b8bfSJed Brown       xs     = grloc[r][0];
744061b8bfSJed Brown       xm     = grloc[r][1];
754061b8bfSJed Brown       ys     = grloc[r][2];
764061b8bfSJed Brown       ym     = grloc[r][3];
774061b8bfSJed Brown       zs     = grloc[r][4];
784061b8bfSJed Brown       zm     = grloc[r][5];
794061b8bfSJed Brown       nnodes = xm*ym*zm;
804061b8bfSJed Brown     }
814061b8bfSJed Brown     maxnnodes = PetscMax(maxnnodes,nnodes);
824061b8bfSJed Brown #if 0
834061b8bfSJed Brown     switch (dim) {
844061b8bfSJed Brown     case 1:
854061b8bfSJed Brown       ierr = PetscFPrintf(comm,fp,"    <Piece Extent=\"%D %D %D %D %D %D\">\n",xs,xs+xm-1,0,0,0,0);CHKERRQ(ierr);
864061b8bfSJed Brown       break;
874061b8bfSJed Brown     case 2:
884061b8bfSJed Brown       ierr = PetscFPrintf(comm,fp,"    <Piece Extent=\"%D %D %D %D %D %D\">\n",xs,xs+xm,ys+ym-1,xs,xs+xm-1,0,0);CHKERRQ(ierr);
894061b8bfSJed Brown       break;
904061b8bfSJed Brown     case 3:
914061b8bfSJed 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);
924061b8bfSJed Brown       break;
93ce94432eSBarry Smith     default: SETERRQ1(PetscObjectComm((PetscObject)da),PETSC_ERR_SUP,"No support for dimension %D",dim);
944061b8bfSJed Brown     }
954061b8bfSJed Brown #endif
964061b8bfSJed 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);
974061b8bfSJed Brown     ierr     = PetscFPrintf(comm,fp,"      <Points>\n");CHKERRQ(ierr);
984061b8bfSJed Brown     ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
994061b8bfSJed Brown     boffset += 3*nnodes*sizeof(PetscScalar) + sizeof(int);
1004061b8bfSJed Brown     ierr     = PetscFPrintf(comm,fp,"      </Points>\n");CHKERRQ(ierr);
1014061b8bfSJed Brown 
1024061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"      <PointData Scalars=\"ScalarPointData\">\n");CHKERRQ(ierr);
1034061b8bfSJed Brown     for (link=vtk->link; link; link=link->next) {
1044061b8bfSJed Brown       Vec        X        = (Vec)link->vec;
1054061b8bfSJed Brown       const char *vecname = "";
1064061b8bfSJed Brown       if (((PetscObject)X)->name || link != vtk->link) { /* If the object is already named, use it. If it is past the first link, name it to disambiguate. */
1074061b8bfSJed Brown         ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
1084061b8bfSJed Brown       }
1094061b8bfSJed Brown       for (i=0; i<bs; i++) {
1104061b8bfSJed Brown         char       buf[256];
1114061b8bfSJed Brown         const char *fieldname;
1124061b8bfSJed Brown         ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
1134061b8bfSJed Brown         if (!fieldname) {
1148caf3d72SBarry Smith           ierr      = PetscSNPrintf(buf,sizeof(buf),"Unnamed%D",i);CHKERRQ(ierr);
1154061b8bfSJed Brown           fieldname = buf;
1164061b8bfSJed Brown         }
1174061b8bfSJed Brown         ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
1184061b8bfSJed Brown         boffset += nnodes*sizeof(PetscScalar) + sizeof(int);
1194061b8bfSJed Brown       }
1204061b8bfSJed Brown     }
1214061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"      </PointData>\n");CHKERRQ(ierr);
1224061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"    </Piece>\n");CHKERRQ(ierr);
1234061b8bfSJed Brown   }
1244061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"  </StructuredGrid>\n");CHKERRQ(ierr);
1254061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
1264061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
1274061b8bfSJed Brown 
1284061b8bfSJed Brown   /* Now write the arrays. */
1294061b8bfSJed Brown   tag  = ((PetscObject)viewer)->tag;
130dcca6d9dSJed Brown   ierr = PetscMalloc2(maxnnodes*PetscMax(3,bs),&array,maxnnodes*3,&array2);CHKERRQ(ierr);
1314061b8bfSJed Brown   for (r=0; r<size; r++) {
1324061b8bfSJed Brown     MPI_Status status;
1334061b8bfSJed Brown     PetscInt   xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
1344061b8bfSJed Brown     if (!rank) {
1354061b8bfSJed Brown       xs     = grloc[r][0];
1364061b8bfSJed Brown       xm     = grloc[r][1];
1374061b8bfSJed Brown       ys     = grloc[r][2];
1384061b8bfSJed Brown       ym     = grloc[r][3];
1394061b8bfSJed Brown       zs     = grloc[r][4];
1404061b8bfSJed Brown       zm     = grloc[r][5];
1414061b8bfSJed Brown       nnodes = xm*ym*zm;
1424061b8bfSJed Brown     } else if (r == rank) {
1434061b8bfSJed Brown       nnodes = info.xm*info.ym*info.zm;
1444061b8bfSJed Brown     }
1454061b8bfSJed Brown 
1462eaa9ef4SLisandro Dalcin     /* Write the coordinates */
1474061b8bfSJed Brown     if (Coords) {
1484061b8bfSJed Brown       const PetscScalar *coords;
1494061b8bfSJed Brown       ierr = VecGetArrayRead(Coords,&coords);CHKERRQ(ierr);
1504061b8bfSJed Brown       if (!rank) {
1514061b8bfSJed Brown         if (r) {
1526622f924SJed Brown           PetscMPIInt nn;
1532eaa9ef4SLisandro Dalcin           ierr = MPI_Recv(array,nnodes*cdim,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr);
1544061b8bfSJed Brown           ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr);
1552eaa9ef4SLisandro Dalcin           if (nn != nnodes*cdim) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
1564061b8bfSJed Brown         } else {
1572eaa9ef4SLisandro Dalcin           ierr = PetscMemcpy(array,coords,nnodes*cdim*sizeof(PetscScalar));CHKERRQ(ierr);
1584061b8bfSJed Brown         }
1594061b8bfSJed Brown         /* Transpose coordinates to VTK (C-style) ordering */
1604061b8bfSJed Brown         for (k=0; k<zm; k++) {
1614061b8bfSJed Brown           for (j=0; j<ym; j++) {
1624061b8bfSJed Brown             for (i=0; i<xm; i++) {
1634061b8bfSJed Brown               PetscInt Iloc = i+xm*(j+ym*k);
1642eaa9ef4SLisandro Dalcin               array2[Iloc*3+0] = array[Iloc*cdim + 0];
1652eaa9ef4SLisandro Dalcin               array2[Iloc*3+1] = cdim > 1 ? array[Iloc*cdim + 1] : 0.0;
1662eaa9ef4SLisandro Dalcin               array2[Iloc*3+2] = cdim > 2 ? array[Iloc*cdim + 2] : 0.0;
1674061b8bfSJed Brown             }
1684061b8bfSJed Brown           }
1694061b8bfSJed Brown         }
1704061b8bfSJed Brown       } else if (r == rank) {
1712eaa9ef4SLisandro Dalcin         ierr = MPI_Send((void*)coords,nnodes*cdim,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr);
1724061b8bfSJed Brown       }
1734061b8bfSJed Brown       ierr = VecRestoreArrayRead(Coords,&coords);CHKERRQ(ierr);
1744061b8bfSJed Brown     } else {       /* Fabricate some coordinates using grid index */
1754061b8bfSJed Brown       for (k=0; k<zm; k++) {
1764061b8bfSJed Brown         for (j=0; j<ym; j++) {
1774061b8bfSJed Brown           for (i=0; i<xm; i++) {
1784061b8bfSJed Brown             PetscInt Iloc = i+xm*(j+ym*k);
1794061b8bfSJed Brown             array2[Iloc*3+0] = xs+i;
1804061b8bfSJed Brown             array2[Iloc*3+1] = ys+j;
1814061b8bfSJed Brown             array2[Iloc*3+2] = zs+k;
1824061b8bfSJed Brown           }
1834061b8bfSJed Brown         }
1844061b8bfSJed Brown       }
1854061b8bfSJed Brown     }
186b263465dSJed Brown     ierr = PetscViewerVTKFWrite(viewer,fp,array2,nnodes*3,PETSC_SCALAR);CHKERRQ(ierr);
1874061b8bfSJed Brown 
1884061b8bfSJed Brown     /* Write each of the objects queued up for this file */
1894061b8bfSJed Brown     for (link=vtk->link; link; link=link->next) {
1904061b8bfSJed Brown       Vec               X = (Vec)link->vec;
1914061b8bfSJed Brown       const PetscScalar *x;
1924061b8bfSJed Brown 
1934061b8bfSJed Brown       ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
1944061b8bfSJed Brown       if (!rank) {
1954061b8bfSJed Brown         if (r) {
1966622f924SJed Brown           PetscMPIInt nn;
1974061b8bfSJed Brown           ierr = MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr);
1984061b8bfSJed Brown           ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr);
1994061b8bfSJed Brown           if (nn != nnodes*bs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %D",r);
2004061b8bfSJed Brown         } else {
2014061b8bfSJed Brown           ierr = PetscMemcpy(array,x,nnodes*bs*sizeof(PetscScalar));CHKERRQ(ierr);
2024061b8bfSJed Brown         }
2034061b8bfSJed Brown         for (f=0; f<bs; f++) {
2044061b8bfSJed Brown           /* Extract and transpose the f'th field */
2054061b8bfSJed Brown           for (k=0; k<zm; k++) {
2064061b8bfSJed Brown             for (j=0; j<ym; j++) {
2074061b8bfSJed Brown               for (i=0; i<xm; i++) {
2084061b8bfSJed Brown                 PetscInt Iloc = i+xm*(j+ym*k);
2094061b8bfSJed Brown                 array2[Iloc] = array[Iloc*bs + f];
2104061b8bfSJed Brown               }
2114061b8bfSJed Brown             }
2124061b8bfSJed Brown           }
213b263465dSJed Brown           ierr = PetscViewerVTKFWrite(viewer,fp,array2,nnodes,PETSC_SCALAR);CHKERRQ(ierr);
2144061b8bfSJed Brown         }
2154061b8bfSJed Brown       } else if (r == rank) {
2164061b8bfSJed Brown         ierr = MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr);
2174061b8bfSJed Brown       }
2184061b8bfSJed Brown       ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
2194061b8bfSJed Brown     }
2204061b8bfSJed Brown   }
2214061b8bfSJed Brown   ierr = PetscFree2(array,array2);CHKERRQ(ierr);
2224061b8bfSJed Brown   ierr = PetscFree(grloc);CHKERRQ(ierr);
2234061b8bfSJed Brown 
2244061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr);
2254061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
2264061b8bfSJed Brown   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
2274061b8bfSJed Brown   PetscFunctionReturn(0);
2284061b8bfSJed Brown }
2294061b8bfSJed Brown 
2304061b8bfSJed Brown 
2314061b8bfSJed Brown #undef __FUNCT__
232a13bc4e3SShao-Ching Huang #define __FUNCT__ "DMDAVTKWriteAll_VTR"
233a13bc4e3SShao-Ching Huang static PetscErrorCode DMDAVTKWriteAll_VTR(DM da,PetscViewer viewer)
234a13bc4e3SShao-Ching Huang {
235a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_REAL_SINGLE)
236a13bc4e3SShao-Ching Huang   const char precision[] = "Float32";
237a13bc4e3SShao-Ching Huang #elif defined(PETSC_USE_REAL_DOUBLE)
238a13bc4e3SShao-Ching Huang   const char precision[] = "Float64";
239a13bc4e3SShao-Ching Huang #else
240a13bc4e3SShao-Ching Huang   const char precision[] = "UnknownPrecision";
241a13bc4e3SShao-Ching Huang #endif
242a13bc4e3SShao-Ching Huang   MPI_Comm                 comm;
243a13bc4e3SShao-Ching Huang   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
244a13bc4e3SShao-Ching Huang   PetscViewerVTKObjectLink link;
245a13bc4e3SShao-Ching Huang   FILE                     *fp;
246a13bc4e3SShao-Ching Huang   PetscMPIInt              rank,size,tag;
247a13bc4e3SShao-Ching Huang   DMDALocalInfo            info;
248a13bc4e3SShao-Ching Huang   PetscInt                 dim,mx,my,mz,bs,boffset,maxnnodes,i,j,k,f,r;
249a13bc4e3SShao-Ching Huang   PetscInt                 rloc[6],(*grloc)[6] = NULL;
250a13bc4e3SShao-Ching Huang   PetscScalar              *array,*array2;
251a13bc4e3SShao-Ching Huang   PetscReal                gmin[3],gmax[3];
252a13bc4e3SShao-Ching Huang   PetscErrorCode           ierr;
253a13bc4e3SShao-Ching Huang 
254a13bc4e3SShao-Ching Huang   PetscFunctionBegin;
255a13bc4e3SShao-Ching Huang   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
256a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_COMPLEX)
257a13bc4e3SShao-Ching Huang   SETERRQ(comm,PETSC_ERR_SUP,"Complex values not supported");
258a13bc4e3SShao-Ching Huang #endif
259a13bc4e3SShao-Ching Huang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
260a13bc4e3SShao-Ching Huang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
261a13bc4e3SShao-Ching Huang   ierr = DMDAGetInfo(da,&dim, &mx,&my,&mz, 0,0,0, &bs,0,0,0,0,0);CHKERRQ(ierr);
262a13bc4e3SShao-Ching Huang   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
263a13bc4e3SShao-Ching Huang   ierr = DMDAGetBoundingBox(da,gmin,gmax);CHKERRQ(ierr);
264a13bc4e3SShao-Ching Huang   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
265a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
266a13bc4e3SShao-Ching Huang #if defined(PETSC_WORDS_BIGENDIAN)
267a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"RectilinearGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr);
268a13bc4e3SShao-Ching Huang #else
269a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"RectilinearGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr);
270a13bc4e3SShao-Ching Huang #endif
271a13bc4e3SShao-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);
272a13bc4e3SShao-Ching Huang 
273785e854fSJed Brown   if (!rank) {ierr = PetscMalloc1(size*6,&grloc);CHKERRQ(ierr);}
274a13bc4e3SShao-Ching Huang   rloc[0] = info.xs;
275a13bc4e3SShao-Ching Huang   rloc[1] = info.xm;
276a13bc4e3SShao-Ching Huang   rloc[2] = info.ys;
277a13bc4e3SShao-Ching Huang   rloc[3] = info.ym;
278a13bc4e3SShao-Ching Huang   rloc[4] = info.zs;
279a13bc4e3SShao-Ching Huang   rloc[5] = info.zm;
280a13bc4e3SShao-Ching Huang   ierr    = MPI_Gather(rloc,6,MPIU_INT,&grloc[0][0],6,MPIU_INT,0,comm);CHKERRQ(ierr);
281a13bc4e3SShao-Ching Huang 
282a13bc4e3SShao-Ching Huang   /* Write XML header */
283a13bc4e3SShao-Ching Huang   maxnnodes = 0;                /* Used for the temporary array size on rank 0 */
284a13bc4e3SShao-Ching Huang   boffset   = 0;                /* Offset into binary file */
285a13bc4e3SShao-Ching Huang   for (r=0; r<size; r++) {
286a13bc4e3SShao-Ching Huang     PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
287a13bc4e3SShao-Ching Huang     if (!rank) {
288a13bc4e3SShao-Ching Huang       xs     = grloc[r][0];
289a13bc4e3SShao-Ching Huang       xm     = grloc[r][1];
290a13bc4e3SShao-Ching Huang       ys     = grloc[r][2];
291a13bc4e3SShao-Ching Huang       ym     = grloc[r][3];
292a13bc4e3SShao-Ching Huang       zs     = grloc[r][4];
293a13bc4e3SShao-Ching Huang       zm     = grloc[r][5];
294a13bc4e3SShao-Ching Huang       nnodes = xm*ym*zm;
295a13bc4e3SShao-Ching Huang     }
296a13bc4e3SShao-Ching Huang     maxnnodes = PetscMax(maxnnodes,nnodes);
297a13bc4e3SShao-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);
298a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"      <Coordinates>\n");CHKERRQ(ierr);
299a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Xcoord\"  format=\"appended\"  offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
300a13bc4e3SShao-Ching Huang     boffset += xm*sizeof(PetscScalar) + sizeof(int);
301a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Ycoord\"  format=\"appended\"  offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
302a13bc4e3SShao-Ching Huang     boffset += ym*sizeof(PetscScalar) + sizeof(int);
303a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Zcoord\"  format=\"appended\"  offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
304a13bc4e3SShao-Ching Huang     boffset += zm*sizeof(PetscScalar) + sizeof(int);
305a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"      </Coordinates>\n");CHKERRQ(ierr);
306a13bc4e3SShao-Ching Huang     ierr = PetscFPrintf(comm,fp,"      <PointData Scalars=\"ScalarPointData\">\n");CHKERRQ(ierr);
307a13bc4e3SShao-Ching Huang     for (link=vtk->link; link; link=link->next) {
308a13bc4e3SShao-Ching Huang       Vec        X        = (Vec)link->vec;
309a13bc4e3SShao-Ching Huang       const char *vecname = "";
310a13bc4e3SShao-Ching Huang       if (((PetscObject)X)->name || link != vtk->link) { /* If the object is already named, use it. If it is past the first link, name it to disambiguate. */
311a13bc4e3SShao-Ching Huang         ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
312a13bc4e3SShao-Ching Huang       }
313a13bc4e3SShao-Ching Huang       for (i=0; i<bs; i++) {
314a13bc4e3SShao-Ching Huang         char       buf[256];
315a13bc4e3SShao-Ching Huang         const char *fieldname;
316a13bc4e3SShao-Ching Huang         ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
317a13bc4e3SShao-Ching Huang         if (!fieldname) {
318a13bc4e3SShao-Ching Huang           ierr      = PetscSNPrintf(buf,sizeof(buf),"Unnamed%D",i);CHKERRQ(ierr);
319a13bc4e3SShao-Ching Huang           fieldname = buf;
320a13bc4e3SShao-Ching Huang         }
321a13bc4e3SShao-Ching Huang         ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
322a13bc4e3SShao-Ching Huang         boffset += nnodes*sizeof(PetscScalar) + sizeof(int);
323a13bc4e3SShao-Ching Huang       }
324a13bc4e3SShao-Ching Huang     }
325a13bc4e3SShao-Ching Huang     ierr = PetscFPrintf(comm,fp,"      </PointData>\n");CHKERRQ(ierr);
326a13bc4e3SShao-Ching Huang     ierr = PetscFPrintf(comm,fp,"    </Piece>\n");CHKERRQ(ierr);
327a13bc4e3SShao-Ching Huang   }
328a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"  </RectilinearGrid>\n");CHKERRQ(ierr);
329a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
330a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
331a13bc4e3SShao-Ching Huang 
332a13bc4e3SShao-Ching Huang   /* Now write the arrays. */
333a13bc4e3SShao-Ching Huang   tag  = ((PetscObject)viewer)->tag;
334dcca6d9dSJed Brown   ierr = PetscMalloc2(PetscMax(maxnnodes*bs,info.xm+info.ym+info.zm),&array,PetscMax(maxnnodes*3,info.xm+info.ym+info.zm),&array2);CHKERRQ(ierr);
335a13bc4e3SShao-Ching Huang   for (r=0; r<size; r++) {
336a13bc4e3SShao-Ching Huang     MPI_Status status;
337a13bc4e3SShao-Ching Huang     PetscInt   xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
338a13bc4e3SShao-Ching Huang     if (!rank) {
339a13bc4e3SShao-Ching Huang       xs     = grloc[r][0];
340a13bc4e3SShao-Ching Huang       xm     = grloc[r][1];
341a13bc4e3SShao-Ching Huang       ys     = grloc[r][2];
342a13bc4e3SShao-Ching Huang       ym     = grloc[r][3];
343a13bc4e3SShao-Ching Huang       zs     = grloc[r][4];
344a13bc4e3SShao-Ching Huang       zm     = grloc[r][5];
345a13bc4e3SShao-Ching Huang       nnodes = xm*ym*zm;
346a13bc4e3SShao-Ching Huang     } else if (r == rank) {
347a13bc4e3SShao-Ching Huang       nnodes = info.xm*info.ym*info.zm;
348a13bc4e3SShao-Ching Huang     }
349a13bc4e3SShao-Ching Huang     {                           /* Write the coordinates */
350a13bc4e3SShao-Ching Huang       Vec Coords;
351a13bc4e3SShao-Ching Huang       ierr = DMGetCoordinates(da,&Coords);CHKERRQ(ierr);
352a13bc4e3SShao-Ching Huang       if (Coords) {
353a13bc4e3SShao-Ching Huang         const PetscScalar *coords;
354a13bc4e3SShao-Ching Huang         ierr = VecGetArrayRead(Coords,&coords);CHKERRQ(ierr);
355a13bc4e3SShao-Ching Huang         if (!rank) {
356a13bc4e3SShao-Ching Huang           if (r) {
357a13bc4e3SShao-Ching Huang             PetscMPIInt nn;
358a13bc4e3SShao-Ching Huang             ierr = MPI_Recv(array,xm+ym+zm,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr);
359a13bc4e3SShao-Ching Huang             ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr);
360a13bc4e3SShao-Ching Huang             if (nn != xm+ym+zm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
361a13bc4e3SShao-Ching Huang           } else {
362a13bc4e3SShao-Ching Huang             /* x coordinates */
363a13bc4e3SShao-Ching Huang             for (j=0, k=0, i=0; i<xm; i++) {
364a13bc4e3SShao-Ching Huang               PetscInt Iloc = i+xm*(j+ym*k);
365a13bc4e3SShao-Ching Huang               array[i] = coords[Iloc*dim + 0];
366a13bc4e3SShao-Ching Huang             }
367a13bc4e3SShao-Ching Huang             /* y coordinates */
368a13bc4e3SShao-Ching Huang             for (i=0, k=0, j=0; j<ym; j++) {
369a13bc4e3SShao-Ching Huang               PetscInt Iloc = i+xm*(j+ym*k);
370a13bc4e3SShao-Ching Huang               array[j+xm] = dim > 1 ? coords[Iloc*dim + 1] : 0;
371a13bc4e3SShao-Ching Huang             }
372a13bc4e3SShao-Ching Huang             /* z coordinates */
373a13bc4e3SShao-Ching Huang             for (i=0, j=0, k=0; k<zm; k++) {
374a13bc4e3SShao-Ching Huang               PetscInt Iloc = i+xm*(j+ym*k);
375a13bc4e3SShao-Ching Huang               array[k+xm+ym] = dim > 2 ? coords[Iloc*dim + 2] : 0;
376a13bc4e3SShao-Ching Huang             }
377a13bc4e3SShao-Ching Huang           }
378a13bc4e3SShao-Ching Huang         } else if (r == rank) {
379a13bc4e3SShao-Ching Huang           xm = info.xm;
380a13bc4e3SShao-Ching Huang           ym = info.ym;
381a13bc4e3SShao-Ching Huang           zm = info.zm;
382a13bc4e3SShao-Ching Huang           for (j=0, k=0, i=0; i<xm; i++) {
383a13bc4e3SShao-Ching Huang             PetscInt Iloc = i+xm*(j+ym*k);
384a13bc4e3SShao-Ching Huang             array2[i] = coords[Iloc*dim + 0];
385a13bc4e3SShao-Ching Huang           }
386a13bc4e3SShao-Ching Huang           for (i=0, k=0, j=0; j<ym; j++) {
387a13bc4e3SShao-Ching Huang             PetscInt Iloc = i+xm*(j+ym*k);
388a13bc4e3SShao-Ching Huang             array2[j+xm] = dim > 1 ? coords[Iloc*dim + 1] : 0;
389a13bc4e3SShao-Ching Huang           }
390a13bc4e3SShao-Ching Huang           for (i=0, j=0, k=0; k<zm; k++) {
391a13bc4e3SShao-Ching Huang             PetscInt Iloc = i+xm*(j+ym*k);
392a13bc4e3SShao-Ching Huang             array2[k+xm+ym] = dim > 2 ? coords[Iloc*dim + 2] : 0;
393a13bc4e3SShao-Ching Huang           }
394a13bc4e3SShao-Ching Huang           ierr = MPI_Send((void*)array2,xm+ym+zm,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr);
395a13bc4e3SShao-Ching Huang         }
396a13bc4e3SShao-Ching Huang         ierr = VecRestoreArrayRead(Coords,&coords);CHKERRQ(ierr);
397a13bc4e3SShao-Ching Huang       } else {       /* Fabricate some coordinates using grid index */
398*3d3eaba7SBarry Smith         for (i=0; i<xm; i++) {
399a13bc4e3SShao-Ching Huang           array[i] = xs+i;
400a13bc4e3SShao-Ching Huang         }
401*3d3eaba7SBarry Smith         for (j=0; j<ym; j++) {
402a13bc4e3SShao-Ching Huang           array[j+xm] = ys+j;
403a13bc4e3SShao-Ching Huang         }
404*3d3eaba7SBarry Smith         for (k=0; k<zm; k++) {
405a13bc4e3SShao-Ching Huang           array[k+xm+ym] = zs+k;
406a13bc4e3SShao-Ching Huang         }
407a13bc4e3SShao-Ching Huang       }
408a13bc4e3SShao-Ching Huang       if (!rank) {
409a13bc4e3SShao-Ching Huang         ierr = PetscViewerVTKFWrite(viewer,fp,&(array[0]),xm,PETSC_SCALAR);CHKERRQ(ierr);
410a13bc4e3SShao-Ching Huang         ierr = PetscViewerVTKFWrite(viewer,fp,&(array[xm]),ym,PETSC_SCALAR);CHKERRQ(ierr);
411a13bc4e3SShao-Ching Huang         ierr = PetscViewerVTKFWrite(viewer,fp,&(array[xm+ym]),zm,PETSC_SCALAR);CHKERRQ(ierr);
412a13bc4e3SShao-Ching Huang       }
413a13bc4e3SShao-Ching Huang     }
414a13bc4e3SShao-Ching Huang 
415a13bc4e3SShao-Ching Huang     /* Write each of the objects queued up for this file */
416a13bc4e3SShao-Ching Huang     for (link=vtk->link; link; link=link->next) {
417a13bc4e3SShao-Ching Huang       Vec               X = (Vec)link->vec;
418a13bc4e3SShao-Ching Huang       const PetscScalar *x;
419a13bc4e3SShao-Ching Huang 
420a13bc4e3SShao-Ching Huang       ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
421a13bc4e3SShao-Ching Huang       if (!rank) {
422a13bc4e3SShao-Ching Huang         if (r) {
423a13bc4e3SShao-Ching Huang           PetscMPIInt nn;
424a13bc4e3SShao-Ching Huang           ierr = MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr);
425a13bc4e3SShao-Ching Huang           ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr);
426a13bc4e3SShao-Ching Huang           if (nn != nnodes*bs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %D",r);
427a13bc4e3SShao-Ching Huang         } else {
428a13bc4e3SShao-Ching Huang           ierr = PetscMemcpy(array,x,nnodes*bs*sizeof(PetscScalar));CHKERRQ(ierr);
429a13bc4e3SShao-Ching Huang         }
430a13bc4e3SShao-Ching Huang         for (f=0; f<bs; f++) {
431a13bc4e3SShao-Ching Huang           /* Extract and transpose the f'th field */
432a13bc4e3SShao-Ching Huang           for (k=0; k<zm; k++) {
433a13bc4e3SShao-Ching Huang             for (j=0; j<ym; j++) {
434a13bc4e3SShao-Ching Huang               for (i=0; i<xm; i++) {
435a13bc4e3SShao-Ching Huang                 PetscInt Iloc = i+xm*(j+ym*k);
436a13bc4e3SShao-Ching Huang                 array2[Iloc] = array[Iloc*bs + f];
437a13bc4e3SShao-Ching Huang               }
438a13bc4e3SShao-Ching Huang             }
439a13bc4e3SShao-Ching Huang           }
440a13bc4e3SShao-Ching Huang           ierr = PetscViewerVTKFWrite(viewer,fp,array2,nnodes,PETSC_SCALAR);CHKERRQ(ierr);
441a13bc4e3SShao-Ching Huang         }
442a13bc4e3SShao-Ching Huang       } else if (r == rank) {
443a13bc4e3SShao-Ching Huang         ierr = MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr);
444a13bc4e3SShao-Ching Huang       }
445a13bc4e3SShao-Ching Huang       ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
446a13bc4e3SShao-Ching Huang     }
447a13bc4e3SShao-Ching Huang   }
448a13bc4e3SShao-Ching Huang   ierr = PetscFree2(array,array2);CHKERRQ(ierr);
449a13bc4e3SShao-Ching Huang   ierr = PetscFree(grloc);CHKERRQ(ierr);
450a13bc4e3SShao-Ching Huang 
451a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr);
452a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
453a13bc4e3SShao-Ching Huang   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
454a13bc4e3SShao-Ching Huang   PetscFunctionReturn(0);
455a13bc4e3SShao-Ching Huang }
456a13bc4e3SShao-Ching Huang 
457a13bc4e3SShao-Ching Huang #undef __FUNCT__
4584061b8bfSJed Brown #define __FUNCT__ "DMDAVTKWriteAll"
4594061b8bfSJed Brown /*@C
4604061b8bfSJed Brown    DMDAVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
4614061b8bfSJed Brown 
4624061b8bfSJed Brown    Collective
4634061b8bfSJed Brown 
4644061b8bfSJed Brown    Input Arguments:
4654061b8bfSJed Brown    odm - DM specifying the grid layout, passed as a PetscObject
4664061b8bfSJed Brown    viewer - viewer of type VTK
4674061b8bfSJed Brown 
4684061b8bfSJed Brown    Level: developer
4694061b8bfSJed Brown 
4704061b8bfSJed Brown    Note:
4714061b8bfSJed Brown    This function is a callback used by the VTK viewer to actually write the file.
4724061b8bfSJed 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.
4734061b8bfSJed Brown    Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
4744061b8bfSJed Brown 
4754061b8bfSJed Brown .seealso: PETSCVIEWERVTK
4764061b8bfSJed Brown @*/
4774061b8bfSJed Brown PetscErrorCode DMDAVTKWriteAll(PetscObject odm,PetscViewer viewer)
4784061b8bfSJed Brown {
4794061b8bfSJed Brown   DM             dm = (DM)odm;
4804061b8bfSJed Brown   PetscBool      isvtk;
4814061b8bfSJed Brown   PetscErrorCode ierr;
4824061b8bfSJed Brown 
4834061b8bfSJed Brown   PetscFunctionBegin;
4844061b8bfSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
4854061b8bfSJed Brown   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
486251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr);
487ce94432eSBarry Smith   if (!isvtk) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_INCOMP,"Cannot use viewer type %s",((PetscObject)viewer)->type_name);
4884061b8bfSJed Brown   switch (viewer->format) {
4894061b8bfSJed Brown   case PETSC_VIEWER_VTK_VTS:
4904061b8bfSJed Brown     ierr = DMDAVTKWriteAll_VTS(dm,viewer);CHKERRQ(ierr);
4914061b8bfSJed Brown     break;
492a13bc4e3SShao-Ching Huang   case PETSC_VIEWER_VTK_VTR:
493a13bc4e3SShao-Ching Huang     ierr = DMDAVTKWriteAll_VTR(dm,viewer);CHKERRQ(ierr);
494a13bc4e3SShao-Ching Huang     break;
495ce94432eSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No support for format '%s'",PetscViewerFormats[viewer->format]);
4964061b8bfSJed Brown   }
4974061b8bfSJed Brown   PetscFunctionReturn(0);
4984061b8bfSJed Brown }
499