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