xref: /libCEED/examples/fluids/src/turb_spanstats.c (revision c91984182f631443696c7d7792303473abfec66c)
151ee423eSJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
251ee423eSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
351ee423eSJames Wright //
451ee423eSJames Wright // SPDX-License-Identifier: BSD-2-Clause
551ee423eSJames Wright //
651ee423eSJames Wright // This file is part of CEED:  http://github.com/ceed
751ee423eSJames Wright /// @file
8a175e481SJames Wright /// Functions for setting up and performing statistics collection
9a175e481SJames Wright 
10a175e481SJames Wright #include "../qfunctions/turb_spanstats.h"
1151ee423eSJames Wright 
1219706a06SJames Wright #include <petscsf.h>
1319706a06SJames Wright 
14a175e481SJames Wright #include "../include/matops.h"
1551ee423eSJames Wright #include "../navierstokes.h"
16a175e481SJames Wright #include "ceed/ceed.h"
17a175e481SJames Wright #include "petscerror.h"
186dcea3beSJames Wright #include "petsclog.h"
19a175e481SJames Wright #include "petscmat.h"
2019706a06SJames Wright #include "petscsys.h"
21a175e481SJames Wright #include "petscvec.h"
228ed52730SJames Wright #include "petscviewer.h"
2351ee423eSJames Wright 
2451ee423eSJames Wright PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree, SimpleBC bc) {
25a175e481SJames Wright   user->spanstats.num_comp_stats = 22;
26a175e481SJames Wright   PetscReal domain_min[3], domain_max[3];
2751ee423eSJames Wright   PetscFunctionBeginUser;
2851ee423eSJames Wright 
29a175e481SJames Wright   // Get spanwise length
30a175e481SJames Wright   PetscCall(DMGetBoundingBox(user->dm, domain_min, domain_max));
31a175e481SJames Wright   user->spanstats.span_width = domain_max[2] - domain_min[1];
32a175e481SJames Wright 
3351ee423eSJames Wright   // Get DM from surface
3451ee423eSJames Wright   {
35*c9198418SJames Wright     PetscSF isoperiodicface;
3651ee423eSJames Wright     DMLabel label;
37*c9198418SJames Wright 
38*c9198418SJames Wright     PetscCall(DMPlexGetIsoperiodicFaceSF(user->dm, &isoperiodicface));
39*c9198418SJames Wright 
40*c9198418SJames Wright     if (isoperiodicface) {
41*c9198418SJames Wright       PetscSF         inv_isoperiodicface;
42*c9198418SJames Wright       PetscInt        nleaves;
43*c9198418SJames Wright       const PetscInt *ilocal;
44*c9198418SJames Wright 
45*c9198418SJames Wright       PetscCall(PetscSFCreateInverseSF(isoperiodicface, &inv_isoperiodicface));
46*c9198418SJames Wright       PetscCall(PetscSFGetGraph(inv_isoperiodicface, NULL, &nleaves, &ilocal, NULL));
47*c9198418SJames Wright       PetscCall(DMCreateLabel(user->dm, "Periodic Face"));
48*c9198418SJames Wright       PetscCall(DMGetLabel(user->dm, "Periodic Face", &label));
49*c9198418SJames Wright       for (PetscInt i = 0; i < nleaves; i++) {
50*c9198418SJames Wright         PetscCall(DMLabelSetValue(label, ilocal[i], 1));
51*c9198418SJames Wright       }
52*c9198418SJames Wright     } else {
5351ee423eSJames Wright       PetscCall(DMGetLabel(user->dm, "Face Sets", &label));
54*c9198418SJames Wright     }
55*c9198418SJames Wright 
5651ee423eSJames Wright     PetscCall(DMPlexLabelComplete(user->dm, label));
5751ee423eSJames Wright     PetscCall(DMPlexFilter(user->dm, label, 1, &user->spanstats.dm));
5851ee423eSJames Wright     PetscCall(DMProjectCoordinates(user->spanstats.dm, NULL));  // Ensure that a coordinate FE exists
5951ee423eSJames Wright   }
6051ee423eSJames Wright 
6151ee423eSJames Wright   PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats"));
6251ee423eSJames Wright   PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "spanstats_"));
6351ee423eSJames Wright   PetscCall(DMSetFromOptions(user->spanstats.dm));
64a175e481SJames Wright   PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view"));  // -spanstats_dm_view
65*c9198418SJames Wright 
6651ee423eSJames Wright   {
6751ee423eSJames Wright     PetscFE fe;
6851ee423eSJames Wright     DMLabel label;
6951ee423eSJames Wright 
7051ee423eSJames Wright     PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim - 1, user->spanstats.num_comp_stats, PETSC_FALSE, degree, PETSC_DECIDE, &fe));
7151ee423eSJames Wright     PetscCall(PetscObjectSetName((PetscObject)fe, "stats"));
7251ee423eSJames Wright     PetscCall(DMAddField(user->spanstats.dm, NULL, (PetscObject)fe));
7351ee423eSJames Wright     PetscCall(DMCreateDS(user->spanstats.dm));
7451ee423eSJames Wright     PetscCall(DMGetLabel(user->spanstats.dm, "Face Sets", &label));
7551ee423eSJames Wright 
7651ee423eSJames Wright     PetscCall(DMPlexSetClosurePermutationTensor(user->spanstats.dm, PETSC_DETERMINE, NULL));
7751ee423eSJames Wright     PetscCall(PetscFEDestroy(&fe));
7851ee423eSJames Wright   }
7951ee423eSJames Wright 
8051ee423eSJames Wright   PetscSection section;
8151ee423eSJames Wright   PetscCall(DMGetLocalSection(user->spanstats.dm, &section));
8251ee423eSJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
83a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 0, "Mean Density"));
84a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 1, "Mean Pressure"));
85a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 2, "Mean Pressure Squared"));
86a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 3, "Mean Pressure Velocity X"));
87a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 4, "Mean Pressure Velocity Y"));
88a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 5, "Mean Pressure Velocity Z"));
89a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 6, "Mean Density Temperature"));
90a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 7, "Mean Density Temperature Flux X"));
91a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 8, "Mean Density Temperature Flux Y"));
92a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 9, "Mean Density Temperature Flux Z"));
93a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 10, "Mean Momentum X"));
94a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 11, "Mean Momentum Y"));
95a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 12, "Mean Momentum Z"));
96a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 13, "Mean Momentum Flux XX"));
97a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 14, "Mean Momentum Flux YY"));
98a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 15, "Mean Momentum Flux ZZ"));
99a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 16, "Mean Momentum Flux YZ"));
100a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 17, "Mean Momentum Flux XZ"));
101a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 18, "Mean Momentum Flux XY"));
102a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 19, "Mean Velocity X"));
103a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 20, "Mean Velocity Y"));
104a175e481SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, 21, "Mean Velocity Z"));
10519706a06SJames Wright 
10619706a06SJames Wright   PetscFunctionReturn(0);
10719706a06SJames Wright }
10819706a06SJames Wright 
1091737222fSJames Wright // Create CeedElemRestriction for collocated data based on associated CeedBasis and CeedElemRestriction
1101737222fSJames Wright // Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction
1111737222fSJames Wright PetscErrorCode CreateElemRestrColloc(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base,
1121737222fSJames Wright                                      CeedElemRestriction *elem_restr_collocated, CeedVector *l_vec, CeedVector *e_vec) {
1131737222fSJames Wright   CeedInt num_elem_qpts, loc_num_elem;
1141737222fSJames Wright   PetscFunctionBeginUser;
1151737222fSJames Wright 
1161737222fSJames Wright   CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts);
1171737222fSJames Wright   CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem);
1181737222fSJames Wright 
1191737222fSJames Wright   const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp};
1201737222fSJames Wright   CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides,
1211737222fSJames Wright                                    elem_restr_collocated);
1221737222fSJames Wright   CeedElemRestrictionCreateVector(*elem_restr_collocated, l_vec, e_vec);
1231737222fSJames Wright   PetscFunctionReturn(0);
1241737222fSJames Wright }
1251737222fSJames Wright 
126a175e481SJames Wright // Get coordinates of quadrature points
12719706a06SJames Wright PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, CeedVector *qx_coords,
12819706a06SJames Wright                                    PetscInt *total_nqpnts) {
12919706a06SJames Wright   CeedQFunction       qf_quad_coords;
13019706a06SJames Wright   CeedOperator        op_quad_coords;
13119706a06SJames Wright   PetscInt            num_comp_x, loc_num_elem, num_elem_qpts;
13219706a06SJames Wright   CeedElemRestriction elem_restr_qx;
13319706a06SJames Wright   PetscFunctionBeginUser;
13419706a06SJames Wright 
13519706a06SJames Wright   // Create Element Restriction and CeedVector for quadrature coordinates
13619706a06SJames Wright   CeedBasisGetNumQuadraturePoints(basis_x, &num_elem_qpts);
13719706a06SJames Wright   CeedElemRestrictionGetNumElements(elem_restr_x, &loc_num_elem);
13819706a06SJames Wright   CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x);
13919706a06SJames Wright   *total_nqpnts = num_elem_qpts * loc_num_elem;
1401737222fSJames Wright   PetscCall(CreateElemRestrColloc(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx, qx_coords, NULL));
14119706a06SJames Wright 
14219706a06SJames Wright   // Create QFunction
14319706a06SJames Wright   CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords);
14419706a06SJames Wright 
14519706a06SJames Wright   // Create Operator
14619706a06SJames Wright   CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords);
14719706a06SJames Wright   CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
14819706a06SJames Wright   CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
14919706a06SJames Wright 
15019706a06SJames Wright   CeedOperatorApply(op_quad_coords, x_coords, *qx_coords, CEED_REQUEST_IMMEDIATE);
15119706a06SJames Wright 
15219706a06SJames Wright   CeedQFunctionDestroy(&qf_quad_coords);
15319706a06SJames Wright   CeedOperatorDestroy(&op_quad_coords);
15419706a06SJames Wright   PetscFunctionReturn(0);
15519706a06SJames Wright }
15619706a06SJames Wright 
157a175e481SJames Wright // Create PetscSF for child-to-parent communication
15819706a06SJames Wright PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, DM parentdm, DM childdm, PetscSF statssf) {
15919706a06SJames Wright   PetscInt   child_num_qpnts, parent_num_qpnts, num_comp_x;
16019706a06SJames Wright   CeedVector child_qx_coords, parent_qx_coords;
16119706a06SJames Wright   PetscReal *child_coords, *parent_coords;
16219706a06SJames Wright   PetscFunctionBeginUser;
16319706a06SJames Wright 
16419706a06SJames Wright   // Assume that child and parent have the same number of components
16519706a06SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x);
16619706a06SJames Wright   const PetscInt num_comp_sf = num_comp_x - 1;  // Number of coord components used in the creation of the SF
16719706a06SJames Wright 
16819706a06SJames Wright   // Get quad_coords for child DM
169a175e481SJames Wright   PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &child_qx_coords, &child_num_qpnts));
17019706a06SJames Wright 
17119706a06SJames Wright   // Get quad_coords for parent DM
17219706a06SJames Wright   PetscCall(GetQuadratureCoords(ceed, parentdm, ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, ceed_data->spanstats.x_coord,
17319706a06SJames Wright                                 &parent_qx_coords, &parent_num_qpnts));
17419706a06SJames Wright 
17519706a06SJames Wright   // Remove z component of coordinates for matching
17619706a06SJames Wright   {
17719706a06SJames Wright     const PetscReal *child_quad_coords, *parent_quad_coords;
17819706a06SJames Wright 
17919706a06SJames Wright     CeedVectorGetArrayRead(child_qx_coords, CEED_MEM_HOST, &child_quad_coords);
18019706a06SJames Wright     CeedVectorGetArrayRead(parent_qx_coords, CEED_MEM_HOST, &parent_quad_coords);
18119706a06SJames Wright 
18219706a06SJames Wright     PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords));
18319706a06SJames Wright     for (int i = 0; i < child_num_qpnts; i++) {
18419706a06SJames Wright       child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x];
18519706a06SJames Wright       child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x];
18619706a06SJames Wright     }
18719706a06SJames Wright     for (int i = 0; i < parent_num_qpnts; i++) {
18819706a06SJames Wright       parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x];
18919706a06SJames Wright       parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x];
19019706a06SJames Wright     }
19119706a06SJames Wright     CeedVectorRestoreArrayRead(child_qx_coords, &child_quad_coords);
19219706a06SJames Wright     CeedVectorRestoreArrayRead(parent_qx_coords, &parent_quad_coords);
19319706a06SJames Wright   }
19419706a06SJames Wright 
19519706a06SJames Wright   PetscCall(PetscSFSetGraphFromCoordinates(statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords));
19619706a06SJames Wright 
197a175e481SJames Wright   PetscCall(PetscSFViewFromOptions(statssf, NULL, "-spanstats_sf_view"));
198a175e481SJames Wright 
19919706a06SJames Wright   PetscCall(PetscFree2(child_coords, parent_coords));
20019706a06SJames Wright   CeedVectorDestroy(&child_qx_coords);
20119706a06SJames Wright   CeedVectorDestroy(&parent_qx_coords);
20219706a06SJames Wright   PetscFunctionReturn(0);
20319706a06SJames Wright }
20419706a06SJames Wright 
205a175e481SJames Wright // Compute mass matrix for statistics projection
206a175e481SJames Wright PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data) {
207a175e481SJames Wright   CeedQFunction qf_mass;
208a175e481SJames Wright   CeedOperator  op_mass;
209a175e481SJames Wright   CeedInt       num_comp_q, q_data_size;
210a175e481SJames Wright   PetscFunctionBeginUser;
211a175e481SJames Wright 
212a175e481SJames Wright   // CEED Restriction
213a175e481SJames Wright   CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_stats, &num_comp_q);
214a175e481SJames Wright   CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_qd, &q_data_size);
215a175e481SJames Wright 
216a175e481SJames Wright   // Create Mass CeedOperator
217a175e481SJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
218a175e481SJames Wright   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
219a175e481SJames Wright   CeedOperatorSetField(op_mass, "q", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE);
220a175e481SJames Wright   CeedOperatorSetField(op_mass, "qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, ceed_data->spanstats.q_data);
221a175e481SJames Wright   CeedOperatorSetField(op_mass, "v", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE);
222a175e481SJames Wright 
223a175e481SJames Wright   // Setup KSP for L^2 projection
224a175e481SJames Wright   {
225a175e481SJames Wright     MatopApplyContext M_ctx;
226a175e481SJames Wright     PetscInt          l_size, g_size;
227a175e481SJames Wright     Mat               mat_mass;
228a175e481SJames Wright     VecType           vec_type;
229a175e481SJames Wright     KSP               ksp;
230a175e481SJames Wright     Vec               ones, M_inv;
231a175e481SJames Wright     CeedVector        x_ceed, y_ceed;
232a175e481SJames Wright 
233a175e481SJames Wright     PetscCall(DMCreateGlobalVector(user->spanstats.dm, &M_inv));
234a175e481SJames Wright     PetscCall(VecGetLocalSize(M_inv, &l_size));
235a175e481SJames Wright     PetscCall(VecGetSize(M_inv, &g_size));
236a175e481SJames Wright     PetscCall(VecGetType(M_inv, &vec_type));
237a175e481SJames Wright 
238a175e481SJames Wright     PetscCall(PetscMalloc1(1, &M_ctx));
239a175e481SJames Wright     PetscCall(MatCreateShell(PETSC_COMM_WORLD, l_size, l_size, g_size, g_size, M_ctx, &mat_mass));
240a175e481SJames Wright     PetscCall(MatShellSetOperation(mat_mass, MATOP_MULT, (void (*)(void))MatMult_Ceed));
241a175e481SJames Wright     PetscCall(MatShellSetOperation(mat_mass, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed));
242a175e481SJames Wright     PetscCall(MatShellSetVecType(mat_mass, vec_type));
243a175e481SJames Wright 
244a175e481SJames Wright     CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &x_ceed, NULL);
245a175e481SJames Wright     CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &y_ceed, NULL);
246a175e481SJames Wright 
247a175e481SJames Wright     PetscCall(SetupMatopApplyCtx(PETSC_COMM_WORLD, user->spanstats.dm, user->ceed, op_mass, x_ceed, y_ceed, NULL, M_ctx));
248a175e481SJames Wright     user->spanstats.M_ctx = M_ctx;
249a175e481SJames Wright 
250a175e481SJames Wright     // Create lumped mass matrix inverse
251a175e481SJames Wright     PetscCall(DMGetGlobalVector(user->spanstats.dm, &ones));
252a175e481SJames Wright     PetscCall(VecZeroEntries(M_inv));
253a175e481SJames Wright     PetscCall(VecSet(ones, 1));
254a175e481SJames Wright     PetscCall(MatMult(mat_mass, ones, M_inv));
255a175e481SJames Wright     PetscCall(VecReciprocal(M_inv));
256a175e481SJames Wright     user->spanstats.M_inv = M_inv;
257a175e481SJames Wright     PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &ones));
258a175e481SJames Wright 
259a175e481SJames Wright     PetscCall(KSPCreate(PETSC_COMM_WORLD, &ksp));
260a175e481SJames Wright     PetscCall(KSPSetOptionsPrefix(ksp, "spanstats_"));
261a175e481SJames Wright     {
262a175e481SJames Wright       PC pc;
263a175e481SJames Wright       PetscCall(KSPGetPC(ksp, &pc));
264a175e481SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
265a175e481SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
266a175e481SJames Wright       PetscCall(KSPSetType(ksp, KSPCG));
267a175e481SJames Wright       PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
268a175e481SJames Wright       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
269a175e481SJames Wright     }
270a175e481SJames Wright     PetscCall(KSPSetOperators(ksp, mat_mass, mat_mass));
271a175e481SJames Wright     PetscCall(KSPSetFromOptions(ksp));
272a175e481SJames Wright     user->spanstats.ksp = ksp;
273a175e481SJames Wright   }
274a175e481SJames Wright 
275a175e481SJames Wright   // Cleanup
276a175e481SJames Wright   CeedQFunctionDestroy(&qf_mass);
277a175e481SJames Wright   PetscFunctionReturn(0);
278a175e481SJames Wright }
279a175e481SJames Wright 
280a175e481SJames Wright // Create CeedOperators and KSP for the statistics collection and processing
281a175e481SJames Wright PetscErrorCode CreateStatisticsOperators(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
282a175e481SJames Wright   CeedInt      num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q;
283a175e481SJames Wright   CeedOperator op_setup_sur;
284a175e481SJames Wright   PetscFunctionBeginUser;
285a175e481SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q);
286a175e481SJames Wright 
287a175e481SJames Wright   // Create Operator for statistics collection
288a175e481SJames Wright   switch (user->phys->state_var) {
289a175e481SJames Wright     case STATEVAR_PRIMITIVE:
290a175e481SJames Wright       CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &ceed_data->spanstats.qf_stats_collect);
291a175e481SJames Wright       break;
292a175e481SJames Wright     case STATEVAR_CONSERVATIVE:
293a175e481SJames Wright       CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &ceed_data->spanstats.qf_stats_collect);
294a175e481SJames Wright       break;
295a175e481SJames Wright   }
296a175e481SJames Wright 
297a175e481SJames Wright   if (user->app_ctx->stats_test) {
298a175e481SJames Wright     CeedQFunctionDestroy(&ceed_data->spanstats.qf_stats_collect);
299a175e481SJames Wright     CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionTest, ChildStatsCollectionTest_loc, &ceed_data->spanstats.qf_stats_collect);
300a175e481SJames Wright   }
301a175e481SJames Wright 
302a175e481SJames Wright   CeedQFunctionSetContext(ceed_data->spanstats.qf_stats_collect, problem->apply_vol_ifunction.qfunction_context);
303a175e481SJames Wright   CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP);
304a175e481SJames Wright   CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE);
305a175e481SJames Wright   CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP);
306a175e481SJames Wright   CeedQFunctionAddOutput(ceed_data->spanstats.qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE);
307a175e481SJames Wright 
308a175e481SJames Wright   CeedOperatorCreate(ceed, ceed_data->spanstats.qf_stats_collect, NULL, NULL, &user->spanstats.op_stats_collect);
309a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
310a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
311a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
312a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_collect, "v", ceed_data->spanstats.elem_restr_child_colloc, CEED_BASIS_COLLOCATED,
313a175e481SJames Wright                        CEED_VECTOR_ACTIVE);
314a175e481SJames Wright 
315823a1283SJames Wright   // Create Operator for L^2 projection of statistics
316a175e481SJames Wright   // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function.
317a175e481SJames Wright   // Therefore, an Identity QF is sufficient
318a175e481SJames Wright   CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &ceed_data->spanstats.qf_stats_proj);
319a175e481SJames Wright 
320a175e481SJames Wright   CeedOperatorCreate(ceed, ceed_data->spanstats.qf_stats_proj, NULL, NULL, &user->spanstats.op_stats_proj);
321a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_proj, "input", ceed_data->spanstats.elem_restr_parent_colloc, CEED_BASIS_COLLOCATED,
322a175e481SJames Wright                        CEED_VECTOR_ACTIVE);
323a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_proj, "output", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats,
324a175e481SJames Wright                        CEED_VECTOR_ACTIVE);
325a175e481SJames Wright 
326a175e481SJames Wright   // Get q_data for lumped mass matrix formation
327a175e481SJames Wright   CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur);
328a175e481SJames Wright   CeedOperatorSetField(op_setup_sur, "dx", ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, CEED_VECTOR_ACTIVE);
329a175e481SJames Wright   CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->spanstats.basis_x, CEED_VECTOR_NONE);
330a175e481SJames Wright   CeedOperatorSetField(op_setup_sur, "surface qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
331a175e481SJames Wright   CeedOperatorApply(op_setup_sur, ceed_data->spanstats.x_coord, ceed_data->spanstats.q_data, CEED_REQUEST_IMMEDIATE);
332a175e481SJames Wright 
333a175e481SJames Wright   CeedOperatorDestroy(&op_setup_sur);
334a175e481SJames Wright   PetscFunctionReturn(0);
335a175e481SJames Wright }
336a175e481SJames Wright 
337823a1283SJames Wright PetscErrorCode SetupErrorTesting(Ceed ceed, User user, CeedData ceed_data) {
338823a1283SJames Wright   CeedInt       num_comp_stats = user->spanstats.num_comp_stats, num_comp_x;
339823a1283SJames Wright   CeedQFunction qf_error;
340823a1283SJames Wright   CeedOperator  op_error;
341823a1283SJames Wright   CeedInt       q_data_size;
342823a1283SJames Wright   CeedVector    x_ceed, y_ceed;
343823a1283SJames Wright   PetscFunctionBeginUser;
344823a1283SJames Wright 
345823a1283SJames Wright   CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_qd, &q_data_size);
346823a1283SJames Wright   CeedBasisGetNumComponents(ceed_data->spanstats.basis_x, &num_comp_x);
347823a1283SJames Wright 
348823a1283SJames Wright   CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionTest_Error, ChildStatsCollectionTest_Error_loc, &qf_error);
349823a1283SJames Wright   CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP);
350823a1283SJames Wright   CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE);
351823a1283SJames Wright   CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP);
352823a1283SJames Wright   CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP);
353823a1283SJames Wright 
354823a1283SJames Wright   CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
355823a1283SJames Wright   CeedOperatorSetField(op_error, "q", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE);
356823a1283SJames Wright   CeedOperatorSetField(op_error, "qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, ceed_data->spanstats.q_data);
357823a1283SJames Wright   CeedOperatorSetField(op_error, "x", ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, ceed_data->spanstats.x_coord);
358823a1283SJames Wright   CeedOperatorSetField(op_error, "v", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE);
359823a1283SJames Wright 
360823a1283SJames Wright   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &x_ceed, NULL);
361823a1283SJames Wright   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &y_ceed, NULL);
362823a1283SJames Wright 
363823a1283SJames Wright   PetscCall(PetscCalloc1(1, &user->spanstats.test_error_ctx));
364823a1283SJames Wright   PetscCall(SetupMatopApplyCtx(PETSC_COMM_WORLD, user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, user->spanstats.test_error_ctx));
365823a1283SJames Wright 
366823a1283SJames Wright   PetscFunctionReturn(0);
367823a1283SJames Wright }
368823a1283SJames Wright 
369a175e481SJames Wright // Setup for statistics collection
37019706a06SJames Wright PetscErrorCode SetupStatsCollection(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
37119706a06SJames Wright   DM                 dm   = user->spanstats.dm;
37219706a06SJames Wright   MPI_Comm           comm = PetscObjectComm((PetscObject)dm);
37319706a06SJames Wright   CeedInt            dim, P, Q, num_comp_x;
37419706a06SJames Wright   Vec                X_loc;
37519706a06SJames Wright   PetscMemType       X_loc_memtype;
37619706a06SJames Wright   const PetscScalar *X_loc_array;
37719706a06SJames Wright   PetscFunctionBeginUser;
37819706a06SJames Wright 
37919706a06SJames Wright   PetscCall(DMGetDimension(dm, &dim));
380a175e481SJames Wright   CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &Q);
381a175e481SJames Wright   CeedBasisGetNumNodes1D(ceed_data->basis_q, &P);
38219706a06SJames Wright 
3831737222fSJames Wright   PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, problem->q_data_size_sur, &ceed_data->spanstats.elem_restr_parent_stats,
38419706a06SJames Wright                                     &ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.elem_restr_parent_qd));
38519706a06SJames Wright   CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_x, &num_comp_x);
38619706a06SJames Wright   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.x_coord, NULL);
387a175e481SJames Wright   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &user->spanstats.rhs_ceed, NULL);
3881737222fSJames Wright   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_qd, &ceed_data->spanstats.q_data, NULL);
38919706a06SJames Wright 
390a175e481SJames Wright   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS, &ceed_data->spanstats.basis_x);
39119706a06SJames Wright   CeedBasisCreateTensorH1Lagrange(ceed, dim, user->spanstats.num_comp_stats, P, Q, CEED_GAUSS, &ceed_data->spanstats.basis_stats);
39219706a06SJames Wright 
3931737222fSJames Wright   PetscCall(CreateElemRestrColloc(ceed, user->spanstats.num_comp_stats, ceed_data->spanstats.basis_stats,
3941737222fSJames Wright                                   ceed_data->spanstats.elem_restr_parent_stats, &ceed_data->spanstats.elem_restr_parent_colloc,
3951737222fSJames Wright                                   &user->spanstats.parent_stats, NULL));
3961737222fSJames Wright   PetscCall(CreateElemRestrColloc(ceed, user->spanstats.num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q,
3971737222fSJames Wright                                   &ceed_data->spanstats.elem_restr_child_colloc, &user->spanstats.child_stats, NULL));
398a175e481SJames Wright   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_child_colloc, &user->spanstats.child_inst_stats, NULL);
399a175e481SJames Wright   CeedVectorSetValue(user->spanstats.child_stats, 0);
4001737222fSJames Wright 
40119706a06SJames Wright   // -- Copy DM coordinates into CeedVector
40219706a06SJames Wright   {
40319706a06SJames Wright     DM cdm;
40419706a06SJames Wright     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
40519706a06SJames Wright     if (cdm) {
40619706a06SJames Wright       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
40719706a06SJames Wright     } else {
40819706a06SJames Wright       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
40919706a06SJames Wright     }
41019706a06SJames Wright   }
41119706a06SJames Wright   PetscCall(VecScale(X_loc, problem->dm_scale));
41219706a06SJames Wright   PetscCall(VecGetArrayReadAndMemType(X_loc, &X_loc_array, &X_loc_memtype));
41319706a06SJames Wright   CeedVectorSetArray(ceed_data->spanstats.x_coord, MemTypeP2C(X_loc_memtype), CEED_COPY_VALUES, (PetscScalar *)X_loc_array);
41419706a06SJames Wright   PetscCall(VecRestoreArrayRead(X_loc, &X_loc_array));
41519706a06SJames Wright 
41619706a06SJames Wright   // Create SF for communicating child data back their respective parents
41719706a06SJames Wright   PetscCall(PetscSFCreate(comm, &user->spanstats.sf));
41819706a06SJames Wright   PetscCall(CreateStatsSF(ceed, ceed_data, user->dm, user->spanstats.dm, user->spanstats.sf));
41951ee423eSJames Wright 
420a175e481SJames Wright   // Create CeedOperators for statistics collection
421a175e481SJames Wright   PetscCall(CreateStatisticsOperators(ceed, user, ceed_data, problem));
422a175e481SJames Wright 
423a175e481SJames Wright   // Setup KSP and Mat for L^2 projection of statistics
424a175e481SJames Wright   PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data));
425a175e481SJames Wright 
426823a1283SJames Wright   if (user->app_ctx->stats_test) {
427823a1283SJames Wright     PetscCall(SetupErrorTesting(ceed, user, ceed_data));
428823a1283SJames Wright   }
429823a1283SJames Wright 
4308ed52730SJames Wright   // Setup stats viewer with prefix
4318ed52730SJames Wright   {
4328ed52730SJames Wright     PetscViewerType viewer_type;
4338ed52730SJames Wright     PetscCall(PetscViewerGetType(user->app_ctx->stats_viewer, &viewer_type));
4348ed52730SJames Wright     PetscCall(PetscOptionsSetValue(NULL, "-stats_viewer_type", viewer_type));
4358ed52730SJames Wright 
4368ed52730SJames Wright     PetscCall(PetscViewerSetOptionsPrefix(user->app_ctx->stats_viewer, "stats_"));
4378ed52730SJames Wright     PetscCall(PetscViewerSetFromOptions(user->app_ctx->stats_viewer));
4388ed52730SJames Wright   }
439a175e481SJames Wright   PetscFunctionReturn(0);
440a175e481SJames Wright }
441a175e481SJames Wright 
442a175e481SJames Wright // Collect statistics based on the solution Q
443a175e481SJames Wright PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) {
444a175e481SJames Wright   PetscMemType       q_mem_type;
445a175e481SJames Wright   const PetscScalar *q_arr;
446a175e481SJames Wright   Vec                Q_loc;
447a175e481SJames Wright   PetscFunctionBeginUser;
448a175e481SJames Wright 
4496dcea3beSJames Wright   PetscLogStage stage_stats_collect;
4506dcea3beSJames Wright   PetscCall(PetscLogStageGetId("Stats Collect", &stage_stats_collect));
4516dcea3beSJames Wright   if (stage_stats_collect == -1) PetscCall(PetscLogStageRegister("Stats Collect", &stage_stats_collect));
4526dcea3beSJames Wright   PetscCall(PetscLogStagePush(stage_stats_collect));
4536dcea3beSJames Wright 
454a175e481SJames Wright   PetscCall(DMGetLocalVector(user->dm, &Q_loc));
455a175e481SJames Wright   PetscCall(VecZeroEntries(Q_loc));
456a175e481SJames Wright   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, Q_loc));
457a175e481SJames Wright 
458a175e481SJames Wright   PetscCall(VecGetArrayReadAndMemType(Q_loc, &q_arr, &q_mem_type));
459a175e481SJames Wright   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q_arr);
460a175e481SJames Wright 
461a175e481SJames Wright   CeedOperatorApply(user->spanstats.op_stats_collect, user->q_ceed, user->spanstats.child_inst_stats, CEED_REQUEST_IMMEDIATE);
462a175e481SJames Wright 
463a175e481SJames Wright   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
464a175e481SJames Wright   PetscCall(VecRestoreArrayReadAndMemType(Q_loc, &q_arr));
465a175e481SJames Wright   PetscCall(DMRestoreLocalVector(user->dm, &Q_loc));
466a175e481SJames Wright 
467a175e481SJames Wright   // Record averaging using left rectangle rule
468a175e481SJames Wright   PetscScalar delta_t           = solution_time - user->spanstats.prev_time;
469a175e481SJames Wright   PetscScalar prev_timeinterval = user->spanstats.prev_time - user->app_ctx->cont_time;
470a175e481SJames Wright   CeedVectorScale(user->spanstats.child_stats, prev_timeinterval / (prev_timeinterval + delta_t));
471a175e481SJames Wright   CeedVectorAXPY(user->spanstats.child_stats, delta_t / (prev_timeinterval + delta_t), user->spanstats.child_inst_stats);
472a175e481SJames Wright   user->spanstats.prev_time = solution_time;
473a175e481SJames Wright 
4746dcea3beSJames Wright   PetscCall(PetscLogStagePop());
475a175e481SJames Wright   PetscFunctionReturn(0);
476a175e481SJames Wright }
477a175e481SJames Wright 
478a175e481SJames Wright // Process the child statistics into parent statistics and project them onto stats
479a175e481SJames Wright PetscErrorCode ProcessStatistics(User user, Vec *stats) {
480a175e481SJames Wright   Span_Stats         user_stats = user->spanstats;
481a175e481SJames Wright   const PetscScalar *child_stats;
482a175e481SJames Wright   PetscScalar       *parent_stats;
483a175e481SJames Wright   MPI_Datatype       unit;
484a175e481SJames Wright   Vec                rhs_loc, rhs;
485a175e481SJames Wright   PetscMemType       rhs_mem_type;
486a175e481SJames Wright   CeedScalar        *rhs_arr;
487a175e481SJames Wright   CeedMemType        ceed_mem_type;
488a175e481SJames Wright   PetscFunctionBeginUser;
489a175e481SJames Wright 
4906dcea3beSJames Wright   PetscLogStage stage_stats_process;
4916dcea3beSJames Wright   PetscCall(PetscLogStageGetId("Stats Process", &stage_stats_process));
4926dcea3beSJames Wright   if (stage_stats_process == -1) PetscCall(PetscLogStageRegister("Stats Process", &stage_stats_process));
4936dcea3beSJames Wright   PetscCall(PetscLogStagePush(stage_stats_process));
4946dcea3beSJames Wright 
495a175e481SJames Wright   CeedGetPreferredMemType(user->ceed, &ceed_mem_type);
496a175e481SJames Wright   CeedVectorSetValue(user_stats.parent_stats, 0);
497a175e481SJames Wright 
498a175e481SJames Wright   CeedVectorGetArrayRead(user_stats.child_stats, ceed_mem_type, &child_stats);
499a175e481SJames Wright   CeedVectorGetArray(user_stats.parent_stats, ceed_mem_type, &parent_stats);
500a175e481SJames Wright 
501a175e481SJames Wright   if (user_stats.num_comp_stats == 1) unit = MPIU_REAL;
502a175e481SJames Wright   else {
503a175e481SJames Wright     PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit));
504a175e481SJames Wright     PetscCallMPI(MPI_Type_commit(&unit));
505a175e481SJames Wright   }
506a175e481SJames Wright 
507a175e481SJames Wright   PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
508a175e481SJames Wright   PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
509a175e481SJames Wright 
510a175e481SJames Wright   CeedVectorRestoreArrayRead(user_stats.child_stats, &child_stats);
511a175e481SJames Wright   CeedVectorRestoreArray(user_stats.parent_stats, &parent_stats);
512a175e481SJames Wright   PetscCallMPI(MPI_Type_free(&unit));
513a175e481SJames Wright 
514a175e481SJames Wright   CeedVectorScale(user_stats.parent_stats, 1 / user_stats.span_width);
515a175e481SJames Wright 
516a175e481SJames Wright   // L^2 projection with the parent_data
517a175e481SJames Wright   PetscCall(DMGetGlobalVector(user_stats.dm, &rhs));
518a175e481SJames Wright   PetscCall(DMGetLocalVector(user_stats.dm, &rhs_loc));
519a175e481SJames Wright   PetscCall(VecZeroEntries(rhs));
520a175e481SJames Wright   PetscCall(VecZeroEntries(rhs_loc));
521a175e481SJames Wright   PetscCall(VecGetArrayWriteAndMemType(rhs_loc, &rhs_arr, &rhs_mem_type));
522a175e481SJames Wright   CeedVectorSetArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), CEED_USE_POINTER, (PetscScalar *)rhs_arr);
523a175e481SJames Wright 
524a175e481SJames Wright   CeedOperatorApply(user_stats.op_stats_proj, user_stats.parent_stats, user_stats.rhs_ceed, CEED_REQUEST_IMMEDIATE);
525a175e481SJames Wright 
526a175e481SJames Wright   CeedVectorTakeArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), &rhs_arr);
527a175e481SJames Wright   PetscCall(VecRestoreArrayAndMemType(rhs_loc, &rhs_arr));
528a175e481SJames Wright   PetscCall(DMLocalToGlobal(user_stats.dm, rhs_loc, ADD_VALUES, rhs));
529a175e481SJames Wright 
530a175e481SJames Wright   PetscCall(VecDuplicate(rhs, stats));
531a175e481SJames Wright   PetscCall(VecPointwiseMult(*stats, rhs, user_stats.M_inv));
532a175e481SJames Wright 
533a175e481SJames Wright   PetscCall(KSPSolve(user_stats.ksp, rhs, *stats));
534a175e481SJames Wright 
5356dcea3beSJames Wright   PetscCall(PetscLogStagePop());
536a175e481SJames Wright   PetscFunctionReturn(0);
537a175e481SJames Wright }
538a175e481SJames Wright 
539a175e481SJames Wright // TSMonitor for the statistics collection and processing
540a175e481SJames Wright PetscErrorCode TSMonitor_Statistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) {
541a175e481SJames Wright   User user = (User)ctx;
542a175e481SJames Wright   Vec  stats;
543a175e481SJames Wright   PetscFunctionBeginUser;
544a175e481SJames Wright 
545a175e481SJames Wright   // Do not collect or process on the first step of the run (ie. on the initial condition)
546522ee345SJames Wright   if (steps == user->app_ctx->cont_steps && !user->spanstats.monitor_final_call) PetscFunctionReturn(0);
547a175e481SJames Wright 
548522ee345SJames Wright   if (steps % user->app_ctx->stats_collect_interval == 0 || user->spanstats.monitor_final_call) {
549a175e481SJames Wright     PetscCall(CollectStatistics(user, solution_time, Q));
550a175e481SJames Wright   }
551a175e481SJames Wright 
5528ed52730SJames Wright   if ((steps % user->app_ctx->stats_viewer_interval == 0 && user->app_ctx->stats_viewer_interval != -1) || user->spanstats.monitor_final_call) {
553a175e481SJames Wright     PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats));
554a175e481SJames Wright     PetscCall(ProcessStatistics(user, &stats));
5558ed52730SJames Wright     PetscCall(DMSetOutputSequenceNumber(user->spanstats.dm, steps, solution_time));
5568ed52730SJames Wright     PetscCall(PetscViewerPushFormat(user->app_ctx->stats_viewer, user->app_ctx->stats_viewer_format));
5578ed52730SJames Wright     PetscCall(VecView(stats, user->app_ctx->stats_viewer));
5588ed52730SJames Wright     PetscCall(PetscViewerPopFormat(user->app_ctx->stats_viewer));
559522ee345SJames Wright     if (user->app_ctx->stats_test && user->spanstats.monitor_final_call) {
560823a1283SJames Wright       Vec error;
561823a1283SJames Wright       PetscCall(VecDuplicate(stats, &error));
562823a1283SJames Wright       PetscCall(ApplyLocal_Ceed(stats, error, user->spanstats.test_error_ctx));
563823a1283SJames Wright       PetscScalar error_sq = 0;
564823a1283SJames Wright       PetscCall(VecSum(error, &error_sq));
565823a1283SJames Wright       PetscScalar l2_error = sqrt(error_sq);
566823a1283SJames Wright       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error));
567823a1283SJames Wright     }
568a175e481SJames Wright     PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats));
569522ee345SJames Wright   }
57051ee423eSJames Wright   PetscFunctionReturn(0);
57151ee423eSJames Wright }
572fd170fd0SJames Wright 
573fd170fd0SJames Wright PetscErrorCode CleanupStats(User user, CeedData ceed_data) {
574fd170fd0SJames Wright   PetscFunctionBeginUser;
575fd170fd0SJames Wright 
576fd170fd0SJames Wright   // -- CeedVectors
577fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.child_stats);
578fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.child_inst_stats);
579fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.parent_stats);
580fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.rhs_ceed);
581fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.x_ceed);
582fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.y_ceed);
583fd170fd0SJames Wright   CeedVectorDestroy(&ceed_data->spanstats.x_coord);
584fd170fd0SJames Wright   CeedVectorDestroy(&ceed_data->spanstats.q_data);
585fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.M_ctx->x_ceed);
586fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.M_ctx->y_ceed);
587fd170fd0SJames Wright   if (user->app_ctx->stats_test) {
588fd170fd0SJames Wright     CeedVectorDestroy(&user->spanstats.test_error_ctx->x_ceed);
589fd170fd0SJames Wright     CeedVectorDestroy(&user->spanstats.test_error_ctx->y_ceed);
590fd170fd0SJames Wright   }
591fd170fd0SJames Wright 
592fd170fd0SJames Wright   // -- QFunctions
593fd170fd0SJames Wright   CeedQFunctionDestroy(&ceed_data->spanstats.qf_stats_collect);
594fd170fd0SJames Wright   CeedQFunctionDestroy(&ceed_data->spanstats.qf_stats_proj);
595fd170fd0SJames Wright 
596fd170fd0SJames Wright   // -- CeedBasis
597fd170fd0SJames Wright   CeedBasisDestroy(&ceed_data->spanstats.basis_stats);
598fd170fd0SJames Wright   CeedBasisDestroy(&ceed_data->spanstats.basis_x);
599fd170fd0SJames Wright 
600fd170fd0SJames Wright   // -- CeedElemRestriction
601fd170fd0SJames Wright   CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_parent_stats);
602fd170fd0SJames Wright   CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_parent_qd);
603fd170fd0SJames Wright   CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_parent_x);
604fd170fd0SJames Wright   CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_parent_colloc);
605fd170fd0SJames Wright   CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_child_colloc);
606fd170fd0SJames Wright 
607fd170fd0SJames Wright   // -- CeedOperators
608fd170fd0SJames Wright   CeedOperatorDestroy(&user->spanstats.op_stats_collect);
609fd170fd0SJames Wright   CeedOperatorDestroy(&user->spanstats.op_stats_proj);
610fd170fd0SJames Wright   CeedOperatorDestroy(&user->spanstats.M_ctx->op);
611fd170fd0SJames Wright   if (user->app_ctx->stats_test) CeedOperatorDestroy(&user->spanstats.test_error_ctx->op);
612fd170fd0SJames Wright 
613fd170fd0SJames Wright   // -- Vec
614fd170fd0SJames Wright   PetscCall(VecDestroy(&user->spanstats.M_inv));
615fd170fd0SJames Wright 
616fd170fd0SJames Wright   // -- KSP
617fd170fd0SJames Wright   PetscCall(KSPDestroy(&user->spanstats.ksp));
618fd170fd0SJames Wright 
619fd170fd0SJames Wright   // -- SF
620fd170fd0SJames Wright   PetscCall(PetscSFDestroy(&user->spanstats.sf));
621fd170fd0SJames Wright 
622fd170fd0SJames Wright   // -- DM
623fd170fd0SJames Wright   PetscCall(DMDestroy(&user->spanstats.dm));
624fd170fd0SJames Wright 
625fd170fd0SJames Wright   PetscFunctionReturn(0);
626fd170fd0SJames Wright }
627