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