Actual source code: swarmpic_view.c
1: #include <petscdmda.h>
2: #include <petsc/private/dmswarmimpl.h>
3: #include "../src/dm/impls/swarm/data_bucket.h"
5: static PetscErrorCode private_PetscViewerCreate_XDMF(MPI_Comm comm, const char filename[], PetscViewer *v)
6: {
7: long int *bytes;
8: PetscContainer container;
9: PetscViewer viewer;
11: PetscFunctionBegin;
12: PetscCall(PetscViewerCreate(comm, &viewer));
13: PetscCall(PetscViewerSetType(viewer, PETSCVIEWERASCII));
14: PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_WRITE));
15: PetscCall(PetscViewerFileSetName(viewer, filename));
17: PetscCall(PetscMalloc1(1, &bytes));
18: bytes[0] = 0;
19: PetscCall(PetscContainerCreate(comm, &container));
20: PetscCall(PetscContainerSetPointer(container, (void *)bytes));
21: PetscCall(PetscObjectCompose((PetscObject)viewer, "XDMFViewerContext", (PetscObject)container));
23: /* write xdmf header */
24: PetscCall(PetscViewerASCIIPrintf(viewer, "<?xml version=\"1.0\" encoding=\"utf-8\"?>\n"));
25: PetscCall(PetscViewerASCIIPrintf(viewer, "<Xdmf xmlns:xi=\"http://www.w3.org/2001/XInclude/\" Version=\"2.99\">\n"));
26: /* write xdmf domain */
27: PetscCall(PetscViewerASCIIPushTab(viewer));
28: PetscCall(PetscViewerASCIIPrintf(viewer, "<Domain>\n"));
29: *v = viewer;
30: PetscFunctionReturn(PETSC_SUCCESS);
31: }
33: static PetscErrorCode private_PetscViewerDestroy_XDMF(PetscViewer *v)
34: {
35: PetscViewer viewer;
36: DM dm = NULL;
37: long int *bytes;
38: PetscContainer container = NULL;
40: PetscFunctionBegin;
41: if (!v) PetscFunctionReturn(PETSC_SUCCESS);
42: viewer = *v;
44: PetscCall(PetscObjectQuery((PetscObject)viewer, "DMSwarm", (PetscObject *)&dm));
45: if (dm) {
46: PetscCall(PetscViewerASCIIPrintf(viewer, "</Grid>\n"));
47: PetscCall(PetscViewerASCIIPopTab(viewer));
48: }
50: /* close xdmf header */
51: PetscCall(PetscViewerASCIIPrintf(viewer, "</Domain>\n"));
52: PetscCall(PetscViewerASCIIPopTab(viewer));
53: PetscCall(PetscViewerASCIIPrintf(viewer, "</Xdmf>\n"));
55: PetscCall(PetscObjectQuery((PetscObject)viewer, "XDMFViewerContext", (PetscObject *)&container));
56: if (container) {
57: PetscCall(PetscContainerGetPointer(container, (void **)&bytes));
58: PetscCall(PetscFree(bytes));
59: PetscCall(PetscContainerDestroy(&container));
60: }
61: PetscCall(PetscViewerDestroy(&viewer));
62: *v = NULL;
63: PetscFunctionReturn(PETSC_SUCCESS);
64: }
66: static PetscErrorCode private_CreateDataFileNameXDMF(const char filename[], char dfilename[])
67: {
68: const char dot_xmf[] = ".xmf";
69: size_t len;
70: char viewername_minus_ext[PETSC_MAX_PATH_LEN];
71: PetscBool flg;
73: PetscFunctionBegin;
74: PetscCall(PetscStrendswith(filename, dot_xmf, &flg));
75: PetscCheck(flg, PETSC_COMM_SELF, PETSC_ERR_SUP, "File extension must be %s", dot_xmf);
76: PetscCall(PetscStrncpy(viewername_minus_ext, filename, sizeof(viewername_minus_ext)));
77: PetscCall(PetscStrlen(filename, &len));
78: len -= sizeof(dot_xmf) - 1;
79: if (sizeof(viewername_minus_ext) > len) viewername_minus_ext[len] = '\0';
80: PetscCall(PetscSNPrintf(dfilename, PETSC_MAX_PATH_LEN - 1, "%s_swarm_fields.pbin", viewername_minus_ext));
81: PetscFunctionReturn(PETSC_SUCCESS);
82: }
84: static PetscErrorCode private_DMSwarmView_XDMF(DM dm, PetscViewer viewer)
85: {
86: PetscBool isswarm = PETSC_FALSE;
87: const char *viewername;
88: char datafile[PETSC_MAX_PATH_LEN];
89: char *datafilename;
90: PetscViewer fviewer;
91: PetscInt k, ng, dim;
92: Vec dvec;
93: long int *bytes = NULL;
94: PetscContainer container = NULL;
95: const char *dmname, *coordname;
97: PetscFunctionBegin;
98: PetscCall(PetscObjectQuery((PetscObject)viewer, "XDMFViewerContext", (PetscObject *)&container));
99: if (container) {
100: PetscCall(PetscContainerGetPointer(container, (void **)&bytes));
101: } else SETERRQ(PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Valid to find attached data XDMFViewerContext");
103: PetscCall(PetscObjectTypeCompare((PetscObject)dm, DMSWARM, &isswarm));
104: PetscCheck(isswarm, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Only valid for DMSwarm");
106: PetscCall(PetscObjectCompose((PetscObject)viewer, "DMSwarm", (PetscObject)dm));
108: PetscCall(PetscViewerASCIIPushTab(viewer));
109: PetscCall(PetscObjectGetName((PetscObject)dm, &dmname));
110: if (!dmname) PetscCall(DMGetOptionsPrefix(dm, &dmname));
111: if (!dmname) {
112: PetscCall(PetscViewerASCIIPrintf(viewer, "<Grid Name=\"DMSwarm\" GridType=\"Uniform\">\n"));
113: } else {
114: PetscCall(PetscViewerASCIIPrintf(viewer, "<Grid Name=\"DMSwarm[%s]\" GridType=\"Uniform\">\n", dmname));
115: }
117: /* create a sub-viewer for topology, geometry and all data fields */
118: /* name is viewer.name + "_swarm_fields.pbin" */
119: PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)viewer), &fviewer));
120: PetscCall(PetscViewerSetType(fviewer, PETSCVIEWERBINARY));
121: PetscCall(PetscViewerBinarySetSkipHeader(fviewer, PETSC_TRUE));
122: PetscCall(PetscViewerBinarySetSkipInfo(fviewer, PETSC_TRUE));
123: PetscCall(PetscViewerFileSetMode(fviewer, FILE_MODE_WRITE));
125: PetscCall(PetscViewerFileGetName(viewer, &viewername));
126: PetscCall(private_CreateDataFileNameXDMF(viewername, datafile));
127: PetscCall(PetscViewerFileSetName(fviewer, datafile));
128: PetscCall(PetscStrrchr(datafile, '/', &datafilename));
130: PetscCall(DMSwarmGetSize(dm, &ng));
132: /* write topology header */
133: PetscCall(PetscViewerASCIIPushTab(viewer));
134: PetscCall(PetscViewerASCIIPrintf(viewer, "<Topology Dimensions=\"%" PetscInt_FMT "\" TopologyType=\"Mixed\">\n", ng));
135: PetscCall(PetscViewerASCIIPushTab(viewer));
136: PetscCall(PetscViewerASCIIPrintf(viewer, "<DataItem Format=\"Binary\" Endian=\"Big\" DataType=\"Int\" Dimensions=\"%" PetscInt_FMT "\" Seek=\"%ld\">\n", ng * 3, bytes[0]));
137: PetscCall(PetscViewerASCIIPushTab(viewer));
138: PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", datafilename));
139: PetscCall(PetscViewerASCIIPopTab(viewer));
140: PetscCall(PetscViewerASCIIPrintf(viewer, "</DataItem>\n"));
141: PetscCall(PetscViewerASCIIPopTab(viewer));
142: PetscCall(PetscViewerASCIIPrintf(viewer, "</Topology>\n"));
143: PetscCall(PetscViewerASCIIPopTab(viewer));
145: /* write topology data */
146: for (k = 0; k < ng; k++) {
147: PetscInt pvertex[3];
149: pvertex[0] = 1;
150: pvertex[1] = 1;
151: pvertex[2] = k;
152: PetscCall(PetscViewerBinaryWrite(fviewer, pvertex, 3, PETSC_INT));
153: }
154: bytes[0] += sizeof(PetscInt) * ng * 3;
156: /* write geometry header */
157: PetscCall(PetscViewerASCIIPushTab(viewer));
158: PetscCall(DMGetDimension(dm, &dim));
159: switch (dim) {
160: case 1:
161: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for 1D");
162: case 2:
163: PetscCall(PetscViewerASCIIPrintf(viewer, "<Geometry Type=\"XY\">\n"));
164: break;
165: case 3:
166: PetscCall(PetscViewerASCIIPrintf(viewer, "<Geometry Type=\"XYZ\">\n"));
167: break;
168: }
169: PetscCall(PetscViewerASCIIPushTab(viewer));
170: PetscCall(PetscViewerASCIIPrintf(viewer, "<DataItem Format=\"Binary\" Endian=\"Big\" DataType=\"Float\" Precision=\"8\" Dimensions=\"%" PetscInt_FMT " %" PetscInt_FMT "\" Seek=\"%ld\">\n", ng, dim, bytes[0]));
171: PetscCall(PetscViewerASCIIPushTab(viewer));
172: PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", datafilename));
173: PetscCall(PetscViewerASCIIPopTab(viewer));
174: PetscCall(PetscViewerASCIIPrintf(viewer, "</DataItem>\n"));
175: PetscCall(PetscViewerASCIIPopTab(viewer));
176: PetscCall(PetscViewerASCIIPrintf(viewer, "</Geometry>\n"));
177: PetscCall(PetscViewerASCIIPopTab(viewer));
179: /* write geometry data */
180: PetscCall(DMSwarmGetCoordinateField(dm, &coordname));
181: PetscCall(DMSwarmCreateGlobalVectorFromField(dm, coordname, &dvec));
182: PetscCall(VecView(dvec, fviewer));
183: PetscCall(DMSwarmDestroyGlobalVectorFromField(dm, coordname, &dvec));
184: bytes[0] += sizeof(PetscReal) * ng * dim;
186: PetscCall(PetscViewerDestroy(&fviewer));
187: PetscFunctionReturn(PETSC_SUCCESS);
188: }
190: static PetscErrorCode private_VecView_Swarm_XDMF(Vec x, PetscViewer viewer)
191: {
192: long int *bytes = NULL;
193: PetscContainer container = NULL;
194: const char *viewername;
195: char datafile[PETSC_MAX_PATH_LEN];
196: char *datafilename;
197: PetscViewer fviewer;
198: PetscInt N, bs;
199: const char *vecname;
200: char fieldname[PETSC_MAX_PATH_LEN];
202: PetscFunctionBegin;
203: PetscCall(PetscObjectQuery((PetscObject)viewer, "XDMFViewerContext", (PetscObject *)&container));
204: PetscCheck(container, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unable to find attached data XDMFViewerContext");
205: PetscCall(PetscContainerGetPointer(container, (void **)&bytes));
206: PetscCall(PetscViewerFileGetName(viewer, &viewername));
207: PetscCall(private_CreateDataFileNameXDMF(viewername, datafile));
209: /* re-open a sub-viewer for all data fields */
210: /* name is viewer.name + "_swarm_fields.pbin" */
211: PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)viewer), &fviewer));
212: PetscCall(PetscViewerSetType(fviewer, PETSCVIEWERBINARY));
213: PetscCall(PetscViewerBinarySetSkipHeader(fviewer, PETSC_TRUE));
214: PetscCall(PetscViewerBinarySetSkipInfo(fviewer, PETSC_TRUE));
215: PetscCall(PetscViewerFileSetMode(fviewer, FILE_MODE_APPEND));
216: PetscCall(PetscViewerFileSetName(fviewer, datafile));
217: PetscCall(PetscStrrchr(datafile, '/', &datafilename));
219: PetscCall(VecGetSize(x, &N));
220: PetscCall(VecGetBlockSize(x, &bs));
221: N = N / bs;
222: PetscCall(PetscObjectGetName((PetscObject)x, &vecname));
223: if (!vecname) {
224: PetscCall(PetscSNPrintf(fieldname, PETSC_MAX_PATH_LEN - 1, "swarmfield_%d", ((PetscObject)x)->tag));
225: } else {
226: PetscCall(PetscSNPrintf(fieldname, PETSC_MAX_PATH_LEN - 1, "%s", vecname));
227: }
229: /* write data header */
230: PetscCall(PetscViewerASCIIPushTab(viewer));
231: PetscCall(PetscViewerASCIIPrintf(viewer, "<Attribute Center=\"Node\" Name=\"%s\" Type=\"None\">\n", fieldname));
232: PetscCall(PetscViewerASCIIPushTab(viewer));
233: if (bs == 1) {
234: PetscCall(PetscViewerASCIIPrintf(viewer, "<DataItem Format=\"Binary\" Endian=\"Big\" DataType=\"Float\" Precision=\"8\" Dimensions=\"%" PetscInt_FMT "\" Seek=\"%ld\">\n", N, bytes[0]));
235: } else {
236: PetscCall(PetscViewerASCIIPrintf(viewer, "<DataItem Format=\"Binary\" Endian=\"Big\" DataType=\"Float\" Precision=\"8\" Dimensions=\"%" PetscInt_FMT " %" PetscInt_FMT "\" Seek=\"%ld\">\n", N, bs, bytes[0]));
237: }
238: PetscCall(PetscViewerASCIIPushTab(viewer));
239: PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", datafilename));
240: PetscCall(PetscViewerASCIIPopTab(viewer));
241: PetscCall(PetscViewerASCIIPrintf(viewer, "</DataItem>\n"));
242: PetscCall(PetscViewerASCIIPopTab(viewer));
243: PetscCall(PetscViewerASCIIPrintf(viewer, "</Attribute>\n"));
244: PetscCall(PetscViewerASCIIPopTab(viewer));
246: /* write data */
247: PetscCall(VecView(x, fviewer));
248: bytes[0] += sizeof(PetscReal) * N * bs;
250: PetscCall(PetscViewerDestroy(&fviewer));
251: PetscFunctionReturn(PETSC_SUCCESS);
252: }
254: static PetscErrorCode private_ISView_Swarm_XDMF(IS is, PetscViewer viewer)
255: {
256: long int *bytes = NULL;
257: PetscContainer container = NULL;
258: const char *viewername;
259: char datafile[PETSC_MAX_PATH_LEN];
260: char *datafilename;
261: PetscViewer fviewer;
262: PetscInt N, bs;
263: const char *vecname;
264: char fieldname[PETSC_MAX_PATH_LEN];
266: PetscFunctionBegin;
267: PetscCall(PetscObjectQuery((PetscObject)viewer, "XDMFViewerContext", (PetscObject *)&container));
268: PetscCheck(container, PetscObjectComm((PetscObject)viewer), PETSC_ERR_SUP, "Unable to find attached data XDMFViewerContext");
269: PetscCall(PetscContainerGetPointer(container, (void **)&bytes));
270: PetscCall(PetscViewerFileGetName(viewer, &viewername));
271: PetscCall(private_CreateDataFileNameXDMF(viewername, datafile));
273: /* re-open a sub-viewer for all data fields */
274: /* name is viewer.name + "_swarm_fields.pbin" */
275: PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)viewer), &fviewer));
276: PetscCall(PetscViewerSetType(fviewer, PETSCVIEWERBINARY));
277: PetscCall(PetscViewerBinarySetSkipHeader(fviewer, PETSC_TRUE));
278: PetscCall(PetscViewerBinarySetSkipInfo(fviewer, PETSC_TRUE));
279: PetscCall(PetscViewerFileSetMode(fviewer, FILE_MODE_APPEND));
280: PetscCall(PetscViewerFileSetName(fviewer, datafile));
281: PetscCall(PetscStrrchr(datafile, '/', &datafilename));
283: PetscCall(ISGetSize(is, &N));
284: PetscCall(ISGetBlockSize(is, &bs));
285: N = N / bs;
286: PetscCall(PetscObjectGetName((PetscObject)is, &vecname));
287: if (!vecname) {
288: PetscCall(PetscSNPrintf(fieldname, PETSC_MAX_PATH_LEN - 1, "swarmfield_%d", ((PetscObject)is)->tag));
289: } else {
290: PetscCall(PetscSNPrintf(fieldname, PETSC_MAX_PATH_LEN - 1, "%s", vecname));
291: }
293: /* write data header */
294: PetscCall(PetscViewerASCIIPushTab(viewer));
295: PetscCall(PetscViewerASCIIPrintf(viewer, "<Attribute Center=\"Node\" Name=\"%s\" Type=\"None\">\n", fieldname));
296: PetscCall(PetscViewerASCIIPushTab(viewer));
297: if (bs == 1) {
298: PetscCall(PetscViewerASCIIPrintf(viewer, "<DataItem Format=\"Binary\" Endian=\"Big\" DataType=\"Int\" Precision=\"4\" Dimensions=\"%" PetscInt_FMT "\" Seek=\"%ld\">\n", N, bytes[0]));
299: } else {
300: PetscCall(PetscViewerASCIIPrintf(viewer, "<DataItem Format=\"Binary\" Endian=\"Big\" DataType=\"Int\" Precision=\"4\" Dimensions=\"%" PetscInt_FMT " %" PetscInt_FMT "\" Seek=\"%ld\">\n", N, bs, bytes[0]));
301: }
302: PetscCall(PetscViewerASCIIPushTab(viewer));
303: PetscCall(PetscViewerASCIIPrintf(viewer, "%s\n", datafilename));
304: PetscCall(PetscViewerASCIIPopTab(viewer));
305: PetscCall(PetscViewerASCIIPrintf(viewer, "</DataItem>\n"));
306: PetscCall(PetscViewerASCIIPopTab(viewer));
307: PetscCall(PetscViewerASCIIPrintf(viewer, "</Attribute>\n"));
308: PetscCall(PetscViewerASCIIPopTab(viewer));
310: /* write data */
311: PetscCall(ISView(is, fviewer));
312: bytes[0] += sizeof(PetscInt) * N * bs;
314: PetscCall(PetscViewerDestroy(&fviewer));
315: PetscFunctionReturn(PETSC_SUCCESS);
316: }
318: /*@C
319: DMSwarmViewFieldsXDMF - Write a selection of DMSwarm fields to an XDMF3 file
321: Collective
323: Input Parameters:
324: + dm - the `DMSWARM`
325: . filename - the file name of the XDMF file (must have the extension .xmf)
326: . nfields - the number of fields to write into the XDMF file
327: - field_name_list - array of length nfields containing the textual name of fields to write
329: Level: beginner
331: Note:
332: Only fields registered with data type `PETSC_DOUBLE` or `PETSC_INT` can be written into the file
334: .seealso: `DM`, `DMSWARM`, `DMSwarmViewXDMF()`
335: @*/
336: PetscErrorCode DMSwarmViewFieldsXDMF(DM dm, const char filename[], PetscInt nfields, const char *field_name_list[])
337: {
338: Vec dvec;
339: PetscInt f, N;
340: PetscViewer viewer;
342: PetscFunctionBegin;
343: PetscCall(private_PetscViewerCreate_XDMF(PetscObjectComm((PetscObject)dm), filename, &viewer));
344: PetscCall(private_DMSwarmView_XDMF(dm, viewer));
345: PetscCall(DMSwarmGetLocalSize(dm, &N));
346: for (f = 0; f < nfields; f++) {
347: void *data;
348: PetscDataType type;
350: PetscCall(DMSwarmGetField(dm, field_name_list[f], NULL, &type, &data));
351: PetscCall(DMSwarmRestoreField(dm, field_name_list[f], NULL, &type, &data));
352: if (type == PETSC_DOUBLE) {
353: PetscCall(DMSwarmCreateGlobalVectorFromField(dm, field_name_list[f], &dvec));
354: PetscCall(PetscObjectSetName((PetscObject)dvec, field_name_list[f]));
355: PetscCall(private_VecView_Swarm_XDMF(dvec, viewer));
356: PetscCall(DMSwarmDestroyGlobalVectorFromField(dm, field_name_list[f], &dvec));
357: } else if (type == PETSC_INT) {
358: IS is;
359: const PetscInt *idx;
361: PetscCall(DMSwarmGetField(dm, field_name_list[f], NULL, &type, &data));
362: idx = (const PetscInt *)data;
364: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), N, idx, PETSC_USE_POINTER, &is));
365: PetscCall(PetscObjectSetName((PetscObject)is, field_name_list[f]));
366: PetscCall(private_ISView_Swarm_XDMF(is, viewer));
367: PetscCall(ISDestroy(&is));
368: PetscCall(DMSwarmRestoreField(dm, field_name_list[f], NULL, &type, &data));
369: } else SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Can only write PETSC_INT and PETSC_DOUBLE");
370: }
371: PetscCall(private_PetscViewerDestroy_XDMF(&viewer));
372: PetscFunctionReturn(PETSC_SUCCESS);
373: }
375: /*@
376: DMSwarmViewXDMF - Write `DMSWARM` fields to an XDMF3 file
378: Collective
380: Input Parameters:
381: + dm - the `DMSWARM`
382: - filename - the file name of the XDMF file (must have the extension .xmf)
384: Level: beginner
386: Note:
387: Only fields user registered with data type `PETSC_DOUBLE` or `PETSC_INT` will be written into the file
389: Developer Note:
390: This should be removed and replaced with the standard use of `PetscViewer`
392: .seealso: `DM`, `DMSWARM`, `DMSwarmViewFieldsXDMF()`
393: @*/
394: PetscErrorCode DMSwarmViewXDMF(DM dm, const char filename[])
395: {
396: DM_Swarm *swarm = (DM_Swarm *)dm->data;
397: Vec dvec;
398: PetscInt f;
399: PetscViewer viewer;
401: PetscFunctionBegin;
402: PetscCall(private_PetscViewerCreate_XDMF(PetscObjectComm((PetscObject)dm), filename, &viewer));
403: PetscCall(private_DMSwarmView_XDMF(dm, viewer));
404: for (f = 4; f < swarm->db->nfields; f++) { /* only examine user defined fields - the first 4 are internally created by DMSwarmPIC */
405: DMSwarmDataField field;
407: /* query field type - accept all those of type PETSC_DOUBLE */
408: field = swarm->db->field[f];
409: if (field->petsc_type == PETSC_DOUBLE) {
410: PetscCall(DMSwarmCreateGlobalVectorFromField(dm, field->name, &dvec));
411: PetscCall(PetscObjectSetName((PetscObject)dvec, field->name));
412: PetscCall(private_VecView_Swarm_XDMF(dvec, viewer));
413: PetscCall(DMSwarmDestroyGlobalVectorFromField(dm, field->name, &dvec));
414: } else if (field->petsc_type == PETSC_INT) {
415: IS is;
416: PetscInt N;
417: const PetscInt *idx;
418: void *data;
420: PetscCall(DMSwarmGetLocalSize(dm, &N));
421: PetscCall(DMSwarmGetField(dm, field->name, NULL, NULL, &data));
422: idx = (const PetscInt *)data;
424: PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), N, idx, PETSC_USE_POINTER, &is));
425: PetscCall(PetscObjectSetName((PetscObject)is, field->name));
426: PetscCall(private_ISView_Swarm_XDMF(is, viewer));
427: PetscCall(ISDestroy(&is));
428: PetscCall(DMSwarmRestoreField(dm, field->name, NULL, NULL, &data));
429: }
430: }
431: PetscCall(private_PetscViewerDestroy_XDMF(&viewer));
432: PetscFunctionReturn(PETSC_SUCCESS);
433: }