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)); 521dca8a05SBarry 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); 631dca8a05SBarry 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"); 174*1baa6e33SBarry Smith } else PetscCall(PetscArraycpy(array,coords,nnodes*cdim)); 1754061b8bfSJed Brown /* Transpose coordinates to VTK (C-style) ordering */ 1764061b8bfSJed Brown for (k=0; k<zm; k++) { 1774061b8bfSJed Brown for (j=0; j<ym; j++) { 1784061b8bfSJed Brown for (i=0; i<xm; i++) { 1794061b8bfSJed Brown PetscInt Iloc = i+xm*(j+ym*k); 1802eaa9ef4SLisandro Dalcin array2[Iloc*3+0] = array[Iloc*cdim + 0]; 1812eaa9ef4SLisandro Dalcin array2[Iloc*3+1] = cdim > 1 ? array[Iloc*cdim + 1] : 0.0; 1822eaa9ef4SLisandro Dalcin array2[Iloc*3+2] = cdim > 2 ? array[Iloc*cdim + 2] : 0.0; 1834061b8bfSJed Brown } 1844061b8bfSJed Brown } 1854061b8bfSJed Brown } 1864061b8bfSJed Brown } else if (r == rank) { 1879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void*)coords,nnodes*cdim,MPIU_SCALAR,0,tag,comm)); 1884061b8bfSJed Brown } 1899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Coords,&coords)); 1904061b8bfSJed Brown } else { /* Fabricate some coordinates using grid index */ 1914061b8bfSJed Brown for (k=0; k<zm; k++) { 1924061b8bfSJed Brown for (j=0; j<ym; j++) { 1934061b8bfSJed Brown for (i=0; i<xm; i++) { 1944061b8bfSJed Brown PetscInt Iloc = i+xm*(j+ym*k); 1954061b8bfSJed Brown array2[Iloc*3+0] = xs+i; 1964061b8bfSJed Brown array2[Iloc*3+1] = ys+j; 1974061b8bfSJed Brown array2[Iloc*3+2] = zs+k; 1984061b8bfSJed Brown } 1994061b8bfSJed Brown } 2004061b8bfSJed Brown } 2014061b8bfSJed Brown } 2029566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer,fp,array2,nnodes*3,MPIU_SCALAR)); 2034061b8bfSJed Brown 2044061b8bfSJed Brown /* Write each of the objects queued up for this file */ 2054061b8bfSJed Brown for (link=vtk->link; link; link=link->next) { 2064061b8bfSJed Brown Vec X = (Vec)link->vec; 2074061b8bfSJed Brown const PetscScalar *x; 208fb5bd1c2SPatrick Sanan PetscInt bs,f; 209ea2d7708SPatrick Sanan DM daCurr; 210fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 2119566063dSJacob Faibussowitsch PetscCall(VecGetDM(X,&daCurr)); 2129566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr,NULL,NULL,NULL,NULL, NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL)); 2139566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 214dd400576SPatrick Sanan if (rank == 0) { 2154061b8bfSJed Brown if (r) { 2166622f924SJed Brown PetscMPIInt nn; 2179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status)); 2189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status,MPIU_SCALAR,&nn)); 21963a3b9bcSJacob Faibussowitsch PetscCheck(nn == nnodes*bs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %" PetscInt_FMT,r); 220*1baa6e33SBarry Smith } else PetscCall(PetscArraycpy(array,x,nnodes*bs)); 221fb5bd1c2SPatrick Sanan 222fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 2239566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr,&fieldsnamed)); 224fb5bd1c2SPatrick Sanan if (fieldsnamed) { 225fb5bd1c2SPatrick Sanan for (f=0; f<bs; f++) { 226fb5bd1c2SPatrick Sanan /* Extract and transpose the f'th field */ 227fb5bd1c2SPatrick Sanan for (k=0; k<zm; k++) { 228fb5bd1c2SPatrick Sanan for (j=0; j<ym; j++) { 229fb5bd1c2SPatrick Sanan for (i=0; i<xm; i++) { 230fb5bd1c2SPatrick Sanan PetscInt Iloc = i+xm*(j+ym*k); 231fb5bd1c2SPatrick Sanan array2[Iloc] = array[Iloc*bs + f]; 232fb5bd1c2SPatrick Sanan } 233fb5bd1c2SPatrick Sanan } 234fb5bd1c2SPatrick Sanan } 2359566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer,fp,array2,nnodes,MPIU_SCALAR)); 236fb5bd1c2SPatrick Sanan } 237*1baa6e33SBarry Smith } else PetscCall(PetscViewerVTKFWrite(viewer,fp,array,bs*nnodes,MPIU_SCALAR)); 238*1baa6e33SBarry Smith } else if (r == rank) PetscCallMPI(MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm)); 2399566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 2404061b8bfSJed Brown } 2414061b8bfSJed Brown } 2429566063dSJacob Faibussowitsch PetscCall(PetscFree2(array,array2)); 2439566063dSJacob Faibussowitsch PetscCall(PetscFree(grloc)); 2444061b8bfSJed Brown 2459566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"\n </AppendedData>\n")); 2469566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"</VTKFile>\n")); 2479566063dSJacob Faibussowitsch PetscCall(PetscFClose(comm,fp)); 2484061b8bfSJed Brown PetscFunctionReturn(0); 2494061b8bfSJed Brown } 2504061b8bfSJed Brown 251a13bc4e3SShao-Ching Huang static PetscErrorCode DMDAVTKWriteAll_VTR(DM da,PetscViewer viewer) 252a13bc4e3SShao-Ching Huang { 25330815ce0SLisandro Dalcin const char *byte_order = PetscBinaryBigEndian() ? "BigEndian" : "LittleEndian"; 254a13bc4e3SShao-Ching Huang #if defined(PETSC_USE_REAL_SINGLE) 255a13bc4e3SShao-Ching Huang const char precision[] = "Float32"; 256a13bc4e3SShao-Ching Huang #elif defined(PETSC_USE_REAL_DOUBLE) 257a13bc4e3SShao-Ching Huang const char precision[] = "Float64"; 258a13bc4e3SShao-Ching Huang #else 259a13bc4e3SShao-Ching Huang const char precision[] = "UnknownPrecision"; 260a13bc4e3SShao-Ching Huang #endif 261a13bc4e3SShao-Ching Huang MPI_Comm comm; 262a13bc4e3SShao-Ching Huang PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data; 263a13bc4e3SShao-Ching Huang PetscViewerVTKObjectLink link; 264a13bc4e3SShao-Ching Huang FILE *fp; 265a13bc4e3SShao-Ching Huang PetscMPIInt rank,size,tag; 266a13bc4e3SShao-Ching Huang DMDALocalInfo info; 267ea2d7708SPatrick Sanan PetscInt dim,mx,my,mz,boffset,maxnnodes,maxbs,i,j,k,r; 268a13bc4e3SShao-Ching Huang PetscInt rloc[6],(*grloc)[6] = NULL; 269a13bc4e3SShao-Ching Huang PetscScalar *array,*array2; 270a13bc4e3SShao-Ching Huang 271a13bc4e3SShao-Ching Huang PetscFunctionBegin; 2729566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)da,&comm)); 2731dca8a05SBarry Smith PetscCheck(!PetscDefined(USE_COMPLEX),PETSC_COMM_SELF,PETSC_ERR_SUP,"Complex values not supported"); 2749566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 2759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 2769566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,&mx,&my,&mz,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 2779566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da,&info)); 2789566063dSJacob Faibussowitsch PetscCall(PetscFOpen(comm,vtk->filename,"wb",&fp)); 2799566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"<?xml version=\"1.0\"?>\n")); 2809566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"<VTKFile type=\"RectilinearGrid\" version=\"0.1\" byte_order=\"%s\">\n",byte_order)); 28163a3b9bcSJacob 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)); 282a13bc4e3SShao-Ching Huang 2839566063dSJacob Faibussowitsch if (rank == 0) PetscCall(PetscMalloc1(size*6,&grloc)); 284a13bc4e3SShao-Ching Huang rloc[0] = info.xs; 285a13bc4e3SShao-Ching Huang rloc[1] = info.xm; 286a13bc4e3SShao-Ching Huang rloc[2] = info.ys; 287a13bc4e3SShao-Ching Huang rloc[3] = info.ym; 288a13bc4e3SShao-Ching Huang rloc[4] = info.zs; 289a13bc4e3SShao-Ching Huang rloc[5] = info.zm; 2909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Gather(rloc,6,MPIU_INT,&grloc[0][0],6,MPIU_INT,0,comm)); 291a13bc4e3SShao-Ching Huang 292a13bc4e3SShao-Ching Huang /* Write XML header */ 293a13bc4e3SShao-Ching Huang maxnnodes = 0; /* Used for the temporary array size on rank 0 */ 294ea2d7708SPatrick Sanan maxbs = 0; /* Used for the temporary array size on rank 0 */ 295a13bc4e3SShao-Ching Huang boffset = 0; /* Offset into binary file */ 296a13bc4e3SShao-Ching Huang for (r=0; r<size; r++) { 297a13bc4e3SShao-Ching Huang PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0; 298dd400576SPatrick Sanan if (rank == 0) { 299a13bc4e3SShao-Ching Huang xs = grloc[r][0]; 300a13bc4e3SShao-Ching Huang xm = grloc[r][1]; 301a13bc4e3SShao-Ching Huang ys = grloc[r][2]; 302a13bc4e3SShao-Ching Huang ym = grloc[r][3]; 303a13bc4e3SShao-Ching Huang zs = grloc[r][4]; 304a13bc4e3SShao-Ching Huang zm = grloc[r][5]; 305a13bc4e3SShao-Ching Huang nnodes = xm*ym*zm; 306a13bc4e3SShao-Ching Huang } 307a13bc4e3SShao-Ching Huang maxnnodes = PetscMax(maxnnodes,nnodes); 30863a3b9bcSJacob 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)); 3099566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <Coordinates>\n")); 31063a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"Xcoord\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,boffset)); 311a13bc4e3SShao-Ching Huang boffset += xm*sizeof(PetscScalar) + sizeof(int); 31263a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"Ycoord\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,boffset)); 313a13bc4e3SShao-Ching Huang boffset += ym*sizeof(PetscScalar) + sizeof(int); 31463a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"Zcoord\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,boffset)); 315a13bc4e3SShao-Ching Huang boffset += zm*sizeof(PetscScalar) + sizeof(int); 3169566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," </Coordinates>\n")); 3179566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <PointData Scalars=\"ScalarPointData\">\n")); 318a13bc4e3SShao-Ching Huang for (link=vtk->link; link; link=link->next) { 319a13bc4e3SShao-Ching Huang Vec X = (Vec)link->vec; 320fb5bd1c2SPatrick Sanan PetscInt bs,f; 321ea2d7708SPatrick Sanan DM daCurr; 322fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 323fb5bd1c2SPatrick Sanan const char *vecname = "Unnamed"; 324ea2d7708SPatrick Sanan 3259566063dSJacob Faibussowitsch PetscCall(VecGetDM(X,&daCurr)); 3269566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL)); 327ea2d7708SPatrick Sanan maxbs = PetscMax(maxbs,bs); 328ea2d7708SPatrick Sanan if (((PetscObject)X)->name || link != vtk->link) { 3299566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)X,&vecname)); 330a13bc4e3SShao-Ching Huang } 331ea2d7708SPatrick Sanan 332fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 3339566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr,&fieldsnamed)); 334fb5bd1c2SPatrick Sanan if (fieldsnamed) { 335fb5bd1c2SPatrick Sanan for (f=0; f<bs; f++) { 336fb5bd1c2SPatrick Sanan char buf[256]; 337fb5bd1c2SPatrick Sanan const char *fieldname; 3389566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldName(daCurr,f,&fieldname)); 339fb5bd1c2SPatrick Sanan if (!fieldname) { 34063a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(buf,sizeof(buf),"%" PetscInt_FMT,f)); 341fb5bd1c2SPatrick Sanan fieldname = buf; 342fb5bd1c2SPatrick Sanan } 34363a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s.%s\" NumberOfComponents=\"1\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,vecname,fieldname,boffset)); 344fb5bd1c2SPatrick Sanan boffset += nnodes*sizeof(PetscScalar) + sizeof(int); 345fb5bd1c2SPatrick Sanan } 346fb5bd1c2SPatrick Sanan } else { 34763a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <DataArray type=\"%s\" Name=\"%s\" NumberOfComponents=\"%" PetscInt_FMT "\" format=\"appended\" offset=\"%" PetscInt_FMT "\" />\n",precision,vecname,bs,boffset)); 348ea2d7708SPatrick Sanan boffset += bs*nnodes*sizeof(PetscScalar) + sizeof(int); 349a13bc4e3SShao-Ching Huang } 350fb5bd1c2SPatrick Sanan } 3519566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," </PointData>\n")); 3529566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," </Piece>\n")); 353a13bc4e3SShao-Ching Huang } 3549566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," </RectilinearGrid>\n")); 3559566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp," <AppendedData encoding=\"raw\">\n")); 3569566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"_")); 357a13bc4e3SShao-Ching Huang 358a13bc4e3SShao-Ching Huang /* Now write the arrays. */ 359a13bc4e3SShao-Ching Huang tag = ((PetscObject)viewer)->tag; 3609566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(PetscMax(maxnnodes*maxbs,info.xm+info.ym+info.zm),&array,PetscMax(maxnnodes*maxbs,info.xm+info.ym+info.zm),&array2)); 361a13bc4e3SShao-Ching Huang for (r=0; r<size; r++) { 362a13bc4e3SShao-Ching Huang MPI_Status status; 363a13bc4e3SShao-Ching Huang PetscInt xs=-1,xm=-1,ys=-1,ym=-1,zs=-1,zm=-1,nnodes = 0; 364dd400576SPatrick Sanan if (rank == 0) { 365a13bc4e3SShao-Ching Huang xs = grloc[r][0]; 366a13bc4e3SShao-Ching Huang xm = grloc[r][1]; 367a13bc4e3SShao-Ching Huang ys = grloc[r][2]; 368a13bc4e3SShao-Ching Huang ym = grloc[r][3]; 369a13bc4e3SShao-Ching Huang zs = grloc[r][4]; 370a13bc4e3SShao-Ching Huang zm = grloc[r][5]; 371a13bc4e3SShao-Ching Huang nnodes = xm*ym*zm; 372a13bc4e3SShao-Ching Huang } else if (r == rank) { 373a13bc4e3SShao-Ching Huang nnodes = info.xm*info.ym*info.zm; 374a13bc4e3SShao-Ching Huang } 375a13bc4e3SShao-Ching Huang { /* Write the coordinates */ 376a13bc4e3SShao-Ching Huang Vec Coords; 377ea2d7708SPatrick Sanan 3789566063dSJacob Faibussowitsch PetscCall(DMGetCoordinates(da,&Coords)); 379a13bc4e3SShao-Ching Huang if (Coords) { 380a13bc4e3SShao-Ching Huang const PetscScalar *coords; 3819566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Coords,&coords)); 382dd400576SPatrick Sanan if (rank == 0) { 383a13bc4e3SShao-Ching Huang if (r) { 384a13bc4e3SShao-Ching Huang PetscMPIInt nn; 3859566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array,xm+ym+zm,MPIU_SCALAR,r,tag,comm,&status)); 3869566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status,MPIU_SCALAR,&nn)); 38708401ef6SPierre Jolivet PetscCheck(nn == xm+ym+zm,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch"); 388a13bc4e3SShao-Ching Huang } else { 389a13bc4e3SShao-Ching Huang /* x coordinates */ 390a13bc4e3SShao-Ching Huang for (j=0, k=0, i=0; i<xm; i++) { 391a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 392a13bc4e3SShao-Ching Huang array[i] = coords[Iloc*dim + 0]; 393a13bc4e3SShao-Ching Huang } 394a13bc4e3SShao-Ching Huang /* y coordinates */ 395a13bc4e3SShao-Ching Huang for (i=0, k=0, j=0; j<ym; j++) { 396a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 397a13bc4e3SShao-Ching Huang array[j+xm] = dim > 1 ? coords[Iloc*dim + 1] : 0; 398a13bc4e3SShao-Ching Huang } 399a13bc4e3SShao-Ching Huang /* z coordinates */ 400a13bc4e3SShao-Ching Huang for (i=0, j=0, k=0; k<zm; k++) { 401a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 402a13bc4e3SShao-Ching Huang array[k+xm+ym] = dim > 2 ? coords[Iloc*dim + 2] : 0; 403a13bc4e3SShao-Ching Huang } 404a13bc4e3SShao-Ching Huang } 405a13bc4e3SShao-Ching Huang } else if (r == rank) { 406a13bc4e3SShao-Ching Huang xm = info.xm; 407a13bc4e3SShao-Ching Huang ym = info.ym; 408a13bc4e3SShao-Ching Huang zm = info.zm; 409a13bc4e3SShao-Ching Huang for (j=0, k=0, i=0; i<xm; i++) { 410a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 411a13bc4e3SShao-Ching Huang array2[i] = coords[Iloc*dim + 0]; 412a13bc4e3SShao-Ching Huang } 413a13bc4e3SShao-Ching Huang for (i=0, k=0, j=0; j<ym; j++) { 414a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 415a13bc4e3SShao-Ching Huang array2[j+xm] = dim > 1 ? coords[Iloc*dim + 1] : 0; 416a13bc4e3SShao-Ching Huang } 417a13bc4e3SShao-Ching Huang for (i=0, j=0, k=0; k<zm; k++) { 418a13bc4e3SShao-Ching Huang PetscInt Iloc = i+xm*(j+ym*k); 419a13bc4e3SShao-Ching Huang array2[k+xm+ym] = dim > 2 ? coords[Iloc*dim + 2] : 0; 420a13bc4e3SShao-Ching Huang } 4219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void*)array2,xm+ym+zm,MPIU_SCALAR,0,tag,comm)); 422a13bc4e3SShao-Ching Huang } 4239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Coords,&coords)); 424a13bc4e3SShao-Ching Huang } else { /* Fabricate some coordinates using grid index */ 4253d3eaba7SBarry Smith for (i=0; i<xm; i++) { 426a13bc4e3SShao-Ching Huang array[i] = xs+i; 427a13bc4e3SShao-Ching Huang } 4283d3eaba7SBarry Smith for (j=0; j<ym; j++) { 429a13bc4e3SShao-Ching Huang array[j+xm] = ys+j; 430a13bc4e3SShao-Ching Huang } 4313d3eaba7SBarry Smith for (k=0; k<zm; k++) { 432a13bc4e3SShao-Ching Huang array[k+xm+ym] = zs+k; 433a13bc4e3SShao-Ching Huang } 434a13bc4e3SShao-Ching Huang } 435dd400576SPatrick Sanan if (rank == 0) { 4369566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer,fp,&(array[0]),xm,MPIU_SCALAR)); 4379566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer,fp,&(array[xm]),ym,MPIU_SCALAR)); 4389566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer,fp,&(array[xm+ym]),zm,MPIU_SCALAR)); 439a13bc4e3SShao-Ching Huang } 440a13bc4e3SShao-Ching Huang } 441a13bc4e3SShao-Ching Huang 442a13bc4e3SShao-Ching Huang /* Write each of the objects queued up for this file */ 443a13bc4e3SShao-Ching Huang for (link=vtk->link; link; link=link->next) { 444a13bc4e3SShao-Ching Huang Vec X = (Vec)link->vec; 445a13bc4e3SShao-Ching Huang const PetscScalar *x; 446fb5bd1c2SPatrick Sanan PetscInt bs,f; 447ea2d7708SPatrick Sanan DM daCurr; 448fb5bd1c2SPatrick Sanan PetscBool fieldsnamed; 4499566063dSJacob Faibussowitsch PetscCall(VecGetDM(X,&daCurr)); 4509566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(daCurr,NULL,NULL,NULL,NULL,NULL,NULL,NULL,&bs,NULL,NULL,NULL,NULL,NULL)); 451a13bc4e3SShao-Ching Huang 4529566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 453dd400576SPatrick Sanan if (rank == 0) { 454a13bc4e3SShao-Ching Huang if (r) { 455a13bc4e3SShao-Ching Huang PetscMPIInt nn; 4569566063dSJacob Faibussowitsch PetscCallMPI(MPI_Recv(array,nnodes*bs,MPIU_SCALAR,r,tag,comm,&status)); 4579566063dSJacob Faibussowitsch PetscCallMPI(MPI_Get_count(&status,MPIU_SCALAR,&nn)); 45863a3b9bcSJacob Faibussowitsch PetscCheck(nn == nnodes*bs,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Array size mismatch receiving from rank %" PetscInt_FMT,r); 459*1baa6e33SBarry Smith } else PetscCall(PetscArraycpy(array,x,nnodes*bs)); 460fb5bd1c2SPatrick Sanan /* If any fields are named, add scalar fields. Otherwise, add a vector field */ 4619566063dSJacob Faibussowitsch PetscCall(DMDAGetFieldsNamed(daCurr,&fieldsnamed)); 462fb5bd1c2SPatrick Sanan if (fieldsnamed) { 463fb5bd1c2SPatrick Sanan for (f=0; f<bs; f++) { 464fb5bd1c2SPatrick Sanan /* Extract and transpose the f'th field */ 465fb5bd1c2SPatrick Sanan for (k=0; k<zm; k++) { 466fb5bd1c2SPatrick Sanan for (j=0; j<ym; j++) { 467fb5bd1c2SPatrick Sanan for (i=0; i<xm; i++) { 468fb5bd1c2SPatrick Sanan PetscInt Iloc = i+xm*(j+ym*k); 469fb5bd1c2SPatrick Sanan array2[Iloc] = array[Iloc*bs + f]; 470fb5bd1c2SPatrick Sanan } 471fb5bd1c2SPatrick Sanan } 472fb5bd1c2SPatrick Sanan } 4739566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer,fp,array2,nnodes,MPIU_SCALAR)); 474fb5bd1c2SPatrick Sanan } 475fb5bd1c2SPatrick Sanan } 4769566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKFWrite(viewer,fp,array,nnodes*bs,MPIU_SCALAR)); 477ea2d7708SPatrick Sanan 478a13bc4e3SShao-Ching Huang } else if (r == rank) { 4799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Send((void*)x,nnodes*bs,MPIU_SCALAR,0,tag,comm)); 480a13bc4e3SShao-Ching Huang } 4819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 482a13bc4e3SShao-Ching Huang } 483a13bc4e3SShao-Ching Huang } 4849566063dSJacob Faibussowitsch PetscCall(PetscFree2(array,array2)); 4859566063dSJacob Faibussowitsch PetscCall(PetscFree(grloc)); 486a13bc4e3SShao-Ching Huang 4879566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"\n </AppendedData>\n")); 4889566063dSJacob Faibussowitsch PetscCall(PetscFPrintf(comm,fp,"</VTKFile>\n")); 4899566063dSJacob Faibussowitsch PetscCall(PetscFClose(comm,fp)); 490a13bc4e3SShao-Ching Huang PetscFunctionReturn(0); 491a13bc4e3SShao-Ching Huang } 492a13bc4e3SShao-Ching Huang 4934061b8bfSJed Brown /*@C 4944061b8bfSJed Brown DMDAVTKWriteAll - Write a file containing all the fields that have been provided to the viewer 4954061b8bfSJed Brown 4964061b8bfSJed Brown Collective 4974061b8bfSJed Brown 4984165533cSJose E. Roman Input Parameters: 4994165533cSJose E. Roman + odm - DM specifying the grid layout, passed as a PetscObject 5004165533cSJose E. Roman - viewer - viewer of type VTK 5014061b8bfSJed Brown 5024061b8bfSJed Brown Level: developer 5034061b8bfSJed Brown 504fb5bd1c2SPatrick Sanan Notes: 5054061b8bfSJed Brown This function is a callback used by the VTK viewer to actually write the file. 5064061b8bfSJed 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. 5074061b8bfSJed Brown Instead, metadata for the entire file needs to be available up-front before you can start writing the file. 5084061b8bfSJed Brown 509fb5bd1c2SPatrick Sanan If any fields have been named (see e.g. DMDASetFieldName()), then individual scalar 510fb5bd1c2SPatrick Sanan fields are written. Otherwise, a single multi-dof (vector) field is written. 511fb5bd1c2SPatrick Sanan 512db781477SPatrick Sanan .seealso: `PETSCVIEWERVTK` 5134061b8bfSJed Brown @*/ 5144061b8bfSJed Brown PetscErrorCode DMDAVTKWriteAll(PetscObject odm,PetscViewer viewer) 5154061b8bfSJed Brown { 5164061b8bfSJed Brown DM dm = (DM)odm; 5174061b8bfSJed Brown PetscBool isvtk; 5184061b8bfSJed Brown 5194061b8bfSJed Brown PetscFunctionBegin; 5204061b8bfSJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 5214061b8bfSJed Brown PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 5229566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk)); 5237a8be351SBarry Smith PetscCheck(isvtk,PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_INCOMP,"Cannot use viewer type %s",((PetscObject)viewer)->type_name); 5244061b8bfSJed Brown switch (viewer->format) { 5254061b8bfSJed Brown case PETSC_VIEWER_VTK_VTS: 5269566063dSJacob Faibussowitsch PetscCall(DMDAVTKWriteAll_VTS(dm,viewer)); 5274061b8bfSJed Brown break; 528a13bc4e3SShao-Ching Huang case PETSC_VIEWER_VTK_VTR: 5299566063dSJacob Faibussowitsch PetscCall(DMDAVTKWriteAll_VTR(dm,viewer)); 530a13bc4e3SShao-Ching Huang break; 53198921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No support for format '%s'",PetscViewerFormats[viewer->format]); 5324061b8bfSJed Brown } 5334061b8bfSJed Brown PetscFunctionReturn(0); 5344061b8bfSJed Brown } 535