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 ierr = DMDAGetDof(da,&bs);CHKERRQ(ierr); 18fb5bd1c2SPatrick Sanan *fieldsnamed = PETSC_FALSE; 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 { 32*30815ce0SLisandro 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 PetscReal gmin[3],gmax[3]; 514061b8bfSJed Brown PetscErrorCode ierr; 524061b8bfSJed Brown 534061b8bfSJed Brown PetscFunctionBegin; 54094921d9SJed Brown ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 554061b8bfSJed Brown #if defined(PETSC_USE_COMPLEX) 563bf036e2SBarry Smith SETERRQ(comm,PETSC_ERR_SUP,"Complex values not supported"); 574061b8bfSJed Brown #endif 584061b8bfSJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 594061b8bfSJed Brown ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 604061b8bfSJed Brown ierr = DMDAGetInfo(da,&dim,&mx,&my,&mz,0,0,0,&bs,0,0,0,0,0);CHKERRQ(ierr); 614061b8bfSJed Brown ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 624061b8bfSJed Brown ierr = DMDAGetBoundingBox(da,gmin,gmax);CHKERRQ(ierr); 632eaa9ef4SLisandro Dalcin ierr = DMGetCoordinates(da,&Coords);CHKERRQ(ierr); 642eaa9ef4SLisandro Dalcin if (Coords) { 652eaa9ef4SLisandro Dalcin PetscInt csize; 662eaa9ef4SLisandro Dalcin ierr = VecGetSize(Coords,&csize);CHKERRQ(ierr); 672eaa9ef4SLisandro Dalcin if (csize % (mx*my*mz)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Coordinate vector size mismatch"); 682eaa9ef4SLisandro Dalcin cdim = csize/(mx*my*mz); 692eaa9ef4SLisandro Dalcin if (cdim < dim || cdim > 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Coordinate vector size mismatch"); 702eaa9ef4SLisandro Dalcin } else { 712eaa9ef4SLisandro Dalcin cdim = dim; 722eaa9ef4SLisandro Dalcin } 734061b8bfSJed Brown 744061b8bfSJed Brown ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr); 754061b8bfSJed Brown ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr); 76*30815ce0SLisandro Dalcin ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"%s\">\n",byte_order);CHKERRQ(ierr); 774061b8bfSJed 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); 784061b8bfSJed Brown 79785e854fSJed Brown if (!rank) {ierr = PetscMalloc1(size*6,&grloc);CHKERRQ(ierr);} 804061b8bfSJed Brown rloc[0] = info.xs; 814061b8bfSJed Brown rloc[1] = info.xm; 824061b8bfSJed Brown rloc[2] = info.ys; 834061b8bfSJed Brown rloc[3] = info.ym; 844061b8bfSJed Brown rloc[4] = info.zs; 854061b8bfSJed Brown rloc[5] = info.zm; 864061b8bfSJed Brown ierr = MPI_Gather(rloc,6,MPIU_INT,&grloc[0][0],6,MPIU_INT,0,comm);CHKERRQ(ierr); 874061b8bfSJed Brown 884061b8bfSJed Brown /* Write XML header */ 894061b8bfSJed Brown maxnnodes = 0; /* Used for the temporary array size on rank 0 */ 90ea2d7708SPatrick Sanan maxbs = 0; /* Used for the temporary array size on rank 0 */ 914061b8bfSJed Brown boffset = 0; /* Offset into binary file */ 924061b8bfSJed Brown for (r=0; r<size; r++) { 934061b8bfSJed Brown PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0; 944061b8bfSJed Brown if (!rank) { 954061b8bfSJed Brown xs = grloc[r][0]; 964061b8bfSJed Brown xm = grloc[r][1]; 974061b8bfSJed Brown ys = grloc[r][2]; 984061b8bfSJed Brown ym = grloc[r][3]; 994061b8bfSJed Brown zs = grloc[r][4]; 1004061b8bfSJed Brown zm = grloc[r][5]; 1014061b8bfSJed Brown nnodes = xm*ym*zm; 1024061b8bfSJed Brown } 1034061b8bfSJed Brown maxnnodes = PetscMax(maxnnodes,nnodes); 1044061b8bfSJed 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); 1054061b8bfSJed Brown ierr = PetscFPrintf(comm,fp," <Points>\n");CHKERRQ(ierr); 1064061b8bfSJed Brown ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr); 1074061b8bfSJed Brown boffset += 3*nnodes*sizeof(PetscScalar) + sizeof(int); 1084061b8bfSJed Brown ierr = PetscFPrintf(comm,fp," </Points>\n");CHKERRQ(ierr); 1094061b8bfSJed Brown 1104061b8bfSJed Brown ierr = PetscFPrintf(comm,fp," <PointData Scalars=\"ScalarPointData\">\n");CHKERRQ(ierr); 1114061b8bfSJed Brown for (link=vtk->link; link; link=link->next) { 1124061b8bfSJed Brown Vec X = (Vec)link->vec; 113fb5bd1c2SPatrick Sanan PetscInt bs,f; 114ea2d7708SPatrick Sanan DM daCurr; 115fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 116fb5bd1c2SPatrick Sanan const char *vecname = "Unnamed"; 117ea2d7708SPatrick Sanan 118ea2d7708SPatrick Sanan ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr); 119ea2d7708SPatrick Sanan ierr = DMDAGetInfo(daCurr,0,0,0,0,0,0,0,&bs,0,0,0,0,0);CHKERRQ(ierr); 120ea2d7708SPatrick Sanan maxbs = PetscMax(maxbs,bs); 121ea2d7708SPatrick Sanan 122ea2d7708SPatrick Sanan if (((PetscObject)X)->name || link != vtk->link) { 1234061b8bfSJed Brown ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr); 1244061b8bfSJed Brown } 125fb5bd1c2SPatrick Sanan 126fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 127fb5bd1c2SPatrick Sanan ierr = DMDAGetFieldsNamed(daCurr,&fieldsnamed);CHKERRQ(ierr); 128fb5bd1c2SPatrick Sanan if (fieldsnamed) { 129fb5bd1c2SPatrick Sanan for (f=0; f<bs; f++) { 130fb5bd1c2SPatrick Sanan char buf[256]; 131fb5bd1c2SPatrick Sanan const char *fieldname; 132fb5bd1c2SPatrick Sanan ierr = DMDAGetFieldName(daCurr,f,&fieldname);CHKERRQ(ierr); 133fb5bd1c2SPatrick Sanan if (!fieldname) { 134fb5bd1c2SPatrick Sanan ierr = PetscSNPrintf(buf,sizeof(buf),"%D",f);CHKERRQ(ierr); 135fb5bd1c2SPatrick Sanan fieldname = buf; 136fb5bd1c2SPatrick Sanan } 137fb5bd1c2SPatrick Sanan ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr); 138fb5bd1c2SPatrick Sanan boffset += nnodes*sizeof(PetscScalar) + sizeof(int); 139fb5bd1c2SPatrick Sanan } 140fb5bd1c2SPatrick Sanan } else { 141ea2d7708SPatrick Sanan ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%D\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,bs,boffset);CHKERRQ(ierr); 142ea2d7708SPatrick Sanan boffset += bs*nnodes*sizeof(PetscScalar) + sizeof(int); 1434061b8bfSJed Brown } 144fb5bd1c2SPatrick Sanan } 1454061b8bfSJed Brown ierr = PetscFPrintf(comm,fp," </PointData>\n");CHKERRQ(ierr); 1464061b8bfSJed Brown ierr = PetscFPrintf(comm,fp," </Piece>\n");CHKERRQ(ierr); 1474061b8bfSJed Brown } 1484061b8bfSJed Brown ierr = PetscFPrintf(comm,fp," </StructuredGrid>\n");CHKERRQ(ierr); 1494061b8bfSJed Brown ierr = PetscFPrintf(comm,fp," <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr); 1504061b8bfSJed Brown ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr); 1514061b8bfSJed Brown 1524061b8bfSJed Brown /* Now write the arrays. */ 1534061b8bfSJed Brown tag = ((PetscObject)viewer)->tag; 154fb5bd1c2SPatrick Sanan ierr = PetscMalloc2(maxnnodes*PetscMax(3,maxbs),&array,maxnnodes*PetscMax(3,maxbs),&array2);CHKERRQ(ierr); 1554061b8bfSJed Brown for (r=0; r<size; r++) { 1564061b8bfSJed Brown MPI_Status status; 1574061b8bfSJed Brown PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0; 1584061b8bfSJed Brown if (!rank) { 1594061b8bfSJed Brown xs = grloc[r][0]; 1604061b8bfSJed Brown xm = grloc[r][1]; 1614061b8bfSJed Brown ys = grloc[r][2]; 1624061b8bfSJed Brown ym = grloc[r][3]; 1634061b8bfSJed Brown zs = grloc[r][4]; 1644061b8bfSJed Brown zm = grloc[r][5]; 1654061b8bfSJed Brown nnodes = xm*ym*zm; 1664061b8bfSJed Brown } else if (r == rank) { 1674061b8bfSJed Brown nnodes = info.xm*info.ym*info.zm; 1684061b8bfSJed Brown } 1694061b8bfSJed Brown 1702eaa9ef4SLisandro Dalcin /* Write the coordinates */ 1714061b8bfSJed Brown if (Coords) { 1724061b8bfSJed Brown const PetscScalar *coords; 1734061b8bfSJed Brown ierr = VecGetArrayRead(Coords,&coords);CHKERRQ(ierr); 1744061b8bfSJed Brown if (!rank) { 1754061b8bfSJed Brown if (r) { 1766622f924SJed Brown PetscMPIInt nn; 1772eaa9ef4SLisandro Dalcin ierr = MPI_Recv(array,nnodes*cdim,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr); 1784061b8bfSJed Brown ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr); 1792eaa9ef4SLisandro Dalcin if (nn != nnodes*cdim) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch"); 1804061b8bfSJed Brown } else { 1812eaa9ef4SLisandro Dalcin ierr = PetscMemcpy(array,coords,nnodes*cdim*sizeof(PetscScalar));CHKERRQ(ierr); 1824061b8bfSJed Brown } 1834061b8bfSJed Brown /* Transpose coordinates to VTK (C-style) ordering */ 1844061b8bfSJed Brown for (k=0; k<zm; k++) { 1854061b8bfSJed Brown for (j=0; j<ym; j++) { 1864061b8bfSJed Brown for (i=0; i<xm; i++) { 1874061b8bfSJed Brown PetscInt Iloc = i+xm*(j+ym*k); 1882eaa9ef4SLisandro Dalcin array2[Iloc*3+0] = array[Iloc*cdim + 0]; 1892eaa9ef4SLisandro Dalcin array2[Iloc*3+1] = cdim > 1 ? array[Iloc*cdim + 1] : 0.0; 1902eaa9ef4SLisandro Dalcin array2[Iloc*3+2] = cdim > 2 ? array[Iloc*cdim + 2] : 0.0; 1914061b8bfSJed Brown } 1924061b8bfSJed Brown } 1934061b8bfSJed Brown } 1944061b8bfSJed Brown } else if (r == rank) { 1952eaa9ef4SLisandro Dalcin ierr = MPI_Send((void*)coords,nnodes*cdim,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr); 1964061b8bfSJed Brown } 1974061b8bfSJed Brown ierr = VecRestoreArrayRead(Coords,&coords);CHKERRQ(ierr); 1984061b8bfSJed Brown } else { /* Fabricate some coordinates using grid index */ 1994061b8bfSJed Brown for (k=0; k<zm; k++) { 2004061b8bfSJed Brown for (j=0; j<ym; j++) { 2014061b8bfSJed Brown for (i=0; i<xm; i++) { 2024061b8bfSJed Brown PetscInt Iloc = i+xm*(j+ym*k); 2034061b8bfSJed Brown array2[Iloc*3+0] = xs+i; 2044061b8bfSJed Brown array2[Iloc*3+1] = ys+j; 2054061b8bfSJed Brown array2[Iloc*3+2] = zs+k; 2064061b8bfSJed Brown } 2074061b8bfSJed Brown } 2084061b8bfSJed Brown } 2094061b8bfSJed Brown } 21094fbd55eSBarry Smith ierr = PetscViewerVTKFWrite(viewer,fp,array2,nnodes*3,MPIU_SCALAR);CHKERRQ(ierr); 2114061b8bfSJed Brown 2124061b8bfSJed Brown /* Write each of the objects queued up for this file */ 2134061b8bfSJed Brown for (link=vtk->link; link; link=link->next) { 2144061b8bfSJed Brown Vec X = (Vec)link->vec; 2154061b8bfSJed Brown const PetscScalar *x; 216fb5bd1c2SPatrick Sanan PetscInt bs,f; 217ea2d7708SPatrick Sanan DM daCurr; 218fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 219ea2d7708SPatrick Sanan ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr); 220ea2d7708SPatrick Sanan ierr = DMDAGetInfo(daCurr,0,0,0,0, 0,0,0,&bs,0,0,0,0,0);CHKERRQ(ierr); 2214061b8bfSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 2224061b8bfSJed Brown if (!rank) { 2234061b8bfSJed Brown if (r) { 2246622f924SJed Brown PetscMPIInt nn; 2254061b8bfSJed Brown ierr = MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr); 2264061b8bfSJed Brown ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr); 2274061b8bfSJed Brown if (nn != nnodes*bs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %D",r); 2284061b8bfSJed Brown } else { 2294061b8bfSJed Brown ierr = PetscMemcpy(array,x,nnodes*bs*sizeof(PetscScalar));CHKERRQ(ierr); 2304061b8bfSJed Brown } 231fb5bd1c2SPatrick Sanan 232fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 233fb5bd1c2SPatrick Sanan ierr = DMDAGetFieldsNamed(daCurr,&fieldsnamed);CHKERRQ(ierr); 234fb5bd1c2SPatrick Sanan if (fieldsnamed) { 235fb5bd1c2SPatrick Sanan for (f=0; f<bs; f++) { 236fb5bd1c2SPatrick Sanan /* Extract and transpose the f'th field */ 237fb5bd1c2SPatrick Sanan for (k=0; k<zm; k++) { 238fb5bd1c2SPatrick Sanan for (j=0; j<ym; j++) { 239fb5bd1c2SPatrick Sanan for (i=0; i<xm; i++) { 240fb5bd1c2SPatrick Sanan PetscInt Iloc = i+xm*(j+ym*k); 241fb5bd1c2SPatrick Sanan array2[Iloc] = array[Iloc*bs + f]; 242fb5bd1c2SPatrick Sanan } 243fb5bd1c2SPatrick Sanan } 244fb5bd1c2SPatrick Sanan } 245fb5bd1c2SPatrick Sanan ierr = PetscViewerVTKFWrite(viewer,fp,array2,nnodes,MPIU_SCALAR);CHKERRQ(ierr); 246fb5bd1c2SPatrick Sanan } 247fb5bd1c2SPatrick Sanan } else { 248ea2d7708SPatrick Sanan ierr = PetscViewerVTKFWrite(viewer,fp,array,bs*nnodes,MPIU_SCALAR);CHKERRQ(ierr); 249fb5bd1c2SPatrick Sanan } 2504061b8bfSJed Brown } else if (r == rank) { 2514061b8bfSJed Brown ierr = MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr); 2524061b8bfSJed Brown } 2534061b8bfSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 2544061b8bfSJed Brown } 2554061b8bfSJed Brown } 2564061b8bfSJed Brown ierr = PetscFree2(array,array2);CHKERRQ(ierr); 2574061b8bfSJed Brown ierr = PetscFree(grloc);CHKERRQ(ierr); 2584061b8bfSJed Brown 2594061b8bfSJed Brown ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr); 2604061b8bfSJed Brown ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr); 2614061b8bfSJed Brown ierr = PetscFClose(comm,fp);CHKERRQ(ierr); 2624061b8bfSJed Brown PetscFunctionReturn(0); 2634061b8bfSJed Brown } 2644061b8bfSJed Brown 2654061b8bfSJed Brown 266a13bc4e3SShao-Ching Huang static PetscErrorCode DMDAVTKWriteAll_VTR(DM da,PetscViewer viewer) 267a13bc4e3SShao-Ching Huang { 268*30815ce0SLisandro Dalcin const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian"; 269a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_REAL_SINGLE) 270a13bc4e3SShao-Ching Huang const char precision[] = "Float32"; 271a13bc4e3SShao-Ching Huang #elif defined(PETSC_USE_REAL_DOUBLE) 272a13bc4e3SShao-Ching Huang const char precision[] = "Float64"; 273a13bc4e3SShao-Ching Huang #else 274a13bc4e3SShao-Ching Huang const char precision[] = "UnknownPrecision"; 275a13bc4e3SShao-Ching Huang #endif 276a13bc4e3SShao-Ching Huang MPI_Comm comm; 277a13bc4e3SShao-Ching Huang PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data; 278a13bc4e3SShao-Ching Huang PetscViewerVTKObjectLink link; 279a13bc4e3SShao-Ching Huang FILE *fp; 280a13bc4e3SShao-Ching Huang PetscMPIInt rank,size,tag; 281a13bc4e3SShao-Ching Huang DMDALocalInfo info; 282ea2d7708SPatrick Sanan PetscInt dim,mx,my,mz,boffset,maxnnodes,maxbs,i,j,k,r; 283a13bc4e3SShao-Ching Huang PetscInt rloc[6],(*grloc)[6] = NULL; 284a13bc4e3SShao-Ching Huang PetscScalar *array,*array2; 285a13bc4e3SShao-Ching Huang PetscReal gmin[3],gmax[3]; 286a13bc4e3SShao-Ching Huang PetscErrorCode ierr; 287a13bc4e3SShao-Ching Huang 288a13bc4e3SShao-Ching Huang PetscFunctionBegin; 289a13bc4e3SShao-Ching Huang ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 290a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_COMPLEX) 291a13bc4e3SShao-Ching Huang SETERRQ(comm,PETSC_ERR_SUP,"Complex values not supported"); 292a13bc4e3SShao-Ching Huang #endif 293a13bc4e3SShao-Ching Huang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 294a13bc4e3SShao-Ching Huang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 295ea2d7708SPatrick Sanan ierr = DMDAGetInfo(da,&dim,&mx,&my,&mz,0,0,0,0,0,0,0,0,0);CHKERRQ(ierr); 296a13bc4e3SShao-Ching Huang ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 297a13bc4e3SShao-Ching Huang ierr = DMDAGetBoundingBox(da,gmin,gmax);CHKERRQ(ierr); 298a13bc4e3SShao-Ching Huang ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr); 299a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr); 300*30815ce0SLisandro Dalcin ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"RectilinearGrid\" version=\"0.1\" byte_order=\"%s\">\n",byte_order);CHKERRQ(ierr); 301a13bc4e3SShao-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); 302a13bc4e3SShao-Ching Huang 303785e854fSJed Brown if (!rank) {ierr = PetscMalloc1(size*6,&grloc);CHKERRQ(ierr);} 304a13bc4e3SShao-Ching Huang rloc[0] = info.xs; 305a13bc4e3SShao-Ching Huang rloc[1] = info.xm; 306a13bc4e3SShao-Ching Huang rloc[2] = info.ys; 307a13bc4e3SShao-Ching Huang rloc[3] = info.ym; 308a13bc4e3SShao-Ching Huang rloc[4] = info.zs; 309a13bc4e3SShao-Ching Huang rloc[5] = info.zm; 310a13bc4e3SShao-Ching Huang ierr = MPI_Gather(rloc,6,MPIU_INT,&grloc[0][0],6,MPIU_INT,0,comm);CHKERRQ(ierr); 311a13bc4e3SShao-Ching Huang 312a13bc4e3SShao-Ching Huang /* Write XML header */ 313a13bc4e3SShao-Ching Huang maxnnodes = 0; /* Used for the temporary array size on rank 0 */ 314ea2d7708SPatrick Sanan maxbs = 0; /* Used for the temporary array size on rank 0 */ 315a13bc4e3SShao-Ching Huang boffset = 0; /* Offset into binary file */ 316a13bc4e3SShao-Ching Huang for (r=0; r<size; r++) { 317a13bc4e3SShao-Ching Huang PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0; 318a13bc4e3SShao-Ching Huang if (!rank) { 319a13bc4e3SShao-Ching Huang xs = grloc[r][0]; 320a13bc4e3SShao-Ching Huang xm = grloc[r][1]; 321a13bc4e3SShao-Ching Huang ys = grloc[r][2]; 322a13bc4e3SShao-Ching Huang ym = grloc[r][3]; 323a13bc4e3SShao-Ching Huang zs = grloc[r][4]; 324a13bc4e3SShao-Ching Huang zm = grloc[r][5]; 325a13bc4e3SShao-Ching Huang nnodes = xm*ym*zm; 326a13bc4e3SShao-Ching Huang } 327a13bc4e3SShao-Ching Huang maxnnodes = PetscMax(maxnnodes,nnodes); 328a13bc4e3SShao-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); 329a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," <Coordinates>\n");CHKERRQ(ierr); 330a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"Xcoord\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr); 331a13bc4e3SShao-Ching Huang boffset += xm*sizeof(PetscScalar) + sizeof(int); 332a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"Ycoord\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr); 333a13bc4e3SShao-Ching Huang boffset += ym*sizeof(PetscScalar) + sizeof(int); 334a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"Zcoord\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr); 335a13bc4e3SShao-Ching Huang boffset += zm*sizeof(PetscScalar) + sizeof(int); 336a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," </Coordinates>\n");CHKERRQ(ierr); 337a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," <PointData Scalars=\"ScalarPointData\">\n");CHKERRQ(ierr); 338a13bc4e3SShao-Ching Huang for (link=vtk->link; link; link=link->next) { 339a13bc4e3SShao-Ching Huang Vec X = (Vec)link->vec; 340fb5bd1c2SPatrick Sanan PetscInt bs,f; 341ea2d7708SPatrick Sanan DM daCurr; 342fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 343fb5bd1c2SPatrick Sanan const char *vecname = "Unnamed"; 344ea2d7708SPatrick Sanan 345ea2d7708SPatrick Sanan ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr); 346ea2d7708SPatrick Sanan ierr = DMDAGetInfo(daCurr,0,0,0,0,0,0,0,&bs,0,0,0,0,0);CHKERRQ(ierr); 347ea2d7708SPatrick Sanan maxbs = PetscMax(maxbs,bs); 348ea2d7708SPatrick Sanan if (((PetscObject)X)->name || link != vtk->link) { 349a13bc4e3SShao-Ching Huang ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr); 350a13bc4e3SShao-Ching Huang } 351ea2d7708SPatrick Sanan 352fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 353fb5bd1c2SPatrick Sanan ierr = DMDAGetFieldsNamed(daCurr,&fieldsnamed);CHKERRQ(ierr); 354fb5bd1c2SPatrick Sanan if (fieldsnamed) { 355fb5bd1c2SPatrick Sanan for (f=0; f<bs; f++) { 356fb5bd1c2SPatrick Sanan char buf[256]; 357fb5bd1c2SPatrick Sanan const char *fieldname; 358fb5bd1c2SPatrick Sanan ierr = DMDAGetFieldName(daCurr,f,&fieldname);CHKERRQ(ierr); 359fb5bd1c2SPatrick Sanan if (!fieldname) { 360fb5bd1c2SPatrick Sanan ierr = PetscSNPrintf(buf,sizeof(buf),"%D",f);CHKERRQ(ierr); 361fb5bd1c2SPatrick Sanan fieldname = buf; 362fb5bd1c2SPatrick Sanan } 363fb5bd1c2SPatrick Sanan ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr); 364fb5bd1c2SPatrick Sanan boffset += nnodes*sizeof(PetscScalar) + sizeof(int); 365fb5bd1c2SPatrick Sanan } 366fb5bd1c2SPatrick Sanan } else { 367ea2d7708SPatrick Sanan ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%D\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,bs,boffset);CHKERRQ(ierr); 368ea2d7708SPatrick Sanan boffset += bs*nnodes*sizeof(PetscScalar) + sizeof(int); 369a13bc4e3SShao-Ching Huang } 370fb5bd1c2SPatrick Sanan } 371a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," </PointData>\n");CHKERRQ(ierr); 372a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," </Piece>\n");CHKERRQ(ierr); 373a13bc4e3SShao-Ching Huang } 374a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," </RectilinearGrid>\n");CHKERRQ(ierr); 375a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr); 376a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr); 377a13bc4e3SShao-Ching Huang 378a13bc4e3SShao-Ching Huang /* Now write the arrays. */ 379a13bc4e3SShao-Ching Huang tag = ((PetscObject)viewer)->tag; 380fb5bd1c2SPatrick 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); 381a13bc4e3SShao-Ching Huang for (r=0; r<size; r++) { 382a13bc4e3SShao-Ching Huang MPI_Status status; 383a13bc4e3SShao-Ching Huang PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0; 384a13bc4e3SShao-Ching Huang if (!rank) { 385a13bc4e3SShao-Ching Huang xs = grloc[r][0]; 386a13bc4e3SShao-Ching Huang xm = grloc[r][1]; 387a13bc4e3SShao-Ching Huang ys = grloc[r][2]; 388a13bc4e3SShao-Ching Huang ym = grloc[r][3]; 389a13bc4e3SShao-Ching Huang zs = grloc[r][4]; 390a13bc4e3SShao-Ching Huang zm = grloc[r][5]; 391a13bc4e3SShao-Ching Huang nnodes = xm*ym*zm; 392a13bc4e3SShao-Ching Huang } else if (r == rank) { 393a13bc4e3SShao-Ching Huang nnodes = info.xm*info.ym*info.zm; 394a13bc4e3SShao-Ching Huang } 395a13bc4e3SShao-Ching Huang { /* Write the coordinates */ 396a13bc4e3SShao-Ching Huang Vec Coords; 397ea2d7708SPatrick Sanan 398a13bc4e3SShao-Ching Huang ierr = DMGetCoordinates(da,&Coords);CHKERRQ(ierr); 399a13bc4e3SShao-Ching Huang if (Coords) { 400a13bc4e3SShao-Ching Huang const PetscScalar *coords; 401a13bc4e3SShao-Ching Huang ierr = VecGetArrayRead(Coords,&coords);CHKERRQ(ierr); 402a13bc4e3SShao-Ching Huang if (!rank) { 403a13bc4e3SShao-Ching Huang if (r) { 404a13bc4e3SShao-Ching Huang PetscMPIInt nn; 405a13bc4e3SShao-Ching Huang ierr = MPI_Recv(array,xm+ym+zm,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 != xm+ym+zm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch"); 408a13bc4e3SShao-Ching Huang } else { 409a13bc4e3SShao-Ching Huang /* x coordinates */ 410a13bc4e3SShao-Ching Huang for (j=0, k=0, i=0; i<xm; i++) { 411a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 412a13bc4e3SShao-Ching Huang array[i] = coords[Iloc*dim + 0]; 413a13bc4e3SShao-Ching Huang } 414a13bc4e3SShao-Ching Huang /* y coordinates */ 415a13bc4e3SShao-Ching Huang for (i=0, k=0, j=0; j<ym; j++) { 416a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 417a13bc4e3SShao-Ching Huang array[j+xm] = dim > 1 ? coords[Iloc*dim + 1] : 0; 418a13bc4e3SShao-Ching Huang } 419a13bc4e3SShao-Ching Huang /* z coordinates */ 420a13bc4e3SShao-Ching Huang for (i=0, j=0, k=0; k<zm; k++) { 421a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 422a13bc4e3SShao-Ching Huang array[k+xm+ym] = dim > 2 ? coords[Iloc*dim + 2] : 0; 423a13bc4e3SShao-Ching Huang } 424a13bc4e3SShao-Ching Huang } 425a13bc4e3SShao-Ching Huang } else if (r == rank) { 426a13bc4e3SShao-Ching Huang xm = info.xm; 427a13bc4e3SShao-Ching Huang ym = info.ym; 428a13bc4e3SShao-Ching Huang zm = info.zm; 429a13bc4e3SShao-Ching Huang for (j=0, k=0, i=0; i<xm; i++) { 430a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 431a13bc4e3SShao-Ching Huang array2[i] = coords[Iloc*dim + 0]; 432a13bc4e3SShao-Ching Huang } 433a13bc4e3SShao-Ching Huang for (i=0, k=0, j=0; j<ym; j++) { 434a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 435a13bc4e3SShao-Ching Huang array2[j+xm] = dim > 1 ? coords[Iloc*dim + 1] : 0; 436a13bc4e3SShao-Ching Huang } 437a13bc4e3SShao-Ching Huang for (i=0, j=0, k=0; k<zm; k++) { 438a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 439a13bc4e3SShao-Ching Huang array2[k+xm+ym] = dim > 2 ? coords[Iloc*dim + 2] : 0; 440a13bc4e3SShao-Ching Huang } 441a13bc4e3SShao-Ching Huang ierr = MPI_Send((void*)array2,xm+ym+zm,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr); 442a13bc4e3SShao-Ching Huang } 443a13bc4e3SShao-Ching Huang ierr = VecRestoreArrayRead(Coords,&coords);CHKERRQ(ierr); 444a13bc4e3SShao-Ching Huang } else { /* Fabricate some coordinates using grid index */ 4453d3eaba7SBarry Smith for (i=0; i<xm; i++) { 446a13bc4e3SShao-Ching Huang array[i] = xs+i; 447a13bc4e3SShao-Ching Huang } 4483d3eaba7SBarry Smith for (j=0; j<ym; j++) { 449a13bc4e3SShao-Ching Huang array[j+xm] = ys+j; 450a13bc4e3SShao-Ching Huang } 4513d3eaba7SBarry Smith for (k=0; k<zm; k++) { 452a13bc4e3SShao-Ching Huang array[k+xm+ym] = zs+k; 453a13bc4e3SShao-Ching Huang } 454a13bc4e3SShao-Ching Huang } 455a13bc4e3SShao-Ching Huang if (!rank) { 45694fbd55eSBarry Smith ierr = PetscViewerVTKFWrite(viewer,fp,&(array[0]),xm,MPIU_SCALAR);CHKERRQ(ierr); 45794fbd55eSBarry Smith ierr = PetscViewerVTKFWrite(viewer,fp,&(array[xm]),ym,MPIU_SCALAR);CHKERRQ(ierr); 45894fbd55eSBarry Smith ierr = PetscViewerVTKFWrite(viewer,fp,&(array[xm+ym]),zm,MPIU_SCALAR);CHKERRQ(ierr); 459a13bc4e3SShao-Ching Huang } 460a13bc4e3SShao-Ching Huang } 461a13bc4e3SShao-Ching Huang 462a13bc4e3SShao-Ching Huang /* Write each of the objects queued up for this file */ 463a13bc4e3SShao-Ching Huang for (link=vtk->link; link; link=link->next) { 464a13bc4e3SShao-Ching Huang Vec X = (Vec)link->vec; 465a13bc4e3SShao-Ching Huang const PetscScalar *x; 466fb5bd1c2SPatrick Sanan PetscInt bs,f; 467ea2d7708SPatrick Sanan DM daCurr; 468fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 469ea2d7708SPatrick Sanan ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr); 470ea2d7708SPatrick Sanan ierr = DMDAGetInfo(daCurr,0,0,0,0,0,0,0,&bs,0,0,0,0,0);CHKERRQ(ierr); 471a13bc4e3SShao-Ching Huang 472a13bc4e3SShao-Ching Huang ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 473a13bc4e3SShao-Ching Huang if (!rank) { 474a13bc4e3SShao-Ching Huang if (r) { 475a13bc4e3SShao-Ching Huang PetscMPIInt nn; 476a13bc4e3SShao-Ching Huang ierr = MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr); 477a13bc4e3SShao-Ching Huang ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr); 478a13bc4e3SShao-Ching Huang if (nn != nnodes*bs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %D",r); 479a13bc4e3SShao-Ching Huang } else { 480a13bc4e3SShao-Ching Huang ierr = PetscMemcpy(array,x,nnodes*bs*sizeof(PetscScalar));CHKERRQ(ierr); 481a13bc4e3SShao-Ching Huang } 482fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 483fb5bd1c2SPatrick Sanan ierr = DMDAGetFieldsNamed(daCurr,&fieldsnamed);CHKERRQ(ierr); 484fb5bd1c2SPatrick Sanan if (fieldsnamed) { 485fb5bd1c2SPatrick Sanan for (f=0; f<bs; f++) { 486fb5bd1c2SPatrick Sanan /* Extract and transpose the f'th field */ 487fb5bd1c2SPatrick Sanan for (k=0; k<zm; k++) { 488fb5bd1c2SPatrick Sanan for (j=0; j<ym; j++) { 489fb5bd1c2SPatrick Sanan for (i=0; i<xm; i++) { 490fb5bd1c2SPatrick Sanan PetscInt Iloc = i+xm*(j+ym*k); 491fb5bd1c2SPatrick Sanan array2[Iloc] = array[Iloc*bs + f]; 492fb5bd1c2SPatrick Sanan } 493fb5bd1c2SPatrick Sanan } 494fb5bd1c2SPatrick Sanan } 495fb5bd1c2SPatrick Sanan ierr = PetscViewerVTKFWrite(viewer,fp,array2,nnodes,MPIU_SCALAR);CHKERRQ(ierr); 496fb5bd1c2SPatrick Sanan } 497fb5bd1c2SPatrick Sanan } 498ea2d7708SPatrick Sanan ierr = PetscViewerVTKFWrite(viewer,fp,array,nnodes*bs,MPIU_SCALAR);CHKERRQ(ierr); 499ea2d7708SPatrick Sanan 500a13bc4e3SShao-Ching Huang } else if (r == rank) { 501a13bc4e3SShao-Ching Huang ierr = MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr); 502a13bc4e3SShao-Ching Huang } 503a13bc4e3SShao-Ching Huang ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 504a13bc4e3SShao-Ching Huang } 505a13bc4e3SShao-Ching Huang } 506a13bc4e3SShao-Ching Huang ierr = PetscFree2(array,array2);CHKERRQ(ierr); 507a13bc4e3SShao-Ching Huang ierr = PetscFree(grloc);CHKERRQ(ierr); 508a13bc4e3SShao-Ching Huang 509a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr); 510a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr); 511a13bc4e3SShao-Ching Huang ierr = PetscFClose(comm,fp);CHKERRQ(ierr); 512a13bc4e3SShao-Ching Huang PetscFunctionReturn(0); 513a13bc4e3SShao-Ching Huang } 514a13bc4e3SShao-Ching Huang 5154061b8bfSJed Brown /*@C 5164061b8bfSJed Brown DMDAVTKWriteAll - Write a file containing all the fields that have been provided to the viewer 5174061b8bfSJed Brown 5184061b8bfSJed Brown Collective 5194061b8bfSJed Brown 5204061b8bfSJed Brown Input Arguments: 5214061b8bfSJed Brown odm - DM specifying the grid layout, passed as a PetscObject 5224061b8bfSJed Brown viewer - viewer of type VTK 5234061b8bfSJed Brown 5244061b8bfSJed Brown Level: developer 5254061b8bfSJed Brown 526fb5bd1c2SPatrick Sanan Notes: 5274061b8bfSJed Brown This function is a callback used by the VTK viewer to actually write the file. 5284061b8bfSJed 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. 5294061b8bfSJed Brown Instead, metadata for the entire file needs to be available up-front before you can start writing the file. 5304061b8bfSJed Brown 531fb5bd1c2SPatrick Sanan If any fields have been named (see e.g. DMDASetFieldName()), then individual scalar 532fb5bd1c2SPatrick Sanan fields are written. Otherwise, a single multi-dof (vector) field is written. 533fb5bd1c2SPatrick Sanan 5344061b8bfSJed Brown .seealso: PETSCVIEWERVTK 5354061b8bfSJed Brown @*/ 5364061b8bfSJed Brown PetscErrorCode DMDAVTKWriteAll(PetscObject odm,PetscViewer viewer) 5374061b8bfSJed Brown { 5384061b8bfSJed Brown DM dm = (DM)odm; 5394061b8bfSJed Brown PetscBool isvtk; 5404061b8bfSJed Brown PetscErrorCode ierr; 5414061b8bfSJed Brown 5424061b8bfSJed Brown PetscFunctionBegin; 5434061b8bfSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 5444061b8bfSJed Brown PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 545251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr); 546ce94432eSBarry Smith if (!isvtk) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_INCOMP,"Cannot use viewer type %s",((PetscObject)viewer)->type_name); 5474061b8bfSJed Brown switch (viewer->format) { 5484061b8bfSJed Brown case PETSC_VIEWER_VTK_VTS: 5494061b8bfSJed Brown ierr = DMDAVTKWriteAll_VTS(dm,viewer);CHKERRQ(ierr); 5504061b8bfSJed Brown break; 551a13bc4e3SShao-Ching Huang case PETSC_VIEWER_VTK_VTR: 552a13bc4e3SShao-Ching Huang ierr = DMDAVTKWriteAll_VTR(dm,viewer);CHKERRQ(ierr); 553a13bc4e3SShao-Ching Huang break; 554ce94432eSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No support for format '%s'",PetscViewerFormats[viewer->format]); 5554061b8bfSJed Brown } 5564061b8bfSJed Brown PetscFunctionReturn(0); 5574061b8bfSJed Brown } 558