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); 544061b8bfSJed Brown #if defined(PETSC_USE_COMPLEX) 553bf036e2SBarry Smith SETERRQ(comm,PETSC_ERR_SUP,"Complex values not supported"); 564061b8bfSJed Brown #endif 574061b8bfSJed Brown ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 584061b8bfSJed Brown ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 59*ea78f98cSLisandro Dalcin ierr = DMDAGetInfo(da,&dim,&mx,&my,&mz,NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 604061b8bfSJed Brown ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 612eaa9ef4SLisandro Dalcin ierr = DMGetCoordinates(da,&Coords);CHKERRQ(ierr); 622eaa9ef4SLisandro Dalcin if (Coords) { 632eaa9ef4SLisandro Dalcin PetscInt csize; 642eaa9ef4SLisandro Dalcin ierr = VecGetSize(Coords,&csize);CHKERRQ(ierr); 652eaa9ef4SLisandro Dalcin if (csize % (mx*my*mz)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Coordinate vector size mismatch"); 662eaa9ef4SLisandro Dalcin cdim = csize/(mx*my*mz); 672eaa9ef4SLisandro Dalcin if (cdim < dim || cdim > 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Coordinate vector size mismatch"); 682eaa9ef4SLisandro Dalcin } else { 692eaa9ef4SLisandro Dalcin cdim = dim; 702eaa9ef4SLisandro Dalcin } 714061b8bfSJed Brown 724061b8bfSJed Brown ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr); 734061b8bfSJed Brown ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr); 7430815ce0SLisandro Dalcin ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"%s\">\n",byte_order);CHKERRQ(ierr); 754061b8bfSJed 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); 764061b8bfSJed Brown 77785e854fSJed Brown if (!rank) {ierr = PetscMalloc1(size*6,&grloc);CHKERRQ(ierr);} 784061b8bfSJed Brown rloc[0] = info.xs; 794061b8bfSJed Brown rloc[1] = info.xm; 804061b8bfSJed Brown rloc[2] = info.ys; 814061b8bfSJed Brown rloc[3] = info.ym; 824061b8bfSJed Brown rloc[4] = info.zs; 834061b8bfSJed Brown rloc[5] = info.zm; 844061b8bfSJed Brown ierr = MPI_Gather(rloc,6,MPIU_INT,&grloc[0][0],6,MPIU_INT,0,comm);CHKERRQ(ierr); 854061b8bfSJed Brown 864061b8bfSJed Brown /* Write XML header */ 874061b8bfSJed Brown maxnnodes = 0; /* Used for the temporary array size on rank 0 */ 88ea2d7708SPatrick Sanan maxbs = 0; /* Used for the temporary array size on rank 0 */ 894061b8bfSJed Brown boffset = 0; /* Offset into binary file */ 904061b8bfSJed Brown for (r=0; r<size; r++) { 914061b8bfSJed Brown PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0; 924061b8bfSJed Brown if (!rank) { 934061b8bfSJed Brown xs = grloc[r][0]; 944061b8bfSJed Brown xm = grloc[r][1]; 954061b8bfSJed Brown ys = grloc[r][2]; 964061b8bfSJed Brown ym = grloc[r][3]; 974061b8bfSJed Brown zs = grloc[r][4]; 984061b8bfSJed Brown zm = grloc[r][5]; 994061b8bfSJed Brown nnodes = xm*ym*zm; 1004061b8bfSJed Brown } 1014061b8bfSJed Brown maxnnodes = PetscMax(maxnnodes,nnodes); 1024061b8bfSJed 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); 1034061b8bfSJed Brown ierr = PetscFPrintf(comm,fp," <Points>\n");CHKERRQ(ierr); 1044061b8bfSJed Brown ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr); 1054061b8bfSJed Brown boffset += 3*nnodes*sizeof(PetscScalar) + sizeof(int); 1064061b8bfSJed Brown ierr = PetscFPrintf(comm,fp," </Points>\n");CHKERRQ(ierr); 1074061b8bfSJed Brown 1084061b8bfSJed Brown ierr = PetscFPrintf(comm,fp," <PointData Scalars=\"ScalarPointData\">\n");CHKERRQ(ierr); 1094061b8bfSJed Brown for (link=vtk->link; link; link=link->next) { 1104061b8bfSJed Brown Vec X = (Vec)link->vec; 111fb5bd1c2SPatrick Sanan PetscInt bs,f; 112ea2d7708SPatrick Sanan DM daCurr; 113fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 114fb5bd1c2SPatrick Sanan const char *vecname = "Unnamed"; 115ea2d7708SPatrick Sanan 116ea2d7708SPatrick Sanan ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr); 117*ea78f98cSLisandro Dalcin ierr = DMDAGetInfo(daCurr,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 118ea2d7708SPatrick Sanan maxbs = PetscMax(maxbs,bs); 119ea2d7708SPatrick Sanan 120ea2d7708SPatrick Sanan if (((PetscObject)X)->name || link != vtk->link) { 1214061b8bfSJed Brown ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr); 1224061b8bfSJed Brown } 123fb5bd1c2SPatrick Sanan 124fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 125fb5bd1c2SPatrick Sanan ierr = DMDAGetFieldsNamed(daCurr,&fieldsnamed);CHKERRQ(ierr); 126fb5bd1c2SPatrick Sanan if (fieldsnamed) { 127fb5bd1c2SPatrick Sanan for (f=0; f<bs; f++) { 128fb5bd1c2SPatrick Sanan char buf[256]; 129fb5bd1c2SPatrick Sanan const char *fieldname; 130fb5bd1c2SPatrick Sanan ierr = DMDAGetFieldName(daCurr,f,&fieldname);CHKERRQ(ierr); 131fb5bd1c2SPatrick Sanan if (!fieldname) { 132fb5bd1c2SPatrick Sanan ierr = PetscSNPrintf(buf,sizeof(buf),"%D",f);CHKERRQ(ierr); 133fb5bd1c2SPatrick Sanan fieldname = buf; 134fb5bd1c2SPatrick Sanan } 135fb5bd1c2SPatrick Sanan ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr); 136fb5bd1c2SPatrick Sanan boffset += nnodes*sizeof(PetscScalar) + sizeof(int); 137fb5bd1c2SPatrick Sanan } 138fb5bd1c2SPatrick Sanan } else { 139ea2d7708SPatrick Sanan ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%D\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,bs,boffset);CHKERRQ(ierr); 140ea2d7708SPatrick Sanan boffset += bs*nnodes*sizeof(PetscScalar) + sizeof(int); 1414061b8bfSJed Brown } 142fb5bd1c2SPatrick Sanan } 1434061b8bfSJed Brown ierr = PetscFPrintf(comm,fp," </PointData>\n");CHKERRQ(ierr); 1444061b8bfSJed Brown ierr = PetscFPrintf(comm,fp," </Piece>\n");CHKERRQ(ierr); 1454061b8bfSJed Brown } 1464061b8bfSJed Brown ierr = PetscFPrintf(comm,fp," </StructuredGrid>\n");CHKERRQ(ierr); 1474061b8bfSJed Brown ierr = PetscFPrintf(comm,fp," <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr); 1484061b8bfSJed Brown ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr); 1494061b8bfSJed Brown 1504061b8bfSJed Brown /* Now write the arrays. */ 1514061b8bfSJed Brown tag = ((PetscObject)viewer)->tag; 152fb5bd1c2SPatrick Sanan ierr = PetscMalloc2(maxnnodes*PetscMax(3,maxbs),&array,maxnnodes*PetscMax(3,maxbs),&array2);CHKERRQ(ierr); 1534061b8bfSJed Brown for (r=0; r<size; r++) { 1544061b8bfSJed Brown MPI_Status status; 1554061b8bfSJed Brown PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0; 1564061b8bfSJed Brown if (!rank) { 1574061b8bfSJed Brown xs = grloc[r][0]; 1584061b8bfSJed Brown xm = grloc[r][1]; 1594061b8bfSJed Brown ys = grloc[r][2]; 1604061b8bfSJed Brown ym = grloc[r][3]; 1614061b8bfSJed Brown zs = grloc[r][4]; 1624061b8bfSJed Brown zm = grloc[r][5]; 1634061b8bfSJed Brown nnodes = xm*ym*zm; 1644061b8bfSJed Brown } else if (r == rank) { 1654061b8bfSJed Brown nnodes = info.xm*info.ym*info.zm; 1664061b8bfSJed Brown } 1674061b8bfSJed Brown 1682eaa9ef4SLisandro Dalcin /* Write the coordinates */ 1694061b8bfSJed Brown if (Coords) { 1704061b8bfSJed Brown const PetscScalar *coords; 1714061b8bfSJed Brown ierr = VecGetArrayRead(Coords,&coords);CHKERRQ(ierr); 1724061b8bfSJed Brown if (!rank) { 1734061b8bfSJed Brown if (r) { 1746622f924SJed Brown PetscMPIInt nn; 1752eaa9ef4SLisandro Dalcin ierr = MPI_Recv(array,nnodes*cdim,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr); 1764061b8bfSJed Brown ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr); 1772eaa9ef4SLisandro Dalcin if (nn != nnodes*cdim) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch"); 1784061b8bfSJed Brown } else { 179580bdb30SBarry Smith ierr = PetscArraycpy(array,coords,nnodes*cdim);CHKERRQ(ierr); 1804061b8bfSJed Brown } 1814061b8bfSJed Brown /* Transpose coordinates to VTK (C-style) ordering */ 1824061b8bfSJed Brown for (k=0; k<zm; k++) { 1834061b8bfSJed Brown for (j=0; j<ym; j++) { 1844061b8bfSJed Brown for (i=0; i<xm; i++) { 1854061b8bfSJed Brown PetscInt Iloc = i+xm*(j+ym*k); 1862eaa9ef4SLisandro Dalcin array2[Iloc*3+0] = array[Iloc*cdim + 0]; 1872eaa9ef4SLisandro Dalcin array2[Iloc*3+1] = cdim > 1 ? array[Iloc*cdim + 1] : 0.0; 1882eaa9ef4SLisandro Dalcin array2[Iloc*3+2] = cdim > 2 ? array[Iloc*cdim + 2] : 0.0; 1894061b8bfSJed Brown } 1904061b8bfSJed Brown } 1914061b8bfSJed Brown } 1924061b8bfSJed Brown } else if (r == rank) { 1932eaa9ef4SLisandro Dalcin ierr = MPI_Send((void*)coords,nnodes*cdim,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr); 1944061b8bfSJed Brown } 1954061b8bfSJed Brown ierr = VecRestoreArrayRead(Coords,&coords);CHKERRQ(ierr); 1964061b8bfSJed Brown } else { /* Fabricate some coordinates using grid index */ 1974061b8bfSJed Brown for (k=0; k<zm; k++) { 1984061b8bfSJed Brown for (j=0; j<ym; j++) { 1994061b8bfSJed Brown for (i=0; i<xm; i++) { 2004061b8bfSJed Brown PetscInt Iloc = i+xm*(j+ym*k); 2014061b8bfSJed Brown array2[Iloc*3+0] = xs+i; 2024061b8bfSJed Brown array2[Iloc*3+1] = ys+j; 2034061b8bfSJed Brown array2[Iloc*3+2] = zs+k; 2044061b8bfSJed Brown } 2054061b8bfSJed Brown } 2064061b8bfSJed Brown } 2074061b8bfSJed Brown } 20894fbd55eSBarry Smith ierr = PetscViewerVTKFWrite(viewer,fp,array2,nnodes*3,MPIU_SCALAR);CHKERRQ(ierr); 2094061b8bfSJed Brown 2104061b8bfSJed Brown /* Write each of the objects queued up for this file */ 2114061b8bfSJed Brown for (link=vtk->link; link; link=link->next) { 2124061b8bfSJed Brown Vec X = (Vec)link->vec; 2134061b8bfSJed Brown const PetscScalar *x; 214fb5bd1c2SPatrick Sanan PetscInt bs,f; 215ea2d7708SPatrick Sanan DM daCurr; 216fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 217ea2d7708SPatrick Sanan ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr); 218*ea78f98cSLisandro Dalcin ierr = DMDAGetInfo(daCurr,NULL,NULL,NULL,NULL, NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 2194061b8bfSJed Brown ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 2204061b8bfSJed Brown if (!rank) { 2214061b8bfSJed Brown if (r) { 2226622f924SJed Brown PetscMPIInt nn; 2234061b8bfSJed Brown ierr = MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr); 2244061b8bfSJed Brown ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr); 2254061b8bfSJed Brown if (nn != nnodes*bs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %D",r); 2264061b8bfSJed Brown } else { 227580bdb30SBarry Smith ierr = PetscArraycpy(array,x,nnodes*bs);CHKERRQ(ierr); 2284061b8bfSJed Brown } 229fb5bd1c2SPatrick Sanan 230fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 231fb5bd1c2SPatrick Sanan ierr = DMDAGetFieldsNamed(daCurr,&fieldsnamed);CHKERRQ(ierr); 232fb5bd1c2SPatrick Sanan if (fieldsnamed) { 233fb5bd1c2SPatrick Sanan for (f=0; f<bs; f++) { 234fb5bd1c2SPatrick Sanan /* Extract and transpose the f'th field */ 235fb5bd1c2SPatrick Sanan for (k=0; k<zm; k++) { 236fb5bd1c2SPatrick Sanan for (j=0; j<ym; j++) { 237fb5bd1c2SPatrick Sanan for (i=0; i<xm; i++) { 238fb5bd1c2SPatrick Sanan PetscInt Iloc = i+xm*(j+ym*k); 239fb5bd1c2SPatrick Sanan array2[Iloc] = array[Iloc*bs + f]; 240fb5bd1c2SPatrick Sanan } 241fb5bd1c2SPatrick Sanan } 242fb5bd1c2SPatrick Sanan } 243fb5bd1c2SPatrick Sanan ierr = PetscViewerVTKFWrite(viewer,fp,array2,nnodes,MPIU_SCALAR);CHKERRQ(ierr); 244fb5bd1c2SPatrick Sanan } 245fb5bd1c2SPatrick Sanan } else { 246ea2d7708SPatrick Sanan ierr = PetscViewerVTKFWrite(viewer,fp,array,bs*nnodes,MPIU_SCALAR);CHKERRQ(ierr); 247fb5bd1c2SPatrick Sanan } 2484061b8bfSJed Brown } else if (r == rank) { 2494061b8bfSJed Brown ierr = MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr); 2504061b8bfSJed Brown } 2514061b8bfSJed Brown ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 2524061b8bfSJed Brown } 2534061b8bfSJed Brown } 2544061b8bfSJed Brown ierr = PetscFree2(array,array2);CHKERRQ(ierr); 2554061b8bfSJed Brown ierr = PetscFree(grloc);CHKERRQ(ierr); 2564061b8bfSJed Brown 2574061b8bfSJed Brown ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr); 2584061b8bfSJed Brown ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr); 2594061b8bfSJed Brown ierr = PetscFClose(comm,fp);CHKERRQ(ierr); 2604061b8bfSJed Brown PetscFunctionReturn(0); 2614061b8bfSJed Brown } 2624061b8bfSJed Brown 2634061b8bfSJed Brown 264a13bc4e3SShao-Ching Huang static PetscErrorCode DMDAVTKWriteAll_VTR(DM da,PetscViewer viewer) 265a13bc4e3SShao-Ching Huang { 26630815ce0SLisandro Dalcin const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian"; 267a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_REAL_SINGLE) 268a13bc4e3SShao-Ching Huang const char precision[] = "Float32"; 269a13bc4e3SShao-Ching Huang #elif defined(PETSC_USE_REAL_DOUBLE) 270a13bc4e3SShao-Ching Huang const char precision[] = "Float64"; 271a13bc4e3SShao-Ching Huang #else 272a13bc4e3SShao-Ching Huang const char precision[] = "UnknownPrecision"; 273a13bc4e3SShao-Ching Huang #endif 274a13bc4e3SShao-Ching Huang MPI_Comm comm; 275a13bc4e3SShao-Ching Huang PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data; 276a13bc4e3SShao-Ching Huang PetscViewerVTKObjectLink link; 277a13bc4e3SShao-Ching Huang FILE *fp; 278a13bc4e3SShao-Ching Huang PetscMPIInt rank,size,tag; 279a13bc4e3SShao-Ching Huang DMDALocalInfo info; 280ea2d7708SPatrick Sanan PetscInt dim,mx,my,mz,boffset,maxnnodes,maxbs,i,j,k,r; 281a13bc4e3SShao-Ching Huang PetscInt rloc[6],(*grloc)[6] = NULL; 282a13bc4e3SShao-Ching Huang PetscScalar *array,*array2; 283a13bc4e3SShao-Ching Huang PetscErrorCode ierr; 284a13bc4e3SShao-Ching Huang 285a13bc4e3SShao-Ching Huang PetscFunctionBegin; 286a13bc4e3SShao-Ching Huang ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr); 287a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_COMPLEX) 288a13bc4e3SShao-Ching Huang SETERRQ(comm,PETSC_ERR_SUP,"Complex values not supported"); 289a13bc4e3SShao-Ching Huang #endif 290a13bc4e3SShao-Ching Huang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 291a13bc4e3SShao-Ching Huang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 292*ea78f98cSLisandro Dalcin ierr = DMDAGetInfo(da,&dim,&mx,&my,&mz,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 293a13bc4e3SShao-Ching Huang ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr); 294a13bc4e3SShao-Ching Huang ierr = PetscFOpen(comm,vtk->filename,"wb",&fp);CHKERRQ(ierr); 295a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n");CHKERRQ(ierr); 29630815ce0SLisandro Dalcin ierr = PetscFPrintf(comm,fp,"<VTKFile type=\"RectilinearGrid\" version=\"0.1\" byte_order=\"%s\">\n",byte_order);CHKERRQ(ierr); 297a13bc4e3SShao-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); 298a13bc4e3SShao-Ching Huang 299785e854fSJed Brown if (!rank) {ierr = PetscMalloc1(size*6,&grloc);CHKERRQ(ierr);} 300a13bc4e3SShao-Ching Huang rloc[0] = info.xs; 301a13bc4e3SShao-Ching Huang rloc[1] = info.xm; 302a13bc4e3SShao-Ching Huang rloc[2] = info.ys; 303a13bc4e3SShao-Ching Huang rloc[3] = info.ym; 304a13bc4e3SShao-Ching Huang rloc[4] = info.zs; 305a13bc4e3SShao-Ching Huang rloc[5] = info.zm; 306a13bc4e3SShao-Ching Huang ierr = MPI_Gather(rloc,6,MPIU_INT,&grloc[0][0],6,MPIU_INT,0,comm);CHKERRQ(ierr); 307a13bc4e3SShao-Ching Huang 308a13bc4e3SShao-Ching Huang /* Write XML header */ 309a13bc4e3SShao-Ching Huang maxnnodes = 0; /* Used for the temporary array size on rank 0 */ 310ea2d7708SPatrick Sanan maxbs = 0; /* Used for the temporary array size on rank 0 */ 311a13bc4e3SShao-Ching Huang boffset = 0; /* Offset into binary file */ 312a13bc4e3SShao-Ching Huang for (r=0; r<size; r++) { 313a13bc4e3SShao-Ching Huang PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0; 314a13bc4e3SShao-Ching Huang if (!rank) { 315a13bc4e3SShao-Ching Huang xs = grloc[r][0]; 316a13bc4e3SShao-Ching Huang xm = grloc[r][1]; 317a13bc4e3SShao-Ching Huang ys = grloc[r][2]; 318a13bc4e3SShao-Ching Huang ym = grloc[r][3]; 319a13bc4e3SShao-Ching Huang zs = grloc[r][4]; 320a13bc4e3SShao-Ching Huang zm = grloc[r][5]; 321a13bc4e3SShao-Ching Huang nnodes = xm*ym*zm; 322a13bc4e3SShao-Ching Huang } 323a13bc4e3SShao-Ching Huang maxnnodes = PetscMax(maxnnodes,nnodes); 324a13bc4e3SShao-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); 325a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," <Coordinates>\n");CHKERRQ(ierr); 326a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"Xcoord\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr); 327a13bc4e3SShao-Ching Huang boffset += xm*sizeof(PetscScalar) + sizeof(int); 328a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"Ycoord\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr); 329a13bc4e3SShao-Ching Huang boffset += ym*sizeof(PetscScalar) + sizeof(int); 330a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"Zcoord\" format=\"appended\" offset=\"%D\" />\n",precision,boffset);CHKERRQ(ierr); 331a13bc4e3SShao-Ching Huang boffset += zm*sizeof(PetscScalar) + sizeof(int); 332a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," </Coordinates>\n");CHKERRQ(ierr); 333a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," <PointData Scalars=\"ScalarPointData\">\n");CHKERRQ(ierr); 334a13bc4e3SShao-Ching Huang for (link=vtk->link; link; link=link->next) { 335a13bc4e3SShao-Ching Huang Vec X = (Vec)link->vec; 336fb5bd1c2SPatrick Sanan PetscInt bs,f; 337ea2d7708SPatrick Sanan DM daCurr; 338fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 339fb5bd1c2SPatrick Sanan const char *vecname = "Unnamed"; 340ea2d7708SPatrick Sanan 341ea2d7708SPatrick Sanan ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr); 342*ea78f98cSLisandro Dalcin ierr = DMDAGetInfo(daCurr,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 343ea2d7708SPatrick Sanan maxbs = PetscMax(maxbs,bs); 344ea2d7708SPatrick Sanan if (((PetscObject)X)->name || link != vtk->link) { 345a13bc4e3SShao-Ching Huang ierr = PetscObjectGetName((PetscObject)X,&vecname);CHKERRQ(ierr); 346a13bc4e3SShao-Ching Huang } 347ea2d7708SPatrick Sanan 348fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 349fb5bd1c2SPatrick Sanan ierr = DMDAGetFieldsNamed(daCurr,&fieldsnamed);CHKERRQ(ierr); 350fb5bd1c2SPatrick Sanan if (fieldsnamed) { 351fb5bd1c2SPatrick Sanan for (f=0; f<bs; f++) { 352fb5bd1c2SPatrick Sanan char buf[256]; 353fb5bd1c2SPatrick Sanan const char *fieldname; 354fb5bd1c2SPatrick Sanan ierr = DMDAGetFieldName(daCurr,f,&fieldname);CHKERRQ(ierr); 355fb5bd1c2SPatrick Sanan if (!fieldname) { 356fb5bd1c2SPatrick Sanan ierr = PetscSNPrintf(buf,sizeof(buf),"%D",f);CHKERRQ(ierr); 357fb5bd1c2SPatrick Sanan fieldname = buf; 358fb5bd1c2SPatrick Sanan } 359fb5bd1c2SPatrick Sanan ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,fieldname,boffset);CHKERRQ(ierr); 360fb5bd1c2SPatrick Sanan boffset += nnodes*sizeof(PetscScalar) + sizeof(int); 361fb5bd1c2SPatrick Sanan } 362fb5bd1c2SPatrick Sanan } else { 363ea2d7708SPatrick Sanan ierr = PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%D\" format=\"appended\" offset=\"%D\" />\n",precision,vecname,bs,boffset);CHKERRQ(ierr); 364ea2d7708SPatrick Sanan boffset += bs*nnodes*sizeof(PetscScalar) + sizeof(int); 365a13bc4e3SShao-Ching Huang } 366fb5bd1c2SPatrick Sanan } 367a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," </PointData>\n");CHKERRQ(ierr); 368a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," </Piece>\n");CHKERRQ(ierr); 369a13bc4e3SShao-Ching Huang } 370a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," </RectilinearGrid>\n");CHKERRQ(ierr); 371a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp," <AppendedData encoding=\"raw\">\n");CHKERRQ(ierr); 372a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp,"_");CHKERRQ(ierr); 373a13bc4e3SShao-Ching Huang 374a13bc4e3SShao-Ching Huang /* Now write the arrays. */ 375a13bc4e3SShao-Ching Huang tag = ((PetscObject)viewer)->tag; 376fb5bd1c2SPatrick 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); 377a13bc4e3SShao-Ching Huang for (r=0; r<size; r++) { 378a13bc4e3SShao-Ching Huang MPI_Status status; 379a13bc4e3SShao-Ching Huang PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0; 380a13bc4e3SShao-Ching Huang if (!rank) { 381a13bc4e3SShao-Ching Huang xs = grloc[r][0]; 382a13bc4e3SShao-Ching Huang xm = grloc[r][1]; 383a13bc4e3SShao-Ching Huang ys = grloc[r][2]; 384a13bc4e3SShao-Ching Huang ym = grloc[r][3]; 385a13bc4e3SShao-Ching Huang zs = grloc[r][4]; 386a13bc4e3SShao-Ching Huang zm = grloc[r][5]; 387a13bc4e3SShao-Ching Huang nnodes = xm*ym*zm; 388a13bc4e3SShao-Ching Huang } else if (r == rank) { 389a13bc4e3SShao-Ching Huang nnodes = info.xm*info.ym*info.zm; 390a13bc4e3SShao-Ching Huang } 391a13bc4e3SShao-Ching Huang { /* Write the coordinates */ 392a13bc4e3SShao-Ching Huang Vec Coords; 393ea2d7708SPatrick Sanan 394a13bc4e3SShao-Ching Huang ierr = DMGetCoordinates(da,&Coords);CHKERRQ(ierr); 395a13bc4e3SShao-Ching Huang if (Coords) { 396a13bc4e3SShao-Ching Huang const PetscScalar *coords; 397a13bc4e3SShao-Ching Huang ierr = VecGetArrayRead(Coords,&coords);CHKERRQ(ierr); 398a13bc4e3SShao-Ching Huang if (!rank) { 399a13bc4e3SShao-Ching Huang if (r) { 400a13bc4e3SShao-Ching Huang PetscMPIInt nn; 401a13bc4e3SShao-Ching Huang ierr = MPI_Recv(array,xm+ym+zm,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr); 402a13bc4e3SShao-Ching Huang ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr); 403a13bc4e3SShao-Ching Huang if (nn != xm+ym+zm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch"); 404a13bc4e3SShao-Ching Huang } else { 405a13bc4e3SShao-Ching Huang /* x coordinates */ 406a13bc4e3SShao-Ching Huang for (j=0, k=0, i=0; i<xm; i++) { 407a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 408a13bc4e3SShao-Ching Huang array[i] = coords[Iloc*dim + 0]; 409a13bc4e3SShao-Ching Huang } 410a13bc4e3SShao-Ching Huang /* y coordinates */ 411a13bc4e3SShao-Ching Huang for (i=0, k=0, j=0; j<ym; j++) { 412a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 413a13bc4e3SShao-Ching Huang array[j+xm] = dim > 1 ? coords[Iloc*dim + 1] : 0; 414a13bc4e3SShao-Ching Huang } 415a13bc4e3SShao-Ching Huang /* z coordinates */ 416a13bc4e3SShao-Ching Huang for (i=0, j=0, k=0; k<zm; k++) { 417a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 418a13bc4e3SShao-Ching Huang array[k+xm+ym] = dim > 2 ? coords[Iloc*dim + 2] : 0; 419a13bc4e3SShao-Ching Huang } 420a13bc4e3SShao-Ching Huang } 421a13bc4e3SShao-Ching Huang } else if (r == rank) { 422a13bc4e3SShao-Ching Huang xm = info.xm; 423a13bc4e3SShao-Ching Huang ym = info.ym; 424a13bc4e3SShao-Ching Huang zm = info.zm; 425a13bc4e3SShao-Ching Huang for (j=0, k=0, i=0; i<xm; i++) { 426a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 427a13bc4e3SShao-Ching Huang array2[i] = coords[Iloc*dim + 0]; 428a13bc4e3SShao-Ching Huang } 429a13bc4e3SShao-Ching Huang for (i=0, k=0, j=0; j<ym; j++) { 430a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 431a13bc4e3SShao-Ching Huang array2[j+xm] = dim > 1 ? coords[Iloc*dim + 1] : 0; 432a13bc4e3SShao-Ching Huang } 433a13bc4e3SShao-Ching Huang for (i=0, j=0, k=0; k<zm; k++) { 434a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 435a13bc4e3SShao-Ching Huang array2[k+xm+ym] = dim > 2 ? coords[Iloc*dim + 2] : 0; 436a13bc4e3SShao-Ching Huang } 437a13bc4e3SShao-Ching Huang ierr = MPI_Send((void*)array2,xm+ym+zm,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr); 438a13bc4e3SShao-Ching Huang } 439a13bc4e3SShao-Ching Huang ierr = VecRestoreArrayRead(Coords,&coords);CHKERRQ(ierr); 440a13bc4e3SShao-Ching Huang } else { /* Fabricate some coordinates using grid index */ 4413d3eaba7SBarry Smith for (i=0; i<xm; i++) { 442a13bc4e3SShao-Ching Huang array[i] = xs+i; 443a13bc4e3SShao-Ching Huang } 4443d3eaba7SBarry Smith for (j=0; j<ym; j++) { 445a13bc4e3SShao-Ching Huang array[j+xm] = ys+j; 446a13bc4e3SShao-Ching Huang } 4473d3eaba7SBarry Smith for (k=0; k<zm; k++) { 448a13bc4e3SShao-Ching Huang array[k+xm+ym] = zs+k; 449a13bc4e3SShao-Ching Huang } 450a13bc4e3SShao-Ching Huang } 451a13bc4e3SShao-Ching Huang if (!rank) { 45294fbd55eSBarry Smith ierr = PetscViewerVTKFWrite(viewer,fp,&(array[0]),xm,MPIU_SCALAR);CHKERRQ(ierr); 45394fbd55eSBarry Smith ierr = PetscViewerVTKFWrite(viewer,fp,&(array[xm]),ym,MPIU_SCALAR);CHKERRQ(ierr); 45494fbd55eSBarry Smith ierr = PetscViewerVTKFWrite(viewer,fp,&(array[xm+ym]),zm,MPIU_SCALAR);CHKERRQ(ierr); 455a13bc4e3SShao-Ching Huang } 456a13bc4e3SShao-Ching Huang } 457a13bc4e3SShao-Ching Huang 458a13bc4e3SShao-Ching Huang /* Write each of the objects queued up for this file */ 459a13bc4e3SShao-Ching Huang for (link=vtk->link; link; link=link->next) { 460a13bc4e3SShao-Ching Huang Vec X = (Vec)link->vec; 461a13bc4e3SShao-Ching Huang const PetscScalar *x; 462fb5bd1c2SPatrick Sanan PetscInt bs,f; 463ea2d7708SPatrick Sanan DM daCurr; 464fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 465ea2d7708SPatrick Sanan ierr = VecGetDM(X,&daCurr);CHKERRQ(ierr); 466*ea78f98cSLisandro Dalcin ierr = DMDAGetInfo(daCurr,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 467a13bc4e3SShao-Ching Huang 468a13bc4e3SShao-Ching Huang ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 469a13bc4e3SShao-Ching Huang if (!rank) { 470a13bc4e3SShao-Ching Huang if (r) { 471a13bc4e3SShao-Ching Huang PetscMPIInt nn; 472a13bc4e3SShao-Ching Huang ierr = MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status);CHKERRQ(ierr); 473a13bc4e3SShao-Ching Huang ierr = MPI_Get_count(&status,MPIU_SCALAR,&nn);CHKERRQ(ierr); 474a13bc4e3SShao-Ching Huang if (nn != nnodes*bs) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %D",r); 475a13bc4e3SShao-Ching Huang } else { 476580bdb30SBarry Smith ierr = PetscArraycpy(array,x,nnodes*bs);CHKERRQ(ierr); 477a13bc4e3SShao-Ching Huang } 478fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 479fb5bd1c2SPatrick Sanan ierr = DMDAGetFieldsNamed(daCurr,&fieldsnamed);CHKERRQ(ierr); 480fb5bd1c2SPatrick Sanan if (fieldsnamed) { 481fb5bd1c2SPatrick Sanan for (f=0; f<bs; f++) { 482fb5bd1c2SPatrick Sanan /* Extract and transpose the f'th field */ 483fb5bd1c2SPatrick Sanan for (k=0; k<zm; k++) { 484fb5bd1c2SPatrick Sanan for (j=0; j<ym; j++) { 485fb5bd1c2SPatrick Sanan for (i=0; i<xm; i++) { 486fb5bd1c2SPatrick Sanan PetscInt Iloc = i+xm*(j+ym*k); 487fb5bd1c2SPatrick Sanan array2[Iloc] = array[Iloc*bs + f]; 488fb5bd1c2SPatrick Sanan } 489fb5bd1c2SPatrick Sanan } 490fb5bd1c2SPatrick Sanan } 491fb5bd1c2SPatrick Sanan ierr = PetscViewerVTKFWrite(viewer,fp,array2,nnodes,MPIU_SCALAR);CHKERRQ(ierr); 492fb5bd1c2SPatrick Sanan } 493fb5bd1c2SPatrick Sanan } 494ea2d7708SPatrick Sanan ierr = PetscViewerVTKFWrite(viewer,fp,array,nnodes*bs,MPIU_SCALAR);CHKERRQ(ierr); 495ea2d7708SPatrick Sanan 496a13bc4e3SShao-Ching Huang } else if (r == rank) { 497a13bc4e3SShao-Ching Huang ierr = MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm);CHKERRQ(ierr); 498a13bc4e3SShao-Ching Huang } 499a13bc4e3SShao-Ching Huang ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 500a13bc4e3SShao-Ching Huang } 501a13bc4e3SShao-Ching Huang } 502a13bc4e3SShao-Ching Huang ierr = PetscFree2(array,array2);CHKERRQ(ierr); 503a13bc4e3SShao-Ching Huang ierr = PetscFree(grloc);CHKERRQ(ierr); 504a13bc4e3SShao-Ching Huang 505a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp,"\n </AppendedData>\n");CHKERRQ(ierr); 506a13bc4e3SShao-Ching Huang ierr = PetscFPrintf(comm,fp,"</VTKFile>\n");CHKERRQ(ierr); 507a13bc4e3SShao-Ching Huang ierr = PetscFClose(comm,fp);CHKERRQ(ierr); 508a13bc4e3SShao-Ching Huang PetscFunctionReturn(0); 509a13bc4e3SShao-Ching Huang } 510a13bc4e3SShao-Ching Huang 5114061b8bfSJed Brown /*@C 5124061b8bfSJed Brown DMDAVTKWriteAll - Write a file containing all the fields that have been provided to the viewer 5134061b8bfSJed Brown 5144061b8bfSJed Brown Collective 5154061b8bfSJed Brown 5164061b8bfSJed Brown Input Arguments: 5174061b8bfSJed Brown odm - DM specifying the grid layout, passed as a PetscObject 5184061b8bfSJed Brown viewer - viewer of type VTK 5194061b8bfSJed Brown 5204061b8bfSJed Brown Level: developer 5214061b8bfSJed Brown 522fb5bd1c2SPatrick Sanan Notes: 5234061b8bfSJed Brown This function is a callback used by the VTK viewer to actually write the file. 5244061b8bfSJed 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. 5254061b8bfSJed Brown Instead, metadata for the entire file needs to be available up-front before you can start writing the file. 5264061b8bfSJed Brown 527fb5bd1c2SPatrick Sanan If any fields have been named (see e.g. DMDASetFieldName()), then individual scalar 528fb5bd1c2SPatrick Sanan fields are written. Otherwise, a single multi-dof (vector) field is written. 529fb5bd1c2SPatrick Sanan 5304061b8bfSJed Brown .seealso: PETSCVIEWERVTK 5314061b8bfSJed Brown @*/ 5324061b8bfSJed Brown PetscErrorCode DMDAVTKWriteAll(PetscObject odm,PetscViewer viewer) 5334061b8bfSJed Brown { 5344061b8bfSJed Brown DM dm = (DM)odm; 5354061b8bfSJed Brown PetscBool isvtk; 5364061b8bfSJed Brown PetscErrorCode ierr; 5374061b8bfSJed Brown 5384061b8bfSJed Brown PetscFunctionBegin; 5394061b8bfSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 5404061b8bfSJed Brown PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 541251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr); 542ce94432eSBarry Smith if (!isvtk) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_INCOMP,"Cannot use viewer type %s",((PetscObject)viewer)->type_name); 5434061b8bfSJed Brown switch (viewer->format) { 5444061b8bfSJed Brown case PETSC_VIEWER_VTK_VTS: 5454061b8bfSJed Brown ierr = DMDAVTKWriteAll_VTS(dm,viewer);CHKERRQ(ierr); 5464061b8bfSJed Brown break; 547a13bc4e3SShao-Ching Huang case PETSC_VIEWER_VTK_VTR: 548a13bc4e3SShao-Ching Huang ierr = DMDAVTKWriteAll_VTR(dm,viewer);CHKERRQ(ierr); 549a13bc4e3SShao-Ching Huang break; 550ce94432eSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No support for format '%s'",PetscViewerFormats[viewer->format]); 5514061b8bfSJed Brown } 5524061b8bfSJed Brown PetscFunctionReturn(0); 5534061b8bfSJed Brown } 554