Actual source code: vtkv.c

  1: #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h>

  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*/

 19: /*@C
 20:    PetscViewerVTKAddField - Add a field to the viewer

 22:    Collective

 24:    Input Parameters:
 25: + viewer - VTK viewer
 26: . dm - DM on which Vec lives
 27: . PetscViewerVTKWriteFunction - function to write this Vec
 28: . fieldnum - which field of the DM to write (PETSC_DEFAULT if the whle vector should be written)
 29: . fieldtype - Either PETSC_VTK_POINT_FIELD or PETSC_VTK_CELL_FIELD
 30: . checkdm - whether to check for identical dm arguments as fields are added
 31: - vec - Vec from which to write

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

 36:    Level: developer

 38: .seealso: PetscViewerVTKOpen(), DMDAVTKWriteAll(), PetscViewerVTKWriteFunction, PetscViewerVTKGetDM()
 39: @*/
 40: PetscErrorCode PetscViewerVTKAddField(PetscViewer viewer,PetscObject dm,PetscErrorCode (*PetscViewerVTKWriteFunction)(PetscObject,PetscViewer),PetscInt fieldnum,PetscViewerVTKFieldType fieldtype,PetscBool checkdm,PetscObject vec)
 41: {
 45:   PetscUseMethod(viewer,"PetscViewerVTKAddField_C",(PetscViewer,PetscObject,PetscErrorCode (*)(PetscObject,PetscViewer),PetscInt,PetscViewerVTKFieldType,PetscBool,PetscObject),(viewer,dm,PetscViewerVTKWriteFunction,fieldnum,fieldtype,checkdm,vec));
 46:   return 0;
 47: }

 49: /*@C
 50:    PetscViewerVTKGetDM - get the DM associated with the viewer

 52:    Collective

 54:    Input Parameters:
 55: + viewer - VTK viewer
 56: - dm - DM associated with the viewer (as PetscObject)

 58:    Level: developer

 60: .seealso: PetscViewerVTKOpen(), DMDAVTKWriteAll(), PetscViewerVTKWriteFunction, PetscViewerVTKAddField()
 61: @*/
 62: PetscErrorCode PetscViewerVTKGetDM(PetscViewer viewer,PetscObject *dm)
 63: {
 65:   PetscUseMethod(viewer,"PetscViewerVTKGetDM_C",(PetscViewer,PetscObject*),(viewer,dm));
 66:   return 0;
 67: }

 69: static PetscErrorCode PetscViewerDestroy_VTK(PetscViewer viewer)
 70: {
 71:   PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;

 73:   PetscFree(vtk->filename);
 74:   PetscFree(vtk);
 75:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);
 76:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);
 77:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);
 78:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetMode_C",NULL);
 79:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerVTKAddField_C",NULL);
 80:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerVTKGetDM_C",NULL);
 81:   return 0;
 82: }

 84: static PetscErrorCode PetscViewerFlush_VTK(PetscViewer viewer)
 85: {
 86:   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
 87:   PetscViewerVTKObjectLink link,next;

 90:   if (vtk->write) (*vtk->write)(vtk->dm,viewer);
 91:   for (link=vtk->link; link; link=next) {
 92:     next = link->next;
 93:     PetscObjectDestroy(&link->vec);
 94:     PetscFree(link);
 95:   }
 96:   PetscObjectDestroy(&vtk->dm);
 97:   vtk->write = NULL;
 98:   vtk->link  = NULL;
 99:   return 0;
100: }

102: PetscErrorCode  PetscViewerFileSetName_VTK(PetscViewer viewer,const char name[])
103: {
104:   PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;
105:   PetscBool       isvtk,isvts,isvtu,isvtr;
106:   size_t          len;

108:   PetscViewerFlush(viewer);
109:   PetscFree(vtk->filename);
110:   PetscStrlen(name,&len);
111:   if (!len) {
112:     isvtk = PETSC_TRUE;
113:   } else {
114:     PetscStrcasecmp(name+len-4,".vtk",&isvtk);
115:     PetscStrcasecmp(name+len-4,".vts",&isvts);
116:     PetscStrcasecmp(name+len-4,".vtu",&isvtu);
117:     PetscStrcasecmp(name+len-4,".vtr",&isvtr);
118:   }
119:   if (isvtk) {
120:     if (viewer->format == PETSC_VIEWER_DEFAULT) viewer->format = PETSC_VIEWER_ASCII_VTK_DEPRECATED;
122:   } else if (isvts) {
123:     if (viewer->format == PETSC_VIEWER_DEFAULT) viewer->format = PETSC_VIEWER_VTK_VTS;
125:   } else if (isvtu) {
126:     if (viewer->format == PETSC_VIEWER_DEFAULT) viewer->format = PETSC_VIEWER_VTK_VTU;
128:   } else if (isvtr) {
129:     if (viewer->format == PETSC_VIEWER_DEFAULT) viewer->format = PETSC_VIEWER_VTK_VTR;
131:   } else SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_UNKNOWN_TYPE,"File '%s' has unrecognized extension",name);
132:   PetscStrallocpy(len ? name : "stdout",&vtk->filename);
133:   return 0;
134: }

136: PetscErrorCode  PetscViewerFileGetName_VTK(PetscViewer viewer,const char **name)
137: {
138:   PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;
139:   *name = vtk->filename;
140:   return 0;
141: }

143: PetscErrorCode  PetscViewerFileSetMode_VTK(PetscViewer viewer,PetscFileMode type)
144: {
145:   PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;

147:   vtk->btype = type;
148:   return 0;
149: }

151: PetscErrorCode  PetscViewerFileGetMode_VTK(PetscViewer viewer,PetscFileMode *type)
152: {
153:   PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;

155:   *type = vtk->btype;
156:   return 0;
157: }

159: PetscErrorCode  PetscViewerVTKAddField_VTK(PetscViewer viewer,PetscObject dm,PetscErrorCode (*PetscViewerVTKWriteFunction)(PetscObject,PetscViewer),PetscInt fieldnum,PetscViewerVTKFieldType fieldtype,PetscBool checkdm,PetscObject vec)
160: {
161:   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
162:   PetscViewerVTKObjectLink link, tail = vtk->link;

164:   if (vtk->dm) {
166:   } else {
167:     PetscObjectReference(dm);
168:     vtk->dm = dm;
169:   }
170:   vtk->write  = PetscViewerVTKWriteFunction;
171:   PetscNew(&link);
172:   link->ft    = fieldtype;
173:   link->vec   = vec;
174:   link->field = fieldnum;
175:   link->next  = NULL;
176:   /* Append to list */
177:   if (tail) {
178:     while (tail->next) tail = tail->next;
179:     tail->next = link;
180:   } else vtk->link = link;
181:   return 0;
182: }

184: PetscErrorCode PetscViewerVTKGetDM_VTK(PetscViewer viewer,PetscObject *dm)
185: {
186:   PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;

188:   *dm = vtk->dm;
189:   return 0;
190: }

192: /*MC
193:    PETSCVIEWERVTK - A viewer that writes to a VTK file

195: .seealso:  PetscViewerVTKOpen(), PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
196:            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
197:            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
198:            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()

200:   Level: beginner
201: M*/

203: PETSC_EXTERN PetscErrorCode PetscViewerCreate_VTK(PetscViewer v)
204: {
205:   PetscViewer_VTK *vtk;

207:   PetscNewLog(v,&vtk);

209:   v->data         = (void*)vtk;
210:   v->ops->destroy = PetscViewerDestroy_VTK;
211:   v->ops->flush   = PetscViewerFlush_VTK;
212:   vtk->btype      = FILE_MODE_UNDEFINED;
213:   vtk->filename   = NULL;

215:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_VTK);
216:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_VTK);
217:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_VTK);
218:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_VTK);
219:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerVTKAddField_C",PetscViewerVTKAddField_VTK);
220:   PetscObjectComposeFunction((PetscObject)v,"PetscViewerVTKGetDM_C",PetscViewerVTKGetDM_VTK);
221:   return 0;
222: }

224: /*@C
225:    PetscViewerVTKOpen - Opens a file for VTK output.

227:    Collective

229:    Input Parameters:
230: +  comm - MPI communicator
231: .  name - name of file
232: -  type - type of file
233: $    FILE_MODE_WRITE - create new file for binary output
234: $    FILE_MODE_READ - open existing file for binary input (not currently supported)
235: $    FILE_MODE_APPEND - open existing file for binary output (not currently supported)

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

240:    Level: beginner

242:    Note:
243:    This PetscViewer should be destroyed with PetscViewerDestroy().

245: .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(),
246:           VecView(), MatView(), VecLoad(), MatLoad(),
247:           PetscFileMode, PetscViewer
248: @*/
249: PetscErrorCode PetscViewerVTKOpen(MPI_Comm comm,const char name[],PetscFileMode type,PetscViewer *vtk)
250: {
251:   PetscViewerCreate(comm,vtk);
252:   PetscViewerSetType(*vtk,PETSCVIEWERVTK);
253:   PetscViewerFileSetMode(*vtk,type);
254:   PetscViewerFileSetName(*vtk,name);
255:   return 0;
256: }

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

261:    Logically collective on PetscViewer

263:    Input Parameters:
264: +  viewer - logically collective viewer, data written from rank 0
265: .  fp - file pointer valid on rank 0
266: .  data - data pointer valid on rank 0
267: .  n - number of data items
268: -  dtype - data type

270:    Level: developer

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

275: .seealso: DMDAVTKWriteAll(), DMPlexVTKWriteAll(), PetscViewerPushFormat(), PetscViewerVTKOpen(), PetscBinaryWrite()
276: @*/
277: PetscErrorCode PetscViewerVTKFWrite(PetscViewer viewer,FILE *fp,const void *data,PetscInt n,MPI_Datatype dtype)
278: {
279:   PetscMPIInt    rank;
280:   MPI_Datatype   vdtype=dtype;
281: #if defined(PETSC_USE_REAL___FLOAT128)
282:   double         *tmp;
283:   PetscInt       i;
284:   PetscReal      *ttmp = (PetscReal*)data;
285: #endif

288:   if (!n) return 0;
289:   MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);
290:   if (rank == 0) {
291:     size_t      count;
292:     PetscMPIInt dsize;
293:     PetscVTKInt bytes;

295: #if defined(PETSC_USE_REAL___FLOAT128)
296:     if (dtype == MPIU___FLOAT128) {
297:       PetscMalloc1(n,&tmp);
298:       for (i=0; i<n; i++) tmp[i] = ttmp[i];
299:       data  = (void*) tmp;
300:       vdtype = MPI_DOUBLE;
301:     }
302: #endif
303:     MPI_Type_size(vdtype,&dsize);
304:     bytes = PetscVTKIntCast(dsize*n);

306:     count = fwrite(&bytes,sizeof(int),1,fp);
308:     count = fwrite(data,dsize,(size_t)n,fp);
310: #if defined(PETSC_USE_REAL___FLOAT128)
311:     if (dtype == MPIU___FLOAT128) {
312:       PetscFree(tmp);
313:     }
314: #endif
315:   }
316:   return 0;
317: }