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