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