xref: /petsc/src/dm/impls/da/grvtk.c (revision 6daa58f615fbca3ea8728881d7825373f9836f4e)
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 
8fb5bd1c2SPatrick Sanan /* Helper function which determines if any DMDA fields are named.  This is used
9fb5bd1c2SPatrick Sanan    as a proxy for the user's intention to use DMDA fields as distinct
10fb5bd1c2SPatrick Sanan    scalar-valued fields as opposed to a single vector-valued field */
11fb5bd1c2SPatrick Sanan static PetscErrorCode DMDAGetFieldsNamed(DM da,PetscBool *fieldsnamed)
12fb5bd1c2SPatrick Sanan {
13fb5bd1c2SPatrick Sanan   PetscErrorCode ierr;
14fb5bd1c2SPatrick Sanan   PetscInt       f,bs;
15fb5bd1c2SPatrick Sanan 
16fb5bd1c2SPatrick Sanan   PetscFunctionBegin;
17fb5bd1c2SPatrick Sanan   *fieldsnamed = PETSC_FALSE;
18277f51e8SBarry Smith   ierr = DMDAGetDof(da,&bs);CHKERRQ(ierr);
19fb5bd1c2SPatrick Sanan   for (f=0; f<bs; ++f) {
20fb5bd1c2SPatrick Sanan     const char * fieldname;
21fb5bd1c2SPatrick Sanan     ierr = DMDAGetFieldName(da,f,&fieldname);CHKERRQ(ierr);
22fb5bd1c2SPatrick Sanan     if (fieldname) {
23fb5bd1c2SPatrick Sanan       *fieldsnamed = PETSC_TRUE;
24fb5bd1c2SPatrick Sanan       break;
25fb5bd1c2SPatrick Sanan     }
26fb5bd1c2SPatrick Sanan   }
27fb5bd1c2SPatrick Sanan   PetscFunctionReturn(0);
28fb5bd1c2SPatrick Sanan }
29fb5bd1c2SPatrick Sanan 
304061b8bfSJed Brown static PetscErrorCode DMDAVTKWriteAll_VTS(DM da,PetscViewer viewer)
314061b8bfSJed Brown {
3230815ce0SLisandro Dalcin   const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian";
334061b8bfSJed Brown #if defined(PETSC_USE_REAL_SINGLE)
344061b8bfSJed Brown   const char precision[] = "Float32";
354061b8bfSJed Brown #elif defined(PETSC_USE_REAL_DOUBLE)
364061b8bfSJed Brown   const char precision[] = "Float64";
374061b8bfSJed Brown #else
384061b8bfSJed Brown   const char precision[] = "UnknownPrecision";
394061b8bfSJed Brown #endif
40ce94432eSBarry Smith   MPI_Comm                 comm;
412eaa9ef4SLisandro Dalcin   Vec                      Coords;
424061b8bfSJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
434061b8bfSJed Brown   PetscViewerVTKObjectLink link;
444061b8bfSJed Brown   FILE                     *fp;
454061b8bfSJed Brown   PetscMPIInt              rank,size,tag;
464061b8bfSJed Brown   DMDALocalInfo            info;
47ea2d7708SPatrick Sanan   PetscInt                 dim,mx,my,mz,cdim,bs,boffset,maxnnodes,maxbs,i,j,k,r;
480298fd71SBarry Smith   PetscInt                 rloc[6],(*grloc)[6] = NULL;
494061b8bfSJed Brown   PetscScalar              *array,*array2;
504061b8bfSJed Brown   PetscErrorCode           ierr;
514061b8bfSJed Brown 
524061b8bfSJed Brown   PetscFunctionBegin;
53094921d9SJed Brown   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
54*6daa58f6SJose E. Roman   if (PetscDefined(USE_COMPLEX)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Complex values not supported");
55ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
56ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
57ea78f98cSLisandro Dalcin   ierr = DMDAGetInfo(da,&dim,&mx,&my,&mz,NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
584061b8bfSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
592eaa9ef4SLisandro Dalcin   ierr = DMGetCoordinates(da,&Coords);CHKERRQ(ierr);
602eaa9ef4SLisandro Dalcin   if (Coords) {
612eaa9ef4SLisandro Dalcin     PetscInt csize;
622eaa9ef4SLisandro Dalcin     ierr = VecGetSize(Coords,&csize);CHKERRQ(ierr);
632eaa9ef4SLisandro Dalcin     if (csize % (mx*my*mz)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Coordinate vector size mismatch");
642eaa9ef4SLisandro Dalcin     cdim = csize/(mx*my*mz);
652eaa9ef4SLisandro Dalcin     if (cdim < dim || cdim > 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Coordinate vector size mismatch");
662eaa9ef4SLisandro Dalcin   } else {
672eaa9ef4SLisandro Dalcin     cdim = dim;
682eaa9ef4SLisandro Dalcin   }
694061b8bfSJed Brown 
704061b8bfSJed Brown   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
714061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
7230815ce0SLisandro Dalcin   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"%s\">\n",byte_order);CHKERRQ(ierr);
734061b8bfSJed 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);
744061b8bfSJed Brown 
75785e854fSJed Brown   if (!rank) {ierr = PetscMalloc1(size*6,&grloc);CHKERRQ(ierr);}
764061b8bfSJed Brown   rloc[0] = info.xs;
774061b8bfSJed Brown   rloc[1] = info.xm;
784061b8bfSJed Brown   rloc[2] = info.ys;
794061b8bfSJed Brown   rloc[3] = info.ym;
804061b8bfSJed Brown   rloc[4] = info.zs;
814061b8bfSJed Brown   rloc[5] = info.zm;
8255b25c41SPierre Jolivet   ierr    = MPI_Gather(rloc,6,MPIU_INT,&grloc[0][0],6,MPIU_INT,0,comm);CHKERRMPI(ierr);
834061b8bfSJed Brown 
844061b8bfSJed Brown   /* Write XML header */
854061b8bfSJed Brown   maxnnodes = 0;                /* Used for the temporary array size on rank 0 */
86ea2d7708SPatrick Sanan   maxbs     = 0;                /* Used for the temporary array size on rank 0 */
874061b8bfSJed Brown   boffset   = 0;                /* Offset into binary file */
884061b8bfSJed Brown   for (r=0; r<size; r++) {
894061b8bfSJed Brown     PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
904061b8bfSJed Brown     if (!rank) {
914061b8bfSJed Brown       xs     = grloc[r][0];
924061b8bfSJed Brown       xm     = grloc[r][1];
934061b8bfSJed Brown       ys     = grloc[r][2];
944061b8bfSJed Brown       ym     = grloc[r][3];
954061b8bfSJed Brown       zs     = grloc[r][4];
964061b8bfSJed Brown       zm     = grloc[r][5];
974061b8bfSJed Brown       nnodes = xm*ym*zm;
984061b8bfSJed Brown     }
994061b8bfSJed Brown     maxnnodes = PetscMax(maxnnodes,nnodes);
1004061b8bfSJed 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);
1014061b8bfSJed Brown     ierr     = PetscFPrintf(comm,fp,"      <Points>\n");CHKERRQ(ierr);
1024061b8bfSJed Brown     ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
1034061b8bfSJed Brown     boffset += 3*nnodes*sizeof(PetscScalar) + sizeof(int);
1044061b8bfSJed Brown     ierr     = PetscFPrintf(comm,fp,"      </Points>\n");CHKERRQ(ierr);
1054061b8bfSJed Brown 
1064061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"      <PointData Scalars=\"ScalarPointData\">\n");CHKERRQ(ierr);
1074061b8bfSJed Brown     for (link=vtk->link; link; link=link->next) {
1084061b8bfSJed Brown       Vec        X = (Vec)link->vec;
109fb5bd1c2SPatrick Sanan       PetscInt   bs,f;
110ea2d7708SPatrick Sanan       DM         daCurr;
111fb5bd1c2SPatrick Sanan       PetscBool  fieldsnamed;
112fb5bd1c2SPatrick Sanan       const char *vecname = "Unnamed";
113ea2d7708SPatrick Sanan 
114ea2d7708SPatrick Sanan       ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr);
115ea78f98cSLisandro Dalcin       ierr = DMDAGetInfo(daCurr,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
116ea2d7708SPatrick Sanan       maxbs = PetscMax(maxbs,bs);
117ea2d7708SPatrick Sanan 
118ea2d7708SPatrick Sanan       if (((PetscObject)X)->name || link != vtk->link) {
1194061b8bfSJed Brown         ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
1204061b8bfSJed Brown       }
121fb5bd1c2SPatrick Sanan 
122fb5bd1c2SPatrick Sanan       /* If any fields are named, add scalar fields. Otherwise, add a vector field */
123fb5bd1c2SPatrick Sanan       ierr = DMDAGetFieldsNamed(daCurr,&fieldsnamed);CHKERRQ(ierr);
124fb5bd1c2SPatrick Sanan       if (fieldsnamed) {
125fb5bd1c2SPatrick Sanan         for (f=0; f<bs; f++) {
126fb5bd1c2SPatrick Sanan           char       buf[256];
127fb5bd1c2SPatrick Sanan           const char *fieldname;
128fb5bd1c2SPatrick Sanan           ierr = DMDAGetFieldName(daCurr,f,&fieldname);CHKERRQ(ierr);
129fb5bd1c2SPatrick Sanan           if (!fieldname) {
130fb5bd1c2SPatrick Sanan             ierr      = PetscSNPrintf(buf,sizeof(buf),"%D",f);CHKERRQ(ierr);
131fb5bd1c2SPatrick Sanan             fieldname = buf;
132fb5bd1c2SPatrick Sanan           }
133fb5bd1c2SPatrick Sanan           ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
134fb5bd1c2SPatrick Sanan           boffset += nnodes*sizeof(PetscScalar) + sizeof(int);
135fb5bd1c2SPatrick Sanan         }
136fb5bd1c2SPatrick Sanan       } else {
137ea2d7708SPatrick Sanan         ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%D\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,bs,boffset);CHKERRQ(ierr);
138ea2d7708SPatrick Sanan         boffset += bs*nnodes*sizeof(PetscScalar) + sizeof(int);
1394061b8bfSJed Brown       }
140fb5bd1c2SPatrick Sanan     }
1414061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"      </PointData>\n");CHKERRQ(ierr);
1424061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"    </Piece>\n");CHKERRQ(ierr);
1434061b8bfSJed Brown   }
1444061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"  </StructuredGrid>\n");CHKERRQ(ierr);
1454061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
1464061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
1474061b8bfSJed Brown 
1484061b8bfSJed Brown   /* Now write the arrays. */
1494061b8bfSJed Brown   tag  = ((PetscObject)viewer)->tag;
150fb5bd1c2SPatrick Sanan   ierr = PetscMalloc2(maxnnodes*PetscMax(3,maxbs),&array,maxnnodes*PetscMax(3,maxbs),&array2);CHKERRQ(ierr);
1514061b8bfSJed Brown   for (r=0; r<size; r++) {
1524061b8bfSJed Brown     MPI_Status status;
1534061b8bfSJed Brown     PetscInt   xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
1544061b8bfSJed Brown     if (!rank) {
1554061b8bfSJed Brown       xs     = grloc[r][0];
1564061b8bfSJed Brown       xm     = grloc[r][1];
1574061b8bfSJed Brown       ys     = grloc[r][2];
1584061b8bfSJed Brown       ym     = grloc[r][3];
1594061b8bfSJed Brown       zs     = grloc[r][4];
1604061b8bfSJed Brown       zm     = grloc[r][5];
1614061b8bfSJed Brown       nnodes = xm*ym*zm;
1624061b8bfSJed Brown     } else if (r == rank) {
1634061b8bfSJed Brown       nnodes = info.xm*info.ym*info.zm;
1644061b8bfSJed Brown     }
1654061b8bfSJed Brown 
1662eaa9ef4SLisandro Dalcin     /* Write the coordinates */
1674061b8bfSJed Brown     if (Coords) {
1684061b8bfSJed Brown       const PetscScalar *coords;
1694061b8bfSJed Brown       ierr = VecGetArrayRead(Coords,&coords);CHKERRQ(ierr);
1704061b8bfSJed Brown       if (!rank) {
1714061b8bfSJed Brown         if (r) {
1726622f924SJed Brown           PetscMPIInt nn;
173ffc4695bSBarry Smith           ierr = MPI_Recv(array,nnodes*cdim,MPIU_SCALAR,r,tag,comm,&status);CHKERRMPI(ierr);
174ffc4695bSBarry Smith           ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRMPI(ierr);
1752eaa9ef4SLisandro Dalcin           if (nn != nnodes*cdim) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
1764061b8bfSJed Brown         } else {
177580bdb30SBarry Smith           ierr = PetscArraycpy(array,coords,nnodes*cdim);CHKERRQ(ierr);
1784061b8bfSJed Brown         }
1794061b8bfSJed Brown         /* Transpose coordinates to VTK (C-style) ordering */
1804061b8bfSJed Brown         for (k=0; k<zm; k++) {
1814061b8bfSJed Brown           for (j=0; j<ym; j++) {
1824061b8bfSJed Brown             for (i=0; i<xm; i++) {
1834061b8bfSJed Brown               PetscInt Iloc = i+xm*(j+ym*k);
1842eaa9ef4SLisandro Dalcin               array2[Iloc*3+0] = array[Iloc*cdim + 0];
1852eaa9ef4SLisandro Dalcin               array2[Iloc*3+1] = cdim > 1 ? array[Iloc*cdim + 1] : 0.0;
1862eaa9ef4SLisandro Dalcin               array2[Iloc*3+2] = cdim > 2 ? array[Iloc*cdim + 2] : 0.0;
1874061b8bfSJed Brown             }
1884061b8bfSJed Brown           }
1894061b8bfSJed Brown         }
1904061b8bfSJed Brown       } else if (r == rank) {
191ffc4695bSBarry Smith         ierr = MPI_Send((void*)coords,nnodes*cdim,MPIU_SCALAR,0,tag,comm);CHKERRMPI(ierr);
1924061b8bfSJed Brown       }
1934061b8bfSJed Brown       ierr = VecRestoreArrayRead(Coords,&coords);CHKERRQ(ierr);
1944061b8bfSJed Brown     } else {       /* Fabricate some coordinates using grid index */
1954061b8bfSJed Brown       for (k=0; k<zm; k++) {
1964061b8bfSJed Brown         for (j=0; j<ym; j++) {
1974061b8bfSJed Brown           for (i=0; i<xm; i++) {
1984061b8bfSJed Brown             PetscInt Iloc = i+xm*(j+ym*k);
1994061b8bfSJed Brown             array2[Iloc*3+0] = xs+i;
2004061b8bfSJed Brown             array2[Iloc*3+1] = ys+j;
2014061b8bfSJed Brown             array2[Iloc*3+2] = zs+k;
2024061b8bfSJed Brown           }
2034061b8bfSJed Brown         }
2044061b8bfSJed Brown       }
2054061b8bfSJed Brown     }
20694fbd55eSBarry Smith     ierr = PetscViewerVTKFWrite(viewer,fp,array2,nnodes*3,MPIU_SCALAR);CHKERRQ(ierr);
2074061b8bfSJed Brown 
2084061b8bfSJed Brown     /* Write each of the objects queued up for this file */
2094061b8bfSJed Brown     for (link=vtk->link; link; link=link->next) {
2104061b8bfSJed Brown       Vec               X = (Vec)link->vec;
2114061b8bfSJed Brown       const PetscScalar *x;
212fb5bd1c2SPatrick Sanan       PetscInt          bs,f;
213ea2d7708SPatrick Sanan       DM                daCurr;
214fb5bd1c2SPatrick Sanan       PetscBool         fieldsnamed;
215ea2d7708SPatrick Sanan       ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr);
216ea78f98cSLisandro Dalcin       ierr = DMDAGetInfo(daCurr,NULL,NULL,NULL,NULL, NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
2174061b8bfSJed Brown       ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
2184061b8bfSJed Brown       if (!rank) {
2194061b8bfSJed Brown         if (r) {
2206622f924SJed Brown           PetscMPIInt nn;
221ffc4695bSBarry Smith           ierr = MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status);CHKERRMPI(ierr);
222ffc4695bSBarry Smith           ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRMPI(ierr);
2234061b8bfSJed Brown           if (nn != nnodes*bs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %D",r);
2244061b8bfSJed Brown         } else {
225580bdb30SBarry Smith           ierr = PetscArraycpy(array,x,nnodes*bs);CHKERRQ(ierr);
2264061b8bfSJed Brown         }
227fb5bd1c2SPatrick Sanan 
228fb5bd1c2SPatrick Sanan         /* If any fields are named, add scalar fields. Otherwise, add a vector field */
229fb5bd1c2SPatrick Sanan         ierr = DMDAGetFieldsNamed(daCurr,&fieldsnamed);CHKERRQ(ierr);
230fb5bd1c2SPatrick Sanan         if (fieldsnamed) {
231fb5bd1c2SPatrick Sanan           for (f=0; f<bs; f++) {
232fb5bd1c2SPatrick Sanan             /* Extract and transpose the f'th field */
233fb5bd1c2SPatrick Sanan             for (k=0; k<zm; k++) {
234fb5bd1c2SPatrick Sanan               for (j=0; j<ym; j++) {
235fb5bd1c2SPatrick Sanan                 for (i=0; i<xm; i++) {
236fb5bd1c2SPatrick Sanan                   PetscInt Iloc = i+xm*(j+ym*k);
237fb5bd1c2SPatrick Sanan                   array2[Iloc] = array[Iloc*bs + f];
238fb5bd1c2SPatrick Sanan                 }
239fb5bd1c2SPatrick Sanan               }
240fb5bd1c2SPatrick Sanan             }
241fb5bd1c2SPatrick Sanan             ierr = PetscViewerVTKFWrite(viewer,fp,array2,nnodes,MPIU_SCALAR);CHKERRQ(ierr);
242fb5bd1c2SPatrick Sanan           }
243fb5bd1c2SPatrick Sanan         } else {
244ea2d7708SPatrick Sanan           ierr = PetscViewerVTKFWrite(viewer,fp,array,bs*nnodes,MPIU_SCALAR);CHKERRQ(ierr);
245fb5bd1c2SPatrick Sanan         }
2464061b8bfSJed Brown       } else if (r == rank) {
247ffc4695bSBarry Smith         ierr = MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm);CHKERRMPI(ierr);
2484061b8bfSJed Brown       }
2494061b8bfSJed Brown       ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
2504061b8bfSJed Brown     }
2514061b8bfSJed Brown   }
2524061b8bfSJed Brown   ierr = PetscFree2(array,array2);CHKERRQ(ierr);
2534061b8bfSJed Brown   ierr = PetscFree(grloc);CHKERRQ(ierr);
2544061b8bfSJed Brown 
2554061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr);
2564061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
2574061b8bfSJed Brown   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
2584061b8bfSJed Brown   PetscFunctionReturn(0);
2594061b8bfSJed Brown }
2604061b8bfSJed Brown 
2614061b8bfSJed Brown 
262a13bc4e3SShao-Ching Huang static PetscErrorCode DMDAVTKWriteAll_VTR(DM da,PetscViewer viewer)
263a13bc4e3SShao-Ching Huang {
26430815ce0SLisandro Dalcin   const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian";
265a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_REAL_SINGLE)
266a13bc4e3SShao-Ching Huang   const char precision[] = "Float32";
267a13bc4e3SShao-Ching Huang #elif defined(PETSC_USE_REAL_DOUBLE)
268a13bc4e3SShao-Ching Huang   const char precision[] = "Float64";
269a13bc4e3SShao-Ching Huang #else
270a13bc4e3SShao-Ching Huang   const char precision[] = "UnknownPrecision";
271a13bc4e3SShao-Ching Huang #endif
272a13bc4e3SShao-Ching Huang   MPI_Comm                 comm;
273a13bc4e3SShao-Ching Huang   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
274a13bc4e3SShao-Ching Huang   PetscViewerVTKObjectLink link;
275a13bc4e3SShao-Ching Huang   FILE                     *fp;
276a13bc4e3SShao-Ching Huang   PetscMPIInt              rank,size,tag;
277a13bc4e3SShao-Ching Huang   DMDALocalInfo            info;
278ea2d7708SPatrick Sanan   PetscInt                 dim,mx,my,mz,boffset,maxnnodes,maxbs,i,j,k,r;
279a13bc4e3SShao-Ching Huang   PetscInt                 rloc[6],(*grloc)[6] = NULL;
280a13bc4e3SShao-Ching Huang   PetscScalar              *array,*array2;
281a13bc4e3SShao-Ching Huang   PetscErrorCode           ierr;
282a13bc4e3SShao-Ching Huang 
283a13bc4e3SShao-Ching Huang   PetscFunctionBegin;
284a13bc4e3SShao-Ching Huang   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
285*6daa58f6SJose E. Roman   if (PetscDefined(USE_COMPLEX)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Complex values not supported");
286ffc4695bSBarry Smith   ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr);
287ffc4695bSBarry Smith   ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr);
288ea78f98cSLisandro Dalcin   ierr = DMDAGetInfo(da,&dim,&mx,&my,&mz,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
289a13bc4e3SShao-Ching Huang   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
290a13bc4e3SShao-Ching Huang   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
291a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
29230815ce0SLisandro Dalcin   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"RectilinearGrid\" version=\"0.1\" byte_order=\"%s\">\n",byte_order);CHKERRQ(ierr);
293a13bc4e3SShao-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);
294a13bc4e3SShao-Ching Huang 
295785e854fSJed Brown   if (!rank) {ierr = PetscMalloc1(size*6,&grloc);CHKERRQ(ierr);}
296a13bc4e3SShao-Ching Huang   rloc[0] = info.xs;
297a13bc4e3SShao-Ching Huang   rloc[1] = info.xm;
298a13bc4e3SShao-Ching Huang   rloc[2] = info.ys;
299a13bc4e3SShao-Ching Huang   rloc[3] = info.ym;
300a13bc4e3SShao-Ching Huang   rloc[4] = info.zs;
301a13bc4e3SShao-Ching Huang   rloc[5] = info.zm;
30255b25c41SPierre Jolivet   ierr    = MPI_Gather(rloc,6,MPIU_INT,&grloc[0][0],6,MPIU_INT,0,comm);CHKERRMPI(ierr);
303a13bc4e3SShao-Ching Huang 
304a13bc4e3SShao-Ching Huang   /* Write XML header */
305a13bc4e3SShao-Ching Huang   maxnnodes = 0;                /* Used for the temporary array size on rank 0 */
306ea2d7708SPatrick Sanan   maxbs     = 0;                /* Used for the temporary array size on rank 0 */
307a13bc4e3SShao-Ching Huang   boffset   = 0;                /* Offset into binary file */
308a13bc4e3SShao-Ching Huang   for (r=0; r<size; r++) {
309a13bc4e3SShao-Ching Huang     PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
310a13bc4e3SShao-Ching Huang     if (!rank) {
311a13bc4e3SShao-Ching Huang       xs     = grloc[r][0];
312a13bc4e3SShao-Ching Huang       xm     = grloc[r][1];
313a13bc4e3SShao-Ching Huang       ys     = grloc[r][2];
314a13bc4e3SShao-Ching Huang       ym     = grloc[r][3];
315a13bc4e3SShao-Ching Huang       zs     = grloc[r][4];
316a13bc4e3SShao-Ching Huang       zm     = grloc[r][5];
317a13bc4e3SShao-Ching Huang       nnodes = xm*ym*zm;
318a13bc4e3SShao-Ching Huang     }
319a13bc4e3SShao-Ching Huang     maxnnodes = PetscMax(maxnnodes,nnodes);
320a13bc4e3SShao-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);
321a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"      <Coordinates>\n");CHKERRQ(ierr);
322a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Xcoord\"  format=\"appended\"  offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
323a13bc4e3SShao-Ching Huang     boffset += xm*sizeof(PetscScalar) + sizeof(int);
324a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Ycoord\"  format=\"appended\"  offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
325a13bc4e3SShao-Ching Huang     boffset += ym*sizeof(PetscScalar) + sizeof(int);
326a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Zcoord\"  format=\"appended\"  offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
327a13bc4e3SShao-Ching Huang     boffset += zm*sizeof(PetscScalar) + sizeof(int);
328a13bc4e3SShao-Ching Huang     ierr     = PetscFPrintf(comm,fp,"      </Coordinates>\n");CHKERRQ(ierr);
329a13bc4e3SShao-Ching Huang     ierr = PetscFPrintf(comm,fp,"      <PointData Scalars=\"ScalarPointData\">\n");CHKERRQ(ierr);
330a13bc4e3SShao-Ching Huang     for (link=vtk->link; link; link=link->next) {
331a13bc4e3SShao-Ching Huang       Vec        X = (Vec)link->vec;
332fb5bd1c2SPatrick Sanan       PetscInt   bs,f;
333ea2d7708SPatrick Sanan       DM         daCurr;
334fb5bd1c2SPatrick Sanan       PetscBool  fieldsnamed;
335fb5bd1c2SPatrick Sanan       const char *vecname = "Unnamed";
336ea2d7708SPatrick Sanan 
337ea2d7708SPatrick Sanan       ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr);
338ea78f98cSLisandro Dalcin       ierr = DMDAGetInfo(daCurr,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
339ea2d7708SPatrick Sanan       maxbs = PetscMax(maxbs,bs);
340ea2d7708SPatrick Sanan       if (((PetscObject)X)->name || link != vtk->link) {
341a13bc4e3SShao-Ching Huang         ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
342a13bc4e3SShao-Ching Huang       }
343ea2d7708SPatrick Sanan 
344fb5bd1c2SPatrick Sanan       /* If any fields are named, add scalar fields. Otherwise, add a vector field */
345fb5bd1c2SPatrick Sanan       ierr = DMDAGetFieldsNamed(daCurr,&fieldsnamed);CHKERRQ(ierr);
346fb5bd1c2SPatrick Sanan       if (fieldsnamed) {
347fb5bd1c2SPatrick Sanan         for (f=0; f<bs; f++) {
348fb5bd1c2SPatrick Sanan           char       buf[256];
349fb5bd1c2SPatrick Sanan           const char *fieldname;
350fb5bd1c2SPatrick Sanan           ierr = DMDAGetFieldName(daCurr,f,&fieldname);CHKERRQ(ierr);
351fb5bd1c2SPatrick Sanan           if (!fieldname) {
352fb5bd1c2SPatrick Sanan             ierr      = PetscSNPrintf(buf,sizeof(buf),"%D",f);CHKERRQ(ierr);
353fb5bd1c2SPatrick Sanan             fieldname = buf;
354fb5bd1c2SPatrick Sanan           }
355fb5bd1c2SPatrick Sanan           ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
356fb5bd1c2SPatrick Sanan           boffset += nnodes*sizeof(PetscScalar) + sizeof(int);
357fb5bd1c2SPatrick Sanan         }
358fb5bd1c2SPatrick Sanan       } else {
359ea2d7708SPatrick Sanan         ierr     = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%D\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,bs,boffset);CHKERRQ(ierr);
360ea2d7708SPatrick Sanan         boffset += bs*nnodes*sizeof(PetscScalar) + sizeof(int);
361a13bc4e3SShao-Ching Huang       }
362fb5bd1c2SPatrick Sanan     }
363a13bc4e3SShao-Ching Huang     ierr = PetscFPrintf(comm,fp,"      </PointData>\n");CHKERRQ(ierr);
364a13bc4e3SShao-Ching Huang     ierr = PetscFPrintf(comm,fp,"    </Piece>\n");CHKERRQ(ierr);
365a13bc4e3SShao-Ching Huang   }
366a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"  </RectilinearGrid>\n");CHKERRQ(ierr);
367a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
368a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
369a13bc4e3SShao-Ching Huang 
370a13bc4e3SShao-Ching Huang   /* Now write the arrays. */
371a13bc4e3SShao-Ching Huang   tag  = ((PetscObject)viewer)->tag;
372fb5bd1c2SPatrick 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);
373a13bc4e3SShao-Ching Huang   for (r=0; r<size; r++) {
374a13bc4e3SShao-Ching Huang     MPI_Status status;
375a13bc4e3SShao-Ching Huang     PetscInt   xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
376a13bc4e3SShao-Ching Huang     if (!rank) {
377a13bc4e3SShao-Ching Huang       xs     = grloc[r][0];
378a13bc4e3SShao-Ching Huang       xm     = grloc[r][1];
379a13bc4e3SShao-Ching Huang       ys     = grloc[r][2];
380a13bc4e3SShao-Ching Huang       ym     = grloc[r][3];
381a13bc4e3SShao-Ching Huang       zs     = grloc[r][4];
382a13bc4e3SShao-Ching Huang       zm     = grloc[r][5];
383a13bc4e3SShao-Ching Huang       nnodes = xm*ym*zm;
384a13bc4e3SShao-Ching Huang     } else if (r == rank) {
385a13bc4e3SShao-Ching Huang       nnodes = info.xm*info.ym*info.zm;
386a13bc4e3SShao-Ching Huang     }
387a13bc4e3SShao-Ching Huang     {                           /* Write the coordinates */
388a13bc4e3SShao-Ching Huang       Vec Coords;
389ea2d7708SPatrick Sanan 
390a13bc4e3SShao-Ching Huang       ierr = DMGetCoordinates(da,&Coords);CHKERRQ(ierr);
391a13bc4e3SShao-Ching Huang       if (Coords) {
392a13bc4e3SShao-Ching Huang         const PetscScalar *coords;
393a13bc4e3SShao-Ching Huang         ierr = VecGetArrayRead(Coords,&coords);CHKERRQ(ierr);
394a13bc4e3SShao-Ching Huang         if (!rank) {
395a13bc4e3SShao-Ching Huang           if (r) {
396a13bc4e3SShao-Ching Huang             PetscMPIInt nn;
397ffc4695bSBarry Smith             ierr = MPI_Recv(array,xm+ym+zm,MPIU_SCALAR,r,tag,comm,&status);CHKERRMPI(ierr);
398ffc4695bSBarry Smith             ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRMPI(ierr);
399a13bc4e3SShao-Ching Huang             if (nn != xm+ym+zm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
400a13bc4e3SShao-Ching Huang           } else {
401a13bc4e3SShao-Ching Huang             /* x coordinates */
402a13bc4e3SShao-Ching Huang             for (j=0, k=0, i=0; i<xm; i++) {
403a13bc4e3SShao-Ching Huang               PetscInt Iloc = i+xm*(j+ym*k);
404a13bc4e3SShao-Ching Huang               array[i] = coords[Iloc*dim + 0];
405a13bc4e3SShao-Ching Huang             }
406a13bc4e3SShao-Ching Huang             /* y coordinates */
407a13bc4e3SShao-Ching Huang             for (i=0, k=0, j=0; j<ym; j++) {
408a13bc4e3SShao-Ching Huang               PetscInt Iloc = i+xm*(j+ym*k);
409a13bc4e3SShao-Ching Huang               array[j+xm] = dim > 1 ? coords[Iloc*dim + 1] : 0;
410a13bc4e3SShao-Ching Huang             }
411a13bc4e3SShao-Ching Huang             /* z coordinates */
412a13bc4e3SShao-Ching Huang             for (i=0, j=0, k=0; k<zm; k++) {
413a13bc4e3SShao-Ching Huang               PetscInt Iloc = i+xm*(j+ym*k);
414a13bc4e3SShao-Ching Huang               array[k+xm+ym] = dim > 2 ? coords[Iloc*dim + 2] : 0;
415a13bc4e3SShao-Ching Huang             }
416a13bc4e3SShao-Ching Huang           }
417a13bc4e3SShao-Ching Huang         } else if (r == rank) {
418a13bc4e3SShao-Ching Huang           xm = info.xm;
419a13bc4e3SShao-Ching Huang           ym = info.ym;
420a13bc4e3SShao-Ching Huang           zm = info.zm;
421a13bc4e3SShao-Ching Huang           for (j=0, k=0, i=0; i<xm; i++) {
422a13bc4e3SShao-Ching Huang             PetscInt Iloc = i+xm*(j+ym*k);
423a13bc4e3SShao-Ching Huang             array2[i] = coords[Iloc*dim + 0];
424a13bc4e3SShao-Ching Huang           }
425a13bc4e3SShao-Ching Huang           for (i=0, k=0, j=0; j<ym; j++) {
426a13bc4e3SShao-Ching Huang             PetscInt Iloc = i+xm*(j+ym*k);
427a13bc4e3SShao-Ching Huang             array2[j+xm] = dim > 1 ? coords[Iloc*dim + 1] : 0;
428a13bc4e3SShao-Ching Huang           }
429a13bc4e3SShao-Ching Huang           for (i=0, j=0, k=0; k<zm; k++) {
430a13bc4e3SShao-Ching Huang             PetscInt Iloc = i+xm*(j+ym*k);
431a13bc4e3SShao-Ching Huang             array2[k+xm+ym] = dim > 2 ? coords[Iloc*dim + 2] : 0;
432a13bc4e3SShao-Ching Huang           }
433ffc4695bSBarry Smith           ierr = MPI_Send((void*)array2,xm+ym+zm,MPIU_SCALAR,0,tag,comm);CHKERRMPI(ierr);
434a13bc4e3SShao-Ching Huang         }
435a13bc4e3SShao-Ching Huang         ierr = VecRestoreArrayRead(Coords,&coords);CHKERRQ(ierr);
436a13bc4e3SShao-Ching Huang       } else {       /* Fabricate some coordinates using grid index */
4373d3eaba7SBarry Smith         for (i=0; i<xm; i++) {
438a13bc4e3SShao-Ching Huang           array[i] = xs+i;
439a13bc4e3SShao-Ching Huang         }
4403d3eaba7SBarry Smith         for (j=0; j<ym; j++) {
441a13bc4e3SShao-Ching Huang           array[j+xm] = ys+j;
442a13bc4e3SShao-Ching Huang         }
4433d3eaba7SBarry Smith         for (k=0; k<zm; k++) {
444a13bc4e3SShao-Ching Huang           array[k+xm+ym] = zs+k;
445a13bc4e3SShao-Ching Huang         }
446a13bc4e3SShao-Ching Huang       }
447a13bc4e3SShao-Ching Huang       if (!rank) {
44894fbd55eSBarry Smith         ierr = PetscViewerVTKFWrite(viewer,fp,&(array[0]),xm,MPIU_SCALAR);CHKERRQ(ierr);
44994fbd55eSBarry Smith         ierr = PetscViewerVTKFWrite(viewer,fp,&(array[xm]),ym,MPIU_SCALAR);CHKERRQ(ierr);
45094fbd55eSBarry Smith         ierr = PetscViewerVTKFWrite(viewer,fp,&(array[xm+ym]),zm,MPIU_SCALAR);CHKERRQ(ierr);
451a13bc4e3SShao-Ching Huang       }
452a13bc4e3SShao-Ching Huang     }
453a13bc4e3SShao-Ching Huang 
454a13bc4e3SShao-Ching Huang     /* Write each of the objects queued up for this file */
455a13bc4e3SShao-Ching Huang     for (link=vtk->link; link; link=link->next) {
456a13bc4e3SShao-Ching Huang       Vec               X = (Vec)link->vec;
457a13bc4e3SShao-Ching Huang       const PetscScalar *x;
458fb5bd1c2SPatrick Sanan       PetscInt          bs,f;
459ea2d7708SPatrick Sanan       DM                daCurr;
460fb5bd1c2SPatrick Sanan       PetscBool         fieldsnamed;
461ea2d7708SPatrick Sanan       ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr);
462ea78f98cSLisandro Dalcin       ierr = DMDAGetInfo(daCurr,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
463a13bc4e3SShao-Ching Huang 
464a13bc4e3SShao-Ching Huang       ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
465a13bc4e3SShao-Ching Huang       if (!rank) {
466a13bc4e3SShao-Ching Huang         if (r) {
467a13bc4e3SShao-Ching Huang           PetscMPIInt nn;
468ffc4695bSBarry Smith           ierr = MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status);CHKERRMPI(ierr);
469ffc4695bSBarry Smith           ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRMPI(ierr);
470a13bc4e3SShao-Ching Huang           if (nn != nnodes*bs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %D",r);
471a13bc4e3SShao-Ching Huang         } else {
472580bdb30SBarry Smith           ierr = PetscArraycpy(array,x,nnodes*bs);CHKERRQ(ierr);
473a13bc4e3SShao-Ching Huang         }
474fb5bd1c2SPatrick Sanan         /* If any fields are named, add scalar fields. Otherwise, add a vector field */
475fb5bd1c2SPatrick Sanan         ierr = DMDAGetFieldsNamed(daCurr,&fieldsnamed);CHKERRQ(ierr);
476fb5bd1c2SPatrick Sanan         if (fieldsnamed) {
477fb5bd1c2SPatrick Sanan           for (f=0; f<bs; f++) {
478fb5bd1c2SPatrick Sanan             /* Extract and transpose the f'th field */
479fb5bd1c2SPatrick Sanan             for (k=0; k<zm; k++) {
480fb5bd1c2SPatrick Sanan               for (j=0; j<ym; j++) {
481fb5bd1c2SPatrick Sanan                 for (i=0; i<xm; i++) {
482fb5bd1c2SPatrick Sanan                   PetscInt Iloc = i+xm*(j+ym*k);
483fb5bd1c2SPatrick Sanan                   array2[Iloc] = array[Iloc*bs + f];
484fb5bd1c2SPatrick Sanan                 }
485fb5bd1c2SPatrick Sanan               }
486fb5bd1c2SPatrick Sanan             }
487fb5bd1c2SPatrick Sanan             ierr = PetscViewerVTKFWrite(viewer,fp,array2,nnodes,MPIU_SCALAR);CHKERRQ(ierr);
488fb5bd1c2SPatrick Sanan           }
489fb5bd1c2SPatrick Sanan         }
490ea2d7708SPatrick Sanan         ierr = PetscViewerVTKFWrite(viewer,fp,array,nnodes*bs,MPIU_SCALAR);CHKERRQ(ierr);
491ea2d7708SPatrick Sanan 
492a13bc4e3SShao-Ching Huang       } else if (r == rank) {
493ffc4695bSBarry Smith         ierr = MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm);CHKERRMPI(ierr);
494a13bc4e3SShao-Ching Huang       }
495a13bc4e3SShao-Ching Huang       ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
496a13bc4e3SShao-Ching Huang     }
497a13bc4e3SShao-Ching Huang   }
498a13bc4e3SShao-Ching Huang   ierr = PetscFree2(array,array2);CHKERRQ(ierr);
499a13bc4e3SShao-Ching Huang   ierr = PetscFree(grloc);CHKERRQ(ierr);
500a13bc4e3SShao-Ching Huang 
501a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr);
502a13bc4e3SShao-Ching Huang   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
503a13bc4e3SShao-Ching Huang   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
504a13bc4e3SShao-Ching Huang   PetscFunctionReturn(0);
505a13bc4e3SShao-Ching Huang }
506a13bc4e3SShao-Ching Huang 
5074061b8bfSJed Brown /*@C
5084061b8bfSJed Brown    DMDAVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
5094061b8bfSJed Brown 
5104061b8bfSJed Brown    Collective
5114061b8bfSJed Brown 
5124061b8bfSJed Brown    Input Arguments:
5134061b8bfSJed Brown    odm - DM specifying the grid layout, passed as a PetscObject
5144061b8bfSJed Brown    viewer - viewer of type VTK
5154061b8bfSJed Brown 
5164061b8bfSJed Brown    Level: developer
5174061b8bfSJed Brown 
518fb5bd1c2SPatrick Sanan    Notes:
5194061b8bfSJed Brown    This function is a callback used by the VTK viewer to actually write the file.
5204061b8bfSJed 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.
5214061b8bfSJed Brown    Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
5224061b8bfSJed Brown 
523fb5bd1c2SPatrick Sanan    If any fields have been named (see e.g. DMDASetFieldName()), then individual scalar
524fb5bd1c2SPatrick Sanan    fields are written. Otherwise, a single multi-dof (vector) field is written.
525fb5bd1c2SPatrick Sanan 
5264061b8bfSJed Brown .seealso: PETSCVIEWERVTK
5274061b8bfSJed Brown @*/
5284061b8bfSJed Brown PetscErrorCode DMDAVTKWriteAll(PetscObject odm,PetscViewer viewer)
5294061b8bfSJed Brown {
5304061b8bfSJed Brown   DM             dm = (DM)odm;
5314061b8bfSJed Brown   PetscBool      isvtk;
5324061b8bfSJed Brown   PetscErrorCode ierr;
5334061b8bfSJed Brown 
5344061b8bfSJed Brown   PetscFunctionBegin;
5354061b8bfSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
5364061b8bfSJed Brown   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
537251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr);
538ce94432eSBarry Smith   if (!isvtk) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_INCOMP,"Cannot use viewer type %s",((PetscObject)viewer)->type_name);
5394061b8bfSJed Brown   switch (viewer->format) {
5404061b8bfSJed Brown   case PETSC_VIEWER_VTK_VTS:
5414061b8bfSJed Brown     ierr = DMDAVTKWriteAll_VTS(dm,viewer);CHKERRQ(ierr);
5424061b8bfSJed Brown     break;
543a13bc4e3SShao-Ching Huang   case PETSC_VIEWER_VTK_VTR:
544a13bc4e3SShao-Ching Huang     ierr = DMDAVTKWriteAll_VTR(dm,viewer);CHKERRQ(ierr);
545a13bc4e3SShao-Ching Huang     break;
546ce94432eSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No support for format '%s'",PetscViewerFormats[viewer->format]);
5474061b8bfSJed Brown   }
5484061b8bfSJed Brown   PetscFunctionReturn(0);
5494061b8bfSJed Brown }
550