xref: /libCEED/examples/fluids/src/turb_spanstats.c (revision 5ebd836c59d60a2e5e1cb67f6731404c7da26f85)
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_qd, elem_restr_parent_colloc, elem_restr_child_colloc;
25   CeedBasis           basis_x, basis_stats;
26   CeedVector          x_coord, q_data;
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[1];
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), PetscObjectComm((PetscObject)dm), 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   PetscInt dim;
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(DMGetDimension(dm, &dim));
196   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &(*stats_data)->elem_restr_parent_stats));
197   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &(*stats_data)->elem_restr_parent_x));
198   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, problem->q_data_size_sur,
199                                                  &(*stats_data)->elem_restr_parent_qd));
200   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents((*stats_data)->elem_restr_parent_x, &num_comp_x));
201   PetscCallCeed(ceed, CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_x, &(*stats_data)->x_coord, NULL));
202   PetscCallCeed(ceed, CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_qd, &(*stats_data)->q_data, NULL));
203 
204   {
205     DM dm_coord;
206     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
207     PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &(*stats_data)->basis_x));
208     PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &(*stats_data)->basis_stats));
209   }
210 
211   PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_stats, (*stats_data)->basis_stats, (*stats_data)->elem_restr_parent_stats,
212                                             &(*stats_data)->elem_restr_parent_colloc));
213   PetscCall(
214       CreateElemRestrColloc_CompMajor(ceed, num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q, &(*stats_data)->elem_restr_child_colloc));
215 
216   {  // -- Copy DM coordinates into CeedVector
217     DM cdm;
218     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
219     if (cdm) {
220       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
221     } else {
222       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
223     }
224   }
225   PetscCall(VecScale(X_loc, problem->dm_scale));
226   PetscCall(VecCopyPetscToCeed(X_loc, (*stats_data)->x_coord));
227   PetscFunctionReturn(PETSC_SUCCESS);
228 }
229 
230 PetscErrorCode SpanStatsSetupDataDestroy(SpanStatsSetupData data) {
231   Ceed ceed;
232 
233   PetscFunctionBeginUser;
234   PetscCall(CeedElemRestrictionGetCeed(data->elem_restr_parent_x, &ceed));
235   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_x));
236   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_stats));
237   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_qd));
238   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_colloc));
239   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_child_colloc));
240 
241   PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_x));
242   PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_stats));
243 
244   PetscCallCeed(ceed, CeedVectorDestroy(&data->x_coord));
245   PetscCallCeed(ceed, CeedVectorDestroy(&data->q_data));
246 
247   PetscCall(PetscFree(data));
248   PetscFunctionReturn(PETSC_SUCCESS);
249 }
250 
251 // Create PetscSF for child-to-parent communication
252 PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, SpanStatsSetupData stats_data, DM parentdm, DM childdm, PetscSF *statssf) {
253   PetscInt child_num_qpnts, parent_num_qpnts;
254   CeedInt  num_comp_x;
255   Vec      Child_qx_coords, Parent_qx_coords;
256 
257   PetscFunctionBeginUser;
258   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)childdm), statssf));
259 
260   // Assume that child and parent have the same number of components
261   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x));
262   const PetscInt num_comp_sf = num_comp_x - 1;  // Number of coord components used in the creation of the SF
263 
264   // Get quad_coords for child and parent DM
265   PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &Child_qx_coords, &child_num_qpnts));
266   PetscCall(GetQuadratureCoords(ceed, parentdm, stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord, &Parent_qx_coords,
267                                 &parent_num_qpnts));
268 
269   {  // Remove z component of coordinates for matching
270     const PetscReal *child_quad_coords, *parent_quad_coords;
271     PetscReal       *child_coords, *parent_coords;
272 
273     PetscCall(VecGetArrayRead(Child_qx_coords, &child_quad_coords));
274     PetscCall(VecGetArrayRead(Parent_qx_coords, &parent_quad_coords));
275 
276     PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords));
277     for (int i = 0; i < child_num_qpnts; i++) {
278       child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x];
279       child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x];
280     }
281     for (int i = 0; i < parent_num_qpnts; i++) {
282       parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x];
283       parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x];
284     }
285     PetscCall(VecRestoreArrayRead(Child_qx_coords, &child_quad_coords));
286     PetscCall(VecRestoreArrayRead(Parent_qx_coords, &parent_quad_coords));
287 
288     PetscCall(PetscSFSetGraphFromCoordinates(*statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords));
289     PetscCall(PetscFree2(child_coords, parent_coords));
290   }
291 
292   PetscCall(PetscSFViewFromOptions(*statssf, NULL, "-spanstats_sf_view"));
293 
294   PetscCall(VecDestroy(&Child_qx_coords));
295   PetscCall(VecDestroy(&Parent_qx_coords));
296   PetscFunctionReturn(PETSC_SUCCESS);
297 }
298 
299 // @brief Setup RHS and LHS for L^2 projection of statistics
300 PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) {
301   CeedOperator  op_mass, op_setup_sur, op_proj_rhs;
302   CeedQFunction qf_mass, qf_stats_proj;
303   CeedInt       q_data_size, num_comp_stats = user->spanstats.num_comp_stats;
304   MPI_Comm      comm = PetscObjectComm((PetscObject)user->spanstats.dm);
305 
306   PetscFunctionBeginUser;
307   // -- Create Operator for RHS of L^2 projection of statistics
308   // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function.
309   // Therefore, an Identity QF is sufficient
310   PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &qf_stats_proj));
311 
312   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_proj, NULL, NULL, &op_proj_rhs));
313   PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "input", stats_data->elem_restr_parent_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
314   PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "output", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
315 
316   PetscCall(OperatorApplyContextCreate(NULL, user->spanstats.dm, ceed, op_proj_rhs, NULL, NULL, NULL, NULL, &user->spanstats.op_proj_rhs_ctx));
317   PetscCall(CeedOperatorCreateLocalVecs(op_proj_rhs, DMReturnVecType(user->spanstats.dm), comm, &user->spanstats.Parent_Stats_loc, NULL));
318 
319   // -- Setup LHS of L^2 projection
320   // Get q_data for mass matrix operator
321   PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur));
322   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", stats_data->elem_restr_parent_x, stats_data->basis_x, CEED_VECTOR_ACTIVE));
323   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, stats_data->basis_x, CEED_VECTOR_NONE));
324   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
325   PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, stats_data->x_coord, stats_data->q_data, CEED_REQUEST_IMMEDIATE));
326 
327   // CEED Restriction
328   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size));
329 
330   // Create Mass CeedOperator
331   PetscCall(CreateMassQFunction(ceed, num_comp_stats, q_data_size, &qf_mass));
332   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
333   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
334   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_NONE, stats_data->q_data));
335   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
336 
337   {  // Setup KSP for L^2 projection
338     Mat mat_mass;
339     KSP ksp;
340 
341     PetscCall(MatCeedCreate(user->spanstats.dm, user->spanstats.dm, op_mass, NULL, &mat_mass));
342 
343     PetscCall(KSPCreate(comm, &ksp));
344     PetscCall(KSPSetOptionsPrefix(ksp, "turbulence_spanstats_"));
345     {
346       PC pc;
347       PetscCall(KSPGetPC(ksp, &pc));
348       PetscCall(PCSetType(pc, PCJACOBI));
349       PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
350       PetscCall(KSPSetType(ksp, KSPCG));
351       PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
352       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
353     }
354     PetscCall(KSPSetFromOptions_WithMatCeed(ksp, mat_mass));
355     user->spanstats.ksp = ksp;
356     PetscCall(MatDestroy(&mat_mass));
357   }
358 
359   // Cleanup
360   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
361   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_proj));
362   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
363   PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur));
364   PetscCallCeed(ceed, CeedOperatorDestroy(&op_proj_rhs));
365   PetscFunctionReturn(PETSC_SUCCESS);
366 }
367 
368 // Create CeedOperator for statistics collection
369 PetscErrorCode CreateStatisticCollectionOperator(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data, ProblemData *problem) {
370   CeedInt                     num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q;
371   Turbulence_SpanStatsContext collect_ctx;
372   NewtonianIdealGasContext    newtonian_ig_ctx;
373   CeedQFunctionContext        collect_context;
374   CeedQFunction               qf_stats_collect;
375   CeedOperator                op_stats_collect;
376 
377   PetscFunctionBeginUser;
378   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q));
379 
380   // Create Operator for statistics collection
381   switch (user->phys->state_var) {
382     case STATEVAR_PRIMITIVE:
383       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &qf_stats_collect));
384       break;
385     case STATEVAR_CONSERVATIVE:
386       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &qf_stats_collect));
387       break;
388     default:
389       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No statisics collection available for chosen state variable");
390   }
391 
392   if (user->spanstats.do_mms_test) {
393     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect));
394     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest, ChildStatsCollectionMMSTest_loc, &qf_stats_collect));
395   }
396 
397   {  // Setup Collection Context
398     PetscCall(PetscNew(&collect_ctx));
399     PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx));
400     collect_ctx->gas = *newtonian_ig_ctx;
401 
402     PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &collect_context));
403     PetscCallCeed(ceed, CeedQFunctionContextSetData(collect_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx));
404     PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(collect_context, CEED_MEM_HOST, FreeContextPetsc));
405 
406     PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(collect_context, "solution time",
407                                                            offsetof(struct Turbulence_SpanStatsContext_, solution_time), 1, "Current solution time"));
408     PetscCallCeed(
409         ceed, CeedQFunctionContextRegisterDouble(collect_context, "previous time", offsetof(struct Turbulence_SpanStatsContext_, previous_time), 1,
410                                                  "Previous time statistics collection was done"));
411 
412     PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx));
413   }
414 
415   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_stats_collect, collect_context));
416   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&collect_context));
417   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP));
418   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE));
419   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP));
420   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE));
421 
422   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_collect, NULL, NULL, &op_stats_collect));
423   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
424   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
425   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
426   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "v", stats_data->elem_restr_child_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
427 
428   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "solution time", &user->spanstats.solution_time_label));
429   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "previous time", &user->spanstats.previous_time_label));
430 
431   PetscCall(OperatorApplyContextCreate(user->dm, user->spanstats.dm, user->ceed, op_stats_collect, user->q_ceed, NULL, NULL, NULL,
432                                        &user->spanstats.op_stats_collect_ctx));
433 
434   PetscCall(CeedOperatorCreateLocalVecs(op_stats_collect, DMReturnVecType(user->spanstats.dm), PetscObjectComm((PetscObject)user->spanstats.dm), NULL,
435                                         &user->spanstats.Child_Stats_loc));
436   PetscCall(VecZeroEntries(user->spanstats.Child_Stats_loc));
437 
438   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect));
439   PetscCallCeed(ceed, CeedOperatorDestroy(&op_stats_collect));
440   PetscFunctionReturn(PETSC_SUCCESS);
441 }
442 
443 // Creates operator for calculating error of method of manufactured solution (MMS) test
444 PetscErrorCode SetupMMSErrorChecking(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) {
445   CeedInt       num_comp_stats = user->spanstats.num_comp_stats, num_comp_x, q_data_size;
446   CeedQFunction qf_error;
447   CeedOperator  op_error;
448   CeedVector    x_ceed, y_ceed;
449 
450   PetscFunctionBeginUser;
451   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size));
452   PetscCallCeed(ceed, CeedBasisGetNumComponents(stats_data->basis_x, &num_comp_x));
453 
454   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest_Error, ChildStatsCollectionMMSTest_Error_loc, &qf_error));
455   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP));
456   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE));
457   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP));
458   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP));
459 
460   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error));
461   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "q", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
462   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_NONE, stats_data->q_data));
463   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "x", stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord));
464   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
465 
466   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &x_ceed, NULL));
467   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &y_ceed, NULL));
468   PetscCall(OperatorApplyContextCreate(user->spanstats.dm, user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, NULL,
469                                        &user->spanstats.mms_error_ctx));
470 
471   PetscCallCeed(ceed, CeedOperatorDestroy(&op_error));
472   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_error));
473   PetscCallCeed(ceed, CeedVectorDestroy(&x_ceed));
474   PetscCallCeed(ceed, CeedVectorDestroy(&y_ceed));
475   PetscFunctionReturn(PETSC_SUCCESS);
476 }
477 
478 // Setup for statistics collection
479 PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
480   SpanStatsSetupData stats_data;
481   PetscLogStage      stage_stats_setup;
482 
483   PetscFunctionBeginUser;
484   PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup));
485   if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup));
486   PetscCall(PetscLogStagePush(stage_stats_setup));
487 
488   // Create parent DM
489   PetscCall(CreateStatsDM(user, problem, user->app_ctx->degree));
490 
491   // Create necessary CeedObjects for setting up statistics
492   PetscCall(SpanStatsSetupDataCreate(ceed, user, ceed_data, problem, &stats_data));
493 
494   // Create SF for communicating child data back their respective parents
495   PetscCall(CreateStatsSF(ceed, ceed_data, stats_data, user->dm, user->spanstats.dm, &user->spanstats.sf));
496 
497   // Create CeedOperators for statistics collection
498   PetscCall(CreateStatisticCollectionOperator(ceed, user, ceed_data, stats_data, problem));
499 
500   // Setup KSP and Mat for L^2 projection of statistics
501   PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data, stats_data));
502 
503   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_turbulence_spanstats_mms", &user->spanstats.do_mms_test, NULL));
504   if (user->spanstats.do_mms_test) {
505     PetscCall(SetupMMSErrorChecking(ceed, user, ceed_data, stats_data));
506   }
507 
508   {  // Setup stats viewer with prefix
509     PetscViewerType viewer_type;
510     PetscCall(PetscViewerGetType(user->app_ctx->turb_spanstats_viewer, &viewer_type));
511     PetscCall(PetscOptionsSetValue(NULL, "-ts_monitor_turbulence_spanstats_viewer_type", viewer_type));
512 
513     PetscCall(PetscViewerSetOptionsPrefix(user->app_ctx->turb_spanstats_viewer, "ts_monitor_turbulence_spanstats_"));
514     PetscCall(PetscViewerSetFromOptions(user->app_ctx->turb_spanstats_viewer));
515   }
516 
517   PetscCall(SpanStatsSetupDataDestroy(stats_data));
518   PetscCall(PetscLogStagePop());
519   PetscFunctionReturn(PETSC_SUCCESS);
520 }
521 
522 // Collect statistics based on the solution Q
523 PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) {
524   SpanStatsData user_stats = user->spanstats;
525 
526   PetscFunctionBeginUser;
527   PetscLogStage stage_stats_collect;
528   PetscCall(PetscLogStageGetId("Stats Collect", &stage_stats_collect));
529   if (stage_stats_collect == -1) PetscCall(PetscLogStageRegister("Stats Collect", &stage_stats_collect));
530   PetscCall(PetscLogStagePush(stage_stats_collect));
531 
532   PetscCall(UpdateBoundaryValues(user, user->Q_loc, solution_time));
533   PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.solution_time_label, &solution_time));
534   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, user->Q_loc));
535   PetscCall(ApplyAddCeedOperatorLocalToLocal(user->Q_loc, user_stats.Child_Stats_loc, user_stats.op_stats_collect_ctx));
536 
537   PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.previous_time_label, &solution_time));
538 
539   PetscCall(PetscLogStagePop());
540   PetscFunctionReturn(PETSC_SUCCESS);
541 }
542 
543 // Process the child statistics into parent statistics and project them onto stats
544 PetscErrorCode ProcessStatistics(User user, Vec stats) {
545   SpanStatsData      user_stats = user->spanstats;
546   const PetscScalar *child_stats;
547   PetscScalar       *parent_stats;
548   MPI_Datatype       unit;
549   Vec                RHS;
550 
551   PetscFunctionBeginUser;
552   PetscLogStage stage_stats_process;
553   PetscCall(PetscLogStageGetId("Stats Process", &stage_stats_process));
554   if (stage_stats_process == -1) PetscCall(PetscLogStageRegister("Stats Process", &stage_stats_process));
555   PetscCall(PetscLogStagePush(stage_stats_process));
556 
557   PetscCall(VecZeroEntries(user_stats.Parent_Stats_loc));
558 
559   PetscCall(VecGetArrayRead(user_stats.Child_Stats_loc, &child_stats));
560   PetscCall(VecGetArray(user_stats.Parent_Stats_loc, &parent_stats));
561 
562   if (user_stats.num_comp_stats == 1) unit = MPIU_REAL;
563   else {
564     PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit));
565     PetscCallMPI(MPI_Type_commit(&unit));
566   }
567 
568   PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
569   PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
570 
571   PetscCall(VecRestoreArrayRead(user_stats.Child_Stats_loc, &child_stats));
572   PetscCall(VecRestoreArray(user_stats.Parent_Stats_loc, &parent_stats));
573   PetscCallMPI(MPI_Type_free(&unit));
574 
575   PetscReal solution_time;
576   PetscCall(DMGetOutputSequenceNumber(user_stats.dm, NULL, &solution_time));
577   PetscReal summing_duration = solution_time - user->app_ctx->cont_time;
578   PetscCall(VecScale(user_stats.Parent_Stats_loc, 1 / (summing_duration * user_stats.span_width)));
579 
580   // L^2 projection with the parent_data
581   PetscCall(DMGetGlobalVector(user_stats.dm, &RHS));
582   PetscCall(ApplyCeedOperatorLocalToGlobal(user_stats.Parent_Stats_loc, RHS, user_stats.op_proj_rhs_ctx));
583 
584   PetscCall(KSPSolve(user_stats.ksp, RHS, stats));
585 
586   PetscCall(DMRestoreGlobalVector(user_stats.dm, &RHS));
587   PetscCall(PetscLogStagePop());
588   PetscFunctionReturn(PETSC_SUCCESS);
589 }
590 
591 // TSMonitor for the statistics collection and processing
592 PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) {
593   User              user = (User)ctx;
594   Vec               stats;
595   TSConvergedReason reason;
596   PetscInt collect_interval = user->app_ctx->turb_spanstats_collect_interval, viewer_interval = user->app_ctx->turb_spanstats_viewer_interval;
597 
598   PetscFunctionBeginUser;
599   PetscCall(TSGetConvergedReason(ts, &reason));
600   // Do not collect or process on the first step of the run (ie. on the initial condition)
601   if (steps == user->app_ctx->cont_steps) PetscFunctionReturn(PETSC_SUCCESS);
602 
603   PetscBool run_processing_and_viewer = (steps % viewer_interval == 0 && viewer_interval != -1) || reason != TS_CONVERGED_ITERATING;
604 
605   if (steps % collect_interval == 0 || run_processing_and_viewer) {
606     PetscCall(CollectStatistics(user, solution_time, Q));
607 
608     if (run_processing_and_viewer) {
609       PetscCall(DMSetOutputSequenceNumber(user->spanstats.dm, steps, solution_time));
610       PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats));
611       PetscCall(ProcessStatistics(user, stats));
612       if (user->app_ctx->test_type == TESTTYPE_NONE) {
613         PetscCall(PetscViewerPushFormat(user->app_ctx->turb_spanstats_viewer, user->app_ctx->turb_spanstats_viewer_format));
614         PetscCall(VecView(stats, user->app_ctx->turb_spanstats_viewer));
615         PetscCall(PetscViewerPopFormat(user->app_ctx->turb_spanstats_viewer));
616       }
617       if (user->app_ctx->test_type == TESTTYPE_TURB_SPANSTATS && reason != TS_CONVERGED_ITERATING) {
618         PetscCall(RegressionTest(user->app_ctx, stats));
619       }
620       if (user->spanstats.do_mms_test && reason != TS_CONVERGED_ITERATING) {
621         Vec error;
622         PetscCall(VecDuplicate(stats, &error));
623         PetscCall(ApplyCeedOperatorGlobalToGlobal(stats, error, user->spanstats.mms_error_ctx));
624         PetscScalar error_sq = 0;
625         PetscCall(VecSum(error, &error_sq));
626         PetscScalar l2_error = sqrt(error_sq);
627         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error));
628       }
629       PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats));
630     }
631   }
632   PetscFunctionReturn(PETSC_SUCCESS);
633 }
634 
635 PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data) {
636   PetscFunctionBeginUser;
637   PetscCall(VecDestroy(&user->spanstats.Child_Stats_loc));
638   PetscCall(VecDestroy(&user->spanstats.Parent_Stats_loc));
639 
640   PetscCall(OperatorApplyContextDestroy(user->spanstats.op_stats_collect_ctx));
641   PetscCall(OperatorApplyContextDestroy(user->spanstats.op_proj_rhs_ctx));
642   PetscCall(OperatorApplyContextDestroy(user->spanstats.mms_error_ctx));
643 
644   PetscCall(KSPDestroy(&user->spanstats.ksp));
645   PetscCall(PetscSFDestroy(&user->spanstats.sf));
646   PetscCall(DMDestroy(&user->spanstats.dm));
647   PetscFunctionReturn(PETSC_SUCCESS);
648 }
649