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