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