xref: /petsc/src/dm/impls/da/grvtk.c (revision ce94432eddcd14845bc7e8083b7f8ea723b9bf7d)
1b45d2f2cSJed Brown #include <petsc-private/daimpl.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
15*ce94432eSBarry 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;
284061b8bfSJed Brown #if defined(PETSC_USE_COMPLEX)
293bf036e2SBarry Smith   SETERRQ(comm,PETSC_ERR_SUP,"Complex values not supported");
304061b8bfSJed Brown #endif
31*ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
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;
82*ce94432eSBarry 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;
1194061b8bfSJed Brown   ierr = PetscMalloc2(maxnnodes*PetscMax(3,bs),PetscScalar,&array,maxnnodes*3,PetscScalar,&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__
2244061b8bfSJed Brown #define __FUNCT__ "DMDAVTKWriteAll"
2254061b8bfSJed Brown /*@C
2264061b8bfSJed Brown    DMDAVTKWriteAll - Write a file containing all the fields that have been provided to the viewer
2274061b8bfSJed Brown 
2284061b8bfSJed Brown    Collective
2294061b8bfSJed Brown 
2304061b8bfSJed Brown    Input Arguments:
2314061b8bfSJed Brown    odm - DM specifying the grid layout, passed as a PetscObject
2324061b8bfSJed Brown    viewer - viewer of type VTK
2334061b8bfSJed Brown 
2344061b8bfSJed Brown    Level: developer
2354061b8bfSJed Brown 
2364061b8bfSJed Brown    Note:
2374061b8bfSJed Brown    This function is a callback used by the VTK viewer to actually write the file.
2384061b8bfSJed 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.
2394061b8bfSJed Brown    Instead, metadata for the entire file needs to be available up-front before you can start writing the file.
2404061b8bfSJed Brown 
2414061b8bfSJed Brown .seealso: PETSCVIEWERVTK
2424061b8bfSJed Brown @*/
2434061b8bfSJed Brown PetscErrorCode DMDAVTKWriteAll(PetscObject odm,PetscViewer viewer)
2444061b8bfSJed Brown {
2454061b8bfSJed Brown   DM             dm = (DM)odm;
2464061b8bfSJed Brown   PetscBool      isvtk;
2474061b8bfSJed Brown   PetscErrorCode ierr;
2484061b8bfSJed Brown 
2494061b8bfSJed Brown   PetscFunctionBegin;
2504061b8bfSJed Brown   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2514061b8bfSJed Brown   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
252251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERVTK,&isvtk);CHKERRQ(ierr);
253*ce94432eSBarry Smith   if (!isvtk) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_INCOMP,"Cannot use viewer type %s",((PetscObject)viewer)->type_name);
2544061b8bfSJed Brown   switch (viewer->format) {
2554061b8bfSJed Brown   case PETSC_VIEWER_VTK_VTS:
2564061b8bfSJed Brown     ierr = DMDAVTKWriteAll_VTS(dm,viewer);CHKERRQ(ierr);
2574061b8bfSJed Brown     break;
258*ce94432eSBarry Smith   default: SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"No support for format '%s'",PetscViewerFormats[viewer->format]);
2594061b8bfSJed Brown   }
2604061b8bfSJed Brown   PetscFunctionReturn(0);
2614061b8bfSJed Brown }
262