xref: /libCEED/examples/fluids/src/turb_spanstats.c (revision 823a1283cec5199b0df409b58df1642984320911)
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 <petscsf.h>
13 
14 #include "../include/matops.h"
15 #include "../navierstokes.h"
16 #include "ceed/ceed.h"
17 #include "petscerror.h"
18 #include "petscmat.h"
19 #include "petscsys.h"
20 #include "petscvec.h"
21 
22 PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree, SimpleBC bc) {
23   user->spanstats.num_comp_stats = 22;
24   PetscReal domain_min[3], domain_max[3];
25   PetscFunctionBeginUser;
26 
27   // Get spanwise length
28   PetscCall(DMGetBoundingBox(user->dm, domain_min, domain_max));
29   user->spanstats.span_width = domain_max[2] - domain_min[1];
30 
31   // Get DM from surface
32   {
33     DMLabel label;
34     PetscCall(DMGetLabel(user->dm, "Face Sets", &label));
35     PetscCall(DMPlexLabelComplete(user->dm, label));
36     PetscCall(DMPlexFilter(user->dm, label, 1, &user->spanstats.dm));
37     PetscCall(DMProjectCoordinates(user->spanstats.dm, NULL));  // Ensure that a coordinate FE exists
38   }
39 
40   PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats"));
41   PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "spanstats_"));
42   PetscCall(DMSetFromOptions(user->spanstats.dm));
43   PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view"));  // -spanstats_dm_view
44   {
45     PetscFE fe;
46     DMLabel label;
47 
48     PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim - 1, user->spanstats.num_comp_stats, PETSC_FALSE, degree, PETSC_DECIDE, &fe));
49     PetscCall(PetscObjectSetName((PetscObject)fe, "stats"));
50     PetscCall(DMAddField(user->spanstats.dm, NULL, (PetscObject)fe));
51     PetscCall(DMCreateDS(user->spanstats.dm));
52     PetscCall(DMGetLabel(user->spanstats.dm, "Face Sets", &label));
53 
54     PetscCall(DMPlexSetClosurePermutationTensor(user->spanstats.dm, PETSC_DETERMINE, NULL));
55     PetscCall(PetscFEDestroy(&fe));
56   }
57 
58   PetscSection section;
59   PetscCall(DMGetLocalSection(user->spanstats.dm, &section));
60   PetscCall(PetscSectionSetFieldName(section, 0, ""));
61   PetscCall(PetscSectionSetComponentName(section, 0, 0, "Mean Density"));
62   PetscCall(PetscSectionSetComponentName(section, 0, 1, "Mean Pressure"));
63   PetscCall(PetscSectionSetComponentName(section, 0, 2, "Mean Pressure Squared"));
64   PetscCall(PetscSectionSetComponentName(section, 0, 3, "Mean Pressure Velocity X"));
65   PetscCall(PetscSectionSetComponentName(section, 0, 4, "Mean Pressure Velocity Y"));
66   PetscCall(PetscSectionSetComponentName(section, 0, 5, "Mean Pressure Velocity Z"));
67   PetscCall(PetscSectionSetComponentName(section, 0, 6, "Mean Density Temperature"));
68   PetscCall(PetscSectionSetComponentName(section, 0, 7, "Mean Density Temperature Flux X"));
69   PetscCall(PetscSectionSetComponentName(section, 0, 8, "Mean Density Temperature Flux Y"));
70   PetscCall(PetscSectionSetComponentName(section, 0, 9, "Mean Density Temperature Flux Z"));
71   PetscCall(PetscSectionSetComponentName(section, 0, 10, "Mean Momentum X"));
72   PetscCall(PetscSectionSetComponentName(section, 0, 11, "Mean Momentum Y"));
73   PetscCall(PetscSectionSetComponentName(section, 0, 12, "Mean Momentum Z"));
74   PetscCall(PetscSectionSetComponentName(section, 0, 13, "Mean Momentum Flux XX"));
75   PetscCall(PetscSectionSetComponentName(section, 0, 14, "Mean Momentum Flux YY"));
76   PetscCall(PetscSectionSetComponentName(section, 0, 15, "Mean Momentum Flux ZZ"));
77   PetscCall(PetscSectionSetComponentName(section, 0, 16, "Mean Momentum Flux YZ"));
78   PetscCall(PetscSectionSetComponentName(section, 0, 17, "Mean Momentum Flux XZ"));
79   PetscCall(PetscSectionSetComponentName(section, 0, 18, "Mean Momentum Flux XY"));
80   PetscCall(PetscSectionSetComponentName(section, 0, 19, "Mean Velocity X"));
81   PetscCall(PetscSectionSetComponentName(section, 0, 20, "Mean Velocity Y"));
82   PetscCall(PetscSectionSetComponentName(section, 0, 21, "Mean Velocity Z"));
83 
84   PetscFunctionReturn(0);
85 }
86 
87 // Create CeedElemRestriction for collocated data based on associated CeedBasis and CeedElemRestriction
88 // Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction
89 PetscErrorCode CreateElemRestrColloc(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base,
90                                      CeedElemRestriction *elem_restr_collocated, CeedVector *l_vec, CeedVector *e_vec) {
91   CeedInt num_elem_qpts, loc_num_elem;
92   PetscFunctionBeginUser;
93 
94   CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts);
95   CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem);
96 
97   const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp};
98   CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides,
99                                    elem_restr_collocated);
100   CeedElemRestrictionCreateVector(*elem_restr_collocated, l_vec, e_vec);
101   PetscFunctionReturn(0);
102 }
103 
104 // Get coordinates of quadrature points
105 PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, CeedVector *qx_coords,
106                                    PetscInt *total_nqpnts) {
107   CeedQFunction       qf_quad_coords;
108   CeedOperator        op_quad_coords;
109   PetscInt            num_comp_x, loc_num_elem, num_elem_qpts;
110   CeedElemRestriction elem_restr_qx;
111   PetscFunctionBeginUser;
112 
113   // Create Element Restriction and CeedVector for quadrature coordinates
114   CeedBasisGetNumQuadraturePoints(basis_x, &num_elem_qpts);
115   CeedElemRestrictionGetNumElements(elem_restr_x, &loc_num_elem);
116   CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x);
117   *total_nqpnts = num_elem_qpts * loc_num_elem;
118   PetscCall(CreateElemRestrColloc(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx, qx_coords, NULL));
119 
120   // Create QFunction
121   CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords);
122 
123   // Create Operator
124   CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords);
125   CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
126   CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
127 
128   CeedOperatorApply(op_quad_coords, x_coords, *qx_coords, CEED_REQUEST_IMMEDIATE);
129 
130   CeedQFunctionDestroy(&qf_quad_coords);
131   CeedOperatorDestroy(&op_quad_coords);
132   PetscFunctionReturn(0);
133 }
134 
135 // Create PetscSF for child-to-parent communication
136 PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, DM parentdm, DM childdm, PetscSF statssf) {
137   PetscInt   child_num_qpnts, parent_num_qpnts, num_comp_x;
138   CeedVector child_qx_coords, parent_qx_coords;
139   PetscReal *child_coords, *parent_coords;
140   PetscFunctionBeginUser;
141 
142   // Assume that child and parent have the same number of components
143   CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x);
144   const PetscInt num_comp_sf = num_comp_x - 1;  // Number of coord components used in the creation of the SF
145 
146   // Get quad_coords for child DM
147   PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &child_qx_coords, &child_num_qpnts));
148 
149   // Get quad_coords for parent DM
150   PetscCall(GetQuadratureCoords(ceed, parentdm, ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, ceed_data->spanstats.x_coord,
151                                 &parent_qx_coords, &parent_num_qpnts));
152 
153   // Remove z component of coordinates for matching
154   {
155     const PetscReal *child_quad_coords, *parent_quad_coords;
156 
157     CeedVectorGetArrayRead(child_qx_coords, CEED_MEM_HOST, &child_quad_coords);
158     CeedVectorGetArrayRead(parent_qx_coords, CEED_MEM_HOST, &parent_quad_coords);
159 
160     PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords));
161     for (int i = 0; i < child_num_qpnts; i++) {
162       child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x];
163       child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x];
164     }
165     for (int i = 0; i < parent_num_qpnts; i++) {
166       parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x];
167       parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x];
168     }
169     CeedVectorRestoreArrayRead(child_qx_coords, &child_quad_coords);
170     CeedVectorRestoreArrayRead(parent_qx_coords, &parent_quad_coords);
171   }
172 
173   PetscCall(PetscSFSetGraphFromCoordinates(statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords));
174 
175   PetscCall(PetscSFViewFromOptions(statssf, NULL, "-spanstats_sf_view"));
176 
177   PetscCall(PetscFree2(child_coords, parent_coords));
178   CeedVectorDestroy(&child_qx_coords);
179   CeedVectorDestroy(&parent_qx_coords);
180   PetscFunctionReturn(0);
181 }
182 
183 // Compute mass matrix for statistics projection
184 PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data) {
185   CeedQFunction qf_mass;
186   CeedOperator  op_mass;
187   CeedInt       num_comp_q, q_data_size;
188   PetscFunctionBeginUser;
189 
190   // CEED Restriction
191   CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_stats, &num_comp_q);
192   CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_qd, &q_data_size);
193 
194   // Create Mass CeedOperator
195   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
196   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
197   CeedOperatorSetField(op_mass, "q", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE);
198   CeedOperatorSetField(op_mass, "qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, ceed_data->spanstats.q_data);
199   CeedOperatorSetField(op_mass, "v", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE);
200 
201   // Setup KSP for L^2 projection
202   {
203     MatopApplyContext M_ctx;
204     PetscInt          l_size, g_size;
205     Mat               mat_mass;
206     VecType           vec_type;
207     KSP               ksp;
208     Vec               ones, M_inv;
209     CeedVector        x_ceed, y_ceed;
210 
211     PetscCall(DMCreateGlobalVector(user->spanstats.dm, &M_inv));
212     PetscCall(VecGetLocalSize(M_inv, &l_size));
213     PetscCall(VecGetSize(M_inv, &g_size));
214     PetscCall(VecGetType(M_inv, &vec_type));
215 
216     PetscCall(PetscMalloc1(1, &M_ctx));
217     PetscCall(MatCreateShell(PETSC_COMM_WORLD, l_size, l_size, g_size, g_size, M_ctx, &mat_mass));
218     PetscCall(MatShellSetOperation(mat_mass, MATOP_MULT, (void (*)(void))MatMult_Ceed));
219     PetscCall(MatShellSetOperation(mat_mass, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed));
220     PetscCall(MatShellSetVecType(mat_mass, vec_type));
221 
222     CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &x_ceed, NULL);
223     CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &y_ceed, NULL);
224 
225     PetscCall(SetupMatopApplyCtx(PETSC_COMM_WORLD, user->spanstats.dm, user->ceed, op_mass, x_ceed, y_ceed, NULL, M_ctx));
226     user->spanstats.M_ctx = M_ctx;
227 
228     // Create lumped mass matrix inverse
229     PetscCall(DMGetGlobalVector(user->spanstats.dm, &ones));
230     PetscCall(VecZeroEntries(M_inv));
231     PetscCall(VecSet(ones, 1));
232     PetscCall(MatMult(mat_mass, ones, M_inv));
233     PetscCall(VecReciprocal(M_inv));
234     user->spanstats.M_inv = M_inv;
235     PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &ones));
236 
237     PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
238     PetscCall(KSPSetOptionsPrefix(ksp, "spanstats_"));
239     {
240       PC pc;
241       PetscCall(KSPGetPC(ksp, &pc));
242       PetscCall(PCSetType(pc, PCJACOBI));
243       PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
244       PetscCall(KSPSetType(ksp, KSPCG));
245       PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
246       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
247     }
248     PetscCall(KSPSetOperators(ksp, mat_mass, mat_mass));
249     PetscCall(KSPSetFromOptions(ksp));
250     user->spanstats.ksp = ksp;
251   }
252 
253   // Cleanup
254   CeedQFunctionDestroy(&qf_mass);
255   PetscFunctionReturn(0);
256 }
257 
258 // Create CeedOperators and KSP for the statistics collection and processing
259 PetscErrorCode CreateStatisticsOperators(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
260   CeedInt      num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q;
261   CeedOperator op_setup_sur;
262   PetscFunctionBeginUser;
263   CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q);
264 
265   // Create Operator for statistics collection
266   switch (user->phys->state_var) {
267     case STATEVAR_PRIMITIVE:
268       CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &ceed_data->spanstats.qf_stats_collect);
269       break;
270     case STATEVAR_CONSERVATIVE:
271       CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &ceed_data->spanstats.qf_stats_collect);
272       break;
273   }
274 
275   if (user->app_ctx->stats_test) {
276     CeedQFunctionDestroy(&ceed_data->spanstats.qf_stats_collect);
277     CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionTest, ChildStatsCollectionTest_loc, &ceed_data->spanstats.qf_stats_collect);
278   }
279 
280   CeedQFunctionSetContext(ceed_data->spanstats.qf_stats_collect, problem->apply_vol_ifunction.qfunction_context);
281   CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP);
282   CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE);
283   CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP);
284   CeedQFunctionAddOutput(ceed_data->spanstats.qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE);
285 
286   CeedOperatorCreate(ceed, ceed_data->spanstats.qf_stats_collect, NULL, NULL, &user->spanstats.op_stats_collect);
287   CeedOperatorSetField(user->spanstats.op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
288   CeedOperatorSetField(user->spanstats.op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
289   CeedOperatorSetField(user->spanstats.op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
290   CeedOperatorSetField(user->spanstats.op_stats_collect, "v", ceed_data->spanstats.elem_restr_child_colloc, CEED_BASIS_COLLOCATED,
291                        CEED_VECTOR_ACTIVE);
292 
293   // Create Operator for L^2 projection of statistics
294   // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function.
295   // Therefore, an Identity QF is sufficient
296   CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &ceed_data->spanstats.qf_stats_proj);
297 
298   CeedOperatorCreate(ceed, ceed_data->spanstats.qf_stats_proj, NULL, NULL, &user->spanstats.op_stats_proj);
299   CeedOperatorSetField(user->spanstats.op_stats_proj, "input", ceed_data->spanstats.elem_restr_parent_colloc, CEED_BASIS_COLLOCATED,
300                        CEED_VECTOR_ACTIVE);
301   CeedOperatorSetField(user->spanstats.op_stats_proj, "output", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats,
302                        CEED_VECTOR_ACTIVE);
303 
304   // Get q_data for lumped mass matrix formation
305   CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur);
306   CeedOperatorSetField(op_setup_sur, "dx", ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, CEED_VECTOR_ACTIVE);
307   CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->spanstats.basis_x, CEED_VECTOR_NONE);
308   CeedOperatorSetField(op_setup_sur, "surface qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
309   CeedOperatorApply(op_setup_sur, ceed_data->spanstats.x_coord, ceed_data->spanstats.q_data, CEED_REQUEST_IMMEDIATE);
310 
311   CeedOperatorDestroy(&op_setup_sur);
312   PetscFunctionReturn(0);
313 }
314 
315 PetscErrorCode SetupErrorTesting(Ceed ceed, User user, CeedData ceed_data) {
316   CeedInt       num_comp_stats = user->spanstats.num_comp_stats, num_comp_x;
317   CeedQFunction qf_error;
318   CeedOperator  op_error;
319   CeedInt       q_data_size;
320   CeedVector    x_ceed, y_ceed;
321   PetscFunctionBeginUser;
322 
323   CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_qd, &q_data_size);
324   CeedBasisGetNumComponents(ceed_data->spanstats.basis_x, &num_comp_x);
325 
326   CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionTest_Error, ChildStatsCollectionTest_Error_loc, &qf_error);
327   CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP);
328   CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE);
329   CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP);
330   CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP);
331 
332   CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
333   CeedOperatorSetField(op_error, "q", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE);
334   CeedOperatorSetField(op_error, "qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, ceed_data->spanstats.q_data);
335   CeedOperatorSetField(op_error, "x", ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, ceed_data->spanstats.x_coord);
336   CeedOperatorSetField(op_error, "v", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE);
337 
338   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &x_ceed, NULL);
339   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &y_ceed, NULL);
340 
341   PetscCall(PetscCalloc1(1, &user->spanstats.test_error_ctx));
342   PetscCall(SetupMatopApplyCtx(PETSC_COMM_WORLD, user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, user->spanstats.test_error_ctx));
343 
344   PetscFunctionReturn(0);
345 }
346 
347 // Setup for statistics collection
348 PetscErrorCode SetupStatsCollection(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
349   DM                 dm   = user->spanstats.dm;
350   MPI_Comm           comm = PetscObjectComm((PetscObject)dm);
351   CeedInt            dim, P, Q, num_comp_x;
352   Vec                X_loc;
353   PetscMemType       X_loc_memtype;
354   const PetscScalar *X_loc_array;
355   PetscFunctionBeginUser;
356 
357   PetscCall(DMGetDimension(dm, &dim));
358   CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &Q);
359   CeedBasisGetNumNodes1D(ceed_data->basis_q, &P);
360 
361   PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, problem->q_data_size_sur, &ceed_data->spanstats.elem_restr_parent_stats,
362                                     &ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.elem_restr_parent_qd));
363   CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_x, &num_comp_x);
364   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.x_coord, NULL);
365   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &user->spanstats.rhs_ceed, NULL);
366   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_qd, &ceed_data->spanstats.q_data, NULL);
367 
368   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS, &ceed_data->spanstats.basis_x);
369   CeedBasisCreateTensorH1Lagrange(ceed, dim, user->spanstats.num_comp_stats, P, Q, CEED_GAUSS, &ceed_data->spanstats.basis_stats);
370 
371   PetscCall(CreateElemRestrColloc(ceed, user->spanstats.num_comp_stats, ceed_data->spanstats.basis_stats,
372                                   ceed_data->spanstats.elem_restr_parent_stats, &ceed_data->spanstats.elem_restr_parent_colloc,
373                                   &user->spanstats.parent_stats, NULL));
374   PetscCall(CreateElemRestrColloc(ceed, user->spanstats.num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q,
375                                   &ceed_data->spanstats.elem_restr_child_colloc, &user->spanstats.child_stats, NULL));
376   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_child_colloc, &user->spanstats.child_inst_stats, NULL);
377   CeedVectorSetValue(user->spanstats.child_stats, 0);
378 
379   // -- Copy DM coordinates into CeedVector
380   {
381     DM cdm;
382     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
383     if (cdm) {
384       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
385     } else {
386       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
387     }
388   }
389   PetscCall(VecScale(X_loc, problem->dm_scale));
390   PetscCall(VecGetArrayReadAndMemType(X_loc, &X_loc_array, &X_loc_memtype));
391   CeedVectorSetArray(ceed_data->spanstats.x_coord, MemTypeP2C(X_loc_memtype), CEED_COPY_VALUES, (PetscScalar *)X_loc_array);
392   PetscCall(VecRestoreArrayRead(X_loc, &X_loc_array));
393 
394   // Create SF for communicating child data back their respective parents
395   PetscCall(PetscSFCreate(comm, &user->spanstats.sf));
396   PetscCall(CreateStatsSF(ceed, ceed_data, user->dm, user->spanstats.dm, user->spanstats.sf));
397 
398   // Create CeedOperators for statistics collection
399   PetscCall(CreateStatisticsOperators(ceed, user, ceed_data, problem));
400 
401   // Setup KSP and Mat for L^2 projection of statistics
402   PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data));
403 
404   if (user->app_ctx->stats_test) {
405     PetscCall(SetupErrorTesting(ceed, user, ceed_data));
406   }
407 
408   PetscFunctionReturn(0);
409 }
410 
411 // Collect statistics based on the solution Q
412 PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) {
413   PetscMemType       q_mem_type;
414   const PetscScalar *q_arr;
415   Vec                Q_loc;
416   PetscFunctionBeginUser;
417 
418   PetscCall(DMGetLocalVector(user->dm, &Q_loc));
419   PetscCall(VecZeroEntries(Q_loc));
420   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
421 
422   PetscCall(VecGetArrayReadAndMemType(Q_loc, &q_arr, &q_mem_type));
423   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q_arr);
424 
425   CeedOperatorApply(user->spanstats.op_stats_collect, user->q_ceed, user->spanstats.child_inst_stats, CEED_REQUEST_IMMEDIATE);
426 
427   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
428   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q_arr));
429   PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
430 
431   // Record averaging using left rectangle rule
432   PetscScalar delta_t           = solution_time - user->spanstats.prev_time;
433   PetscScalar prev_timeinterval = user->spanstats.prev_time - user->app_ctx->cont_time;
434   CeedVectorScale(user->spanstats.child_stats, prev_timeinterval / (prev_timeinterval + delta_t));
435   CeedVectorAXPY(user->spanstats.child_stats, delta_t / (prev_timeinterval + delta_t), user->spanstats.child_inst_stats);
436   user->spanstats.prev_time = solution_time;
437 
438   PetscFunctionReturn(0);
439 }
440 
441 // Process the child statistics into parent statistics and project them onto stats
442 PetscErrorCode ProcessStatistics(User user, Vec *stats) {
443   Span_Stats         user_stats = user->spanstats;
444   const PetscScalar *child_stats;
445   PetscScalar       *parent_stats;
446   MPI_Datatype       unit;
447   Vec                rhs_loc, rhs;
448   PetscMemType       rhs_mem_type;
449   CeedScalar        *rhs_arr;
450   CeedMemType        ceed_mem_type;
451   PetscFunctionBeginUser;
452 
453   CeedGetPreferredMemType(user->ceed, &ceed_mem_type);
454   CeedVectorSetValue(user_stats.parent_stats, 0);
455 
456   CeedVectorGetArrayRead(user_stats.child_stats, ceed_mem_type, &child_stats);
457   CeedVectorGetArray(user_stats.parent_stats, ceed_mem_type, &parent_stats);
458 
459   if (user_stats.num_comp_stats == 1) unit = MPIU_REAL;
460   else {
461     PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit));
462     PetscCallMPI(MPI_Type_commit(&unit));
463   }
464 
465   PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
466   PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
467 
468   CeedVectorRestoreArrayRead(user_stats.child_stats, &child_stats);
469   CeedVectorRestoreArray(user_stats.parent_stats, &parent_stats);
470   PetscCallMPI(MPI_Type_free(&unit));
471 
472   CeedVectorScale(user_stats.parent_stats, 1 / user_stats.span_width);
473 
474   // L^2 projection with the parent_data
475   PetscCall(DMGetGlobalVector(user_stats.dm, &rhs));
476   PetscCall(DMGetLocalVector(user_stats.dm, &rhs_loc));
477   PetscCall(VecZeroEntries(rhs));
478   PetscCall(VecZeroEntries(rhs_loc));
479   PetscCall(VecGetArrayWriteAndMemType(rhs_loc, &rhs_arr, &rhs_mem_type));
480   CeedVectorSetArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), CEED_USE_POINTER, (PetscScalar *)rhs_arr);
481 
482   CeedOperatorApply(user_stats.op_stats_proj, user_stats.parent_stats, user_stats.rhs_ceed, CEED_REQUEST_IMMEDIATE);
483 
484   CeedVectorTakeArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), &rhs_arr);
485   PetscCall(VecRestoreArrayAndMemType(rhs_loc, &rhs_arr));
486   PetscCall(DMLocalToGlobal(user_stats.dm, rhs_loc, ADD_VALUES, rhs));
487 
488   PetscCall(VecDuplicate(rhs, stats));
489   PetscCall(VecPointwiseMult(*stats, rhs, user_stats.M_inv));
490 
491   PetscCall(KSPSolve(user_stats.ksp, rhs, *stats));
492 
493   PetscFunctionReturn(0);
494 }
495 
496 // TSMonitor for the statistics collection and processing
497 PetscErrorCode TSMonitor_Statistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) {
498   User user = (User)ctx;
499   Vec  stats;
500   PetscFunctionBeginUser;
501 
502   // Do not collect or process on the first step of the run (ie. on the initial condition)
503   if (steps == user->app_ctx->cont_steps) PetscFunctionReturn(0);
504 
505   if (steps % user->app_ctx->stats_collect_interval == 0) {
506     PetscCall(CollectStatistics(user, solution_time, Q));
507   }
508 
509   if (steps % user->app_ctx->stats_write_interval == 0 && user->app_ctx->stats_write_interval != -1) {
510     PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats));
511     PetscCall(ProcessStatistics(user, &stats));
512     PetscCall(VecViewFromOptions(stats, NULL, "-stats_write_view"));
513     PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats));
514   }
515   PetscFunctionReturn(0);
516 }
517 
518 // Function to be called at the end of a simulation
519 PetscErrorCode StatsCollectFinalCall(User user, PetscReal solution_time, Vec Q) {
520   Vec stats;
521   PetscFunctionBeginUser;
522 
523   PetscCall(CollectStatistics(user, solution_time, Q));
524 
525   PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats));
526   PetscCall(ProcessStatistics(user, &stats));
527   PetscCall(VecViewFromOptions(stats, NULL, "-stats_write_view"));
528 
529   if (user->app_ctx->stats_test) {
530     Vec error;
531     PetscCall(VecDuplicate(stats, &error));
532     PetscCall(ApplyLocal_Ceed(stats, error, user->spanstats.test_error_ctx));
533     PetscScalar error_sq = 0;
534     PetscCall(VecSum(error, &error_sq));
535     PetscScalar l2_error = sqrt(error_sq);
536     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error));
537   }
538   PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats));
539 
540   PetscFunctionReturn(0);
541 }
542