xref: /honee/src/spanstats/spanstats.c (revision cf8f12c89a11a3279f97d0071ccf240f989bc25c)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 /// @file
4 /// Functions for setting up and performing spanwise-statistics collection
5 ///
6 /// "Parent" refers to the 2D plane on which statistics are collected *onto*.
7 /// "Child" refers to the 3D domain where statistics are gathered *from*.
8 /// Each quadrature point on the parent plane has several children in the child domain that it performs spanwise averaging with.
9 
10 #include <ceed.h>
11 #include <petscdmplex.h>
12 #include <petscsf.h>
13 #include <spanstats.h>
14 
15 #include <navierstokes.h>
16 #include <petsc_ops.h>
17 
18 PetscErrorCode SpanStatsCtxDestroy(SpanStatsCtx *spanstats) {
19   SpanStatsCtx spanstats_ = *spanstats;
20 
21   PetscFunctionBeginUser;
22   if (spanstats_ == NULL) PetscFunctionReturn(PETSC_SUCCESS);
23   PetscCall(VecDestroy(&spanstats_->Child_Stats_loc));
24   PetscCall(VecDestroy(&spanstats_->Parent_Stats_loc));
25 
26   PetscCall(OperatorApplyContextDestroy(spanstats_->op_stats_collect_ctx));
27   PetscCall(OperatorApplyContextDestroy(spanstats_->op_proj_rhs_ctx));
28   PetscCall(OperatorApplyContextDestroy(spanstats_->mms_error_ctx));
29 
30   PetscCall(KSPDestroy(&spanstats_->ksp));
31   PetscCall(PetscSFDestroy(&spanstats_->sf));
32   PetscCall(DMDestroy(&spanstats_->dm));
33   PetscCall(PetscFree(spanstats_->prefix));
34   PetscCall(PetscFree(spanstats_));
35   *spanstats = NULL;
36   PetscFunctionReturn(PETSC_SUCCESS);
37 }
38 
39 static PetscErrorCode SpanwiseStatisticssCreateDM(Honee honee, SpanStatsCtx spanstats, PetscInt degree, PetscInt num_comps) {
40   PetscReal     domain_min[3], domain_max[3];
41   PetscLogStage stage_stats_setup;
42   MPI_Comm      comm = PetscObjectComm((PetscObject)honee->dm);
43 
44   PetscFunctionBeginUser;
45   PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup));
46   if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup));
47   PetscCall(PetscLogStagePush(stage_stats_setup));
48 
49   spanstats->num_comp_stats = num_comps;
50 
51   // Get spanwise length
52   PetscCall(DMGetBoundingBox(honee->dm, domain_min, domain_max));
53   spanstats->span_width = domain_max[2] - domain_min[2];
54 
55   {  // Get DM from surface
56     DM             parent_distributed_dm;
57     const PetscSF *isoperiodicface;
58     PetscInt       num_isoperiodicface;
59     DMLabel        label;
60     PetscMPIInt    size;
61 
62     PetscCall(DMPlexGetIsoperiodicFaceSF(honee->dm, &num_isoperiodicface, &isoperiodicface));
63 
64     if (isoperiodicface) {
65       PetscSF         inv_isoperiodicface;
66       PetscInt        nleaves, isoperiodicface_index = -1;
67       const PetscInt *ilocal;
68       char            isoperiodicface_key[1024];
69 
70       PetscCall(PetscSNPrintf(isoperiodicface_key, sizeof isoperiodicface_key, "-%sisoperiodic_facesf", spanstats->prefix));
71       PetscCall(PetscOptionsGetInt(NULL, NULL, isoperiodicface_key, &isoperiodicface_index, NULL));
72       isoperiodicface_index = isoperiodicface_index == -1 ? num_isoperiodicface - 1 : isoperiodicface_index;
73       PetscCall(PetscSFCreateInverseSF(isoperiodicface[isoperiodicface_index], &inv_isoperiodicface));
74       PetscCall(PetscSFGetGraph(inv_isoperiodicface, NULL, &nleaves, &ilocal, NULL));
75       PetscCall(DMCreateLabel(honee->dm, "Periodic Face"));
76       PetscCall(DMGetLabel(honee->dm, "Periodic Face", &label));
77       for (PetscInt i = 0; i < nleaves; i++) {
78         PetscCall(DMLabelSetValue(label, ilocal[i], 1));
79       }
80       PetscCall(PetscSFDestroy(&inv_isoperiodicface));
81     } else {
82       PetscCall(DMGetLabel(honee->dm, "Face Sets", &label));
83     }
84 
85     PetscCall(DMPlexLabelComplete(honee->dm, label));
86     PetscCall(DMPlexFilter(honee->dm, label, 1, PETSC_FALSE, PETSC_FALSE, NULL, &spanstats->dm));
87     PetscCall(DMSetCoordinateDisc(spanstats->dm, NULL, PETSC_FALSE, PETSC_TRUE));  // Ensure that a coordinate FE exists
88 
89     PetscCall(DMPlexDistribute(spanstats->dm, 0, NULL, &parent_distributed_dm));
90     PetscCallMPI(MPI_Comm_size(comm, &size));
91     if (parent_distributed_dm) {
92       PetscCall(DMDestroy(&spanstats->dm));
93       spanstats->dm = parent_distributed_dm;
94     } else if (size > 1) {
95       PetscCall(PetscPrintf(comm, "WARNING: Spanwise statistics: parent DM could not be distributed accross %d ranks.\n", size));
96     }
97   }
98   {
99     PetscBool is_simplex = PETSC_FALSE;
100     PetscCall(DMPlexIsSimplex(spanstats->dm, &is_simplex));
101     PetscCheck(is_simplex != PETSC_TRUE, comm, PETSC_ERR_ARG_WRONGSTATE, "Spanwise statistics is not implemented for non-tensor product grids");
102   }
103 
104   PetscCall(PetscObjectSetName((PetscObject)spanstats->dm, "Spanwise_Stats"));
105   PetscCall(DMSetOptionsPrefix(spanstats->dm, spanstats->prefix));
106   PetscCall(DMSetFromOptions(spanstats->dm));
107   PetscCall(DMViewFromOptions(spanstats->dm, NULL, "-dm_view"));
108 
109   // Create FE space for parent DM
110   PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, honee->app_ctx->degree, 1, honee->app_ctx->q_extra, 1, &spanstats->num_comp_stats,
111                                spanstats->dm));
112 
113   PetscCall(PetscLogStagePop());
114   PetscFunctionReturn(PETSC_SUCCESS);
115 }
116 
117 /** @brief Create CeedElemRestriction for collocated data in component-major order.
118 a. Sets the strides of the restriction to component-major order
119  Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction.
120 */
121 static PetscErrorCode CreateElemRestrColloc_CompMajor(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base,
122                                                       CeedElemRestriction *elem_restr_collocated) {
123   CeedInt num_elem_qpts, loc_num_elem;
124 
125   PetscFunctionBeginUser;
126   PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts));
127   PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem));
128 
129   const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp};
130   PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides,
131                                                        elem_restr_collocated));
132   PetscFunctionReturn(PETSC_SUCCESS);
133 }
134 
135 // Get coordinates of quadrature points
136 static PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, Vec *Qx_coords,
137                                           PetscInt *total_nqpnts) {
138   CeedElemRestriction  elem_restr_qx;
139   CeedQFunction        qf_quad_coords;
140   CeedOperator         op_quad_coords;
141   CeedInt              num_comp_x;
142   CeedSize             l_vec_size;
143   OperatorApplyContext op_quad_coords_ctx;
144 
145   PetscFunctionBeginUser;
146   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
147   PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx));
148   PetscCallCeed(ceed, CeedElemRestrictionGetLVectorSize(elem_restr_qx, &l_vec_size));
149   *total_nqpnts = l_vec_size / num_comp_x;
150 
151   // Create QFunction
152   PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords));
153 
154   // Create Operator
155   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords));
156   PetscCallCeed(ceed, CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, x_coords));
157   PetscCallCeed(ceed, CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
158 
159   PetscCall(CeedOperatorCreateLocalVecs(op_quad_coords, DMReturnVecType(dm), PETSC_COMM_SELF, NULL, Qx_coords));
160   PetscCall(OperatorApplyContextCreate(NULL, NULL, ceed, op_quad_coords, CEED_VECTOR_NONE, NULL, NULL, NULL, &op_quad_coords_ctx));
161 
162   PetscCall(ApplyCeedOperatorLocalToLocal(NULL, *Qx_coords, op_quad_coords_ctx));
163 
164   PetscCall(OperatorApplyContextDestroy(op_quad_coords_ctx));
165   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qx));
166   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_quad_coords));
167   PetscCallCeed(ceed, CeedOperatorDestroy(&op_quad_coords));
168   PetscFunctionReturn(PETSC_SUCCESS);
169 }
170 
171 static PetscErrorCode SpanwiseStatisticsSetupDataCreate(Honee honee, SpanStatsCtx spanstats, SpanStatsSetupData *stats_setup_data) {
172   Ceed                ceed = honee->ceed;
173   DM                  dm   = spanstats->dm;
174   Vec                 X_loc;
175   CeedInt             num_comp_x, num_comp_stats = spanstats->num_comp_stats;
176   PetscInt            height = 0, dm_field = 0;
177   CeedElemRestriction elem_restr_q;
178   CeedBasis           basis_q;
179 
180   PetscFunctionBeginUser;
181   PetscCall(PetscNew(stats_setup_data));
182 
183   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field,
184                                             &(*stats_setup_data)->elem_restr_parent_stats));
185   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height,
186                                                       &(*stats_setup_data)->elem_restr_parent_x));
187   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents((*stats_setup_data)->elem_restr_parent_x, &num_comp_x));
188   PetscCallCeed(ceed, CeedElemRestrictionCreateVector((*stats_setup_data)->elem_restr_parent_x, &(*stats_setup_data)->x_coord, NULL));
189   PetscCall(DMPlexCeedBasisCoordinateCreate(ceed, dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, &(*stats_setup_data)->basis_x));
190   PetscCall(DMPlexCeedBasisCreate(ceed, dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, height, dm_field, &(*stats_setup_data)->basis_stats));
191 
192   PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_stats, (*stats_setup_data)->basis_stats, (*stats_setup_data)->elem_restr_parent_stats,
193                                             &(*stats_setup_data)->elem_restr_parent_colloc));
194   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &elem_restr_q));
195   PetscCall(DMPlexCeedBasisCreate(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, 0, 0, &basis_q));
196   PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_stats, basis_q, elem_restr_q, &(*stats_setup_data)->elem_restr_child_colloc));
197   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
198   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
199 
200   {  // -- Copy DM coordinates into CeedVector
201     DM cdm;
202     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
203     if (cdm) {
204       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
205     } else {
206       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
207     }
208   }
209   // Do not need to scale the coordinates. The DMClone operation has already copied the scaled coordinates
210   PetscCall(VecCopyPetscToCeed(X_loc, (*stats_setup_data)->x_coord));
211   PetscFunctionReturn(PETSC_SUCCESS);
212 }
213 
214 static PetscErrorCode SpanwiseStatisticsSetupDataDestroy(SpanStatsSetupData data) {
215   Ceed ceed;
216 
217   PetscFunctionBeginUser;
218   PetscCall(CeedElemRestrictionGetCeed(data->elem_restr_parent_x, &ceed));
219   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_x));
220   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_stats));
221   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_colloc));
222   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_child_colloc));
223   PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_x));
224   PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_stats));
225   PetscCallCeed(ceed, CeedVectorDestroy(&data->x_coord));
226   PetscCall(PetscFree(data));
227   PetscCheck(CeedDestroy(&ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "Destroying Ceed object failed");
228   PetscFunctionReturn(PETSC_SUCCESS);
229 }
230 
231 // Create PetscSF for child-to-parent communication
232 static PetscErrorCode SpanwiseStatisticsCreateSF(Honee honee, SpanStatsCtx spanstats, SpanStatsSetupData stats_setup_data, DM parentdm, DM childdm,
233                                                  PetscSF *statssf) {
234   Ceed     ceed = honee->ceed;
235   PetscInt child_num_qpnts, parent_num_qpnts;
236   CeedInt  num_comp_x;
237   Vec      Child_qx_coords, Parent_qx_coords;
238 
239   PetscFunctionBeginUser;
240   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)childdm), statssf));
241 
242   // Assume that child and parent have the same number of components
243   PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_x, &num_comp_x));
244   const PetscInt num_comp_sf = num_comp_x - 1;  // Number of coord components used in the creation of the SF
245 
246   // Get quad_coords for child and parent DM
247   PetscCall(GetQuadratureCoords(ceed, childdm, honee->elem_restr_x, honee->basis_x, honee->x_coord, &Child_qx_coords, &child_num_qpnts));
248   PetscCall(GetQuadratureCoords(ceed, parentdm, stats_setup_data->elem_restr_parent_x, stats_setup_data->basis_x, stats_setup_data->x_coord,
249                                 &Parent_qx_coords, &parent_num_qpnts));
250 
251   {  // Remove z component of coordinates for matching
252     const PetscReal *child_quad_coords, *parent_quad_coords;
253     PetscReal       *child_coords, *parent_coords;
254 
255     PetscCall(VecGetArrayRead(Child_qx_coords, &child_quad_coords));
256     PetscCall(VecGetArrayRead(Parent_qx_coords, &parent_quad_coords));
257 
258     PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords));
259     for (int i = 0; i < child_num_qpnts; i++) {
260       child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x];
261       child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x];
262     }
263     for (int i = 0; i < parent_num_qpnts; i++) {
264       parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x];
265       parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x];
266     }
267     PetscCall(VecRestoreArrayRead(Child_qx_coords, &child_quad_coords));
268     PetscCall(VecRestoreArrayRead(Parent_qx_coords, &parent_quad_coords));
269 
270     PetscCall(PetscSFSetGraphFromCoordinates(*statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords));
271     PetscCall(PetscFree2(child_coords, parent_coords));
272   }
273 
274   {
275     char spanstats_sf_key[1024];
276     PetscCall(PetscSNPrintf(spanstats_sf_key, sizeof spanstats_sf_key, "-%ssf_view", spanstats->prefix));
277     PetscCall(PetscObjectSetName((PetscObject)*statssf, spanstats->prefix));
278     PetscCall(PetscSFViewFromOptions(*statssf, NULL, spanstats_sf_key));
279   }
280 
281   PetscCall(VecDestroy(&Child_qx_coords));
282   PetscCall(VecDestroy(&Parent_qx_coords));
283   PetscFunctionReturn(PETSC_SUCCESS);
284 }
285 
286 // @brief Setup RHS and LHS for L^2 projection of statistics
287 static PetscErrorCode SpanwiseStatisticsSetupL2Projection(Honee honee, SpanStatsCtx spanstats, SpanStatsSetupData stats_setup_data) {
288   Ceed                ceed = honee->ceed;
289   CeedOperator        op_mass, op_proj_rhs;
290   CeedQFunction       qf_mass, qf_stats_proj;
291   CeedInt             q_data_size, num_comp_stats = spanstats->num_comp_stats;
292   CeedElemRestriction elem_restr_qd;
293   CeedVector          q_data;
294 
295   PetscFunctionBeginUser;
296   // -- Create Operator for RHS of L^2 projection of statistics
297   // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function.
298   // Therefore, an Identity QF is sufficient
299   PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &qf_stats_proj));
300 
301   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_proj, NULL, NULL, &op_proj_rhs));
302   PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "input", stats_setup_data->elem_restr_parent_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
303   PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "output", stats_setup_data->elem_restr_parent_stats, stats_setup_data->basis_stats,
304                                            CEED_VECTOR_ACTIVE));
305 
306   PetscCall(OperatorApplyContextCreate(NULL, spanstats->dm, ceed, op_proj_rhs, NULL, NULL, NULL, NULL, &spanstats->op_proj_rhs_ctx));
307   PetscCall(CeedOperatorCreateLocalVecs(op_proj_rhs, DMReturnVecType(spanstats->dm), PETSC_COMM_SELF, &spanstats->Parent_Stats_loc, NULL));
308   PetscCall(QDataGet(ceed, spanstats->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, stats_setup_data->elem_restr_parent_x, stats_setup_data->basis_x,
309                      stats_setup_data->x_coord, &elem_restr_qd, &q_data, &q_data_size));
310 
311   // Create Mass CeedOperator
312   PetscCall(HoneeMassQFunctionCreate(ceed, num_comp_stats, q_data_size, &qf_mass));
313   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
314   PetscCallCeed(ceed,
315                 CeedOperatorSetField(op_mass, "u", stats_setup_data->elem_restr_parent_stats, stats_setup_data->basis_stats, CEED_VECTOR_ACTIVE));
316   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
317   PetscCallCeed(ceed,
318                 CeedOperatorSetField(op_mass, "v", stats_setup_data->elem_restr_parent_stats, stats_setup_data->basis_stats, CEED_VECTOR_ACTIVE));
319 
320   {  // Setup KSP for L^2 projection
321     Mat mat_mass;
322     KSP ksp;
323 
324     PetscCall(MatCreateCeed(spanstats->dm, spanstats->dm, op_mass, NULL, &mat_mass));
325 
326     PetscCall(KSPCreate(PetscObjectComm((PetscObject)spanstats->dm), &ksp));
327     PetscCall(KSPSetOptionsPrefix(ksp, spanstats->prefix));
328     {
329       PC pc;
330       PetscCall(KSPGetPC(ksp, &pc));
331       PetscCall(PCSetType(pc, PCJACOBI));
332       PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
333       PetscCall(KSPSetType(ksp, KSPCG));
334       PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
335       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
336     }
337     PetscCall(KSPSetFromOptions_WithMatCeed(ksp, mat_mass));
338     spanstats->ksp = ksp;
339     PetscCall(MatDestroy(&mat_mass));
340   }
341 
342   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
343   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
344   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
345   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_proj));
346   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
347   PetscCallCeed(ceed, CeedOperatorDestroy(&op_proj_rhs));
348   PetscFunctionReturn(PETSC_SUCCESS);
349 }
350 
351 PetscErrorCode SpanwiseStatisticsSetupInitialize(Honee honee, PetscInt degree, PetscInt num_comps, const char *prefix,
352                                                  SpanStatsSetupData *stats_setup_data, SpanStatsCtx *spanstats) {
353   PetscLogStage stage_stats_setup;
354   SpanStatsCtx  spanstats_;
355 
356   PetscFunctionBeginUser;
357   PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup));
358   if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup));
359   PetscCall(PetscLogStagePush(stage_stats_setup));
360 
361   PetscCall(PetscNew(&spanstats_));
362   PetscCall(PetscStrallocpy(prefix, &spanstats_->prefix));
363 
364   spanstats_->collect_interval = 1;
365   PetscCall(PetscOptionsGetInt(NULL, prefix, "-collect_interval", &spanstats_->collect_interval, NULL));
366 
367   // Create parent DM
368   PetscCall(SpanwiseStatisticssCreateDM(honee, spanstats_, degree, num_comps));
369 
370   // Create necessary CeedObjects for setting up statistics
371   PetscCall(SpanwiseStatisticsSetupDataCreate(honee, spanstats_, stats_setup_data));
372   //
373   // Create SF for communicating child data back their respective parents
374   PetscCall(SpanwiseStatisticsCreateSF(honee, spanstats_, *stats_setup_data, honee->dm, spanstats_->dm, &spanstats_->sf));
375 
376   *spanstats = spanstats_;
377   PetscFunctionReturn(PETSC_SUCCESS);
378 }
379 
380 PetscErrorCode SpanwiseStatisticsSetupFinalize(TS ts, Honee honee, SpanStatsCtx spanstats, PetscViewerAndFormat *ctx,
381                                                SpanStatsSetupData *stats_setup_data) {
382   PetscFunctionBeginUser;
383   // Setup KSP and Mat for L^2 projection of statistics
384   PetscCall(SpanwiseStatisticsSetupL2Projection(honee, spanstats, *stats_setup_data));
385 
386   PetscCall(PetscViewerSetOptionsPrefix(ctx->viewer, spanstats->prefix));
387   PetscCall(PetscViewerSetFromOptions(ctx->viewer));
388 
389   PetscCall(TSGetTime(ts, &spanstats->initial_solution_time));
390   PetscCall(TSGetStepNumber(ts, &spanstats->initial_solution_step));
391   CeedScalar initial_solution_time = spanstats->initial_solution_time;  // done for type conversion
392   PetscCallCeed(honee->ceed,
393                 CeedOperatorSetContextDouble(spanstats->op_stats_collect_ctx->op, spanstats->previous_time_label, &initial_solution_time));
394 
395   PetscCall(SpanwiseStatisticsSetupDataDestroy(*stats_setup_data));
396   PetscCall(PetscLogStagePop());
397   PetscFunctionReturn(PETSC_SUCCESS);
398 }
399 
400 // Collect statistics based on the solution Q
401 PetscErrorCode SpanwiseStatisticsCollect(Honee honee, SpanStatsCtx spanstats, PetscScalar solution_time, Vec Q) {
402   PetscFunctionBeginUser;
403   PetscLogStage stage_stats_collect;
404   PetscCall(PetscLogStageGetId("Stats Collect", &stage_stats_collect));
405   if (stage_stats_collect == -1) PetscCall(PetscLogStageRegister("Stats Collect", &stage_stats_collect));
406   PetscCall(PetscLogStagePush(stage_stats_collect));
407 
408   PetscCall(UpdateBoundaryValues(honee, honee->Q_loc, solution_time));
409   PetscCallCeed(honee->ceed, CeedOperatorSetContextDouble(spanstats->op_stats_collect_ctx->op, spanstats->solution_time_label, &solution_time));
410   PetscCall(DMGlobalToLocal(honee->dm, Q, INSERT_VALUES, honee->Q_loc));
411   PetscCall(ApplyAddCeedOperatorLocalToLocal(honee->Q_loc, spanstats->Child_Stats_loc, spanstats->op_stats_collect_ctx));
412 
413   PetscCallCeed(honee->ceed, CeedOperatorSetContextDouble(spanstats->op_stats_collect_ctx->op, spanstats->previous_time_label, &solution_time));
414 
415   PetscCall(PetscLogStagePop());
416   PetscFunctionReturn(PETSC_SUCCESS);
417 }
418 
419 // Process the child statistics into parent statistics and project them onto stats
420 PetscErrorCode SpanwiseStatisticsProcess(Honee honee, SpanStatsCtx spanstats, Vec stats) {
421   const PetscScalar *child_stats;
422   PetscScalar       *parent_stats;
423   MPI_Datatype       unit;
424   Vec                RHS;
425 
426   PetscFunctionBeginUser;
427   PetscLogStage stage_stats_process;
428   PetscCall(PetscLogStageGetId("Stats Process", &stage_stats_process));
429   if (stage_stats_process == -1) PetscCall(PetscLogStageRegister("Stats Process", &stage_stats_process));
430   PetscCall(PetscLogStagePush(stage_stats_process));
431 
432   PetscCall(VecZeroEntries(spanstats->Parent_Stats_loc));
433 
434   PetscCall(VecGetArrayRead(spanstats->Child_Stats_loc, &child_stats));
435   PetscCall(VecGetArray(spanstats->Parent_Stats_loc, &parent_stats));
436 
437   if (spanstats->num_comp_stats == 1) unit = MPIU_REAL;
438   else {
439     PetscCallMPI(MPI_Type_contiguous(spanstats->num_comp_stats, MPIU_REAL, &unit));
440     PetscCallMPI(MPI_Type_commit(&unit));
441   }
442 
443   PetscCall(PetscSFReduceBegin(spanstats->sf, unit, child_stats, parent_stats, MPI_SUM));
444   PetscCall(PetscSFReduceEnd(spanstats->sf, unit, child_stats, parent_stats, MPI_SUM));
445 
446   PetscCall(VecRestoreArrayRead(spanstats->Child_Stats_loc, &child_stats));
447   PetscCall(VecRestoreArray(spanstats->Parent_Stats_loc, &parent_stats));
448   PetscCallMPI(MPI_Type_free(&unit));
449 
450   PetscReal solution_time;
451   PetscCall(DMGetOutputSequenceNumber(spanstats->dm, NULL, &solution_time));
452   PetscReal summing_duration = solution_time - honee->app_ctx->cont_time;
453   PetscCall(VecScale(spanstats->Parent_Stats_loc, 1 / (summing_duration * spanstats->span_width)));
454 
455   // L^2 projection with the parent_data
456   PetscCall(DMGetGlobalVector(spanstats->dm, &RHS));
457   PetscCall(ApplyCeedOperatorLocalToGlobal(spanstats->Parent_Stats_loc, RHS, spanstats->op_proj_rhs_ctx));
458 
459   PetscCall(KSPSolve(spanstats->ksp, RHS, stats));
460 
461   PetscCall(DMRestoreGlobalVector(spanstats->dm, &RHS));
462   PetscCall(PetscLogStagePop());
463   PetscFunctionReturn(PETSC_SUCCESS);
464 }
465