xref: /petsc/src/dm/impls/da/grvtk.c (revision 4061b8bf102ddefbaac656687f58545bb9370df1)
1*4061b8bfSJed Brown #include <private/daimpl.h>
2*4061b8bfSJed Brown #include <../src/sys/viewer/impls/vtk/vtkvimpl.h>
3*4061b8bfSJed Brown 
4*4061b8bfSJed Brown #if defined(PETSC_HAVE_STDINT_H) /* The VTK format requires a 32-bit integer */
5*4061b8bfSJed Brown typedef int32_t PetscVTKInt;
6*4061b8bfSJed Brown #else
7*4061b8bfSJed Brown typedef int PetscVTKInt;
8*4061b8bfSJed Brown #endif
9*4061b8bfSJed Brown #define PETSC_VTK_INT_MAX  2147483647
10*4061b8bfSJed Brown #define PETSC_VTK_INT_MIN -2147483647
11*4061b8bfSJed Brown #if defined(PETSC_USE_64BIT_INDICES)
12*4061b8bfSJed Brown #  define PetscVTKIntCheck(a)  if ((a) > PETSC_VTK_INT_MAX) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Array too long for 32-bit VTK binary format")
13*4061b8bfSJed Brown #  define PetscVTKIntCast(a) (PetscVTKInt)(a);PetscVTKIntCheck(a)
14*4061b8bfSJed Brown #else
15*4061b8bfSJed Brown #  define PetscVTKIntCheck(a)
16*4061b8bfSJed Brown #  define PetscVTKIntCast(a) a
17*4061b8bfSJed Brown #endif
18*4061b8bfSJed Brown 
19*4061b8bfSJed Brown #undef __FUNCT__
20*4061b8bfSJed Brown #define __FUNCT__ "PetscFWrite_VTK"
21*4061b8bfSJed Brown /* Write binary data preceded by 32-bit int length (in bytes), does not do byte swapping. */
22*4061b8bfSJed Brown static PetscErrorCode PetscFWrite_VTK(MPI_Comm comm,FILE *fp,void *data,PetscInt n,PetscDataType dtype)
23*4061b8bfSJed Brown {
24*4061b8bfSJed Brown   PetscErrorCode ierr;
25*4061b8bfSJed Brown   PetscMPIInt rank;
26*4061b8bfSJed Brown 
27*4061b8bfSJed Brown   PetscFunctionBegin;
28*4061b8bfSJed Brown   if (n < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Trying to write a negative amount of data %D",n);
29*4061b8bfSJed Brown   if (!n) PetscFunctionReturn(0);
30*4061b8bfSJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
31*4061b8bfSJed Brown   if (!rank) {
32*4061b8bfSJed Brown     size_t count;
33*4061b8bfSJed Brown     PetscInt size;
34*4061b8bfSJed Brown     PetscVTKInt bytes;
35*4061b8bfSJed Brown     switch (dtype) {
36*4061b8bfSJed Brown     case PETSC_DOUBLE:
37*4061b8bfSJed Brown       size = sizeof(double);
38*4061b8bfSJed Brown       break;
39*4061b8bfSJed Brown     case PETSC_FLOAT:
40*4061b8bfSJed Brown       size = sizeof(float);
41*4061b8bfSJed Brown       break;
42*4061b8bfSJed Brown     case PETSC_INT:
43*4061b8bfSJed Brown       size = sizeof(PetscInt);
44*4061b8bfSJed Brown       break;
45*4061b8bfSJed Brown     case PETSC_CHAR:
46*4061b8bfSJed Brown       size = sizeof(char);
47*4061b8bfSJed Brown       break;
48*4061b8bfSJed Brown     default: SETERRQ(comm,PETSC_ERR_SUP,"Data type not supported");
49*4061b8bfSJed Brown     }
50*4061b8bfSJed Brown     bytes = PetscVTKIntCast(size*n);
51*4061b8bfSJed Brown 
52*4061b8bfSJed Brown     count = fwrite(&bytes,sizeof(int),1,fp);
53*4061b8bfSJed Brown     if (count != 1) {
54*4061b8bfSJed Brown       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_WRITE,"Error writing byte count");
55*4061b8bfSJed Brown     }
56*4061b8bfSJed Brown     count = fwrite(data,size,(size_t)n,fp);
57*4061b8bfSJed Brown     if ((PetscInt)count != n) {
58*4061b8bfSJed Brown       SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_FILE_WRITE,"Wrote %D/%D array members of size %D",(PetscInt)count,n,(PetscInt)size);
59*4061b8bfSJed Brown     }
60*4061b8bfSJed Brown   }
61*4061b8bfSJed Brown   PetscFunctionReturn(0);
62*4061b8bfSJed Brown }
63*4061b8bfSJed Brown 
64*4061b8bfSJed Brown #undef __FUNCT__
65*4061b8bfSJed Brown #define __FUNCT__ "DMDAVTKWriteAll_VTS"
66*4061b8bfSJed Brown static PetscErrorCode DMDAVTKWriteAll_VTS(DM da,PetscViewer viewer)
67*4061b8bfSJed Brown {
68*4061b8bfSJed Brown #if defined(PETSC_USE_REAL_SINGLE)
69*4061b8bfSJed Brown   const char precision[]  = "Float32";
70*4061b8bfSJed Brown #elif defined(PETSC_USE_REAL_DOUBLE)
71*4061b8bfSJed Brown   const char precision[]  = "Float64";
72*4061b8bfSJed Brown #else
73*4061b8bfSJed Brown   const char precision[]  = "UnknownPrecision";
74*4061b8bfSJed Brown #endif
75*4061b8bfSJed Brown   MPI_Comm                 comm = ((PetscObject)da)->comm;
76*4061b8bfSJed Brown   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
77*4061b8bfSJed Brown   PetscViewerVTKObjectLink link;
78*4061b8bfSJed Brown   FILE                     *fp;
79*4061b8bfSJed Brown   PetscMPIInt              rank,size,tag;
80*4061b8bfSJed Brown   DMDALocalInfo            info;
81*4061b8bfSJed Brown   PetscInt                 dim,mx,my,mz,bs,boffset,maxnnodes,i,j,k,f,r;
82*4061b8bfSJed Brown   PetscInt                 rloc[6],(*grloc)[6] = PETSC_NULL;
83*4061b8bfSJed Brown   PetscScalar              *array,*array2;
84*4061b8bfSJed Brown   PetscReal                gmin[3],gmax[3];
85*4061b8bfSJed Brown   PetscErrorCode           ierr;
86*4061b8bfSJed Brown 
87*4061b8bfSJed Brown   PetscFunctionBegin;
88*4061b8bfSJed Brown #if defined(PETSC_USE_COMPLEX)
89*4061b8bfSJed Brown   SETERRQ(((PetscObject)da)->comm,PETSC_ERR_SUP,"Complex values not supported");
90*4061b8bfSJed Brown #endif
91*4061b8bfSJed Brown 
92*4061b8bfSJed Brown   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
93*4061b8bfSJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
94*4061b8bfSJed Brown   ierr = DMDAGetInfo(da,&dim, &mx,&my,&mz, 0,0,0, &bs,0,0,0,0,0);CHKERRQ(ierr);
95*4061b8bfSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
96*4061b8bfSJed Brown   ierr = DMDAGetBoundingBox(da,gmin,gmax);CHKERRQ(ierr);
97*4061b8bfSJed Brown 
98*4061b8bfSJed Brown   ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr);
99*4061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr);
100*4061b8bfSJed Brown #ifdef PETSC_WORDS_BIGENDIAN
101*4061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"BigEndian\">\n");CHKERRQ(ierr);
102*4061b8bfSJed Brown #else
103*4061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"LittleEndian\">\n");CHKERRQ(ierr);
104*4061b8bfSJed Brown #endif
105*4061b8bfSJed 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);
106*4061b8bfSJed Brown 
107*4061b8bfSJed Brown   if (!rank) {ierr = PetscMalloc(size*6*sizeof(PetscInt),&grloc);CHKERRQ(ierr);}
108*4061b8bfSJed Brown   rloc[0] = info.xs;
109*4061b8bfSJed Brown   rloc[1] = info.xm;
110*4061b8bfSJed Brown   rloc[2] = info.ys;
111*4061b8bfSJed Brown   rloc[3] = info.ym;
112*4061b8bfSJed Brown   rloc[4] = info.zs;
113*4061b8bfSJed Brown   rloc[5] = info.zm;
114*4061b8bfSJed Brown   ierr = MPI_Gather(rloc,6,MPIU_INT,&grloc[0][0],6,MPIU_INT,0,comm);CHKERRQ(ierr);
115*4061b8bfSJed Brown 
116*4061b8bfSJed Brown   /* Write XML header */
117*4061b8bfSJed Brown   maxnnodes = 0;                /* Used for the temporary array size on rank 0 */
118*4061b8bfSJed Brown   boffset = 0;                  /* Offset into binary file */
119*4061b8bfSJed Brown   for (r=0; r<size; r++) {
120*4061b8bfSJed Brown     PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
121*4061b8bfSJed Brown     if (!rank) {
122*4061b8bfSJed Brown       xs = grloc[r][0];
123*4061b8bfSJed Brown       xm = grloc[r][1];
124*4061b8bfSJed Brown       ys = grloc[r][2];
125*4061b8bfSJed Brown       ym = grloc[r][3];
126*4061b8bfSJed Brown       zs = grloc[r][4];
127*4061b8bfSJed Brown       zm = grloc[r][5];
128*4061b8bfSJed Brown       nnodes = xm*ym*zm;
129*4061b8bfSJed Brown     }
130*4061b8bfSJed Brown     maxnnodes = PetscMax(maxnnodes,nnodes);
131*4061b8bfSJed Brown #if 0
132*4061b8bfSJed Brown     switch (dim) {
133*4061b8bfSJed Brown     case 1:
134*4061b8bfSJed Brown       ierr = PetscFPrintf(comm,fp,"    <Piece Extent=\"%D %D %D %D %D %D\">\n",xs,xs+xm-1,0,0,0,0);CHKERRQ(ierr);
135*4061b8bfSJed Brown       break;
136*4061b8bfSJed Brown     case 2:
137*4061b8bfSJed 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);
138*4061b8bfSJed Brown       break;
139*4061b8bfSJed Brown     case 3:
140*4061b8bfSJed 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);
141*4061b8bfSJed Brown       break;
142*4061b8bfSJed Brown     default: SETERRQ1(((PetscObject)da)->comm,PETSC_ERR_SUP,"No support for dimension %D",dim);
143*4061b8bfSJed Brown     }
144*4061b8bfSJed Brown #endif
145*4061b8bfSJed 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);
146*4061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"      <Points>\n");CHKERRQ(ierr);
147*4061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr);
148*4061b8bfSJed Brown     boffset += 3*nnodes*sizeof(PetscScalar) + sizeof(int);
149*4061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"      </Points>\n");CHKERRQ(ierr);
150*4061b8bfSJed Brown 
151*4061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"      <PointData Scalars=\"ScalarPointData\">\n");CHKERRQ(ierr);
152*4061b8bfSJed Brown     for (link=vtk->link; link; link=link->next) {
153*4061b8bfSJed Brown       Vec X = (Vec)link->vec;
154*4061b8bfSJed Brown       const char *vecname = "";
155*4061b8bfSJed 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. */
156*4061b8bfSJed Brown         ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr);
157*4061b8bfSJed Brown       }
158*4061b8bfSJed Brown       for (i=0; i<bs; i++) {
159*4061b8bfSJed Brown         char buf[256];
160*4061b8bfSJed Brown         const char *fieldname;
161*4061b8bfSJed Brown         ierr = DMDAGetFieldName(da,i,&fieldname);CHKERRQ(ierr);
162*4061b8bfSJed Brown         if (!fieldname) {
163*4061b8bfSJed Brown           ierr = PetscSNPrintf(buf,sizeof buf,"Unnamed%D",i);CHKERRQ(ierr);
164*4061b8bfSJed Brown           fieldname = buf;
165*4061b8bfSJed Brown         }
166*4061b8bfSJed Brown         ierr = PetscFPrintf(comm,fp,"        <DataArray type=\"%s\" Name=\"%s%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr);
167*4061b8bfSJed Brown         boffset += nnodes*sizeof(PetscScalar) + sizeof(int);
168*4061b8bfSJed Brown       }
169*4061b8bfSJed Brown     }
170*4061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"      </PointData>\n");CHKERRQ(ierr);
171*4061b8bfSJed Brown     ierr = PetscFPrintf(comm,fp,"    </Piece>\n");CHKERRQ(ierr);
172*4061b8bfSJed Brown   }
173*4061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"  </StructuredGrid>\n");CHKERRQ(ierr);
174*4061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"  <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr);
175*4061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr);
176*4061b8bfSJed Brown 
177*4061b8bfSJed Brown   /* Now write the arrays. */
178*4061b8bfSJed Brown   tag  = ((PetscObject)viewer)->tag;
179*4061b8bfSJed Brown   ierr = PetscMalloc2(maxnnodes*PetscMax(3,bs),PetscScalar,&array,maxnnodes*3,PetscScalar,&array2);CHKERRQ(ierr);
180*4061b8bfSJed Brown   for (r=0; r<size; r++) {
181*4061b8bfSJed Brown     MPI_Status status;
182*4061b8bfSJed Brown     PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0;
183*4061b8bfSJed Brown     if (!rank) {
184*4061b8bfSJed Brown       xs = grloc[r][0];
185*4061b8bfSJed Brown       xm = grloc[r][1];
186*4061b8bfSJed Brown       ys = grloc[r][2];
187*4061b8bfSJed Brown       ym = grloc[r][3];
188*4061b8bfSJed Brown       zs = grloc[r][4];
189*4061b8bfSJed Brown       zm = grloc[r][5];
190*4061b8bfSJed Brown       nnodes = xm*ym*zm;
191*4061b8bfSJed Brown     } else if (r == rank) {
192*4061b8bfSJed Brown       nnodes = info.xm*info.ym*info.zm;
193*4061b8bfSJed Brown     }
194*4061b8bfSJed Brown 
195*4061b8bfSJed Brown     {                           /* Write the coordinates */
196*4061b8bfSJed Brown       Vec Coords;
197*4061b8bfSJed Brown       ierr = DMDAGetCoordinates(da,&Coords);CHKERRQ(ierr);
198*4061b8bfSJed Brown       if (Coords) {
199*4061b8bfSJed Brown         const PetscScalar *coords;
200*4061b8bfSJed Brown         ierr = VecGetArrayRead(Coords,&coords);CHKERRQ(ierr);
201*4061b8bfSJed Brown         if (!rank) {
202*4061b8bfSJed Brown           if (r) {
203*4061b8bfSJed Brown             PetscInt nn;
204*4061b8bfSJed Brown             ierr = MPI_Recv(array,nnodes*3,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr);
205*4061b8bfSJed Brown             ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr);
206*4061b8bfSJed Brown             if (nn != nnodes*3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch");
207*4061b8bfSJed Brown           } else {
208*4061b8bfSJed Brown             ierr = PetscMemcpy(array,coords,nnodes*3*sizeof(PetscScalar));CHKERRQ(ierr);
209*4061b8bfSJed Brown           }
210*4061b8bfSJed Brown           /* Transpose coordinates to VTK (C-style) ordering */
211*4061b8bfSJed Brown           for (k=0; k<zm; k++) {
212*4061b8bfSJed Brown             for (j=0; j<ym; j++) {
213*4061b8bfSJed Brown               for (i=0; i<xm; i++) {
214*4061b8bfSJed Brown                 PetscInt Iloc = i+xm*(j+ym*k);
215*4061b8bfSJed Brown                 array2[Iloc*3+0] = array[Iloc*dim + 0];
216*4061b8bfSJed Brown                 array2[Iloc*3+1] = dim > 1 ? array[Iloc*dim + 1] : 0;
217*4061b8bfSJed Brown                 array2[Iloc*3+2] = dim > 2 ? array[Iloc*dim + 2] : 0;
218*4061b8bfSJed Brown               }
219*4061b8bfSJed Brown             }
220*4061b8bfSJed Brown           }
221*4061b8bfSJed Brown         } else if (r == rank) {
222*4061b8bfSJed Brown           ierr = MPI_Send((void*)coords,nnodes*3,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr);
223*4061b8bfSJed Brown         }
224*4061b8bfSJed Brown         ierr = VecRestoreArrayRead(Coords,&coords);CHKERRQ(ierr);
225*4061b8bfSJed Brown       } else {       /* Fabricate some coordinates using grid index */
226*4061b8bfSJed Brown         for (k=0; k<zm; k++) {
227*4061b8bfSJed Brown           for (j=0; j<ym; j++) {
228*4061b8bfSJed Brown             for (i=0; i<xm; i++) {
229*4061b8bfSJed Brown               PetscInt Iloc = i+xm*(j+ym*k);
230*4061b8bfSJed Brown               array2[Iloc*3+0] = xs+i;
231*4061b8bfSJed Brown               array2[Iloc*3+1] = ys+j;
232*4061b8bfSJed Brown               array2[Iloc*3+2] = zs+k;
233*4061b8bfSJed Brown             }
234*4061b8bfSJed Brown           }
235*4061b8bfSJed Brown         }
236*4061b8bfSJed Brown       }
237*4061b8bfSJed Brown       ierr = PetscFWrite_VTK(comm,fp,array2,nnodes*3,PETSC_SCALAR);CHKERRQ(ierr);
238*4061b8bfSJed Brown     }
239*4061b8bfSJed Brown 
240*4061b8bfSJed Brown     /* Write each of the objects queued up for this file */
241*4061b8bfSJed Brown     for (link=vtk->link; link; link=link->next) {
242*4061b8bfSJed Brown       Vec X = (Vec)link->vec;
243*4061b8bfSJed Brown       const PetscScalar *x;
244*4061b8bfSJed Brown 
245*4061b8bfSJed Brown       ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
246*4061b8bfSJed Brown       if (!rank) {
247*4061b8bfSJed Brown         if (r) {
248*4061b8bfSJed Brown           PetscInt nn;
249*4061b8bfSJed Brown           ierr = MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr);
250*4061b8bfSJed Brown           ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr);
251*4061b8bfSJed Brown           if (nn != nnodes*bs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %D",r);
252*4061b8bfSJed Brown         } else {
253*4061b8bfSJed Brown           ierr = PetscMemcpy(array,x,nnodes*bs*sizeof(PetscScalar));CHKERRQ(ierr);
254*4061b8bfSJed Brown         }
255*4061b8bfSJed Brown         for (f=0; f<bs; f++) {
256*4061b8bfSJed Brown           /* Extract and transpose the f'th field */
257*4061b8bfSJed Brown           for (k=0; k<zm; k++) {
258*4061b8bfSJed Brown             for (j=0; j<ym; j++) {
259*4061b8bfSJed Brown               for (i=0; i<xm; i++) {
260*4061b8bfSJed Brown                 PetscInt Iloc = i+xm*(j+ym*k);
261*4061b8bfSJed Brown                 array2[Iloc] = array[Iloc*bs + f];
262*4061b8bfSJed Brown               }
263*4061b8bfSJed Brown             }
264*4061b8bfSJed Brown           }
265*4061b8bfSJed Brown           ierr = PetscFWrite_VTK(comm,fp,array2,nnodes,PETSC_SCALAR);CHKERRQ(ierr);
266*4061b8bfSJed Brown         }
267*4061b8bfSJed Brown       } else if (r == rank) {
268*4061b8bfSJed Brown         ierr = MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr);
269*4061b8bfSJed Brown       }
270*4061b8bfSJed Brown       ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
271*4061b8bfSJed Brown     }
272*4061b8bfSJed Brown   }
273*4061b8bfSJed Brown   ierr = PetscFree2(array,array2);CHKERRQ(ierr);
274*4061b8bfSJed Brown   ierr = PetscFree(grloc);CHKERRQ(ierr);
275*4061b8bfSJed Brown 
276*4061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr);
277*4061b8bfSJed Brown   ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr);
278*4061b8bfSJed Brown   ierr = PetscFClose(comm,fp);CHKERRQ(ierr);
279*4061b8bfSJed Brown   PetscFunctionReturn(0);
280*4061b8bfSJed Brown }
281*4061b8bfSJed Brown 
282*4061b8bfSJed Brown 
283*4061b8bfSJed Brown #undef __FUNCT__
284*4061b8bfSJed Brown #define __FUNCT__ "DMDAVTKWriteAll"
285*4061b8bfSJed Brown /*@C
286*4061b8bfSJed Brown    DMDAVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
287*4061b8bfSJed Brown 
288*4061b8bfSJed Brown    Collective
289*4061b8bfSJed Brown 
290*4061b8bfSJed Brown    Input Arguments:
291*4061b8bfSJed Brown    odm - DM specifying the grid layout, passed as a PetscObject
292*4061b8bfSJed Brown    viewer - viewer of type VTK
293*4061b8bfSJed Brown 
294*4061b8bfSJed Brown    Level: developer
295*4061b8bfSJed Brown 
296*4061b8bfSJed Brown    Note:
297*4061b8bfSJed Brown    This function is a callback used by the VTK viewer to actually write the file.
298*4061b8bfSJed 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.
299*4061b8bfSJed Brown    Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
300*4061b8bfSJed Brown 
301*4061b8bfSJed Brown .seealso: PETSCVIEWERVTK
302*4061b8bfSJed Brown @*/
303*4061b8bfSJed Brown PetscErrorCode DMDAVTKWriteAll(PetscObject odm,PetscViewer viewer)
304*4061b8bfSJed Brown {
305*4061b8bfSJed Brown   DM              dm   = (DM)odm;
306*4061b8bfSJed Brown   PetscBool       isvtk;
307*4061b8bfSJed Brown   PetscErrorCode  ierr;
308*4061b8bfSJed Brown 
309*4061b8bfSJed Brown   PetscFunctionBegin;
310*4061b8bfSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
311*4061b8bfSJed Brown   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
312*4061b8bfSJed Brown   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr);
313*4061b8bfSJed Brown   if (!isvtk) SETERRQ1(((PetscObject)viewer)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use viewer type %s",((PetscObject)viewer)->type_name);
314*4061b8bfSJed Brown   switch (viewer->format) {
315*4061b8bfSJed Brown   case PETSC_VIEWER_VTK_VTS:
316*4061b8bfSJed Brown     ierr = DMDAVTKWriteAll_VTS(dm,viewer);CHKERRQ(ierr);
317*4061b8bfSJed Brown     break;
318*4061b8bfSJed Brown   default: SETERRQ1(((PetscObject)dm)->comm,PETSC_ERR_SUP,"No support for format '%s'",PetscViewerFormats[viewer->format]);
319*4061b8bfSJed Brown   }
320*4061b8bfSJed Brown   PetscFunctionReturn(0);
321*4061b8bfSJed Brown }
322