Actual source code: vtkv.c
petsc-3.6.4 2016-04-12
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) {PetscViewerSetFormat(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) {PetscViewerSetFormat(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) {PetscViewerSetFormat(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) {PetscViewerSetFormat(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: }
147: vtk->dm = dm;
148: vtk->write = PetscViewerVTKWriteFunction;
149: PetscMalloc(sizeof(struct _n_PetscViewerVTKObjectLink),&link);
150: link->ft = fieldtype;
151: link->vec = vec;
152: link->next = NULL;
153: /* Append to list */
154: if (tail) {
155: while (tail->next) tail = tail->next;
156: tail->next = link;
157: } else vtk->link = link;
158: return(0);
159: }
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: }
185: /*@C
186: PetscViewerVTKOpen - Opens a file for VTK output.
188: Collective on MPI_Comm
190: Input Parameters:
191: + comm - MPI communicator
192: . name - name of file
193: - type - type of file
194: $ FILE_MODE_WRITE - create new file for binary output
195: $ FILE_MODE_READ - open existing file for binary input (not currently supported)
196: $ FILE_MODE_APPEND - open existing file for binary output (not currently supported)
198: Output Parameter:
199: . vtk - PetscViewer for VTK input/output to use with the specified file
201: Level: beginner
203: Note:
204: This PetscViewer should be destroyed with PetscViewerDestroy().
206: Concepts: VTK files
207: Concepts: PetscViewer^creating
209: .seealso: PetscViewerASCIIOpen(), PetscViewerSetFormat(), PetscViewerDestroy(),
210: VecView(), MatView(), VecLoad(), MatLoad(),
211: PetscFileMode, PetscViewer
212: @*/
213: PetscErrorCode PetscViewerVTKOpen(MPI_Comm comm,const char name[],PetscFileMode type,PetscViewer *vtk)
214: {
218: PetscViewerCreate(comm,vtk);
219: PetscViewerSetType(*vtk,PETSCVIEWERVTK);
220: PetscViewerFileSetMode(*vtk,type);
221: PetscViewerFileSetName(*vtk,name);
222: return(0);
223: }
227: /*@C
228: PetscViewerVTKFWrite - write binary data preceded by 32-bit int length (in bytes), does not do byte swapping.
230: Logically collective on PetscViewer
232: Input Parameters:
233: + viewer - logically collective viewer, data written from rank 0
234: . fp - file pointer valid on rank 0
235: . data - data pointer valid on rank 0
236: . n - number of data items
237: - dtype - data type
239: Level: developer
241: Notes: If PetscScalar is __float128 then the binary files are written in double precision
243: Concepts: VTK files
244: Concepts: PetscViewer^creating
246: .seealso: DMDAVTKWriteAll(), DMComplexVTKWriteAll(), PetscViewerSetFormat(), PetscViewerVTKOpen(), PetscBinaryWrite()
247: @*/
248: PetscErrorCode PetscViewerVTKFWrite(PetscViewer viewer,FILE *fp,const void *data,PetscInt n,PetscDataType dtype)
249: {
251: PetscMPIInt rank;
252: #if defined(PETSC_USE_REAL___FLOAT128)
253: PetscInt i;
254: double *tmp = NULL;
255: PetscReal *ttmp = (PetscReal*)data;
256: #endif
259: if (n < 0) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_OUTOFRANGE,"Trying to write a negative amount of data %D",n);
260: if (!n) return(0);
261: MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);
262: if (!rank) {
263: size_t count;
264: PetscInt size;
265: PetscVTKInt bytes;
266: switch (dtype) {
267: case PETSC_DOUBLE:
268: size = sizeof(double);
269: break;
270: case PETSC_FLOAT:
271: size = sizeof(float);
272: break;
273: #if defined(PETSC_USE_REAL___FLOAT128)
274: case PETSC___FLOAT128:
275: size = sizeof(double);
276: PetscMalloc1(n,&tmp);
277: for (i=0; i<n; i++) tmp[i] = ttmp[i];
278: data = (void*) tmp;
279: break;
280: #endif
281: case PETSC_INT:
282: size = sizeof(PetscInt);
283: break;
284: case PETSC_ENUM:
285: size = sizeof(PetscEnum);
286: break;
287: case PETSC_CHAR:
288: size = sizeof(char);
289: break;
290: default: SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Data type not supported");
291: }
292: bytes = PetscVTKIntCast(size*n);
294: count = fwrite(&bytes,sizeof(int),1,fp);
295: if (count != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_WRITE,"Error writing byte count");
296: count = fwrite(data,size,(size_t)n,fp);
297: 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);
298: #if defined(PETSC_USE_REAL___FLOAT128)
299: PetscFree(tmp);
300: #endif
301: }
302: return(0);
303: }