xref: /petsc/src/dm/impls/swarm/swarmpic_view.c (revision 3fc23ce8277e00a65c4b1be6e2b4535b56ba7978)
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(0);
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(0);
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(0);
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(0);
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) {
114     PetscCall(DMGetOptionsPrefix(dm,&dmname));
115   }
116   if (!dmname) {
117     PetscCall(PetscViewerASCIIPrintf(viewer,"<Grid Name=\"DMSwarm\" GridType=\"Uniform\">\n"));
118   } else {
119     PetscCall(PetscViewerASCIIPrintf(viewer,"<Grid Name=\"DMSwarm[%s]\" GridType=\"Uniform\">\n",dmname));
120   }
121 
122   /* create a sub-viewer for topology, geometry and all data fields */
123   /* name is viewer.name + "_swarm_fields.pbin" */
124   PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)viewer),&fviewer));
125   PetscCall(PetscViewerSetType(fviewer,PETSCVIEWERBINARY));
126   PetscCall(PetscViewerBinarySetSkipHeader(fviewer,PETSC_TRUE));
127   PetscCall(PetscViewerBinarySetSkipInfo(fviewer,PETSC_TRUE));
128   PetscCall(PetscViewerFileSetMode(fviewer,FILE_MODE_WRITE));
129 
130   PetscCall(PetscViewerFileGetName(viewer,&viewername));
131   PetscCall(private_CreateDataFileNameXDMF(viewername,datafile));
132   PetscCall(PetscViewerFileSetName(fviewer,datafile));
133   PetscCall(PetscStrrchr(datafile,'/',&datafilename));
134 
135   PetscCall(DMSwarmGetSize(dm,&ng));
136 
137   /* write topology header */
138   PetscCall(PetscViewerASCIIPushTab(viewer));
139   PetscCall(PetscViewerASCIIPrintf(viewer,"<Topology Dimensions=\"%" PetscInt_FMT "\" TopologyType=\"Mixed\">\n",ng));
140   PetscCall(PetscViewerASCIIPushTab(viewer));
141   PetscCall(PetscViewerASCIIPrintf(viewer,"<DataItem Format=\"Binary\" Endian=\"Big\" DataType=\"Int\" Dimensions=\"%" PetscInt_FMT "\" Seek=\"%ld\">\n",ng*3,bytes[0]));
142   PetscCall(PetscViewerASCIIPushTab(viewer));
143   PetscCall(PetscViewerASCIIPrintf(viewer,"%s\n",datafilename));
144   PetscCall(PetscViewerASCIIPopTab(viewer));
145   PetscCall(PetscViewerASCIIPrintf(viewer,"</DataItem>\n"));
146   PetscCall(PetscViewerASCIIPopTab(viewer));
147   PetscCall(PetscViewerASCIIPrintf(viewer,"</Topology>\n"));
148   PetscCall(PetscViewerASCIIPopTab(viewer));
149 
150   /* write topology data */
151   for (k=0; k<ng; k++) {
152     PetscInt pvertex[3];
153 
154     pvertex[0] = 1;
155     pvertex[1] = 1;
156     pvertex[2] = k;
157     PetscCall(PetscViewerBinaryWrite(fviewer,pvertex,3,PETSC_INT));
158   }
159   bytes[0] += sizeof(PetscInt) * ng * 3;
160 
161   /* write geometry header */
162   PetscCall(PetscViewerASCIIPushTab(viewer));
163   PetscCall(DMGetDimension(dm,&dim));
164   switch (dim) {
165     case 1:
166       SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for 1D");
167     case 2:
168       PetscCall(PetscViewerASCIIPrintf(viewer,"<Geometry Type=\"XY\">\n"));
169       break;
170     case 3:
171       PetscCall(PetscViewerASCIIPrintf(viewer,"<Geometry Type=\"XYZ\">\n"));
172       break;
173   }
174   PetscCall(PetscViewerASCIIPushTab(viewer));
175   PetscCall(PetscViewerASCIIPrintf(viewer,"<DataItem Format=\"Binary\" Endian=\"Big\" DataType=\"Float\" Precision=\"8\" Dimensions=\"%" PetscInt_FMT " %" PetscInt_FMT "\" Seek=\"%ld\">\n",ng,dim,bytes[0]));
176   PetscCall(PetscViewerASCIIPushTab(viewer));
177   PetscCall(PetscViewerASCIIPrintf(viewer,"%s\n",datafilename));
178   PetscCall(PetscViewerASCIIPopTab(viewer));
179   PetscCall(PetscViewerASCIIPrintf(viewer,"</DataItem>\n"));
180   PetscCall(PetscViewerASCIIPopTab(viewer));
181   PetscCall(PetscViewerASCIIPrintf(viewer,"</Geometry>\n"));
182   PetscCall(PetscViewerASCIIPopTab(viewer));
183 
184   /* write geometry data */
185   PetscCall(DMSwarmCreateGlobalVectorFromField(dm,DMSwarmPICField_coor,&dvec));
186   PetscCall(VecView(dvec,fviewer));
187   PetscCall(DMSwarmDestroyGlobalVectorFromField(dm,DMSwarmPICField_coor,&dvec));
188   bytes[0] += sizeof(PetscReal) * ng * dim;
189 
190   PetscCall(PetscViewerDestroy(&fviewer));
191   PetscFunctionReturn(0);
192 }
193 
194 PetscErrorCode private_VecView_Swarm_XDMF(Vec x,PetscViewer viewer)
195 {
196   long int       *bytes = NULL;
197   PetscContainer container = NULL;
198   const char     *viewername;
199   char           datafile[PETSC_MAX_PATH_LEN];
200   char           *datafilename;
201   PetscViewer    fviewer;
202   PetscInt       N,bs;
203   const char     *vecname;
204   char           fieldname[PETSC_MAX_PATH_LEN];
205 
206   PetscFunctionBegin;
207   PetscCall(PetscObjectQuery((PetscObject)viewer,"XDMFViewerContext",(PetscObject*)&container));
208   PetscCheck(container,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Unable to find attached data XDMFViewerContext");
209   PetscCall(PetscContainerGetPointer(container,(void**)&bytes));
210   PetscCall(PetscViewerFileGetName(viewer,&viewername));
211   PetscCall(private_CreateDataFileNameXDMF(viewername,datafile));
212 
213   /* re-open a sub-viewer for all data fields */
214   /* name is viewer.name + "_swarm_fields.pbin" */
215   PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)viewer),&fviewer));
216   PetscCall(PetscViewerSetType(fviewer,PETSCVIEWERBINARY));
217   PetscCall(PetscViewerBinarySetSkipHeader(fviewer,PETSC_TRUE));
218   PetscCall(PetscViewerBinarySetSkipInfo(fviewer,PETSC_TRUE));
219   PetscCall(PetscViewerFileSetMode(fviewer,FILE_MODE_APPEND));
220   PetscCall(PetscViewerFileSetName(fviewer,datafile));
221   PetscCall(PetscStrrchr(datafile,'/',&datafilename));
222 
223   PetscCall(VecGetSize(x,&N));
224   PetscCall(VecGetBlockSize(x,&bs));
225   N = N/bs;
226   PetscCall(PetscObjectGetName((PetscObject)x,&vecname));
227   if (!vecname) {
228     PetscCall(PetscSNPrintf(fieldname,PETSC_MAX_PATH_LEN-1,"swarmfield_%d",((PetscObject)x)->tag));
229   } else {
230     PetscCall(PetscSNPrintf(fieldname,PETSC_MAX_PATH_LEN-1,"%s",vecname));
231   }
232 
233   /* write data header */
234   PetscCall(PetscViewerASCIIPushTab(viewer));
235   PetscCall(PetscViewerASCIIPrintf(viewer,"<Attribute Center=\"Node\" Name=\"%s\" Type=\"None\">\n",fieldname));
236   PetscCall(PetscViewerASCIIPushTab(viewer));
237   if (bs == 1) {
238     PetscCall(PetscViewerASCIIPrintf(viewer,"<DataItem Format=\"Binary\" Endian=\"Big\" DataType=\"Float\" Precision=\"8\" Dimensions=\"%" PetscInt_FMT "\" Seek=\"%ld\">\n",N,bytes[0]));
239   } else {
240     PetscCall(PetscViewerASCIIPrintf(viewer,"<DataItem Format=\"Binary\" Endian=\"Big\" DataType=\"Float\" Precision=\"8\" Dimensions=\"%" PetscInt_FMT " %" PetscInt_FMT "\" Seek=\"%ld\">\n",N,bs,bytes[0]));
241   }
242   PetscCall(PetscViewerASCIIPushTab(viewer));
243   PetscCall(PetscViewerASCIIPrintf(viewer,"%s\n",datafilename));
244   PetscCall(PetscViewerASCIIPopTab(viewer));
245   PetscCall(PetscViewerASCIIPrintf(viewer,"</DataItem>\n"));
246   PetscCall(PetscViewerASCIIPopTab(viewer));
247   PetscCall(PetscViewerASCIIPrintf(viewer,"</Attribute>\n"));
248   PetscCall(PetscViewerASCIIPopTab(viewer));
249 
250   /* write data */
251   PetscCall(VecView(x,fviewer));
252   bytes[0] += sizeof(PetscReal) * N * bs;
253 
254   PetscCall(PetscViewerDestroy(&fviewer));
255   PetscFunctionReturn(0);
256 }
257 
258 PetscErrorCode private_ISView_Swarm_XDMF(IS is,PetscViewer viewer)
259 {
260   long int       *bytes = NULL;
261   PetscContainer container = NULL;
262   const char     *viewername;
263   char           datafile[PETSC_MAX_PATH_LEN];
264   char           *datafilename;
265   PetscViewer    fviewer;
266   PetscInt       N,bs;
267   const char     *vecname;
268   char           fieldname[PETSC_MAX_PATH_LEN];
269 
270   PetscFunctionBegin;
271   PetscCall(PetscObjectQuery((PetscObject)viewer,"XDMFViewerContext",(PetscObject*)&container));
272   PetscCheck(container,PetscObjectComm((PetscObject)viewer),PETSC_ERR_SUP,"Unable to find attached data XDMFViewerContext");
273   PetscCall(PetscContainerGetPointer(container,(void**)&bytes));
274   PetscCall(PetscViewerFileGetName(viewer,&viewername));
275   PetscCall(private_CreateDataFileNameXDMF(viewername,datafile));
276 
277   /* re-open a sub-viewer for all data fields */
278   /* name is viewer.name + "_swarm_fields.pbin" */
279   PetscCall(PetscViewerCreate(PetscObjectComm((PetscObject)viewer),&fviewer));
280   PetscCall(PetscViewerSetType(fviewer,PETSCVIEWERBINARY));
281   PetscCall(PetscViewerBinarySetSkipHeader(fviewer,PETSC_TRUE));
282   PetscCall(PetscViewerBinarySetSkipInfo(fviewer,PETSC_TRUE));
283   PetscCall(PetscViewerFileSetMode(fviewer,FILE_MODE_APPEND));
284   PetscCall(PetscViewerFileSetName(fviewer,datafile));
285   PetscCall(PetscStrrchr(datafile,'/',&datafilename));
286 
287   PetscCall(ISGetSize(is,&N));
288   PetscCall(ISGetBlockSize(is,&bs));
289   N = N/bs;
290   PetscCall(PetscObjectGetName((PetscObject)is,&vecname));
291   if (!vecname) {
292     PetscCall(PetscSNPrintf(fieldname,PETSC_MAX_PATH_LEN-1,"swarmfield_%d",((PetscObject)is)->tag));
293   } else {
294     PetscCall(PetscSNPrintf(fieldname,PETSC_MAX_PATH_LEN-1,"%s",vecname));
295   }
296 
297   /* write data header */
298   PetscCall(PetscViewerASCIIPushTab(viewer));
299   PetscCall(PetscViewerASCIIPrintf(viewer,"<Attribute Center=\"Node\" Name=\"%s\" Type=\"None\">\n",fieldname));
300   PetscCall(PetscViewerASCIIPushTab(viewer));
301   if (bs == 1) {
302     PetscCall(PetscViewerASCIIPrintf(viewer,"<DataItem Format=\"Binary\" Endian=\"Big\" DataType=\"Int\" Precision=\"4\" Dimensions=\"%" PetscInt_FMT "\" Seek=\"%ld\">\n",N,bytes[0]));
303   } else {
304     PetscCall(PetscViewerASCIIPrintf(viewer,"<DataItem Format=\"Binary\" Endian=\"Big\" DataType=\"Int\" Precision=\"4\" Dimensions=\"%" PetscInt_FMT " %" PetscInt_FMT "\" Seek=\"%ld\">\n",N,bs,bytes[0]));
305   }
306   PetscCall(PetscViewerASCIIPushTab(viewer));
307   PetscCall(PetscViewerASCIIPrintf(viewer,"%s\n",datafilename));
308   PetscCall(PetscViewerASCIIPopTab(viewer));
309   PetscCall(PetscViewerASCIIPrintf(viewer,"</DataItem>\n"));
310   PetscCall(PetscViewerASCIIPopTab(viewer));
311   PetscCall(PetscViewerASCIIPrintf(viewer,"</Attribute>\n"));
312   PetscCall(PetscViewerASCIIPopTab(viewer));
313 
314   /* write data */
315   PetscCall(ISView(is,fviewer));
316   bytes[0] += sizeof(PetscInt) * N * bs;
317 
318   PetscCall(PetscViewerDestroy(&fviewer));
319   PetscFunctionReturn(0);
320 }
321 
322 /*@C
323    DMSwarmViewFieldsXDMF - Write a selection of DMSwarm fields to an XDMF3 file
324 
325    Collective on dm
326 
327    Input parameters:
328 +  dm - the DMSwarm
329 .  filename - the file name of the XDMF file (must have the extension .xmf)
330 .  nfields - the number of fields to write into the XDMF file
331 -  field_name_list - array of length nfields containing the textual name of fields to write
332 
333    Level: beginner
334 
335    Notes:
336    Only fields registered with data type PETSC_DOUBLE or PETSC_INT can be written into the file
337 
338 .seealso: DMSwarmViewXDMF()
339 @*/
340 PETSC_EXTERN PetscErrorCode DMSwarmViewFieldsXDMF(DM dm,const char filename[],PetscInt nfields,const char *field_name_list[])
341 {
342   Vec            dvec;
343   PetscInt       f,N;
344   PetscViewer    viewer;
345 
346   PetscFunctionBegin;
347   PetscCall(private_PetscViewerCreate_XDMF(PetscObjectComm((PetscObject)dm),filename,&viewer));
348   PetscCall(private_DMSwarmView_XDMF(dm,viewer));
349   PetscCall(DMSwarmGetLocalSize(dm,&N));
350   for (f=0; f<nfields; f++) {
351     void          *data;
352     PetscDataType type;
353 
354     PetscCall(DMSwarmGetField(dm,field_name_list[f],NULL,&type,&data));
355     PetscCall(DMSwarmRestoreField(dm,field_name_list[f],NULL,&type,&data));
356     if (type == PETSC_DOUBLE) {
357       PetscCall(DMSwarmCreateGlobalVectorFromField(dm,field_name_list[f],&dvec));
358       PetscCall(PetscObjectSetName((PetscObject)dvec,field_name_list[f]));
359       PetscCall(private_VecView_Swarm_XDMF(dvec,viewer));
360       PetscCall(DMSwarmDestroyGlobalVectorFromField(dm,field_name_list[f],&dvec));
361     } else if (type == PETSC_INT) {
362       IS is;
363       const PetscInt *idx;
364 
365       PetscCall(DMSwarmGetField(dm,field_name_list[f],NULL,&type,&data));
366       idx = (const PetscInt*)data;
367 
368       PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm),N,idx,PETSC_USE_POINTER,&is));
369       PetscCall(PetscObjectSetName((PetscObject)is,field_name_list[f]));
370       PetscCall(private_ISView_Swarm_XDMF(is,viewer));
371       PetscCall(ISDestroy(&is));
372       PetscCall(DMSwarmRestoreField(dm,field_name_list[f],NULL,&type,&data));
373     } else SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Can only write PETSC_INT and PETSC_DOUBLE");
374 
375   }
376   PetscCall(private_PetscViewerDestroy_XDMF(&viewer));
377   PetscFunctionReturn(0);
378 }
379 
380 /*@C
381    DMSwarmViewXDMF - Write DMSwarm fields to an XDMF3 file
382 
383    Collective on dm
384 
385    Input parameters:
386 +  dm - the DMSwarm
387 -  filename - the file name of the XDMF file (must have the extension .xmf)
388 
389    Level: beginner
390 
391    Notes:
392      Only fields user registered with data type PETSC_DOUBLE or PETSC_INT will be written into the file
393 
394    Developer Notes:
395      This should be removed and replaced with the standard use of PetscViewer
396 
397 .seealso: DMSwarmViewFieldsXDMF()
398 @*/
399 PETSC_EXTERN PetscErrorCode DMSwarmViewXDMF(DM dm,const char filename[])
400 {
401   DM_Swarm       *swarm = (DM_Swarm*)dm->data;
402   Vec            dvec;
403   PetscInt       f;
404   PetscViewer    viewer;
405 
406   PetscFunctionBegin;
407   PetscCall(private_PetscViewerCreate_XDMF(PetscObjectComm((PetscObject)dm),filename,&viewer));
408   PetscCall(private_DMSwarmView_XDMF(dm,viewer));
409   for (f=4; f<swarm->db->nfields; f++) { /* only examine user defined fields - the first 4 are internally created by DMSwarmPIC */
410     DMSwarmDataField field;
411 
412     /* query field type - accept all those of type PETSC_DOUBLE */
413     field = swarm->db->field[f];
414     if (field->petsc_type == PETSC_DOUBLE) {
415       PetscCall(DMSwarmCreateGlobalVectorFromField(dm,field->name,&dvec));
416       PetscCall(PetscObjectSetName((PetscObject)dvec,field->name));
417       PetscCall(private_VecView_Swarm_XDMF(dvec,viewer));
418       PetscCall(DMSwarmDestroyGlobalVectorFromField(dm,field->name,&dvec));
419     } else if (field->petsc_type == PETSC_INT) {
420       IS             is;
421       PetscInt       N;
422       const PetscInt *idx;
423       void           *data;
424 
425       PetscCall(DMSwarmGetLocalSize(dm,&N));
426       PetscCall(DMSwarmGetField(dm,field->name,NULL,NULL,&data));
427       idx = (const PetscInt*)data;
428 
429       PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm),N,idx,PETSC_USE_POINTER,&is));
430       PetscCall(PetscObjectSetName((PetscObject)is,field->name));
431       PetscCall(private_ISView_Swarm_XDMF(is,viewer));
432       PetscCall(ISDestroy(&is));
433       PetscCall(DMSwarmRestoreField(dm,field->name,NULL,NULL,&data));
434     }
435   }
436   PetscCall(private_PetscViewerDestroy_XDMF(&viewer));
437   PetscFunctionReturn(0);
438 }
439