Actual source code: vtkv.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  1: #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h> /*I "petscviewer.h" I*/

  3: /*MC
  4:     PetscViewerVTKWriteFunction - functional form used to provide writer to the PetscViewerVTK

  6:      Synopsis:
  7:      #include <petscviewer.h>
  8:      PetscViewerVTKWriteFunction(PetscObject object,PetscViewer viewer)

 10:      Input Parameters:
 11: +      object - the PETSc object to be written
 12: -      viewer - viewer it is to be written to

 14:    Level: developer

 16: .seealso:   PetscViewerVTKAddField()
 17: M*/

 21: /*@C
 22:    PetscViewerVTKAddField - Add a field to the viewer

 24:    Collective

 26:    Input Arguments:
 27: + viewer - VTK viewer
 28: . dm - DM on which Vec lives
 29: . PetscViewerVTKWriteFunction - function to write this Vec
 30: . fieldtype - Either PETSC_VTK_POINT_FIELD or PETSC_VTK_CELL_FIELD
 31: - vec - Vec to write

 33:    Level: developer

 35:    Note:
 36:    This routine keeps exclusive ownership of the Vec. The caller should not use or destroy the Vec after adding it.

 38: .seealso: PetscViewerVTKOpen(), DMDAVTKWriteAll(), PetscViewerVTKWriteFunction
 39: @*/
 40: PetscErrorCode PetscViewerVTKAddField(PetscViewer viewer,PetscObject dm,PetscErrorCode (*PetscViewerVTKWriteFunction)(PetscObject,PetscViewer),PetscViewerVTKFieldType fieldtype,PetscObject vec)
 41: {

 48:   PetscUseMethod(viewer,"PetscViewerVTKAddField_C",(PetscViewer,PetscObject,PetscErrorCode (*)(PetscObject,PetscViewer),PetscViewerVTKFieldType,PetscObject),(viewer,dm,PetscViewerVTKWriteFunction,fieldtype,vec));
 49:   return(0);
 50: }

 54: static PetscErrorCode PetscViewerDestroy_VTK(PetscViewer viewer)
 55: {
 56:   PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;
 57:   PetscErrorCode  ierr;

 60:   PetscFree(vtk->filename);
 61:   PetscFree(vtk);
 62:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);
 63:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);
 64:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerVTKAddField_C",NULL);
 65:   return(0);
 66: }

 70: static PetscErrorCode PetscViewerFlush_VTK(PetscViewer viewer)
 71: {
 72:   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
 73:   PetscErrorCode           ierr;
 74:   PetscViewerVTKObjectLink link,next;

 77:   if (vtk->link && (!vtk->dm || !vtk->write)) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_WRONGSTATE,"No fields or no grid");
 78:   if (vtk->write) {(*vtk->write)(vtk->dm,viewer);}
 79:   for (link=vtk->link; link; link=next) {
 80:     next = link->next;
 81:     PetscObjectDestroy(&link->vec);
 82:     PetscFree(link);
 83:   }
 84:   PetscObjectDestroy(&vtk->dm);
 85:   vtk->write = NULL;
 86:   return(0);
 87: }

 91: PetscErrorCode  PetscViewerFileSetName_VTK(PetscViewer viewer,const char name[])
 92: {
 93:   PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;
 94:   PetscErrorCode  ierr;
 95:   PetscBool       isvtk,isvts,isvtu,isvtr;
 96:   size_t          len;

 99:   PetscViewerFlush(viewer);
100:   PetscFree(vtk->filename);
101:   PetscStrlen(name,&len);
102:   PetscStrcasecmp(name+len-4,".vtk",&isvtk);
103:   PetscStrcasecmp(name+len-4,".vts",&isvts);
104:   PetscStrcasecmp(name+len-4,".vtu",&isvtu);
105:   PetscStrcasecmp(name+len-4,".vtr",&isvtr);
106:   if (isvtk) {
107:     if (viewer->format == PETSC_VIEWER_DEFAULT) {PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_VTK);}
108:     if (viewer->format != PETSC_VIEWER_ASCII_VTK) SETERRQ2(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_INCOMP,"Cannot use file '%s' with format %s, should have '.vtk' extension",name,PetscViewerFormats[viewer->format]);
109:   } else if (isvts) {
110:     if (viewer->format == PETSC_VIEWER_DEFAULT) {PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTS);}
111:     if (viewer->format != PETSC_VIEWER_VTK_VTS) SETERRQ2(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_INCOMP,"Cannot use file '%s' with format %s, should have '.vts' extension",name,PetscViewerFormats[viewer->format]);
112:   } else if (isvtu) {
113:     if (viewer->format == PETSC_VIEWER_DEFAULT) {PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU);}
114:     if (viewer->format != PETSC_VIEWER_VTK_VTU) SETERRQ2(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_INCOMP,"Cannot use file '%s' with format %s, should have '.vtu' extension",name,PetscViewerFormats[viewer->format]);
115:   } else if (isvtr) {
116:     if (viewer->format == PETSC_VIEWER_DEFAULT) {PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTR);}
117:     if (viewer->format != PETSC_VIEWER_VTK_VTR) SETERRQ2(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_INCOMP,"Cannot use file '%s' with format %s, should have '.vtr' extension",name,PetscViewerFormats[viewer->format]);
118:   } else SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_UNKNOWN_TYPE,"File '%s' has unrecognized extension",name);
119:   PetscStrallocpy(name,&vtk->filename);
120:   return(0);
121: }

125: PetscErrorCode  PetscViewerFileSetMode_VTK(PetscViewer viewer,PetscFileMode type)
126: {
127:   PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;

131:   vtk->btype = type;
132:   return(0);
133: }

137: PetscErrorCode  PetscViewerVTKAddField_VTK(PetscViewer viewer,PetscObject dm,PetscErrorCode (*PetscViewerVTKWriteFunction)(PetscObject,PetscViewer),PetscViewerVTKFieldType fieldtype,PetscObject vec)
138: {
139:   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
140:   PetscViewerVTKObjectLink link, tail = vtk->link;
141:   PetscErrorCode           ierr;

144:   if (vtk->dm) {
145:     if (dm != vtk->dm) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_INCOMP,"Cannot write a field from more than one grid to the same VTK file");
146:   } else {
147:     PetscObjectReference(dm);
148:     vtk->dm = dm;
149:   }
150:   vtk->write = PetscViewerVTKWriteFunction;
151:   PetscMalloc(sizeof(struct _n_PetscViewerVTKObjectLink),&link);
152:   link->ft   = fieldtype;
153:   link->vec  = vec;
154:   link->next = NULL;
155:   /* Append to list */
156:   if (tail) {
157:     while (tail->next) tail = tail->next;
158:     tail->next = link;
159:   } else vtk->link = link;
160:   return(0);
161: }

165: PETSC_EXTERN PetscErrorCode PetscViewerCreate_VTK(PetscViewer v)
166: {
167:   PetscViewer_VTK *vtk;
168:   PetscErrorCode  ierr;

171:   PetscNewLog(v,&vtk);

173:   v->data         = (void*)vtk;
174:   v->ops->destroy = PetscViewerDestroy_VTK;
175:   v->ops->flush   = PetscViewerFlush_VTK;
176:   vtk->btype      = (PetscFileMode) -1;
177:   vtk->filename   = 0;

179:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_VTK);
180:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_VTK);
181:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerVTKAddField_C",PetscViewerVTKAddField_VTK);
182:   return(0);
183: }

187: /*@C
188:    PetscViewerVTKOpen - Opens a file for VTK output.

190:    Collective on MPI_Comm

192:    Input Parameters:
193: +  comm - MPI communicator
194: .  name - name of file
195: -  type - type of file
196: $    FILE_MODE_WRITE - create new file for binary output
197: $    FILE_MODE_READ - open existing file for binary input (not currently supported)
198: $    FILE_MODE_APPEND - open existing file for binary output (not currently supported)

200:    Output Parameter:
201: .  vtk - PetscViewer for VTK input/output to use with the specified file

203:    Level: beginner

205:    Note:
206:    This PetscViewer should be destroyed with PetscViewerDestroy().

208:    Concepts: VTK files
209:    Concepts: PetscViewer^creating

211: .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(),
212:           VecView(), MatView(), VecLoad(), MatLoad(),
213:           PetscFileMode, PetscViewer
214: @*/
215: PetscErrorCode PetscViewerVTKOpen(MPI_Comm comm,const char name[],PetscFileMode type,PetscViewer *vtk)
216: {

220:   PetscViewerCreate(comm,vtk);
221:   PetscViewerSetType(*vtk,PETSCVIEWERVTK);
222:   PetscViewerFileSetMode(*vtk,type);
223:   PetscViewerFileSetName(*vtk,name);
224:   return(0);
225: }

229: /*@C
230:    PetscViewerVTKFWrite - write binary data preceded by 32-bit int length (in bytes), does not do byte swapping.

232:    Logically collective on PetscViewer

234:    Input Parameters:
235: +  viewer - logically collective viewer, data written from rank 0
236: .  fp - file pointer valid on rank 0
237: .  data - data pointer valid on rank 0
238: .  n - number of data items
239: -  dtype - data type

241:    Level: developer

243:    Notes: If PetscScalar is __float128 then the binary files are written in double precision

245:    Concepts: VTK files
246:    Concepts: PetscViewer^creating

248: .seealso: DMDAVTKWriteAll(), DMComplexVTKWriteAll(), PetscViewerPushFormat(), PetscViewerVTKOpen(), PetscBinaryWrite()
249: @*/
250: PetscErrorCode PetscViewerVTKFWrite(PetscViewer viewer,FILE *fp,const void *data,PetscInt n,PetscDataType dtype)
251: {
253:   PetscMPIInt    rank;
254: #if defined(PETSC_USE_REAL___FLOAT128)
255:   PetscInt       i;
256:   double         *tmp = NULL;
257:   PetscReal      *ttmp = (PetscReal*)data;
258: #endif

261:   if (n < 0) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_OUTOFRANGE,"Trying to write a negative amount of data %D",n);
262:   if (!n) return(0);
263:   MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);
264:   if (!rank) {
265:     size_t      count;
266:     PetscInt    size;
267:     PetscVTKInt bytes;
268:     switch (dtype) {
269:     case PETSC_DOUBLE:
270:       size = sizeof(double);
271:       break;
272:     case PETSC_FLOAT:
273:       size = sizeof(float);
274:       break;
275: #if defined(PETSC_USE_REAL___FLOAT128)
276:     case PETSC___FLOAT128:
277:       size = sizeof(double);
278:       PetscMalloc1(n,&tmp);
279:       for (i=0; i<n; i++) tmp[i] = ttmp[i];
280:       data = (void*) tmp;
281:       break;
282: #endif
283:     case PETSC_INT:
284:       size = sizeof(PetscInt);
285:       break;
286:     case PETSC_ENUM:
287:       size = sizeof(PetscEnum);
288:       break;
289:     case PETSC_CHAR:
290:       size = sizeof(char);
291:       break;
292:     default: SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Data type not supported");
293:     }
294:     bytes = PetscVTKIntCast(size*n);

296:     count = fwrite(&bytes,sizeof(int),1,fp);
297:     if (count != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_WRITE,"Error writing byte count");
298:     count = fwrite(data,size,(size_t)n,fp);
299:     if ((PetscInt)count != n) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_FILE_WRITE,"Wrote %D/%D array members of size %D",(PetscInt)count,n,(PetscInt)size);
300: #if defined(PETSC_USE_REAL___FLOAT128)
301:     PetscFree(tmp);
302: #endif
303:   }
304:   return(0);
305: }