Actual source code: vtkv.c

petsc-3.8.4 2018-03-24
Report Typos and Errors
  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 Arguments:
 25: + viewer - VTK viewer
 26: . dm - DM on which Vec lives
 27: . PetscViewerVTKWriteFunction - function to write this Vec
 28: . fieldtype - Either PETSC_VTK_POINT_FIELD or PETSC_VTK_CELL_FIELD
 29: - vec - Vec to write

 31:    Level: developer

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

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

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

 50: static PetscErrorCode PetscViewerDestroy_VTK(PetscViewer viewer)
 51: {
 52:   PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;
 53:   PetscErrorCode  ierr;

 56:   PetscFree(vtk->filename);
 57:   PetscFree(vtk);
 58:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);
 59:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);
 60:   PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerVTKAddField_C",NULL);
 61:   return(0);
 62: }

 64: static PetscErrorCode PetscViewerFlush_VTK(PetscViewer viewer)
 65: {
 66:   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
 67:   PetscErrorCode           ierr;
 68:   PetscViewerVTKObjectLink link,next;

 71:   if (vtk->link && (!vtk->dm || !vtk->write)) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_WRONGSTATE,"No fields or no grid");
 72:   if (vtk->write) {(*vtk->write)(vtk->dm,viewer);}
 73:   for (link=vtk->link; link; link=next) {
 74:     next = link->next;
 75:     PetscObjectDestroy(&link->vec);
 76:     PetscFree(link);
 77:   }
 78:   PetscObjectDestroy(&vtk->dm);
 79:   vtk->write = NULL;
 80:   return(0);
 81: }

 83: PetscErrorCode  PetscViewerFileSetName_VTK(PetscViewer viewer,const char name[])
 84: {
 85:   PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;
 86:   PetscErrorCode  ierr;
 87:   PetscBool       isvtk,isvts,isvtu,isvtr;
 88:   size_t          len;

 91:   PetscViewerFlush(viewer);
 92:   PetscFree(vtk->filename);
 93:   PetscStrlen(name,&len);
 94:   PetscStrcasecmp(name+len-4,".vtk",&isvtk);
 95:   PetscStrcasecmp(name+len-4,".vts",&isvts);
 96:   PetscStrcasecmp(name+len-4,".vtu",&isvtu);
 97:   PetscStrcasecmp(name+len-4,".vtr",&isvtr);
 98:   if (isvtk) {
 99:     if (viewer->format == PETSC_VIEWER_DEFAULT) {PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_VTK);}
100:     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]);
101:   } else if (isvts) {
102:     if (viewer->format == PETSC_VIEWER_DEFAULT) {PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTS);}
103:     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]);
104:   } else if (isvtu) {
105:     if (viewer->format == PETSC_VIEWER_DEFAULT) {PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTU);}
106:     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]);
107:   } else if (isvtr) {
108:     if (viewer->format == PETSC_VIEWER_DEFAULT) {PetscViewerPushFormat(viewer,PETSC_VIEWER_VTK_VTR);}
109:     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]);
110:   } else SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_UNKNOWN_TYPE,"File '%s' has unrecognized extension",name);
111:   PetscStrallocpy(name,&vtk->filename);
112:   return(0);
113: }

115: PetscErrorCode  PetscViewerFileSetMode_VTK(PetscViewer viewer,PetscFileMode type)
116: {
117:   PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;

121:   vtk->btype = type;
122:   return(0);
123: }

125: PetscErrorCode  PetscViewerVTKAddField_VTK(PetscViewer viewer,PetscObject dm,PetscErrorCode (*PetscViewerVTKWriteFunction)(PetscObject,PetscViewer),PetscViewerVTKFieldType fieldtype,PetscObject vec)
126: {
127:   PetscViewer_VTK          *vtk = (PetscViewer_VTK*)viewer->data;
128:   PetscViewerVTKObjectLink link, tail = vtk->link;
129:   PetscErrorCode           ierr;

132:   if (vtk->dm) {
133:     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");
134:   } else {
135:     PetscObjectReference(dm);
136:     vtk->dm = dm;
137:   }
138:   vtk->write = PetscViewerVTKWriteFunction;
139:   PetscNew(&link);
140:   link->ft   = fieldtype;
141:   link->vec  = vec;
142:   link->next = NULL;
143:   /* Append to list */
144:   if (tail) {
145:     while (tail->next) tail = tail->next;
146:     tail->next = link;
147:   } else vtk->link = link;
148:   return(0);
149: }

151: /*MC
152:    PETSCVIEWERVTK - A viewer that writes to an VTK file


155: .seealso:  PetscViewerVTKOpen(), PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
156:            PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
157:            PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
158:            PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()

160:   Level: beginner
161: M*/

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

169:   PetscNewLog(v,&vtk);

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

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

183: /*@C
184:    PetscViewerVTKOpen - Opens a file for VTK output.

186:    Collective on MPI_Comm

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

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

199:    Level: beginner

201:    Note:
202:    This PetscViewer should be destroyed with PetscViewerDestroy().

204:    Concepts: VTK files
205:    Concepts: PetscViewer^creating

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

216:   PetscViewerCreate(comm,vtk);
217:   PetscViewerSetType(*vtk,PETSCVIEWERVTK);
218:   PetscViewerFileSetMode(*vtk,type);
219:   PetscViewerFileSetName(*vtk,name);
220:   return(0);
221: }

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

226:    Logically collective on PetscViewer

228:    Input Parameters:
229: +  viewer - logically collective viewer, data written from rank 0
230: .  fp - file pointer valid on rank 0
231: .  data - data pointer valid on rank 0
232: .  n - number of data items
233: -  dtype - data type

235:    Level: developer

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

239:    Concepts: VTK files
240:    Concepts: PetscViewer^creating

242: .seealso: DMDAVTKWriteAll(), DMComplexVTKWriteAll(), PetscViewerPushFormat(), PetscViewerVTKOpen(), PetscBinaryWrite()
243: @*/
244: PetscErrorCode PetscViewerVTKFWrite(PetscViewer viewer,FILE *fp,const void *data,PetscInt n,PetscDataType dtype)
245: {
247:   PetscMPIInt    rank;
248: #if defined(PETSC_USE_REAL___FLOAT128)
249:   PetscInt       i;
250:   double         *tmp = NULL;
251:   PetscReal      *ttmp = (PetscReal*)data;
252: #endif

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

290:     count = fwrite(&bytes,sizeof(int),1,fp);
291:     if (count != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_WRITE,"Error writing byte count");
292:     count = fwrite(data,size,(size_t)n,fp);
293:     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);
294: #if defined(PETSC_USE_REAL___FLOAT128)
295:     PetscFree(tmp);
296: #endif
297:   }
298:   return(0);
299: }