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