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