Actual source code: vtkv.c
1: #include <../src/sys/classes/viewer/impls/vtk/vtkvimpl.h>
3: /*MC
4: PetscViewerVTKWriteFunction - functional form used to provide a 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: [](sec_viewers), `PETSCVIEWERVTK`, `PetscViewerVTKAddField()`
17: M*/
19: /*@C
20: PetscViewerVTKAddField - Add a field to the viewer
22: Collective
24: Input Parameters:
25: + viewer - `PETSCVIEWERVTK`
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 whole 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: Level: developer
35: Note:
36: This routine keeps exclusive ownership of the `Vec`. The caller should not use or destroy the `Vec` after calling it.
38: .seealso: [](sec_viewers), `PETSCVIEWERVTK`, `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: {
42: PetscFunctionBegin;
46: PetscUseMethod(viewer, "PetscViewerVTKAddField_C", (PetscViewer, PetscObject, PetscErrorCode(*)(PetscObject, PetscViewer), PetscInt, PetscViewerVTKFieldType, PetscBool, PetscObject), (viewer, dm, PetscViewerVTKWriteFunction, fieldnum, fieldtype, checkdm, vec));
47: PetscFunctionReturn(PETSC_SUCCESS);
48: }
50: /*@C
51: PetscViewerVTKGetDM - get the `DM` associated with the `PETSCVIEWERVTK` viewer
53: Collective
55: Input Parameters:
56: + viewer - `PETSCVIEWERVTK` viewer
57: - dm - `DM` associated with the viewer (as a `PetscObject`)
59: Level: developer
61: .seealso: [](sec_viewers), `PETSCVIEWERVTK`, `PetscViewerVTKOpen()`, `DMDAVTKWriteAll()`, `PetscViewerVTKWriteFunction`, `PetscViewerVTKAddField()`
62: @*/
63: PetscErrorCode PetscViewerVTKGetDM(PetscViewer viewer, PetscObject *dm)
64: {
65: PetscFunctionBegin;
67: PetscUseMethod(viewer, "PetscViewerVTKGetDM_C", (PetscViewer, PetscObject *), (viewer, dm));
68: PetscFunctionReturn(PETSC_SUCCESS);
69: }
71: static PetscErrorCode PetscViewerDestroy_VTK(PetscViewer viewer)
72: {
73: PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data;
75: PetscFunctionBegin;
76: PetscCall(PetscFree(vtk->filename));
77: PetscCall(PetscFree(vtk));
78: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetName_C", NULL));
79: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetName_C", NULL));
80: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileSetMode_C", NULL));
81: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerFileGetMode_C", NULL));
82: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerVTKAddField_C", NULL));
83: PetscCall(PetscObjectComposeFunction((PetscObject)viewer, "PetscViewerVTKGetDM_C", NULL));
84: PetscFunctionReturn(PETSC_SUCCESS);
85: }
87: static PetscErrorCode PetscViewerFlush_VTK(PetscViewer viewer)
88: {
89: PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data;
90: PetscViewerVTKObjectLink link, next;
92: PetscFunctionBegin;
93: PetscCheck(!vtk->link || !(!vtk->dm || !vtk->write), PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONGSTATE, "No fields or no grid");
94: if (vtk->write) PetscCall((*vtk->write)(vtk->dm, viewer));
95: for (link = vtk->link; link; link = next) {
96: next = link->next;
97: PetscCall(PetscObjectDestroy(&link->vec));
98: PetscCall(PetscFree(link));
99: }
100: PetscCall(PetscObjectDestroy(&vtk->dm));
101: vtk->write = NULL;
102: vtk->link = NULL;
103: PetscFunctionReturn(PETSC_SUCCESS);
104: }
106: static PetscErrorCode PetscViewerFileSetName_VTK(PetscViewer viewer, const char name[])
107: {
108: PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data;
109: PetscBool isvtk, isvts, isvtu, isvtr;
110: size_t len;
112: PetscFunctionBegin;
113: PetscCall(PetscViewerFlush(viewer));
114: PetscCall(PetscFree(vtk->filename));
115: PetscCall(PetscStrlen(name, &len));
116: if (!len) {
117: isvtk = PETSC_TRUE;
118: } else {
119: PetscCall(PetscStrcasecmp(name + len - 4, ".vtk", &isvtk));
120: PetscCall(PetscStrcasecmp(name + len - 4, ".vts", &isvts));
121: PetscCall(PetscStrcasecmp(name + len - 4, ".vtu", &isvtu));
122: PetscCall(PetscStrcasecmp(name + len - 4, ".vtr", &isvtr));
123: }
124: if (isvtk) {
125: if (viewer->format == PETSC_VIEWER_DEFAULT) viewer->format = PETSC_VIEWER_ASCII_VTK_DEPRECATED;
126: PetscCheck(viewer->format == PETSC_VIEWER_ASCII_VTK_DEPRECATED, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use file '%s' with format %s, should have '.vtk' extension", name, PetscViewerFormats[viewer->format]);
127: } else if (isvts) {
128: if (viewer->format == PETSC_VIEWER_DEFAULT) viewer->format = PETSC_VIEWER_VTK_VTS;
129: PetscCheck(viewer->format == PETSC_VIEWER_VTK_VTS, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use file '%s' with format %s, should have '.vts' extension", name, PetscViewerFormats[viewer->format]);
130: } else if (isvtu) {
131: if (viewer->format == PETSC_VIEWER_DEFAULT) viewer->format = PETSC_VIEWER_VTK_VTU;
132: PetscCheck(viewer->format == PETSC_VIEWER_VTK_VTU, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use file '%s' with format %s, should have '.vtu' extension", name, PetscViewerFormats[viewer->format]);
133: } else if (isvtr) {
134: if (viewer->format == PETSC_VIEWER_DEFAULT) viewer->format = PETSC_VIEWER_VTK_VTR;
135: PetscCheck(viewer->format == PETSC_VIEWER_VTK_VTR, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_INCOMP, "Cannot use file '%s' with format %s, should have '.vtr' extension", name, PetscViewerFormats[viewer->format]);
136: } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_UNKNOWN_TYPE, "File '%s' has unrecognized extension", name);
137: PetscCall(PetscStrallocpy(len ? name : "stdout", &vtk->filename));
138: PetscFunctionReturn(PETSC_SUCCESS);
139: }
141: static PetscErrorCode PetscViewerFileGetName_VTK(PetscViewer viewer, const char **name)
142: {
143: PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data;
145: PetscFunctionBegin;
146: *name = vtk->filename;
147: PetscFunctionReturn(PETSC_SUCCESS);
148: }
150: static PetscErrorCode PetscViewerFileSetMode_VTK(PetscViewer viewer, PetscFileMode type)
151: {
152: PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data;
154: PetscFunctionBegin;
155: vtk->btype = type;
156: PetscFunctionReturn(PETSC_SUCCESS);
157: }
159: static PetscErrorCode PetscViewerFileGetMode_VTK(PetscViewer viewer, PetscFileMode *type)
160: {
161: PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data;
163: PetscFunctionBegin;
164: *type = vtk->btype;
165: PetscFunctionReturn(PETSC_SUCCESS);
166: }
168: static PetscErrorCode PetscViewerVTKAddField_VTK(PetscViewer viewer, PetscObject dm, PetscErrorCode (*PetscViewerVTKWriteFunction)(PetscObject, PetscViewer), PetscInt fieldnum, PetscViewerVTKFieldType fieldtype, PetscBool checkdm, PetscObject vec)
169: {
170: PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data;
171: PetscViewerVTKObjectLink link, tail = vtk->link;
173: PetscFunctionBegin;
174: if (vtk->dm) {
175: PetscCheck(!checkdm || dm == vtk->dm, 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.");
176: } else {
177: PetscCall(PetscObjectReference(dm));
178: vtk->dm = dm;
179: }
180: vtk->write = PetscViewerVTKWriteFunction;
181: PetscCall(PetscNew(&link));
182: link->ft = fieldtype;
183: link->vec = vec;
184: link->field = fieldnum;
185: link->next = NULL;
186: /* Append to list */
187: if (tail) {
188: while (tail->next) tail = tail->next;
189: tail->next = link;
190: } else vtk->link = link;
191: PetscFunctionReturn(PETSC_SUCCESS);
192: }
194: static PetscErrorCode PetscViewerVTKGetDM_VTK(PetscViewer viewer, PetscObject *dm)
195: {
196: PetscViewer_VTK *vtk = (PetscViewer_VTK *)viewer->data;
198: PetscFunctionBegin;
199: *dm = vtk->dm;
200: PetscFunctionReturn(PETSC_SUCCESS);
201: }
203: /*MC
204: PETSCVIEWERVTK - A viewer that writes to a VTK file
206: Level: beginner
208: .seealso: [](sec_viewers), `PetscViewerVTKOpen()`, `PetscViewerHDF5Open()`, `PetscViewerStringSPrintf()`, `PetscViewerSocketOpen()`, `PetscViewerDrawOpen()`, `PETSCVIEWERSOCKET`,
209: `PetscViewerCreate()`, `PetscViewerASCIIOpen()`, `PetscViewerBinaryOpen()`, `PETSCVIEWERBINARY`, `PETSCVIEWERDRAW`, `PETSCVIEWERSTRING`,
210: `PetscViewerMatlabOpen()`, `VecView()`, `DMView()`, `PetscViewerMatlabPutArray()`, `PETSCVIEWERASCII`, `PETSCVIEWERMATLAB`,
211: `PetscViewerFileSetName()`, `PetscViewerFileSetMode()`, `PetscViewerFormat`, `PetscViewerType`, `PetscViewerSetType()`
212: M*/
214: PETSC_EXTERN PetscErrorCode PetscViewerCreate_VTK(PetscViewer v)
215: {
216: PetscViewer_VTK *vtk;
218: PetscFunctionBegin;
219: PetscCall(PetscNew(&vtk));
221: v->data = (void *)vtk;
222: v->ops->destroy = PetscViewerDestroy_VTK;
223: v->ops->flush = PetscViewerFlush_VTK;
224: vtk->btype = FILE_MODE_UNDEFINED;
225: vtk->filename = NULL;
227: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetName_C", PetscViewerFileSetName_VTK));
228: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetName_C", PetscViewerFileGetName_VTK));
229: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileSetMode_C", PetscViewerFileSetMode_VTK));
230: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerFileGetMode_C", PetscViewerFileGetMode_VTK));
231: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerVTKAddField_C", PetscViewerVTKAddField_VTK));
232: PetscCall(PetscObjectComposeFunction((PetscObject)v, "PetscViewerVTKGetDM_C", PetscViewerVTKGetDM_VTK));
233: PetscFunctionReturn(PETSC_SUCCESS);
234: }
236: /*@C
237: PetscViewerVTKOpen - Opens a `PETSCVIEWERVTK` viewer file.
239: Collective
241: Input Parameters:
242: + comm - MPI communicator
243: . name - name of file
244: - type - type of file
245: .vb
246: FILE_MODE_WRITE - create new file for binary output
247: FILE_MODE_READ - open existing file for binary input (not currently supported)
248: FILE_MODE_APPEND - open existing file for binary output (not currently supported)
249: .ve
251: Output Parameter:
252: . vtk - `PetscViewer` for VTK input/output to use with the specified file
254: Level: beginner
256: .seealso: [](sec_viewers), `PETSCVIEWERVTK`, `PetscViewerASCIIOpen()`, `PetscViewerPushFormat()`, `PetscViewerDestroy()`,
257: `VecView()`, `MatView()`, `VecLoad()`, `MatLoad()`,
258: `PetscFileMode`, `PetscViewer`
259: @*/
260: PetscErrorCode PetscViewerVTKOpen(MPI_Comm comm, const char name[], PetscFileMode type, PetscViewer *vtk)
261: {
262: PetscFunctionBegin;
263: PetscCall(PetscViewerCreate(comm, vtk));
264: PetscCall(PetscViewerSetType(*vtk, PETSCVIEWERVTK));
265: PetscCall(PetscViewerFileSetMode(*vtk, type));
266: PetscCall(PetscViewerFileSetName(*vtk, name));
267: PetscFunctionReturn(PETSC_SUCCESS);
268: }
270: /*@C
271: PetscViewerVTKFWrite - write binary data preceded by 32-bit int length (in bytes), does not do byte swapping.
273: Logically Collective
275: Input Parameters:
276: + viewer - logically collective viewer, data written from rank 0
277: . fp - file pointer valid on rank 0
278: . data - data pointer valid on rank 0
279: . n - number of data items
280: - dtype - data type
282: Level: developer
284: Note:
285: If `PetscScalar` is `__float128` then the binary files are written in double precision
287: .seealso: [](sec_viewers), `PETSCVIEWERVTK`, `DMDAVTKWriteAll()`, `DMPlexVTKWriteAll()`, `PetscViewerPushFormat()`, `PetscViewerVTKOpen()`, `PetscBinaryWrite()`
288: @*/
289: PetscErrorCode PetscViewerVTKFWrite(PetscViewer viewer, FILE *fp, const void *data, PetscInt n, MPI_Datatype dtype)
290: {
291: PetscMPIInt rank;
292: MPI_Datatype vdtype = dtype;
293: #if defined(PETSC_USE_REAL___FLOAT128)
294: double *tmp;
295: PetscInt i;
296: PetscReal *ttmp = (PetscReal *)data;
297: #endif
299: PetscFunctionBegin;
300: PetscCheck(n >= 0, PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_OUTOFRANGE, "Trying to write a negative amount of data %" PetscInt_FMT, n);
301: if (!n) PetscFunctionReturn(PETSC_SUCCESS);
302: PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank));
303: if (rank == 0) {
304: size_t count;
305: PetscMPIInt dsize;
306: PetscInt64 bytes;
308: #if defined(PETSC_USE_REAL___FLOAT128)
309: if (dtype == MPIU___FLOAT128) {
310: PetscCall(PetscMalloc1(n, &tmp));
311: for (i = 0; i < n; i++) tmp[i] = ttmp[i];
312: data = (void *)tmp;
313: vdtype = MPI_DOUBLE;
314: }
315: #endif
316: PetscCallMPI(MPI_Type_size(vdtype, &dsize));
317: bytes = (PetscInt64)dsize * n;
319: count = fwrite(&bytes, sizeof(bytes), 1, fp);
320: PetscCheck(count == 1, PETSC_COMM_SELF, PETSC_ERR_FILE_WRITE, "Error writing byte count");
321: count = fwrite(data, dsize, (size_t)n, fp);
322: PetscCheck((PetscInt)count == n, PETSC_COMM_SELF, PETSC_ERR_FILE_WRITE, "Wrote %" PetscInt_FMT "/%" PetscInt_FMT " array members of size %d", (PetscInt)count, n, dsize);
323: #if defined(PETSC_USE_REAL___FLOAT128)
324: if (dtype == MPIU___FLOAT128) PetscCall(PetscFree(tmp));
325: #endif
326: }
327: PetscFunctionReturn(PETSC_SUCCESS);
328: }