Actual source code: vtkv.c
petsc-3.10.5 2019-03-28
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: . checkdm - whether to check for identical dm arguments as fields are added
30: - vec - Vec to write
32: Note:
33: This routine keeps exclusive ownership of the Vec. The caller should not use or destroy the Vec after adding it.
35: Level: developer
37: .seealso: PetscViewerVTKOpen(), DMDAVTKWriteAll(), PetscViewerVTKWriteFunction, PetscViewerVTKGetDM()
38: @*/
39: PetscErrorCode PetscViewerVTKAddField(PetscViewer viewer,PetscObject dm,PetscErrorCode (*PetscViewerVTKWriteFunction)(PetscObject,PetscViewer),PetscViewerVTKFieldType fieldtype,PetscBool checkdm,PetscObject vec)
40: {
47: PetscUseMethod(viewer,"PetscViewerVTKAddField_C",(PetscViewer,PetscObject,PetscErrorCode (*)(PetscObject,PetscViewer),PetscViewerVTKFieldType,PetscBool,PetscObject),(viewer,dm,PetscViewerVTKWriteFunction,fieldtype,checkdm,vec));
48: return(0);
49: }
51: /*@C
52: PetscViewerVTKGetDM - get the DM associated with the viewer
54: Collective
56: Input Arguments:
57: + viewer - VTK viewer
58: - dm - DM associated with the viewer (as PetscObject)
60: Level: developer
62: .seealso: PetscViewerVTKOpen(), DMDAVTKWriteAll(), PetscViewerVTKWriteFunction, PetscViewerVTKAddField()
63: @*/
64: PetscErrorCode PetscViewerVTKGetDM(PetscViewer viewer,PetscObject *dm)
65: {
70: PetscUseMethod(viewer,"PetscViewerVTKGetDM_C",(PetscViewer,PetscObject*),(viewer,dm));
71: return(0);
72: }
74: static PetscErrorCode PetscViewerDestroy_VTK(PetscViewer viewer)
75: {
76: PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;
77: PetscErrorCode ierr;
80: PetscFree(vtk->filename);
81: PetscFree(vtk);
82: PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetName_C",NULL);
83: PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetName_C",NULL);
84: PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileSetMode_C",NULL);
85: PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerFileGetMode_C",NULL);
86: PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerVTKAddField_C",NULL);
87: PetscObjectComposeFunction((PetscObject)viewer,"PetscViewerVTKGetDM_C",NULL);
88: return(0);
89: }
91: static PetscErrorCode PetscViewerFlush_VTK(PetscViewer viewer)
92: {
93: PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;
94: PetscErrorCode ierr;
95: PetscViewerVTKObjectLink link,next;
98: if (vtk->link && (!vtk->dm || !vtk->write)) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_WRONGSTATE,"No fields or no grid");
99: if (vtk->write) {(*vtk->write)(vtk->dm,viewer);}
100: for (link=vtk->link; link; link=next) {
101: next = link->next;
102: PetscObjectDestroy(&link->vec);
103: PetscFree(link);
104: }
105: PetscObjectDestroy(&vtk->dm);
106: vtk->write = NULL;
107: return(0);
108: }
110: PetscErrorCode PetscViewerFileSetName_VTK(PetscViewer viewer,const char name[])
111: {
112: PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;
113: PetscErrorCode ierr;
114: PetscBool isvtk,isvts,isvtu,isvtr;
115: size_t len;
118: PetscViewerFlush(viewer);
119: PetscFree(vtk->filename);
120: PetscStrlen(name,&len);
121: if (!len) {
122: isvtk = PETSC_TRUE;
123: } else {
124: PetscStrcasecmp(name+len-4,".vtk",&isvtk);
125: PetscStrcasecmp(name+len-4,".vts",&isvts);
126: PetscStrcasecmp(name+len-4,".vtu",&isvtu);
127: PetscStrcasecmp(name+len-4,".vtr",&isvtr);
128: }
129: if (isvtk) {
130: if (viewer->format == PETSC_VIEWER_DEFAULT) viewer->format = PETSC_VIEWER_ASCII_VTK;
131: 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]);
132: } else if (isvts) {
133: if (viewer->format == PETSC_VIEWER_DEFAULT) viewer->format = PETSC_VIEWER_VTK_VTS;
134: 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]);
135: } else if (isvtu) {
136: if (viewer->format == PETSC_VIEWER_DEFAULT) viewer->format = PETSC_VIEWER_VTK_VTU;
137: 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]);
138: } else if (isvtr) {
139: if (viewer->format == PETSC_VIEWER_DEFAULT) viewer->format = PETSC_VIEWER_VTK_VTR;
140: 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]);
141: } else SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_UNKNOWN_TYPE,"File '%s' has unrecognized extension",name);
142: PetscStrallocpy(len ? name : "stdout",&vtk->filename);
143: return(0);
144: }
146: PetscErrorCode PetscViewerFileGetName_VTK(PetscViewer viewer,const char **name)
147: {
148: PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;
150: *name = vtk->filename;
151: return(0);
152: }
154: PetscErrorCode PetscViewerFileSetMode_VTK(PetscViewer viewer,PetscFileMode type)
155: {
156: PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;
159: vtk->btype = type;
160: return(0);
161: }
163: PetscErrorCode PetscViewerFileGetMode_VTK(PetscViewer viewer,PetscFileMode *type)
164: {
165: PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;
168: *type = vtk->btype;
169: return(0);
170: }
172: PetscErrorCode PetscViewerVTKAddField_VTK(PetscViewer viewer,PetscObject dm,PetscErrorCode (*PetscViewerVTKWriteFunction)(PetscObject,PetscViewer),PetscViewerVTKFieldType fieldtype,PetscBool checkdm,PetscObject vec)
173: {
174: PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;
175: PetscViewerVTKObjectLink link, tail = vtk->link;
176: PetscErrorCode ierr;
179: if (vtk->dm) {
180: if (checkdm && dm != vtk->dm) SETERRQ(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_INCOMP,"Refusing to write a field from more than one grid to the same VTK file. Set checkdm = PETSC_FALSE to skip this check.");
181: } else {
182: PetscObjectReference(dm);
183: vtk->dm = dm;
184: }
185: vtk->write = PetscViewerVTKWriteFunction;
186: PetscNew(&link);
187: link->ft = fieldtype;
188: link->vec = vec;
189: link->next = NULL;
190: /* Append to list */
191: if (tail) {
192: while (tail->next) tail = tail->next;
193: tail->next = link;
194: } else vtk->link = link;
195: return(0);
196: }
198: PetscErrorCode PetscViewerVTKGetDM_VTK(PetscViewer viewer,PetscObject *dm)
199: {
200: PetscViewer_VTK *vtk = (PetscViewer_VTK*)viewer->data;
203: *dm = vtk->dm;
204: return(0);
205: }
207: /*MC
208: PETSCVIEWERVTK - A viewer that writes to a VTK file
211: .seealso: PetscViewerVTKOpen(), PetscViewerHDF5Open(), PetscViewerStringSPrintf(), PetscViewerSocketOpen(), PetscViewerDrawOpen(), PETSCVIEWERSOCKET,
212: PetscViewerCreate(), PetscViewerASCIIOpen(), PetscViewerBinaryOpen(), PETSCVIEWERBINARY, PETSCVIEWERDRAW, PETSCVIEWERSTRING,
213: PetscViewerMatlabOpen(), VecView(), DMView(), PetscViewerMatlabPutArray(), PETSCVIEWERASCII, PETSCVIEWERMATLAB,
214: PetscViewerFileSetName(), PetscViewerFileSetMode(), PetscViewerFormat, PetscViewerType, PetscViewerSetType()
216: Level: beginner
217: M*/
219: PETSC_EXTERN PetscErrorCode PetscViewerCreate_VTK(PetscViewer v)
220: {
221: PetscViewer_VTK *vtk;
222: PetscErrorCode ierr;
225: PetscNewLog(v,&vtk);
227: v->data = (void*)vtk;
228: v->ops->destroy = PetscViewerDestroy_VTK;
229: v->ops->flush = PetscViewerFlush_VTK;
230: vtk->btype = (PetscFileMode) -1;
231: vtk->filename = 0;
233: PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetName_C",PetscViewerFileSetName_VTK);
234: PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetName_C",PetscViewerFileGetName_VTK);
235: PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileSetMode_C",PetscViewerFileSetMode_VTK);
236: PetscObjectComposeFunction((PetscObject)v,"PetscViewerFileGetMode_C",PetscViewerFileGetMode_VTK);
237: PetscObjectComposeFunction((PetscObject)v,"PetscViewerVTKAddField_C",PetscViewerVTKAddField_VTK);
238: PetscObjectComposeFunction((PetscObject)v,"PetscViewerVTKGetDM_C",PetscViewerVTKGetDM_VTK);
239: return(0);
240: }
242: /*@C
243: PetscViewerVTKOpen - Opens a file for VTK output.
245: Collective on MPI_Comm
247: Input Parameters:
248: + comm - MPI communicator
249: . name - name of file
250: - type - type of file
251: $ FILE_MODE_WRITE - create new file for binary output
252: $ FILE_MODE_READ - open existing file for binary input (not currently supported)
253: $ FILE_MODE_APPEND - open existing file for binary output (not currently supported)
255: Output Parameter:
256: . vtk - PetscViewer for VTK input/output to use with the specified file
258: Level: beginner
260: Note:
261: This PetscViewer should be destroyed with PetscViewerDestroy().
263: Concepts: VTK files
264: Concepts: PetscViewer^creating
266: .seealso: PetscViewerASCIIOpen(), PetscViewerPushFormat(), PetscViewerDestroy(),
267: VecView(), MatView(), VecLoad(), MatLoad(),
268: PetscFileMode, PetscViewer
269: @*/
270: PetscErrorCode PetscViewerVTKOpen(MPI_Comm comm,const char name[],PetscFileMode type,PetscViewer *vtk)
271: {
275: PetscViewerCreate(comm,vtk);
276: PetscViewerSetType(*vtk,PETSCVIEWERVTK);
277: PetscViewerFileSetMode(*vtk,type);
278: PetscViewerFileSetName(*vtk,name);
279: return(0);
280: }
282: /*@C
283: PetscViewerVTKFWrite - write binary data preceded by 32-bit int length (in bytes), does not do byte swapping.
285: Logically collective on PetscViewer
287: Input Parameters:
288: + viewer - logically collective viewer, data written from rank 0
289: . fp - file pointer valid on rank 0
290: . data - data pointer valid on rank 0
291: . n - number of data items
292: - dtype - data type
294: Level: developer
296: Notes:
297: If PetscScalar is __float128 then the binary files are written in double precision
299: Concepts: VTK files
300: Concepts: PetscViewer^creating
302: .seealso: DMDAVTKWriteAll(), DMComplexVTKWriteAll(), PetscViewerPushFormat(), PetscViewerVTKOpen(), PetscBinaryWrite()
303: @*/
304: PetscErrorCode PetscViewerVTKFWrite(PetscViewer viewer,FILE *fp,const void *data,PetscInt n,MPI_Datatype dtype)
305: {
307: PetscMPIInt rank;
308: MPI_Datatype vdtype=dtype;
309: #if defined(PETSC_USE_REAL___FLOAT128)
310: double *tmp;
311: PetscInt i;
312: PetscReal *ttmp = (PetscReal*)data;
313: #endif
316: if (n < 0) SETERRQ1(PetscObjectComm((PetscObject)viewer),PETSC_ERR_ARG_OUTOFRANGE,"Trying to write a negative amount of data %D",n);
317: if (!n) return(0);
318: MPI_Comm_rank(PetscObjectComm((PetscObject)viewer),&rank);
319: if (!rank) {
320: size_t count;
321: PetscMPIInt dsize;
322: PetscVTKInt bytes;
324: #if defined(PETSC_USE_REAL___FLOAT128)
325: if (dtype == MPIU___FLOAT128) {
326: PetscMalloc1(n,&tmp);
327: for (i=0; i<n; i++) tmp[i] = ttmp[i];
328: data = (void*) tmp;
329: vdtype = MPI_DOUBLE;
330: }
331: #endif
332: MPI_Type_size(vdtype,&dsize);
333: bytes = PetscVTKIntCast(dsize*n);
335: count = fwrite(&bytes,sizeof(int),1,fp);
336: if (count != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FILE_WRITE,"Error writing byte count");
337: count = fwrite(data,dsize,(size_t)n,fp);
338: if ((PetscInt)count != n) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_FILE_WRITE,"Wrote %D/%D array members of size %d",(PetscInt)count,n,dsize);
339: #if defined(PETSC_USE_REAL___FLOAT128)
340: if (dtype == MPIU___FLOAT128) {
341: PetscFree(tmp);
342: }
343: #endif
344: }
345: return(0);
346: }