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 PetscInt f,bs; 14fb5bd1c2SPatrick Sanan 15fb5bd1c2SPatrick Sanan PetscFunctionBegin; 16fb5bd1c2SPatrick Sanan *fieldsnamed = PETSC_FALSE; 179566063dSJacob Faibussowitsch PetscCall(DMDAGetDof(da,&bs)); 18fb5bd1c2SPatrick Sanan for (f=0; f<bs; ++f) { 19fb5bd1c2SPatrick Sanan const char * fieldname; 209566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(da,f,&fieldname)); 21fb5bd1c2SPatrick Sanan if (fieldname) { 22fb5bd1c2SPatrick Sanan *fieldsnamed = PETSC_TRUE; 23fb5bd1c2SPatrick Sanan break; 24fb5bd1c2SPatrick Sanan } 25fb5bd1c2SPatrick Sanan } 26fb5bd1c2SPatrick Sanan PetscFunctionReturn(0); 27fb5bd1c2SPatrick Sanan } 28fb5bd1c2SPatrick Sanan 294061b8bfSJed Brown static PetscErrorCode DMDAVTKWriteAll_VTS(DM da,PetscViewer viewer) 304061b8bfSJed Brown { 3130815ce0SLisandro Dalcin const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian"; 324061b8bfSJed Brown #if defined(PETSC_USE_REAL_SINGLE) 334061b8bfSJed Brown const char precision[] = "Float32"; 344061b8bfSJed Brown #elif defined(PETSC_USE_REAL_DOUBLE) 354061b8bfSJed Brown const char precision[] = "Float64"; 364061b8bfSJed Brown #else 374061b8bfSJed Brown const char precision[] = "UnknownPrecision"; 384061b8bfSJed Brown #endif 39ce94432eSBarry Smith MPI_Comm comm; 402eaa9ef4SLisandro Dalcin Vec Coords; 414061b8bfSJed Brown PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data; 424061b8bfSJed Brown PetscViewerVTKObjectLink link; 434061b8bfSJed Brown FILE *fp; 444061b8bfSJed Brown PetscMPIInt rank,size,tag; 454061b8bfSJed Brown DMDALocalInfo info; 46ea2d7708SPatrick Sanan PetscInt dim,mx,my,mz,cdim,bs,boffset,maxnnodes,maxbs,i,j,k,r; 470298fd71SBarry Smith PetscInt rloc[6],(*grloc)[6] = NULL; 484061b8bfSJed Brown PetscScalar *array,*array2; 494061b8bfSJed Brown 504061b8bfSJed Brown PetscFunctionBegin; 519566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da,&comm)); 52*1dca8a05SBarry Smith PetscCheck(!PetscDefined(USE_COMPLEX),PETSC_COMM_SELF,PETSC_ERR_SUP,"Complex values not supported"); 539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 559566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,&mx,&my,&mz,NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL)); 569566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da,&info)); 579566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da,&Coords)); 582eaa9ef4SLisandro Dalcin if (Coords) { 592eaa9ef4SLisandro Dalcin PetscInt csize; 609566063dSJacob Faibussowitsch PetscCall(VecGetSize(Coords,&csize)); 6108401ef6SPierre Jolivet PetscCheck(csize % (mx*my*mz) == 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Coordinate vector size mismatch"); 622eaa9ef4SLisandro Dalcin cdim = csize/(mx*my*mz); 63*1dca8a05SBarry Smith PetscCheck(cdim >= dim && cdim <= 3,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Coordinate vector size mismatch"); 642eaa9ef4SLisandro Dalcin } else { 652eaa9ef4SLisandro Dalcin cdim = dim; 662eaa9ef4SLisandro Dalcin } 674061b8bfSJed Brown 689566063dSJacob Faibussowitsch PetscCall(PetscFOpen(comm,vtk->filename,"wb",&fp)); 699566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n")); 709566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"<VTKFile type=\"StructuredGrid\" version=\"0.1\" byte_order=\"%s\">\n",byte_order)); 7163a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <StructuredGrid WholeExtent=\"%d %" PetscInt_FMT " %d %" PetscInt_FMT " %d %" PetscInt_FMT "\">\n",0,mx-1,0,my-1,0,mz-1)); 724061b8bfSJed Brown 739566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc1(size*6,&grloc)); 744061b8bfSJed Brown rloc[0] = info.xs; 754061b8bfSJed Brown rloc[1] = info.xm; 764061b8bfSJed Brown rloc[2] = info.ys; 774061b8bfSJed Brown rloc[3] = info.ym; 784061b8bfSJed Brown rloc[4] = info.zs; 794061b8bfSJed Brown rloc[5] = info.zm; 809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gather(rloc,6,MPIU_INT,&grloc[0][0],6,MPIU_INT,0,comm)); 814061b8bfSJed Brown 824061b8bfSJed Brown /* Write XML header */ 834061b8bfSJed Brown maxnnodes = 0; /* Used for the temporary array size on rank 0 */ 84ea2d7708SPatrick Sanan maxbs = 0; /* Used for the temporary array size on rank 0 */ 854061b8bfSJed Brown boffset = 0; /* Offset into binary file */ 864061b8bfSJed Brown for (r=0; r<size; r++) { 874061b8bfSJed Brown PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0; 88dd400576SPatrick Sanan if (rank == 0) { 894061b8bfSJed Brown xs = grloc[r][0]; 904061b8bfSJed Brown xm = grloc[r][1]; 914061b8bfSJed Brown ys = grloc[r][2]; 924061b8bfSJed Brown ym = grloc[r][3]; 934061b8bfSJed Brown zs = grloc[r][4]; 944061b8bfSJed Brown zm = grloc[r][5]; 954061b8bfSJed Brown nnodes = xm*ym*zm; 964061b8bfSJed Brown } 974061b8bfSJed Brown maxnnodes = PetscMax(maxnnodes,nnodes); 9863a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <Piece Extent=\"%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\">\n",xs,xs+xm-1,ys,ys+ym-1,zs,zs+zm-1)); 999566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <Points>\n")); 10063a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"Position\" NumberOfComponents=\"3\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,boffset)); 1014061b8bfSJed Brown boffset += 3*nnodes*sizeof(PetscScalar) + sizeof(int); 1029566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," </Points>\n")); 1034061b8bfSJed Brown 1049566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <PointData Scalars=\"ScalarPointData\">\n")); 1054061b8bfSJed Brown for (link=vtk->link; link; link=link->next) { 1064061b8bfSJed Brown Vec X = (Vec)link->vec; 107fb5bd1c2SPatrick Sanan PetscInt bs,f; 108ea2d7708SPatrick Sanan DM daCurr; 109fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 110fb5bd1c2SPatrick Sanan const char *vecname = "Unnamed"; 111ea2d7708SPatrick Sanan 1129566063dSJacob Faibussowitsch PetscCall(VecGetDM(X,&daCurr)); 1139566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL)); 114ea2d7708SPatrick Sanan maxbs = PetscMax(maxbs,bs); 115ea2d7708SPatrick Sanan 116ea2d7708SPatrick Sanan if (((PetscObject)X)->name || link != vtk->link) { 1179566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)X,&vecname)); 1184061b8bfSJed Brown } 119fb5bd1c2SPatrick Sanan 120fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 1219566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr,&fieldsnamed)); 122fb5bd1c2SPatrick Sanan if (fieldsnamed) { 123fb5bd1c2SPatrick Sanan for (f=0; f<bs; f++) { 124fb5bd1c2SPatrick Sanan char buf[256]; 125fb5bd1c2SPatrick Sanan const char *fieldname; 1269566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(daCurr,f,&fieldname)); 127fb5bd1c2SPatrick Sanan if (!fieldname) { 12863a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf,sizeof(buf),"%" PetscInt_FMT,f)); 129fb5bd1c2SPatrick Sanan fieldname = buf; 130fb5bd1c2SPatrick Sanan } 13163a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,vecname,fieldname,boffset)); 132fb5bd1c2SPatrick Sanan boffset += nnodes*sizeof(PetscScalar) + sizeof(int); 133fb5bd1c2SPatrick Sanan } 134fb5bd1c2SPatrick Sanan } else { 13563a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%" PetscInt_FMT "\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,vecname,bs,boffset)); 136ea2d7708SPatrick Sanan boffset += bs*nnodes*sizeof(PetscScalar) + sizeof(int); 1374061b8bfSJed Brown } 138fb5bd1c2SPatrick Sanan } 1399566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," </PointData>\n")); 1409566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," </Piece>\n")); 1414061b8bfSJed Brown } 1429566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," </StructuredGrid>\n")); 1439566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <AppendedData encoding=\"raw\">\n")); 1449566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"_")); 1454061b8bfSJed Brown 1464061b8bfSJed Brown /* Now write the arrays. */ 1474061b8bfSJed Brown tag = ((PetscObject)viewer)->tag; 1489566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(maxnnodes*PetscMax(3,maxbs),&array,maxnnodes*PetscMax(3,maxbs),&array2)); 1494061b8bfSJed Brown for (r=0; r<size; r++) { 1504061b8bfSJed Brown MPI_Status status; 1514061b8bfSJed Brown PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0; 152dd400576SPatrick Sanan if (rank == 0) { 1534061b8bfSJed Brown xs = grloc[r][0]; 1544061b8bfSJed Brown xm = grloc[r][1]; 1554061b8bfSJed Brown ys = grloc[r][2]; 1564061b8bfSJed Brown ym = grloc[r][3]; 1574061b8bfSJed Brown zs = grloc[r][4]; 1584061b8bfSJed Brown zm = grloc[r][5]; 1594061b8bfSJed Brown nnodes = xm*ym*zm; 1604061b8bfSJed Brown } else if (r == rank) { 1614061b8bfSJed Brown nnodes = info.xm*info.ym*info.zm; 1624061b8bfSJed Brown } 1634061b8bfSJed Brown 1642eaa9ef4SLisandro Dalcin /* Write the coordinates */ 1654061b8bfSJed Brown if (Coords) { 1664061b8bfSJed Brown const PetscScalar *coords; 1679566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Coords,&coords)); 168dd400576SPatrick Sanan if (rank == 0) { 1694061b8bfSJed Brown if (r) { 1706622f924SJed Brown PetscMPIInt nn; 1719566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array,nnodes*cdim,MPIU_SCALAR,r,tag,comm,&status)); 1729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status,MPIU_SCALAR,&nn)); 17308401ef6SPierre Jolivet PetscCheck(nn == nnodes*cdim,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch"); 1744061b8bfSJed Brown } else { 1759566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(array,coords,nnodes*cdim)); 1764061b8bfSJed Brown } 1774061b8bfSJed Brown /* Transpose coordinates to VTK (C-style) ordering */ 1784061b8bfSJed Brown for (k=0; k<zm; k++) { 1794061b8bfSJed Brown for (j=0; j<ym; j++) { 1804061b8bfSJed Brown for (i=0; i<xm; i++) { 1814061b8bfSJed Brown PetscInt Iloc = i+xm*(j+ym*k); 1822eaa9ef4SLisandro Dalcin array2[Iloc*3+0] = array[Iloc*cdim + 0]; 1832eaa9ef4SLisandro Dalcin array2[Iloc*3+1] = cdim > 1 ? array[Iloc*cdim + 1] : 0.0; 1842eaa9ef4SLisandro Dalcin array2[Iloc*3+2] = cdim > 2 ? array[Iloc*cdim + 2] : 0.0; 1854061b8bfSJed Brown } 1864061b8bfSJed Brown } 1874061b8bfSJed Brown } 1884061b8bfSJed Brown } else if (r == rank) { 1899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void*)coords,nnodes*cdim,MPIU_SCALAR,0,tag,comm)); 1904061b8bfSJed Brown } 1919566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Coords,&coords)); 1924061b8bfSJed Brown } else { /* Fabricate some coordinates using grid index */ 1934061b8bfSJed Brown for (k=0; k<zm; k++) { 1944061b8bfSJed Brown for (j=0; j<ym; j++) { 1954061b8bfSJed Brown for (i=0; i<xm; i++) { 1964061b8bfSJed Brown PetscInt Iloc = i+xm*(j+ym*k); 1974061b8bfSJed Brown array2[Iloc*3+0] = xs+i; 1984061b8bfSJed Brown array2[Iloc*3+1] = ys+j; 1994061b8bfSJed Brown array2[Iloc*3+2] = zs+k; 2004061b8bfSJed Brown } 2014061b8bfSJed Brown } 2024061b8bfSJed Brown } 2034061b8bfSJed Brown } 2049566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer,fp,array2,nnodes*3,MPIU_SCALAR)); 2054061b8bfSJed Brown 2064061b8bfSJed Brown /* Write each of the objects queued up for this file */ 2074061b8bfSJed Brown for (link=vtk->link; link; link=link->next) { 2084061b8bfSJed Brown Vec X = (Vec)link->vec; 2094061b8bfSJed Brown const PetscScalar *x; 210fb5bd1c2SPatrick Sanan PetscInt bs,f; 211ea2d7708SPatrick Sanan DM daCurr; 212fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 2139566063dSJacob Faibussowitsch PetscCall(VecGetDM(X,&daCurr)); 2149566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr,NULL,NULL,NULL,NULL, NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL)); 2159566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 216dd400576SPatrick Sanan if (rank == 0) { 2174061b8bfSJed Brown if (r) { 2186622f924SJed Brown PetscMPIInt nn; 2199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status)); 2209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status,MPIU_SCALAR,&nn)); 22163a3b9bcSJacob Faibussowitsch PetscCheck(nn == nnodes*bs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %" PetscInt_FMT,r); 2224061b8bfSJed Brown } else { 2239566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(array,x,nnodes*bs)); 2244061b8bfSJed Brown } 225fb5bd1c2SPatrick Sanan 226fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 2279566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr,&fieldsnamed)); 228fb5bd1c2SPatrick Sanan if (fieldsnamed) { 229fb5bd1c2SPatrick Sanan for (f=0; f<bs; f++) { 230fb5bd1c2SPatrick Sanan /* Extract and transpose the f'th field */ 231fb5bd1c2SPatrick Sanan for (k=0; k<zm; k++) { 232fb5bd1c2SPatrick Sanan for (j=0; j<ym; j++) { 233fb5bd1c2SPatrick Sanan for (i=0; i<xm; i++) { 234fb5bd1c2SPatrick Sanan PetscInt Iloc = i+xm*(j+ym*k); 235fb5bd1c2SPatrick Sanan array2[Iloc] = array[Iloc*bs + f]; 236fb5bd1c2SPatrick Sanan } 237fb5bd1c2SPatrick Sanan } 238fb5bd1c2SPatrick Sanan } 2399566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer,fp,array2,nnodes,MPIU_SCALAR)); 240fb5bd1c2SPatrick Sanan } 241fb5bd1c2SPatrick Sanan } else { 2429566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer,fp,array,bs*nnodes,MPIU_SCALAR)); 243fb5bd1c2SPatrick Sanan } 2444061b8bfSJed Brown } else if (r == rank) { 2459566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm)); 2464061b8bfSJed Brown } 2479566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 2484061b8bfSJed Brown } 2494061b8bfSJed Brown } 2509566063dSJacob Faibussowitsch PetscCall(PetscFree2(array,array2)); 2519566063dSJacob Faibussowitsch PetscCall(PetscFree(grloc)); 2524061b8bfSJed Brown 2539566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"\n </AppendedData>\n")); 2549566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"</VTKFile>\n")); 2559566063dSJacob Faibussowitsch PetscCall(PetscFClose(comm,fp)); 2564061b8bfSJed Brown PetscFunctionReturn(0); 2574061b8bfSJed Brown } 2584061b8bfSJed Brown 259a13bc4e3SShao-Ching Huang static PetscErrorCode DMDAVTKWriteAll_VTR(DM da,PetscViewer viewer) 260a13bc4e3SShao-Ching Huang { 26130815ce0SLisandro Dalcin const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian"; 262a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_REAL_SINGLE) 263a13bc4e3SShao-Ching Huang const char precision[] = "Float32"; 264a13bc4e3SShao-Ching Huang #elif defined(PETSC_USE_REAL_DOUBLE) 265a13bc4e3SShao-Ching Huang const char precision[] = "Float64"; 266a13bc4e3SShao-Ching Huang #else 267a13bc4e3SShao-Ching Huang const char precision[] = "UnknownPrecision"; 268a13bc4e3SShao-Ching Huang #endif 269a13bc4e3SShao-Ching Huang MPI_Comm comm; 270a13bc4e3SShao-Ching Huang PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data; 271a13bc4e3SShao-Ching Huang PetscViewerVTKObjectLink link; 272a13bc4e3SShao-Ching Huang FILE *fp; 273a13bc4e3SShao-Ching Huang PetscMPIInt rank,size,tag; 274a13bc4e3SShao-Ching Huang DMDALocalInfo info; 275ea2d7708SPatrick Sanan PetscInt dim,mx,my,mz,boffset,maxnnodes,maxbs,i,j,k,r; 276a13bc4e3SShao-Ching Huang PetscInt rloc[6],(*grloc)[6] = NULL; 277a13bc4e3SShao-Ching Huang PetscScalar *array,*array2; 278a13bc4e3SShao-Ching Huang 279a13bc4e3SShao-Ching Huang PetscFunctionBegin; 2809566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da,&comm)); 281*1dca8a05SBarry Smith PetscCheck(!PetscDefined(USE_COMPLEX),PETSC_COMM_SELF,PETSC_ERR_SUP,"Complex values not supported"); 2829566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 2839566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 2849566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,&mx,&my,&mz,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 2859566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da,&info)); 2869566063dSJacob Faibussowitsch PetscCall(PetscFOpen(comm,vtk->filename,"wb",&fp)); 2879566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n")); 2889566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"<VTKFile type=\"RectilinearGrid\" version=\"0.1\" byte_order=\"%s\">\n",byte_order)); 28963a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <RectilinearGrid WholeExtent=\"%d %" PetscInt_FMT " %d %" PetscInt_FMT " %d %" PetscInt_FMT "\">\n",0,mx-1,0,my-1,0,mz-1)); 290a13bc4e3SShao-Ching Huang 2919566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc1(size*6,&grloc)); 292a13bc4e3SShao-Ching Huang rloc[0] = info.xs; 293a13bc4e3SShao-Ching Huang rloc[1] = info.xm; 294a13bc4e3SShao-Ching Huang rloc[2] = info.ys; 295a13bc4e3SShao-Ching Huang rloc[3] = info.ym; 296a13bc4e3SShao-Ching Huang rloc[4] = info.zs; 297a13bc4e3SShao-Ching Huang rloc[5] = info.zm; 2989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gather(rloc,6,MPIU_INT,&grloc[0][0],6,MPIU_INT,0,comm)); 299a13bc4e3SShao-Ching Huang 300a13bc4e3SShao-Ching Huang /* Write XML header */ 301a13bc4e3SShao-Ching Huang maxnnodes = 0; /* Used for the temporary array size on rank 0 */ 302ea2d7708SPatrick Sanan maxbs = 0; /* Used for the temporary array size on rank 0 */ 303a13bc4e3SShao-Ching Huang boffset = 0; /* Offset into binary file */ 304a13bc4e3SShao-Ching Huang for (r=0; r<size; r++) { 305a13bc4e3SShao-Ching Huang PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0; 306dd400576SPatrick Sanan if (rank == 0) { 307a13bc4e3SShao-Ching Huang xs = grloc[r][0]; 308a13bc4e3SShao-Ching Huang xm = grloc[r][1]; 309a13bc4e3SShao-Ching Huang ys = grloc[r][2]; 310a13bc4e3SShao-Ching Huang ym = grloc[r][3]; 311a13bc4e3SShao-Ching Huang zs = grloc[r][4]; 312a13bc4e3SShao-Ching Huang zm = grloc[r][5]; 313a13bc4e3SShao-Ching Huang nnodes = xm*ym*zm; 314a13bc4e3SShao-Ching Huang } 315a13bc4e3SShao-Ching Huang maxnnodes = PetscMax(maxnnodes,nnodes); 31663a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <Piece Extent=\"%" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\">\n",xs,xs+xm-1,ys,ys+ym-1,zs,zs+zm-1)); 3179566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <Coordinates>\n")); 31863a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"Xcoord\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,boffset)); 319a13bc4e3SShao-Ching Huang boffset += xm*sizeof(PetscScalar) + sizeof(int); 32063a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"Ycoord\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,boffset)); 321a13bc4e3SShao-Ching Huang boffset += ym*sizeof(PetscScalar) + sizeof(int); 32263a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"Zcoord\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,boffset)); 323a13bc4e3SShao-Ching Huang boffset += zm*sizeof(PetscScalar) + sizeof(int); 3249566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," </Coordinates>\n")); 3259566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <PointData Scalars=\"ScalarPointData\">\n")); 326a13bc4e3SShao-Ching Huang for (link=vtk->link; link; link=link->next) { 327a13bc4e3SShao-Ching Huang Vec X = (Vec)link->vec; 328fb5bd1c2SPatrick Sanan PetscInt bs,f; 329ea2d7708SPatrick Sanan DM daCurr; 330fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 331fb5bd1c2SPatrick Sanan const char *vecname = "Unnamed"; 332ea2d7708SPatrick Sanan 3339566063dSJacob Faibussowitsch PetscCall(VecGetDM(X,&daCurr)); 3349566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL)); 335ea2d7708SPatrick Sanan maxbs = PetscMax(maxbs,bs); 336ea2d7708SPatrick Sanan if (((PetscObject)X)->name || link != vtk->link) { 3379566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)X,&vecname)); 338a13bc4e3SShao-Ching Huang } 339ea2d7708SPatrick Sanan 340fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 3419566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr,&fieldsnamed)); 342fb5bd1c2SPatrick Sanan if (fieldsnamed) { 343fb5bd1c2SPatrick Sanan for (f=0; f<bs; f++) { 344fb5bd1c2SPatrick Sanan char buf[256]; 345fb5bd1c2SPatrick Sanan const char *fieldname; 3469566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(daCurr,f,&fieldname)); 347fb5bd1c2SPatrick Sanan if (!fieldname) { 34863a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf,sizeof(buf),"%" PetscInt_FMT,f)); 349fb5bd1c2SPatrick Sanan fieldname = buf; 350fb5bd1c2SPatrick Sanan } 35163a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,vecname,fieldname,boffset)); 352fb5bd1c2SPatrick Sanan boffset += nnodes*sizeof(PetscScalar) + sizeof(int); 353fb5bd1c2SPatrick Sanan } 354fb5bd1c2SPatrick Sanan } else { 35563a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%" PetscInt_FMT "\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,vecname,bs,boffset)); 356ea2d7708SPatrick Sanan boffset += bs*nnodes*sizeof(PetscScalar) + sizeof(int); 357a13bc4e3SShao-Ching Huang } 358fb5bd1c2SPatrick Sanan } 3599566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," </PointData>\n")); 3609566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," </Piece>\n")); 361a13bc4e3SShao-Ching Huang } 3629566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," </RectilinearGrid>\n")); 3639566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <AppendedData encoding=\"raw\">\n")); 3649566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"_")); 365a13bc4e3SShao-Ching Huang 366a13bc4e3SShao-Ching Huang /* Now write the arrays. */ 367a13bc4e3SShao-Ching Huang tag = ((PetscObject)viewer)->tag; 3689566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(PetscMax(maxnnodes*maxbs,info.xm+info.ym+info.zm),&array,PetscMax(maxnnodes*maxbs,info.xm+info.ym+info.zm),&array2)); 369a13bc4e3SShao-Ching Huang for (r=0; r<size; r++) { 370a13bc4e3SShao-Ching Huang MPI_Status status; 371a13bc4e3SShao-Ching Huang PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0; 372dd400576SPatrick Sanan if (rank == 0) { 373a13bc4e3SShao-Ching Huang xs = grloc[r][0]; 374a13bc4e3SShao-Ching Huang xm = grloc[r][1]; 375a13bc4e3SShao-Ching Huang ys = grloc[r][2]; 376a13bc4e3SShao-Ching Huang ym = grloc[r][3]; 377a13bc4e3SShao-Ching Huang zs = grloc[r][4]; 378a13bc4e3SShao-Ching Huang zm = grloc[r][5]; 379a13bc4e3SShao-Ching Huang nnodes = xm*ym*zm; 380a13bc4e3SShao-Ching Huang } else if (r == rank) { 381a13bc4e3SShao-Ching Huang nnodes = info.xm*info.ym*info.zm; 382a13bc4e3SShao-Ching Huang } 383a13bc4e3SShao-Ching Huang { /* Write the coordinates */ 384a13bc4e3SShao-Ching Huang Vec Coords; 385ea2d7708SPatrick Sanan 3869566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da,&Coords)); 387a13bc4e3SShao-Ching Huang if (Coords) { 388a13bc4e3SShao-Ching Huang const PetscScalar *coords; 3899566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Coords,&coords)); 390dd400576SPatrick Sanan if (rank == 0) { 391a13bc4e3SShao-Ching Huang if (r) { 392a13bc4e3SShao-Ching Huang PetscMPIInt nn; 3939566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array,xm+ym+zm,MPIU_SCALAR,r,tag,comm,&status)); 3949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status,MPIU_SCALAR,&nn)); 39508401ef6SPierre Jolivet PetscCheck(nn == xm+ym+zm,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch"); 396a13bc4e3SShao-Ching Huang } else { 397a13bc4e3SShao-Ching Huang /* x coordinates */ 398a13bc4e3SShao-Ching Huang for (j=0, k=0, i=0; i<xm; i++) { 399a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 400a13bc4e3SShao-Ching Huang array[i] = coords[Iloc*dim + 0]; 401a13bc4e3SShao-Ching Huang } 402a13bc4e3SShao-Ching Huang /* y coordinates */ 403a13bc4e3SShao-Ching Huang for (i=0, k=0, j=0; j<ym; j++) { 404a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 405a13bc4e3SShao-Ching Huang array[j+xm] = dim > 1 ? coords[Iloc*dim + 1] : 0; 406a13bc4e3SShao-Ching Huang } 407a13bc4e3SShao-Ching Huang /* z coordinates */ 408a13bc4e3SShao-Ching Huang for (i=0, j=0, k=0; k<zm; k++) { 409a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 410a13bc4e3SShao-Ching Huang array[k+xm+ym] = dim > 2 ? coords[Iloc*dim + 2] : 0; 411a13bc4e3SShao-Ching Huang } 412a13bc4e3SShao-Ching Huang } 413a13bc4e3SShao-Ching Huang } else if (r == rank) { 414a13bc4e3SShao-Ching Huang xm = info.xm; 415a13bc4e3SShao-Ching Huang ym = info.ym; 416a13bc4e3SShao-Ching Huang zm = info.zm; 417a13bc4e3SShao-Ching Huang for (j=0, k=0, i=0; i<xm; i++) { 418a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 419a13bc4e3SShao-Ching Huang array2[i] = coords[Iloc*dim + 0]; 420a13bc4e3SShao-Ching Huang } 421a13bc4e3SShao-Ching Huang for (i=0, k=0, j=0; j<ym; j++) { 422a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 423a13bc4e3SShao-Ching Huang array2[j+xm] = dim > 1 ? coords[Iloc*dim + 1] : 0; 424a13bc4e3SShao-Ching Huang } 425a13bc4e3SShao-Ching Huang for (i=0, j=0, k=0; k<zm; k++) { 426a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 427a13bc4e3SShao-Ching Huang array2[k+xm+ym] = dim > 2 ? coords[Iloc*dim + 2] : 0; 428a13bc4e3SShao-Ching Huang } 4299566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void*)array2,xm+ym+zm,MPIU_SCALAR,0,tag,comm)); 430a13bc4e3SShao-Ching Huang } 4319566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Coords,&coords)); 432a13bc4e3SShao-Ching Huang } else { /* Fabricate some coordinates using grid index */ 4333d3eaba7SBarry Smith for (i=0; i<xm; i++) { 434a13bc4e3SShao-Ching Huang array[i] = xs+i; 435a13bc4e3SShao-Ching Huang } 4363d3eaba7SBarry Smith for (j=0; j<ym; j++) { 437a13bc4e3SShao-Ching Huang array[j+xm] = ys+j; 438a13bc4e3SShao-Ching Huang } 4393d3eaba7SBarry Smith for (k=0; k<zm; k++) { 440a13bc4e3SShao-Ching Huang array[k+xm+ym] = zs+k; 441a13bc4e3SShao-Ching Huang } 442a13bc4e3SShao-Ching Huang } 443dd400576SPatrick Sanan if (rank == 0) { 4449566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer,fp,&(array[0]),xm,MPIU_SCALAR)); 4459566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer,fp,&(array[xm]),ym,MPIU_SCALAR)); 4469566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer,fp,&(array[xm+ym]),zm,MPIU_SCALAR)); 447a13bc4e3SShao-Ching Huang } 448a13bc4e3SShao-Ching Huang } 449a13bc4e3SShao-Ching Huang 450a13bc4e3SShao-Ching Huang /* Write each of the objects queued up for this file */ 451a13bc4e3SShao-Ching Huang for (link=vtk->link; link; link=link->next) { 452a13bc4e3SShao-Ching Huang Vec X = (Vec)link->vec; 453a13bc4e3SShao-Ching Huang const PetscScalar *x; 454fb5bd1c2SPatrick Sanan PetscInt bs,f; 455ea2d7708SPatrick Sanan DM daCurr; 456fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 4579566063dSJacob Faibussowitsch PetscCall(VecGetDM(X,&daCurr)); 4589566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL)); 459a13bc4e3SShao-Ching Huang 4609566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 461dd400576SPatrick Sanan if (rank == 0) { 462a13bc4e3SShao-Ching Huang if (r) { 463a13bc4e3SShao-Ching Huang PetscMPIInt nn; 4649566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status)); 4659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status,MPIU_SCALAR,&nn)); 46663a3b9bcSJacob Faibussowitsch PetscCheck(nn == nnodes*bs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %" PetscInt_FMT,r); 467a13bc4e3SShao-Ching Huang } else { 4689566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(array,x,nnodes*bs)); 469a13bc4e3SShao-Ching Huang } 470fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 4719566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr,&fieldsnamed)); 472fb5bd1c2SPatrick Sanan if (fieldsnamed) { 473fb5bd1c2SPatrick Sanan for (f=0; f<bs; f++) { 474fb5bd1c2SPatrick Sanan /* Extract and transpose the f'th field */ 475fb5bd1c2SPatrick Sanan for (k=0; k<zm; k++) { 476fb5bd1c2SPatrick Sanan for (j=0; j<ym; j++) { 477fb5bd1c2SPatrick Sanan for (i=0; i<xm; i++) { 478fb5bd1c2SPatrick Sanan PetscInt Iloc = i+xm*(j+ym*k); 479fb5bd1c2SPatrick Sanan array2[Iloc] = array[Iloc*bs + f]; 480fb5bd1c2SPatrick Sanan } 481fb5bd1c2SPatrick Sanan } 482fb5bd1c2SPatrick Sanan } 4839566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer,fp,array2,nnodes,MPIU_SCALAR)); 484fb5bd1c2SPatrick Sanan } 485fb5bd1c2SPatrick Sanan } 4869566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer,fp,array,nnodes*bs,MPIU_SCALAR)); 487ea2d7708SPatrick Sanan 488a13bc4e3SShao-Ching Huang } else if (r == rank) { 4899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm)); 490a13bc4e3SShao-Ching Huang } 4919566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 492a13bc4e3SShao-Ching Huang } 493a13bc4e3SShao-Ching Huang } 4949566063dSJacob Faibussowitsch PetscCall(PetscFree2(array,array2)); 4959566063dSJacob Faibussowitsch PetscCall(PetscFree(grloc)); 496a13bc4e3SShao-Ching Huang 4979566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"\n </AppendedData>\n")); 4989566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"</VTKFile>\n")); 4999566063dSJacob Faibussowitsch PetscCall(PetscFClose(comm,fp)); 500a13bc4e3SShao-Ching Huang PetscFunctionReturn(0); 501a13bc4e3SShao-Ching Huang } 502a13bc4e3SShao-Ching Huang 5034061b8bfSJed Brown /*@C 5044061b8bfSJed Brown DMDAVTKWriteAll - Write a file containing all the fields that have been provided to the viewer 5054061b8bfSJed Brown 5064061b8bfSJed Brown Collective 5074061b8bfSJed Brown 5084165533cSJose E. Roman Input Parameters: 5094165533cSJose E. Roman + odm - DM specifying the grid layout, passed as a PetscObject 5104165533cSJose E. Roman - viewer - viewer of type VTK 5114061b8bfSJed Brown 5124061b8bfSJed Brown Level: developer 5134061b8bfSJed Brown 514fb5bd1c2SPatrick Sanan Notes: 5154061b8bfSJed Brown This function is a callback used by the VTK viewer to actually write the file. 5164061b8bfSJed 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. 5174061b8bfSJed Brown Instead, metadata for the entire file needs to be available up-front before you can start writing the file. 5184061b8bfSJed Brown 519fb5bd1c2SPatrick Sanan If any fields have been named (see e.g. DMDASetFieldName()), then individual scalar 520fb5bd1c2SPatrick Sanan fields are written. Otherwise, a single multi-dof (vector) field is written. 521fb5bd1c2SPatrick Sanan 5224061b8bfSJed Brown .seealso: PETSCVIEWERVTK 5234061b8bfSJed Brown @*/ 5244061b8bfSJed Brown PetscErrorCode DMDAVTKWriteAll(PetscObject odm,PetscViewer viewer) 5254061b8bfSJed Brown { 5264061b8bfSJed Brown DM dm = (DM)odm; 5274061b8bfSJed Brown PetscBool isvtk; 5284061b8bfSJed Brown 5294061b8bfSJed Brown PetscFunctionBegin; 5304061b8bfSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 5314061b8bfSJed Brown PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 5329566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk)); 5337a8be351SBarry Smith PetscCheck(isvtk,PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_INCOMP,"Cannot use viewer type %s",((PetscObject)viewer)->type_name); 5344061b8bfSJed Brown switch (viewer->format) { 5354061b8bfSJed Brown case PETSC_VIEWER_VTK_VTS: 5369566063dSJacob Faibussowitsch PetscCall(DMDAVTKWriteAll_VTS(dm,viewer)); 5374061b8bfSJed Brown break; 538a13bc4e3SShao-Ching Huang case PETSC_VIEWER_VTK_VTR: 5399566063dSJacob Faibussowitsch PetscCall(DMDAVTKWriteAll_VTR(dm,viewer)); 540a13bc4e3SShao-Ching Huang break; 54198921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No support for format '%s'",PetscViewerFormats[viewer->format]); 5424061b8bfSJed Brown } 5434061b8bfSJed Brown PetscFunctionReturn(0); 5444061b8bfSJed Brown } 545