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