xref: /libCEED/examples/fluids/src/turb_spanstats.c (revision e64bb3f3ed2986a0c10dec3b47522d734c6e367d)
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 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     PetscSF     isoperiodicface;
48     DMLabel     label;
49     PetscMPIInt size;
50 
51     PetscCall(DMPlexGetIsoperiodicFaceSF(user->dm, &isoperiodicface));
52 
53     if (isoperiodicface) {
54       PetscSF         inv_isoperiodicface;
55       PetscInt        nleaves;
56       const PetscInt *ilocal;
57 
58       PetscCall(PetscSFCreateInverseSF(isoperiodicface, &inv_isoperiodicface));
59       PetscCall(PetscSFGetGraph(inv_isoperiodicface, NULL, &nleaves, &ilocal, NULL));
60       PetscCall(DMCreateLabel(user->dm, "Periodic Face"));
61       PetscCall(DMGetLabel(user->dm, "Periodic Face", &label));
62       for (PetscInt i = 0; i < nleaves; i++) {
63         PetscCall(DMLabelSetValue(label, ilocal[i], 1));
64       }
65     } else {
66       PetscCall(DMGetLabel(user->dm, "Face Sets", &label));
67     }
68 
69     PetscCall(DMPlexLabelComplete(user->dm, label));
70     PetscCall(DMPlexFilter(user->dm, label, 1, &user->spanstats.dm));
71     PetscCall(DMProjectCoordinates(user->spanstats.dm, NULL));  // Ensure that a coordinate FE exists
72 
73     PetscCall(DMPlexDistribute(user->spanstats.dm, 0, NULL, &parent_distributed_dm));
74     PetscCallMPI(MPI_Comm_size(comm, &size));
75     if (parent_distributed_dm) {
76       PetscCall(DMDestroy(&user->spanstats.dm));
77       user->spanstats.dm = parent_distributed_dm;
78     } else if (size > 1) {
79       PetscCall(PetscPrintf(comm, "WARNING: Turbulent spanwise statistics: parent DM could not be distributed accross %d ranks.\n", size));
80     }
81   }
82   {
83     PetscBool is_simplex = PETSC_FALSE;
84     PetscCall(DMPlexIsSimplex(user->spanstats.dm, &is_simplex));
85     PetscCheck(is_simplex != PETSC_TRUE, comm, PETSC_ERR_ARG_WRONGSTATE, "Spanwise Statistics is not implemented for non-tensor product grids");
86   }
87 
88   PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats"));
89   PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "turbulence_spanstats_"));
90   PetscCall(DMSetFromOptions(user->spanstats.dm));
91   PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view"));
92 
93   // Create FE space for parent DM
94   PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, user->app_ctx->degree, 1, user->app_ctx->q_extra, 1, &user->spanstats.num_comp_stats,
95                                user->spanstats.dm));
96 
97   // Create Section for data
98   PetscCall(DMGetLocalSection(user->spanstats.dm, &section));
99   PetscCall(PetscSectionSetFieldName(section, 0, ""));
100   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY, "MeanDensity"));
101   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE, "MeanPressure"));
102   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_SQUARED, "MeanPressureSquared"));
103   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_X, "MeanPressureVelocityX"));
104   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Y, "MeanPressureVelocityY"));
105   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Z, "MeanPressureVelocityZ"));
106   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE, "MeanDensityTemperature"));
107   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_X, "MeanDensityTemperatureFluxX"));
108   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Y, "MeanDensityTemperatureFluxY"));
109   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Z, "MeanDensityTemperatureFluxZ"));
110   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_X, "MeanMomentumX"));
111   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Y, "MeanMomentumY"));
112   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Z, "MeanMomentumZ"));
113   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XX, "MeanMomentumFluxXX"));
114   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YY, "MeanMomentumFluxYY"));
115   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_ZZ, "MeanMomentumFluxZZ"));
116   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YZ, "MeanMomentumFluxYZ"));
117   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XZ, "MeanMomentumFluxXZ"));
118   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XY, "MeanMomentumFluxXY"));
119   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_X, "MeanVelocityX"));
120   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Y, "MeanVelocityY"));
121   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Z, "MeanVelocityZ"));
122 
123   PetscCall(PetscLogStagePop());
124   PetscFunctionReturn(PETSC_SUCCESS);
125 }
126 
127 /** @brief Create CeedElemRestriction for collocated data in component-major order.
128 a. Sets the strides of the restriction to component-major order
129  Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction.
130 */
131 static PetscErrorCode CreateElemRestrColloc_CompMajor(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base,
132                                                       CeedElemRestriction *elem_restr_collocated) {
133   CeedInt num_elem_qpts, loc_num_elem;
134 
135   PetscFunctionBeginUser;
136   PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts));
137   PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem));
138 
139   const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp};
140   PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides,
141                                                        elem_restr_collocated));
142   PetscFunctionReturn(PETSC_SUCCESS);
143 }
144 
145 // Get coordinates of quadrature points
146 PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, Vec *Qx_coords,
147                                    PetscInt *total_nqpnts) {
148   CeedElemRestriction  elem_restr_qx;
149   CeedQFunction        qf_quad_coords;
150   CeedOperator         op_quad_coords;
151   CeedInt              num_comp_x;
152   CeedSize             l_vec_size;
153   OperatorApplyContext op_quad_coords_ctx;
154 
155   PetscFunctionBeginUser;
156   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
157   PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx));
158   PetscCallCeed(ceed, CeedElemRestrictionGetLVectorSize(elem_restr_qx, &l_vec_size));
159   *total_nqpnts = l_vec_size / num_comp_x;
160 
161   // Create QFunction
162   PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords));
163 
164   // Create Operator
165   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords));
166   PetscCallCeed(ceed, CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, x_coords));
167   PetscCallCeed(ceed, CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
168 
169   PetscCall(CeedOperatorCreateLocalVecs(op_quad_coords, DMReturnVecType(dm), PetscObjectComm((PetscObject)dm), NULL, Qx_coords));
170   PetscCall(OperatorApplyContextCreate(NULL, NULL, ceed, op_quad_coords, CEED_VECTOR_NONE, NULL, NULL, NULL, &op_quad_coords_ctx));
171 
172   PetscCall(ApplyCeedOperatorLocalToLocal(NULL, *Qx_coords, op_quad_coords_ctx));
173 
174   PetscCall(OperatorApplyContextDestroy(op_quad_coords_ctx));
175   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qx));
176   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_quad_coords));
177   PetscCallCeed(ceed, CeedOperatorDestroy(&op_quad_coords));
178   PetscFunctionReturn(PETSC_SUCCESS);
179 }
180 
181 PetscErrorCode SpanStatsSetupDataCreate(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem, SpanStatsSetupData *stats_data) {
182   DM       dm = user->spanstats.dm;
183   PetscInt dim;
184   CeedInt  num_comp_x, num_comp_stats = user->spanstats.num_comp_stats;
185   Vec      X_loc;
186   DMLabel  domain_label = NULL;
187   PetscInt label_value = 0, height = 0, dm_field = 0;
188 
189   PetscFunctionBeginUser;
190   PetscCall(PetscNew(stats_data));
191 
192   PetscCall(DMGetDimension(dm, &dim));
193   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &(*stats_data)->elem_restr_parent_stats));
194   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &(*stats_data)->elem_restr_parent_x));
195   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, problem->q_data_size_sur,
196                                                  &(*stats_data)->elem_restr_parent_qd));
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   PetscCallCeed(ceed, CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_qd, &(*stats_data)->q_data, NULL));
200 
201   {
202     DM dm_coord;
203     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
204     PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &(*stats_data)->basis_x));
205     PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &(*stats_data)->basis_stats));
206   }
207 
208   PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_stats, (*stats_data)->basis_stats, (*stats_data)->elem_restr_parent_stats,
209                                             &(*stats_data)->elem_restr_parent_colloc));
210   PetscCall(
211       CreateElemRestrColloc_CompMajor(ceed, num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q, &(*stats_data)->elem_restr_child_colloc));
212 
213   {  // -- Copy DM coordinates into CeedVector
214     DM cdm;
215     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
216     if (cdm) {
217       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
218     } else {
219       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
220     }
221   }
222   PetscCall(VecScale(X_loc, problem->dm_scale));
223   PetscCall(VecCopyP2C(X_loc, (*stats_data)->x_coord));
224   PetscFunctionReturn(PETSC_SUCCESS);
225 }
226 
227 PetscErrorCode SpanStatsSetupDataDestroy(SpanStatsSetupData data) {
228   Ceed ceed;
229 
230   PetscFunctionBeginUser;
231   PetscCall(CeedElemRestrictionGetCeed(data->elem_restr_parent_x, &ceed));
232   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_x));
233   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_stats));
234   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_qd));
235   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_colloc));
236   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_child_colloc));
237 
238   PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_x));
239   PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_stats));
240 
241   PetscCallCeed(ceed, CeedVectorDestroy(&data->x_coord));
242   PetscCallCeed(ceed, CeedVectorDestroy(&data->q_data));
243 
244   PetscCall(PetscFree(data));
245   PetscFunctionReturn(PETSC_SUCCESS);
246 }
247 
248 // Create PetscSF for child-to-parent communication
249 PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, SpanStatsSetupData stats_data, DM parentdm, DM childdm, PetscSF *statssf) {
250   PetscInt child_num_qpnts, parent_num_qpnts;
251   CeedInt  num_comp_x;
252   Vec      Child_qx_coords, Parent_qx_coords;
253 
254   PetscFunctionBeginUser;
255   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)childdm), statssf));
256 
257   // Assume that child and parent have the same number of components
258   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x));
259   const PetscInt num_comp_sf = num_comp_x - 1;  // Number of coord components used in the creation of the SF
260 
261   // Get quad_coords for child and parent DM
262   PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &Child_qx_coords, &child_num_qpnts));
263   PetscCall(GetQuadratureCoords(ceed, parentdm, stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord, &Parent_qx_coords,
264                                 &parent_num_qpnts));
265 
266   {  // Remove z component of coordinates for matching
267     const PetscReal *child_quad_coords, *parent_quad_coords;
268     PetscReal       *child_coords, *parent_coords;
269 
270     PetscCall(VecGetArrayRead(Child_qx_coords, &child_quad_coords));
271     PetscCall(VecGetArrayRead(Parent_qx_coords, &parent_quad_coords));
272 
273     PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords));
274     for (int i = 0; i < child_num_qpnts; i++) {
275       child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x];
276       child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x];
277     }
278     for (int i = 0; i < parent_num_qpnts; i++) {
279       parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x];
280       parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x];
281     }
282     PetscCall(VecRestoreArrayRead(Child_qx_coords, &child_quad_coords));
283     PetscCall(VecRestoreArrayRead(Parent_qx_coords, &parent_quad_coords));
284 
285     PetscCall(PetscSFSetGraphFromCoordinates(*statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords));
286     PetscCall(PetscFree2(child_coords, parent_coords));
287   }
288 
289   PetscCall(PetscSFViewFromOptions(*statssf, NULL, "-spanstats_sf_view"));
290 
291   PetscCall(VecDestroy(&Child_qx_coords));
292   PetscCall(VecDestroy(&Parent_qx_coords));
293   PetscFunctionReturn(PETSC_SUCCESS);
294 }
295 
296 // @brief Setup RHS and LHS for L^2 projection of statistics
297 PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) {
298   CeedOperator  op_mass, op_setup_sur, op_proj_rhs;
299   CeedQFunction qf_mass, qf_stats_proj;
300   CeedInt       q_data_size, num_comp_stats = user->spanstats.num_comp_stats;
301   MPI_Comm      comm = PetscObjectComm((PetscObject)user->spanstats.dm);
302 
303   PetscFunctionBeginUser;
304   // -- Create Operator for RHS of L^2 projection of statistics
305   // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function.
306   // Therefore, an Identity QF is sufficient
307   PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &qf_stats_proj));
308 
309   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_proj, NULL, NULL, &op_proj_rhs));
310   PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "input", stats_data->elem_restr_parent_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
311   PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "output", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
312 
313   PetscCall(OperatorApplyContextCreate(NULL, user->spanstats.dm, ceed, op_proj_rhs, NULL, NULL, NULL, NULL, &user->spanstats.op_proj_rhs_ctx));
314   PetscCall(CeedOperatorCreateLocalVecs(op_proj_rhs, DMReturnVecType(user->spanstats.dm), comm, &user->spanstats.Parent_Stats_loc, NULL));
315 
316   // -- Setup LHS of L^2 projection
317   // Get q_data for mass matrix operator
318   PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur));
319   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", stats_data->elem_restr_parent_x, stats_data->basis_x, CEED_VECTOR_ACTIVE));
320   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, stats_data->basis_x, CEED_VECTOR_NONE));
321   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
322   PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, stats_data->x_coord, stats_data->q_data, CEED_REQUEST_IMMEDIATE));
323 
324   // CEED Restriction
325   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size));
326 
327   // Create Mass CeedOperator
328   PetscCall(CreateMassQFunction(ceed, num_comp_stats, q_data_size, &qf_mass));
329   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
330   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
331   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_NONE, stats_data->q_data));
332   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
333 
334   {  // Setup KSP for L^2 projection
335     OperatorApplyContext M_ctx;
336     Mat                  mat_mass;
337     KSP                  ksp;
338 
339     PetscCall(OperatorApplyContextCreate(user->spanstats.dm, user->spanstats.dm, user->ceed, op_mass, NULL, NULL, NULL, NULL, &M_ctx));
340     PetscCall(CreateMatShell_Ceed(M_ctx, &mat_mass));
341 
342     PetscCall(KSPCreate(comm, &ksp));
343     PetscCall(KSPSetOptionsPrefix(ksp, "turbulence_spanstats_"));
344     {
345       PC pc;
346       PetscCall(KSPGetPC(ksp, &pc));
347       PetscCall(PCSetType(pc, PCJACOBI));
348       PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
349       PetscCall(KSPSetType(ksp, KSPCG));
350       PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
351       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
352     }
353     PetscCall(KSPSetOperators(ksp, mat_mass, mat_mass));
354     PetscCall(KSPSetFromOptions(ksp));
355     user->spanstats.ksp = ksp;
356   }
357 
358   // Cleanup
359   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
360   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_proj));
361   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
362   PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur));
363   PetscCallCeed(ceed, CeedOperatorDestroy(&op_proj_rhs));
364   PetscFunctionReturn(PETSC_SUCCESS);
365 }
366 
367 // Create CeedOperator for statistics collection
368 PetscErrorCode CreateStatisticCollectionOperator(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data, ProblemData *problem) {
369   CeedInt                     num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q;
370   Turbulence_SpanStatsContext collect_ctx;
371   NewtonianIdealGasContext    newtonian_ig_ctx;
372   CeedQFunctionContext        collect_context;
373   CeedQFunction               qf_stats_collect;
374   CeedOperator                op_stats_collect;
375 
376   PetscFunctionBeginUser;
377   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q));
378 
379   // Create Operator for statistics collection
380   switch (user->phys->state_var) {
381     case STATEVAR_PRIMITIVE:
382       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &qf_stats_collect));
383       break;
384     case STATEVAR_CONSERVATIVE:
385       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &qf_stats_collect));
386       break;
387     default:
388       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No statisics collection available for chosen state variable");
389   }
390 
391   if (user->spanstats.do_mms_test) {
392     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect));
393     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest, ChildStatsCollectionMMSTest_loc, &qf_stats_collect));
394   }
395 
396   {  // Setup Collection Context
397     PetscCall(PetscNew(&collect_ctx));
398     PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx));
399     collect_ctx->gas = *newtonian_ig_ctx;
400 
401     PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &collect_context));
402     PetscCallCeed(ceed, CeedQFunctionContextSetData(collect_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx));
403     PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(collect_context, CEED_MEM_HOST, FreeContextPetsc));
404 
405     PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(collect_context, "solution time",
406                                                            offsetof(struct Turbulence_SpanStatsContext_, solution_time), 1, "Current solution time"));
407     PetscCallCeed(
408         ceed, CeedQFunctionContextRegisterDouble(collect_context, "previous time", offsetof(struct Turbulence_SpanStatsContext_, previous_time), 1,
409                                                  "Previous time statistics collection was done"));
410 
411     PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx));
412   }
413 
414   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_stats_collect, collect_context));
415   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&collect_context));
416   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP));
417   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE));
418   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP));
419   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE));
420 
421   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_collect, NULL, NULL, &op_stats_collect));
422   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
423   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
424   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
425   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "v", stats_data->elem_restr_child_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
426 
427   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "solution time", &user->spanstats.solution_time_label));
428   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "previous time", &user->spanstats.previous_time_label));
429 
430   PetscCall(OperatorApplyContextCreate(user->dm, user->spanstats.dm, user->ceed, op_stats_collect, user->q_ceed, NULL, NULL, NULL,
431                                        &user->spanstats.op_stats_collect_ctx));
432 
433   PetscCall(CeedOperatorCreateLocalVecs(op_stats_collect, DMReturnVecType(user->spanstats.dm), PetscObjectComm((PetscObject)user->spanstats.dm), NULL,
434                                         &user->spanstats.Child_Stats_loc));
435   PetscCall(VecZeroEntries(user->spanstats.Child_Stats_loc));
436 
437   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect));
438   PetscCallCeed(ceed, CeedOperatorDestroy(&op_stats_collect));
439   PetscFunctionReturn(PETSC_SUCCESS);
440 }
441 
442 // Creates operator for calculating error of method of manufactured solution (MMS) test
443 PetscErrorCode SetupMMSErrorChecking(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) {
444   CeedInt       num_comp_stats = user->spanstats.num_comp_stats, num_comp_x, q_data_size;
445   CeedQFunction qf_error;
446   CeedOperator  op_error;
447   CeedVector    x_ceed, y_ceed;
448 
449   PetscFunctionBeginUser;
450   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size));
451   PetscCallCeed(ceed, CeedBasisGetNumComponents(stats_data->basis_x, &num_comp_x));
452 
453   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest_Error, ChildStatsCollectionMMSTest_Error_loc, &qf_error));
454   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP));
455   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE));
456   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP));
457   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP));
458 
459   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error));
460   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "q", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
461   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_NONE, stats_data->q_data));
462   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "x", stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord));
463   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
464 
465   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &x_ceed, NULL));
466   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &y_ceed, NULL));
467   PetscCall(OperatorApplyContextCreate(user->spanstats.dm, user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, NULL,
468                                        &user->spanstats.mms_error_ctx));
469 
470   PetscCallCeed(ceed, CeedOperatorDestroy(&op_error));
471   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_error));
472   PetscCallCeed(ceed, CeedVectorDestroy(&x_ceed));
473   PetscCallCeed(ceed, CeedVectorDestroy(&y_ceed));
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 && reason == TS_CONVERGED_ITERATING) 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