xref: /libCEED/examples/fluids/src/turb_spanstats.c (revision 3196072fa7d47cffcf474f8b421f80aad90fa9c5)
1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 /// @file
8 /// Functions for setting up and performing spanwise-statistics collection
9 ///
10 /// "Parent" refers to the 2D plane on which statistics are collected *onto*.
11 /// "Child" refers to the 3D domain where statistics are gathered *from*.
12 /// Each quadrature point on the parent plane has several children in the child domain that it performs spanwise averaging with.
13 
14 #include "../qfunctions/turb_spanstats.h"
15 
16 #include <ceed.h>
17 #include <petscdmplex.h>
18 #include <petscsf.h>
19 
20 #include "../include/petsc_ops.h"
21 #include "../navierstokes.h"
22 
23 typedef struct {
24   CeedElemRestriction elem_restr_parent_x, elem_restr_parent_stats, elem_restr_parent_colloc, elem_restr_child_colloc;
25   CeedBasis           basis_x, basis_stats;
26   CeedVector          x_coord;
27 } *SpanStatsSetupData;
28 
29 PetscErrorCode CreateStatsDM(User user, ProblemData problem, PetscInt degree) {
30   user->spanstats.num_comp_stats = TURB_NUM_COMPONENTS;
31   PetscReal     domain_min[3], domain_max[3];
32   PetscSection  section;
33   PetscLogStage stage_stats_setup;
34   MPI_Comm      comm = PetscObjectComm((PetscObject)user->dm);
35 
36   PetscFunctionBeginUser;
37   PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup));
38   if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup));
39   PetscCall(PetscLogStagePush(stage_stats_setup));
40 
41   // Get spanwise length
42   PetscCall(DMGetBoundingBox(user->dm, domain_min, domain_max));
43   user->spanstats.span_width = domain_max[2] - domain_min[2];
44 
45   {  // Get DM from surface
46     DM             parent_distributed_dm;
47     const PetscSF *isoperiodicface;
48     PetscInt       num_isoperiodicface;
49     DMLabel        label;
50     PetscMPIInt    size;
51 
52     PetscCall(DMPlexGetIsoperiodicFaceSF(user->dm, &num_isoperiodicface, &isoperiodicface));
53 
54     if (isoperiodicface) {
55       PetscSF         inv_isoperiodicface;
56       PetscInt        nleaves, isoperiodicface_index = -1;
57       const PetscInt *ilocal;
58 
59       PetscCall(PetscOptionsGetInt(NULL, NULL, "-ts_monitor_turbulence_spanstats_isoperiodic_facesf", &isoperiodicface_index, NULL));
60       isoperiodicface_index = isoperiodicface_index == -1 ? num_isoperiodicface - 1 : isoperiodicface_index;
61       PetscCall(PetscSFCreateInverseSF(isoperiodicface[isoperiodicface_index], &inv_isoperiodicface));
62       PetscCall(PetscSFGetGraph(inv_isoperiodicface, NULL, &nleaves, &ilocal, NULL));
63       PetscCall(DMCreateLabel(user->dm, "Periodic Face"));
64       PetscCall(DMGetLabel(user->dm, "Periodic Face", &label));
65       for (PetscInt i = 0; i < nleaves; i++) {
66         PetscCall(DMLabelSetValue(label, ilocal[i], 1));
67       }
68     } else {
69       PetscCall(DMGetLabel(user->dm, "Face Sets", &label));
70     }
71 
72     PetscCall(DMPlexLabelComplete(user->dm, label));
73     PetscCall(DMPlexFilter(user->dm, label, 1, PETSC_FALSE, PETSC_FALSE, NULL, &user->spanstats.dm));
74     PetscCall(DMSetCoordinateDisc(user->spanstats.dm, NULL, PETSC_TRUE));  // Ensure that a coordinate FE exists
75 
76     PetscCall(DMPlexDistribute(user->spanstats.dm, 0, NULL, &parent_distributed_dm));
77     PetscCallMPI(MPI_Comm_size(comm, &size));
78     if (parent_distributed_dm) {
79       PetscCall(DMDestroy(&user->spanstats.dm));
80       user->spanstats.dm = parent_distributed_dm;
81     } else if (size > 1) {
82       PetscCall(PetscPrintf(comm, "WARNING: Turbulent spanwise statistics: parent DM could not be distributed accross %d ranks.\n", size));
83     }
84   }
85   {
86     PetscBool is_simplex = PETSC_FALSE;
87     PetscCall(DMPlexIsSimplex(user->spanstats.dm, &is_simplex));
88     PetscCheck(is_simplex != PETSC_TRUE, comm, PETSC_ERR_ARG_WRONGSTATE, "Spanwise Statistics is not implemented for non-tensor product grids");
89   }
90 
91   PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats"));
92   PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "turbulence_spanstats_"));
93   PetscCall(DMSetFromOptions(user->spanstats.dm));
94   PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view"));
95 
96   // Create FE space for parent DM
97   PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, user->app_ctx->degree, 1, user->app_ctx->q_extra, 1, &user->spanstats.num_comp_stats,
98                                user->spanstats.dm));
99 
100   // Create Section for data
101   PetscCall(DMGetLocalSection(user->spanstats.dm, &section));
102   PetscCall(PetscSectionSetFieldName(section, 0, ""));
103   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY, "MeanDensity"));
104   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE, "MeanPressure"));
105   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_SQUARED, "MeanPressureSquared"));
106   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_X, "MeanPressureVelocityX"));
107   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Y, "MeanPressureVelocityY"));
108   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Z, "MeanPressureVelocityZ"));
109   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE, "MeanDensityTemperature"));
110   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_X, "MeanDensityTemperatureFluxX"));
111   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Y, "MeanDensityTemperatureFluxY"));
112   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Z, "MeanDensityTemperatureFluxZ"));
113   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_X, "MeanMomentumX"));
114   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Y, "MeanMomentumY"));
115   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Z, "MeanMomentumZ"));
116   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XX, "MeanMomentumFluxXX"));
117   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YY, "MeanMomentumFluxYY"));
118   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_ZZ, "MeanMomentumFluxZZ"));
119   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YZ, "MeanMomentumFluxYZ"));
120   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XZ, "MeanMomentumFluxXZ"));
121   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XY, "MeanMomentumFluxXY"));
122   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_X, "MeanVelocityX"));
123   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Y, "MeanVelocityY"));
124   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Z, "MeanVelocityZ"));
125 
126   PetscCall(PetscLogStagePop());
127   PetscFunctionReturn(PETSC_SUCCESS);
128 }
129 
130 /** @brief Create CeedElemRestriction for collocated data in component-major order.
131 a. Sets the strides of the restriction to component-major order
132  Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction.
133 */
134 static PetscErrorCode CreateElemRestrColloc_CompMajor(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base,
135                                                       CeedElemRestriction *elem_restr_collocated) {
136   CeedInt num_elem_qpts, loc_num_elem;
137 
138   PetscFunctionBeginUser;
139   PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts));
140   PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem));
141 
142   const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp};
143   PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides,
144                                                        elem_restr_collocated));
145   PetscFunctionReturn(PETSC_SUCCESS);
146 }
147 
148 // Get coordinates of quadrature points
149 PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, Vec *Qx_coords,
150                                    PetscInt *total_nqpnts) {
151   CeedElemRestriction  elem_restr_qx;
152   CeedQFunction        qf_quad_coords;
153   CeedOperator         op_quad_coords;
154   CeedInt              num_comp_x;
155   CeedSize             l_vec_size;
156   OperatorApplyContext op_quad_coords_ctx;
157 
158   PetscFunctionBeginUser;
159   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
160   PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx));
161   PetscCallCeed(ceed, CeedElemRestrictionGetLVectorSize(elem_restr_qx, &l_vec_size));
162   *total_nqpnts = l_vec_size / num_comp_x;
163 
164   // Create QFunction
165   PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords));
166 
167   // Create Operator
168   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords));
169   PetscCallCeed(ceed, CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, x_coords));
170   PetscCallCeed(ceed, CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
171 
172   PetscCall(CeedOperatorCreateLocalVecs(op_quad_coords, DMReturnVecType(dm), PETSC_COMM_SELF, NULL, Qx_coords));
173   PetscCall(OperatorApplyContextCreate(NULL, NULL, ceed, op_quad_coords, CEED_VECTOR_NONE, NULL, NULL, NULL, &op_quad_coords_ctx));
174 
175   PetscCall(ApplyCeedOperatorLocalToLocal(NULL, *Qx_coords, op_quad_coords_ctx));
176 
177   PetscCall(OperatorApplyContextDestroy(op_quad_coords_ctx));
178   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qx));
179   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_quad_coords));
180   PetscCallCeed(ceed, CeedOperatorDestroy(&op_quad_coords));
181   PetscFunctionReturn(PETSC_SUCCESS);
182 }
183 
184 PetscErrorCode SpanStatsSetupDataCreate(Ceed ceed, User user, CeedData ceed_data, ProblemData problem, SpanStatsSetupData *stats_data) {
185   DM       dm = user->spanstats.dm;
186   CeedInt  num_comp_x, num_comp_stats = user->spanstats.num_comp_stats;
187   Vec      X_loc;
188   DMLabel  domain_label = NULL;
189   PetscInt label_value = 0, height = 0, dm_field = 0;
190 
191   PetscFunctionBeginUser;
192   PetscCall(PetscNew(stats_data));
193 
194   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &(*stats_data)->elem_restr_parent_stats));
195   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &(*stats_data)->elem_restr_parent_x));
196   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents((*stats_data)->elem_restr_parent_x, &num_comp_x));
197   PetscCallCeed(ceed, CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_x, &(*stats_data)->x_coord, NULL));
198 
199   {
200     DM dm_coord;
201     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
202     PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &(*stats_data)->basis_x));
203     PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &(*stats_data)->basis_stats));
204   }
205 
206   PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_stats, (*stats_data)->basis_stats, (*stats_data)->elem_restr_parent_stats,
207                                             &(*stats_data)->elem_restr_parent_colloc));
208   PetscCall(
209       CreateElemRestrColloc_CompMajor(ceed, num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q, &(*stats_data)->elem_restr_child_colloc));
210 
211   {  // -- Copy DM coordinates into CeedVector
212     DM cdm;
213     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
214     if (cdm) {
215       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
216     } else {
217       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
218     }
219   }
220   PetscCall(VecScale(X_loc, problem->dm_scale));
221   PetscCall(VecCopyPetscToCeed(X_loc, (*stats_data)->x_coord));
222   PetscFunctionReturn(PETSC_SUCCESS);
223 }
224 
225 PetscErrorCode SpanStatsSetupDataDestroy(SpanStatsSetupData data) {
226   Ceed ceed;
227 
228   PetscFunctionBeginUser;
229   PetscCall(CeedElemRestrictionGetCeed(data->elem_restr_parent_x, &ceed));
230   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_x));
231   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_stats));
232   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_colloc));
233   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_child_colloc));
234 
235   PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_x));
236   PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_stats));
237 
238   PetscCallCeed(ceed, CeedVectorDestroy(&data->x_coord));
239 
240   PetscCall(PetscFree(data));
241   PetscFunctionReturn(PETSC_SUCCESS);
242 }
243 
244 // Create PetscSF for child-to-parent communication
245 PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, SpanStatsSetupData stats_data, DM parentdm, DM childdm, PetscSF *statssf) {
246   PetscInt child_num_qpnts, parent_num_qpnts;
247   CeedInt  num_comp_x;
248   Vec      Child_qx_coords, Parent_qx_coords;
249 
250   PetscFunctionBeginUser;
251   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)childdm), statssf));
252 
253   // Assume that child and parent have the same number of components
254   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x));
255   const PetscInt num_comp_sf = num_comp_x - 1;  // Number of coord components used in the creation of the SF
256 
257   // Get quad_coords for child and parent DM
258   PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &Child_qx_coords, &child_num_qpnts));
259   PetscCall(GetQuadratureCoords(ceed, parentdm, stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord, &Parent_qx_coords,
260                                 &parent_num_qpnts));
261 
262   {  // Remove z component of coordinates for matching
263     const PetscReal *child_quad_coords, *parent_quad_coords;
264     PetscReal       *child_coords, *parent_coords;
265 
266     PetscCall(VecGetArrayRead(Child_qx_coords, &child_quad_coords));
267     PetscCall(VecGetArrayRead(Parent_qx_coords, &parent_quad_coords));
268 
269     PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords));
270     for (int i = 0; i < child_num_qpnts; i++) {
271       child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x];
272       child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x];
273     }
274     for (int i = 0; i < parent_num_qpnts; i++) {
275       parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x];
276       parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x];
277     }
278     PetscCall(VecRestoreArrayRead(Child_qx_coords, &child_quad_coords));
279     PetscCall(VecRestoreArrayRead(Parent_qx_coords, &parent_quad_coords));
280 
281     PetscCall(PetscSFSetGraphFromCoordinates(*statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords));
282     PetscCall(PetscFree2(child_coords, parent_coords));
283   }
284 
285   PetscCall(PetscSFViewFromOptions(*statssf, NULL, "-spanstats_sf_view"));
286 
287   PetscCall(VecDestroy(&Child_qx_coords));
288   PetscCall(VecDestroy(&Parent_qx_coords));
289   PetscFunctionReturn(PETSC_SUCCESS);
290 }
291 
292 // @brief Setup RHS and LHS for L^2 projection of statistics
293 PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) {
294   CeedOperator        op_mass, op_proj_rhs;
295   CeedQFunction       qf_mass, qf_stats_proj;
296   CeedInt             q_data_size, num_comp_stats = user->spanstats.num_comp_stats;
297   CeedElemRestriction elem_restr_qd;
298   CeedVector          q_data;
299   DMLabel             domain_label = NULL;
300   PetscInt            label_value  = 0;
301 
302   PetscFunctionBeginUser;
303   // -- Create Operator for RHS of L^2 projection of statistics
304   // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function.
305   // Therefore, an Identity QF is sufficient
306   PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &qf_stats_proj));
307 
308   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_proj, NULL, NULL, &op_proj_rhs));
309   PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "input", stats_data->elem_restr_parent_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
310   PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "output", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
311 
312   PetscCall(OperatorApplyContextCreate(NULL, user->spanstats.dm, ceed, op_proj_rhs, NULL, NULL, NULL, NULL, &user->spanstats.op_proj_rhs_ctx));
313   PetscCall(CeedOperatorCreateLocalVecs(op_proj_rhs, DMReturnVecType(user->spanstats.dm), PETSC_COMM_SELF, &user->spanstats.Parent_Stats_loc, NULL));
314   PetscCall(QDataGet(ceed, user->spanstats.dm, domain_label, label_value, stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord,
315                      &elem_restr_qd, &q_data, &q_data_size));
316 
317   // Create Mass CeedOperator
318   PetscCall(CreateMassQFunction(ceed, num_comp_stats, q_data_size, &qf_mass));
319   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
320   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
321   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
322   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
323 
324   {  // Setup KSP for L^2 projection
325     Mat mat_mass;
326     KSP ksp;
327 
328     PetscCall(MatCeedCreate(user->spanstats.dm, user->spanstats.dm, op_mass, NULL, &mat_mass));
329 
330     PetscCall(KSPCreate(PetscObjectComm((PetscObject)user->spanstats.dm), &ksp));
331     PetscCall(KSPSetOptionsPrefix(ksp, "turbulence_spanstats_"));
332     {
333       PC pc;
334       PetscCall(KSPGetPC(ksp, &pc));
335       PetscCall(PCSetType(pc, PCJACOBI));
336       PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
337       PetscCall(KSPSetType(ksp, KSPCG));
338       PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
339       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
340     }
341     PetscCall(KSPSetFromOptions_WithMatCeed(ksp, mat_mass));
342     user->spanstats.ksp = ksp;
343     PetscCall(MatDestroy(&mat_mass));
344   }
345 
346   // Cleanup
347   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
348   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
349   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
350   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_proj));
351   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
352   PetscCallCeed(ceed, CeedOperatorDestroy(&op_proj_rhs));
353   PetscFunctionReturn(PETSC_SUCCESS);
354 }
355 
356 // Create CeedOperator for statistics collection
357 PetscErrorCode CreateStatisticCollectionOperator(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data, ProblemData problem) {
358   CeedInt                     num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q;
359   Turbulence_SpanStatsContext collect_ctx;
360   NewtonianIdealGasContext    newtonian_ig_ctx;
361   CeedQFunctionContext        collect_context;
362   CeedQFunction               qf_stats_collect;
363   CeedOperator                op_stats_collect;
364 
365   PetscFunctionBeginUser;
366   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q));
367 
368   // Create Operator for statistics collection
369   switch (user->phys->state_var) {
370     case STATEVAR_PRIMITIVE:
371       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &qf_stats_collect));
372       break;
373     case STATEVAR_CONSERVATIVE:
374       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &qf_stats_collect));
375       break;
376     case STATEVAR_ENTROPY:
377       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Entropy, ChildStatsCollection_Entropy_loc, &qf_stats_collect));
378       break;
379   }
380 
381   if (user->spanstats.do_mms_test) {
382     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect));
383     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest, ChildStatsCollectionMMSTest_loc, &qf_stats_collect));
384   }
385 
386   {  // Setup Collection Context
387     PetscCall(PetscNew(&collect_ctx));
388     PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx));
389     collect_ctx->gas = *newtonian_ig_ctx;
390 
391     PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &collect_context));
392     PetscCallCeed(ceed, CeedQFunctionContextSetData(collect_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx));
393     PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(collect_context, CEED_MEM_HOST, FreeContextPetsc));
394 
395     PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(collect_context, "solution time",
396                                                            offsetof(struct Turbulence_SpanStatsContext_, solution_time), 1, "Current solution time"));
397     PetscCallCeed(
398         ceed, CeedQFunctionContextRegisterDouble(collect_context, "previous time", offsetof(struct Turbulence_SpanStatsContext_, previous_time), 1,
399                                                  "Previous time statistics collection was done"));
400 
401     PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx));
402   }
403 
404   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_stats_collect, collect_context));
405   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&collect_context));
406   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP));
407   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE));
408   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP));
409   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE));
410 
411   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_collect, NULL, NULL, &op_stats_collect));
412   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
413   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
414   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
415   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "v", stats_data->elem_restr_child_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
416 
417   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "solution time", &user->spanstats.solution_time_label));
418   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "previous time", &user->spanstats.previous_time_label));
419 
420   PetscCall(OperatorApplyContextCreate(user->dm, user->spanstats.dm, user->ceed, op_stats_collect, user->q_ceed, NULL, NULL, NULL,
421                                        &user->spanstats.op_stats_collect_ctx));
422 
423   PetscCall(
424       CeedOperatorCreateLocalVecs(op_stats_collect, DMReturnVecType(user->spanstats.dm), PETSC_COMM_SELF, NULL, &user->spanstats.Child_Stats_loc));
425   PetscCall(VecZeroEntries(user->spanstats.Child_Stats_loc));
426 
427   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect));
428   PetscCallCeed(ceed, CeedOperatorDestroy(&op_stats_collect));
429   PetscFunctionReturn(PETSC_SUCCESS);
430 }
431 
432 // Creates operator for calculating error of method of manufactured solution (MMS) test
433 PetscErrorCode SetupMMSErrorChecking(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) {
434   CeedInt             num_comp_stats = user->spanstats.num_comp_stats, num_comp_x, q_data_size;
435   CeedQFunction       qf_error;
436   CeedOperator        op_error;
437   CeedVector          x_ceed, y_ceed;
438   DMLabel             domain_label = NULL;
439   PetscInt            label_value  = 0;
440   CeedVector          q_data;
441   CeedElemRestriction elem_restr_parent_qd;
442 
443   PetscFunctionBeginUser;
444   PetscCall(QDataGet(ceed, user->spanstats.dm, domain_label, label_value, stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord,
445                      &elem_restr_parent_qd, &q_data, &q_data_size));
446   PetscCallCeed(ceed, CeedBasisGetNumComponents(stats_data->basis_x, &num_comp_x));
447 
448   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest_Error, ChildStatsCollectionMMSTest_Error_loc, &qf_error));
449   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP));
450   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE));
451   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP));
452   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP));
453 
454   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error));
455   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "q", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
456   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "qdata", elem_restr_parent_qd, CEED_BASIS_NONE, q_data));
457   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "x", stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord));
458   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
459 
460   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &x_ceed, NULL));
461   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &y_ceed, NULL));
462   PetscCall(OperatorApplyContextCreate(user->spanstats.dm, user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, NULL,
463                                        &user->spanstats.mms_error_ctx));
464 
465   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
466   PetscCallCeed(ceed, CeedVectorDestroy(&x_ceed));
467   PetscCallCeed(ceed, CeedVectorDestroy(&y_ceed));
468   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_parent_qd));
469   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_error));
470   PetscCallCeed(ceed, CeedOperatorDestroy(&op_error));
471   PetscFunctionReturn(PETSC_SUCCESS);
472 }
473 
474 // Setup for statistics collection
475 PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) {
476   SpanStatsSetupData stats_data;
477   PetscLogStage      stage_stats_setup;
478 
479   PetscFunctionBeginUser;
480   PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup));
481   if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup));
482   PetscCall(PetscLogStagePush(stage_stats_setup));
483 
484   // Create parent DM
485   PetscCall(CreateStatsDM(user, problem, user->app_ctx->degree));
486 
487   // Create necessary CeedObjects for setting up statistics
488   PetscCall(SpanStatsSetupDataCreate(ceed, user, ceed_data, problem, &stats_data));
489 
490   // Create SF for communicating child data back their respective parents
491   PetscCall(CreateStatsSF(ceed, ceed_data, stats_data, user->dm, user->spanstats.dm, &user->spanstats.sf));
492 
493   // Create CeedOperators for statistics collection
494   PetscCall(CreateStatisticCollectionOperator(ceed, user, ceed_data, stats_data, problem));
495 
496   // Setup KSP and Mat for L^2 projection of statistics
497   PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data, stats_data));
498 
499   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_turbulence_spanstats_mms", &user->spanstats.do_mms_test, NULL));
500   if (user->spanstats.do_mms_test) {
501     PetscCall(SetupMMSErrorChecking(ceed, user, ceed_data, stats_data));
502   }
503 
504   {  // Setup stats viewer with prefix
505     PetscViewerType viewer_type;
506     PetscCall(PetscViewerGetType(user->app_ctx->turb_spanstats_viewer, &viewer_type));
507     PetscCall(PetscOptionsSetValue(NULL, "-ts_monitor_turbulence_spanstats_viewer_type", viewer_type));
508 
509     PetscCall(PetscViewerSetOptionsPrefix(user->app_ctx->turb_spanstats_viewer, "ts_monitor_turbulence_spanstats_"));
510     PetscCall(PetscViewerSetFromOptions(user->app_ctx->turb_spanstats_viewer));
511   }
512 
513   PetscCall(SpanStatsSetupDataDestroy(stats_data));
514   PetscCall(PetscLogStagePop());
515   PetscFunctionReturn(PETSC_SUCCESS);
516 }
517 
518 // Collect statistics based on the solution Q
519 PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) {
520   SpanStatsData user_stats = user->spanstats;
521 
522   PetscFunctionBeginUser;
523   PetscLogStage stage_stats_collect;
524   PetscCall(PetscLogStageGetId("Stats Collect", &stage_stats_collect));
525   if (stage_stats_collect == -1) PetscCall(PetscLogStageRegister("Stats Collect", &stage_stats_collect));
526   PetscCall(PetscLogStagePush(stage_stats_collect));
527 
528   PetscCall(UpdateBoundaryValues(user, user->Q_loc, solution_time));
529   PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.solution_time_label, &solution_time));
530   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, user->Q_loc));
531   PetscCall(ApplyAddCeedOperatorLocalToLocal(user->Q_loc, user_stats.Child_Stats_loc, user_stats.op_stats_collect_ctx));
532 
533   PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.previous_time_label, &solution_time));
534 
535   PetscCall(PetscLogStagePop());
536   PetscFunctionReturn(PETSC_SUCCESS);
537 }
538 
539 // Process the child statistics into parent statistics and project them onto stats
540 PetscErrorCode ProcessStatistics(User user, Vec stats) {
541   SpanStatsData      user_stats = user->spanstats;
542   const PetscScalar *child_stats;
543   PetscScalar       *parent_stats;
544   MPI_Datatype       unit;
545   Vec                RHS;
546 
547   PetscFunctionBeginUser;
548   PetscLogStage stage_stats_process;
549   PetscCall(PetscLogStageGetId("Stats Process", &stage_stats_process));
550   if (stage_stats_process == -1) PetscCall(PetscLogStageRegister("Stats Process", &stage_stats_process));
551   PetscCall(PetscLogStagePush(stage_stats_process));
552 
553   PetscCall(VecZeroEntries(user_stats.Parent_Stats_loc));
554 
555   PetscCall(VecGetArrayRead(user_stats.Child_Stats_loc, &child_stats));
556   PetscCall(VecGetArray(user_stats.Parent_Stats_loc, &parent_stats));
557 
558   if (user_stats.num_comp_stats == 1) unit = MPIU_REAL;
559   else {
560     PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit));
561     PetscCallMPI(MPI_Type_commit(&unit));
562   }
563 
564   PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
565   PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
566 
567   PetscCall(VecRestoreArrayRead(user_stats.Child_Stats_loc, &child_stats));
568   PetscCall(VecRestoreArray(user_stats.Parent_Stats_loc, &parent_stats));
569   PetscCallMPI(MPI_Type_free(&unit));
570 
571   PetscReal solution_time;
572   PetscCall(DMGetOutputSequenceNumber(user_stats.dm, NULL, &solution_time));
573   PetscReal summing_duration = solution_time - user->app_ctx->cont_time;
574   PetscCall(VecScale(user_stats.Parent_Stats_loc, 1 / (summing_duration * user_stats.span_width)));
575 
576   // L^2 projection with the parent_data
577   PetscCall(DMGetGlobalVector(user_stats.dm, &RHS));
578   PetscCall(ApplyCeedOperatorLocalToGlobal(user_stats.Parent_Stats_loc, RHS, user_stats.op_proj_rhs_ctx));
579 
580   PetscCall(KSPSolve(user_stats.ksp, RHS, stats));
581 
582   PetscCall(DMRestoreGlobalVector(user_stats.dm, &RHS));
583   PetscCall(PetscLogStagePop());
584   PetscFunctionReturn(PETSC_SUCCESS);
585 }
586 
587 // TSMonitor for the statistics collection and processing
588 PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) {
589   User              user = (User)ctx;
590   Vec               stats;
591   TSConvergedReason reason;
592   PetscInt collect_interval = user->app_ctx->turb_spanstats_collect_interval, viewer_interval = user->app_ctx->turb_spanstats_viewer_interval;
593 
594   PetscFunctionBeginUser;
595   PetscCall(TSGetConvergedReason(ts, &reason));
596   // Do not collect or process on the first step of the run (ie. on the initial condition)
597   if (steps == user->app_ctx->cont_steps) PetscFunctionReturn(PETSC_SUCCESS);
598 
599   PetscBool run_processing_and_viewer = (steps % viewer_interval == 0 && viewer_interval != -1) || reason != TS_CONVERGED_ITERATING;
600 
601   if (steps % collect_interval == 0 || run_processing_and_viewer) {
602     PetscCall(CollectStatistics(user, solution_time, Q));
603 
604     if (run_processing_and_viewer) {
605       PetscCall(DMSetOutputSequenceNumber(user->spanstats.dm, steps, solution_time));
606       PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats));
607       PetscCall(ProcessStatistics(user, stats));
608       if (user->app_ctx->test_type == TESTTYPE_NONE) {
609         PetscCall(PetscViewerPushFormat(user->app_ctx->turb_spanstats_viewer, user->app_ctx->turb_spanstats_viewer_format));
610         PetscCall(VecView(stats, user->app_ctx->turb_spanstats_viewer));
611         PetscCall(PetscViewerPopFormat(user->app_ctx->turb_spanstats_viewer));
612       }
613       if (user->app_ctx->test_type == TESTTYPE_TURB_SPANSTATS && reason != TS_CONVERGED_ITERATING) {
614         PetscCall(RegressionTest(user->app_ctx, stats));
615       }
616       if (user->spanstats.do_mms_test && reason != TS_CONVERGED_ITERATING) {
617         Vec error;
618         PetscCall(VecDuplicate(stats, &error));
619         PetscCall(ApplyCeedOperatorGlobalToGlobal(stats, error, user->spanstats.mms_error_ctx));
620         PetscScalar error_sq = 0;
621         PetscCall(VecSum(error, &error_sq));
622         PetscScalar l2_error = sqrt(error_sq);
623         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error));
624       }
625       PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats));
626     }
627   }
628   PetscFunctionReturn(PETSC_SUCCESS);
629 }
630 
631 PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data) {
632   PetscFunctionBeginUser;
633   PetscCall(VecDestroy(&user->spanstats.Child_Stats_loc));
634   PetscCall(VecDestroy(&user->spanstats.Parent_Stats_loc));
635 
636   PetscCall(OperatorApplyContextDestroy(user->spanstats.op_stats_collect_ctx));
637   PetscCall(OperatorApplyContextDestroy(user->spanstats.op_proj_rhs_ctx));
638   PetscCall(OperatorApplyContextDestroy(user->spanstats.mms_error_ctx));
639 
640   PetscCall(KSPDestroy(&user->spanstats.ksp));
641   PetscCall(PetscSFDestroy(&user->spanstats.sf));
642   PetscCall(DMDestroy(&user->spanstats.dm));
643   PetscFunctionReturn(PETSC_SUCCESS);
644 }
645