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