xref: /libCEED/examples/fluids/src/turb_spanstats.c (revision f967ad79f48738b8a2e570b74f942e62880c38e7)
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) {
253a4208e6SJames Wright   user->spanstats.num_comp_stats = TURB_NUM_COMPONENTS;
26a175e481SJames Wright   PetscReal     domain_min[3], domain_max[3];
2778837792SJames Wright   PetscFE       fe;
2878837792SJames Wright   PetscSection  section;
294eed8630SJames Wright   PetscLogStage stage_stats_setup;
3051ee423eSJames Wright   PetscFunctionBeginUser;
3151ee423eSJames Wright 
324eed8630SJames Wright   PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup));
334eed8630SJames Wright   if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup));
344eed8630SJames Wright   PetscCall(PetscLogStagePush(stage_stats_setup));
354eed8630SJames Wright 
36a175e481SJames Wright   // Get spanwise length
37a175e481SJames Wright   PetscCall(DMGetBoundingBox(user->dm, domain_min, domain_max));
38a175e481SJames Wright   user->spanstats.span_width = domain_max[2] - domain_min[1];
39a175e481SJames Wright 
40*f967ad79SJames Wright   {  // Get DM from surface
41c9198418SJames Wright     PetscSF isoperiodicface;
4251ee423eSJames Wright     DMLabel label;
43c9198418SJames Wright 
44c9198418SJames Wright     PetscCall(DMPlexGetIsoperiodicFaceSF(user->dm, &isoperiodicface));
45c9198418SJames Wright 
46c9198418SJames Wright     if (isoperiodicface) {
47c9198418SJames Wright       PetscSF         inv_isoperiodicface;
48c9198418SJames Wright       PetscInt        nleaves;
49c9198418SJames Wright       const PetscInt *ilocal;
50c9198418SJames Wright 
51c9198418SJames Wright       PetscCall(PetscSFCreateInverseSF(isoperiodicface, &inv_isoperiodicface));
52c9198418SJames Wright       PetscCall(PetscSFGetGraph(inv_isoperiodicface, NULL, &nleaves, &ilocal, NULL));
53c9198418SJames Wright       PetscCall(DMCreateLabel(user->dm, "Periodic Face"));
54c9198418SJames Wright       PetscCall(DMGetLabel(user->dm, "Periodic Face", &label));
55c9198418SJames Wright       for (PetscInt i = 0; i < nleaves; i++) {
56c9198418SJames Wright         PetscCall(DMLabelSetValue(label, ilocal[i], 1));
57c9198418SJames Wright       }
58c9198418SJames Wright     } else {
5951ee423eSJames Wright       PetscCall(DMGetLabel(user->dm, "Face Sets", &label));
60c9198418SJames Wright     }
61c9198418SJames Wright 
6251ee423eSJames Wright     PetscCall(DMPlexLabelComplete(user->dm, label));
6351ee423eSJames Wright     PetscCall(DMPlexFilter(user->dm, label, 1, &user->spanstats.dm));
6451ee423eSJames Wright     PetscCall(DMProjectCoordinates(user->spanstats.dm, NULL));  // Ensure that a coordinate FE exists
6551ee423eSJames Wright   }
6651ee423eSJames Wright 
6751ee423eSJames Wright   PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats"));
68b7d66439SJames Wright   PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "turbulence_spanstats_"));
6951ee423eSJames Wright   PetscCall(DMSetFromOptions(user->spanstats.dm));
70a175e481SJames Wright   PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view"));  // -spanstats_dm_view
71c9198418SJames Wright 
7278837792SJames Wright   // Create FE space for parent DM
7351ee423eSJames Wright   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim - 1, user->spanstats.num_comp_stats, PETSC_FALSE, degree, PETSC_DECIDE, &fe));
7451ee423eSJames Wright   PetscCall(PetscObjectSetName((PetscObject)fe, "stats"));
7551ee423eSJames Wright   PetscCall(DMAddField(user->spanstats.dm, NULL, (PetscObject)fe));
7651ee423eSJames Wright   PetscCall(DMCreateDS(user->spanstats.dm));
7751ee423eSJames Wright   PetscCall(DMPlexSetClosurePermutationTensor(user->spanstats.dm, PETSC_DETERMINE, NULL));
7851ee423eSJames Wright 
7978837792SJames Wright   // Create Section for data
8051ee423eSJames Wright   PetscCall(DMGetLocalSection(user->spanstats.dm, &section));
8151ee423eSJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
823a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY, "MeanDensity"));
833a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE, "MeanPressure"));
843a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_SQUARED, "MeanPressureSquared"));
853a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_X, "MeanPressureVelocityX"));
863a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Y, "MeanPressureVelocityY"));
873a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Z, "MeanPressureVelocityZ"));
883a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE, "MeanDensityTemperature"));
893a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_X, "MeanDensityTemperatureFluxX"));
903a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Y, "MeanDensityTemperatureFluxY"));
913a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Z, "MeanDensityTemperatureFluxZ"));
923a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_X, "MeanMomentumX"));
933a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Y, "MeanMomentumY"));
943a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Z, "MeanMomentumZ"));
953a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XX, "MeanMomentumFluxXX"));
963a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YY, "MeanMomentumFluxYY"));
973a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_ZZ, "MeanMomentumFluxZZ"));
983a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YZ, "MeanMomentumFluxYZ"));
993a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XZ, "MeanMomentumFluxXZ"));
1003a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XY, "MeanMomentumFluxXY"));
1013a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_X, "MeanVelocityX"));
1023a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Y, "MeanVelocityY"));
1033a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Z, "MeanVelocityZ"));
10419706a06SJames Wright 
10578837792SJames Wright   // Cleanup
10678837792SJames Wright   PetscCall(PetscFEDestroy(&fe));
10778837792SJames Wright 
1084eed8630SJames Wright   PetscCall(PetscLogStagePop());
10919706a06SJames Wright   PetscFunctionReturn(0);
11019706a06SJames Wright }
11119706a06SJames Wright 
1121737222fSJames Wright // Create CeedElemRestriction for collocated data based on associated CeedBasis and CeedElemRestriction
1131737222fSJames Wright // Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction
1141737222fSJames Wright PetscErrorCode CreateElemRestrColloc(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base,
1151737222fSJames Wright                                      CeedElemRestriction *elem_restr_collocated, CeedVector *l_vec, CeedVector *e_vec) {
1161737222fSJames Wright   CeedInt num_elem_qpts, loc_num_elem;
1171737222fSJames Wright   PetscFunctionBeginUser;
1181737222fSJames Wright 
1191737222fSJames Wright   CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts);
1201737222fSJames Wright   CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem);
1211737222fSJames Wright 
1221737222fSJames Wright   const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp};
1231737222fSJames Wright   CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides,
1241737222fSJames Wright                                    elem_restr_collocated);
1251737222fSJames Wright   CeedElemRestrictionCreateVector(*elem_restr_collocated, l_vec, e_vec);
1261737222fSJames Wright   PetscFunctionReturn(0);
1271737222fSJames Wright }
1281737222fSJames Wright 
129a175e481SJames Wright // Get coordinates of quadrature points
13019706a06SJames Wright PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, CeedVector *qx_coords,
13119706a06SJames Wright                                    PetscInt *total_nqpnts) {
13219706a06SJames Wright   CeedQFunction       qf_quad_coords;
13319706a06SJames Wright   CeedOperator        op_quad_coords;
13419706a06SJames Wright   PetscInt            num_comp_x, loc_num_elem, num_elem_qpts;
13519706a06SJames Wright   CeedElemRestriction elem_restr_qx;
13619706a06SJames Wright   PetscFunctionBeginUser;
13719706a06SJames Wright 
13819706a06SJames Wright   // Create Element Restriction and CeedVector for quadrature coordinates
13919706a06SJames Wright   CeedBasisGetNumQuadraturePoints(basis_x, &num_elem_qpts);
14019706a06SJames Wright   CeedElemRestrictionGetNumElements(elem_restr_x, &loc_num_elem);
14119706a06SJames Wright   CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x);
14219706a06SJames Wright   *total_nqpnts = num_elem_qpts * loc_num_elem;
1431737222fSJames Wright   PetscCall(CreateElemRestrColloc(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx, qx_coords, NULL));
14419706a06SJames Wright 
14519706a06SJames Wright   // Create QFunction
14619706a06SJames Wright   CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords);
14719706a06SJames Wright 
14819706a06SJames Wright   // Create Operator
14919706a06SJames Wright   CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords);
15019706a06SJames Wright   CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
15119706a06SJames Wright   CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
15219706a06SJames Wright 
15319706a06SJames Wright   CeedOperatorApply(op_quad_coords, x_coords, *qx_coords, CEED_REQUEST_IMMEDIATE);
15419706a06SJames Wright 
15519706a06SJames Wright   CeedQFunctionDestroy(&qf_quad_coords);
15619706a06SJames Wright   CeedOperatorDestroy(&op_quad_coords);
15719706a06SJames Wright   PetscFunctionReturn(0);
15819706a06SJames Wright }
15919706a06SJames Wright 
160a175e481SJames Wright // Create PetscSF for child-to-parent communication
16119706a06SJames Wright PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, DM parentdm, DM childdm, PetscSF statssf) {
16219706a06SJames Wright   PetscInt   child_num_qpnts, parent_num_qpnts, num_comp_x;
16319706a06SJames Wright   CeedVector child_qx_coords, parent_qx_coords;
16419706a06SJames Wright   PetscReal *child_coords, *parent_coords;
16519706a06SJames Wright   PetscFunctionBeginUser;
16619706a06SJames Wright 
16719706a06SJames Wright   // Assume that child and parent have the same number of components
16819706a06SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x);
16919706a06SJames Wright   const PetscInt num_comp_sf = num_comp_x - 1;  // Number of coord components used in the creation of the SF
17019706a06SJames Wright 
17119706a06SJames Wright   // Get quad_coords for child DM
172a175e481SJames Wright   PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &child_qx_coords, &child_num_qpnts));
17319706a06SJames Wright 
17419706a06SJames Wright   // Get quad_coords for parent DM
17519706a06SJames Wright   PetscCall(GetQuadratureCoords(ceed, parentdm, ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, ceed_data->spanstats.x_coord,
17619706a06SJames Wright                                 &parent_qx_coords, &parent_num_qpnts));
17719706a06SJames Wright 
17819706a06SJames Wright   // Remove z component of coordinates for matching
17919706a06SJames Wright   {
18019706a06SJames Wright     const PetscReal *child_quad_coords, *parent_quad_coords;
18119706a06SJames Wright 
18219706a06SJames Wright     CeedVectorGetArrayRead(child_qx_coords, CEED_MEM_HOST, &child_quad_coords);
18319706a06SJames Wright     CeedVectorGetArrayRead(parent_qx_coords, CEED_MEM_HOST, &parent_quad_coords);
18419706a06SJames Wright 
18519706a06SJames Wright     PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords));
18619706a06SJames Wright     for (int i = 0; i < child_num_qpnts; i++) {
18719706a06SJames Wright       child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x];
18819706a06SJames Wright       child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x];
18919706a06SJames Wright     }
19019706a06SJames Wright     for (int i = 0; i < parent_num_qpnts; i++) {
19119706a06SJames Wright       parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x];
19219706a06SJames Wright       parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x];
19319706a06SJames Wright     }
19419706a06SJames Wright     CeedVectorRestoreArrayRead(child_qx_coords, &child_quad_coords);
19519706a06SJames Wright     CeedVectorRestoreArrayRead(parent_qx_coords, &parent_quad_coords);
19619706a06SJames Wright   }
19719706a06SJames Wright 
19819706a06SJames Wright   PetscCall(PetscSFSetGraphFromCoordinates(statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords));
19919706a06SJames Wright 
200a175e481SJames Wright   PetscCall(PetscSFViewFromOptions(statssf, NULL, "-spanstats_sf_view"));
201a175e481SJames Wright 
20219706a06SJames Wright   PetscCall(PetscFree2(child_coords, parent_coords));
20319706a06SJames Wright   CeedVectorDestroy(&child_qx_coords);
20419706a06SJames Wright   CeedVectorDestroy(&parent_qx_coords);
20519706a06SJames Wright   PetscFunctionReturn(0);
20619706a06SJames Wright }
20719706a06SJames Wright 
208a175e481SJames Wright // Compute mass matrix for statistics projection
209a175e481SJames Wright PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data) {
210a175e481SJames Wright   CeedQFunction qf_mass;
211a175e481SJames Wright   CeedOperator  op_mass;
212a175e481SJames Wright   CeedInt       num_comp_q, q_data_size;
213*f967ad79SJames Wright   MPI_Comm      comm = PetscObjectComm((PetscObject)user->spanstats.dm);
214a175e481SJames Wright   PetscFunctionBeginUser;
215a175e481SJames Wright 
216a175e481SJames Wright   // CEED Restriction
217a175e481SJames Wright   CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_stats, &num_comp_q);
218a175e481SJames Wright   CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_qd, &q_data_size);
219a175e481SJames Wright 
220a175e481SJames Wright   // Create Mass CeedOperator
221a175e481SJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
222a175e481SJames Wright   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
223a175e481SJames Wright   CeedOperatorSetField(op_mass, "q", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE);
224a175e481SJames Wright   CeedOperatorSetField(op_mass, "qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, ceed_data->spanstats.q_data);
225a175e481SJames Wright   CeedOperatorSetField(op_mass, "v", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE);
226a175e481SJames Wright 
227*f967ad79SJames Wright   {  // Setup KSP for L^2 projection
228a175e481SJames Wright     MatopApplyContext M_ctx;
229a175e481SJames Wright     PetscInt          l_size, g_size;
230a175e481SJames Wright     Mat               mat_mass;
231a175e481SJames Wright     VecType           vec_type;
232a175e481SJames Wright     KSP               ksp;
233a175e481SJames Wright     Vec               ones, M_inv;
234a175e481SJames Wright     CeedVector        x_ceed, y_ceed;
235a175e481SJames Wright 
236a175e481SJames Wright     PetscCall(DMCreateGlobalVector(user->spanstats.dm, &M_inv));
237a175e481SJames Wright     PetscCall(VecGetLocalSize(M_inv, &l_size));
238a175e481SJames Wright     PetscCall(VecGetSize(M_inv, &g_size));
239a175e481SJames Wright     PetscCall(VecGetType(M_inv, &vec_type));
240a175e481SJames Wright 
241a175e481SJames Wright     PetscCall(PetscMalloc1(1, &M_ctx));
242*f967ad79SJames Wright     PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, M_ctx, &mat_mass));
243a175e481SJames Wright     PetscCall(MatShellSetOperation(mat_mass, MATOP_MULT, (void (*)(void))MatMult_Ceed));
244a175e481SJames Wright     PetscCall(MatShellSetOperation(mat_mass, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed));
245a175e481SJames Wright     PetscCall(MatShellSetVecType(mat_mass, vec_type));
246a175e481SJames Wright 
247a175e481SJames Wright     CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &x_ceed, NULL);
248a175e481SJames Wright     CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &y_ceed, NULL);
249a175e481SJames Wright 
250*f967ad79SJames Wright     PetscCall(SetupMatopApplyCtx(comm, user->spanstats.dm, user->ceed, op_mass, x_ceed, y_ceed, NULL, M_ctx));
251a175e481SJames Wright     user->spanstats.M_ctx = M_ctx;
252a175e481SJames Wright 
253a175e481SJames Wright     // Create lumped mass matrix inverse
254a175e481SJames Wright     PetscCall(DMGetGlobalVector(user->spanstats.dm, &ones));
255a175e481SJames Wright     PetscCall(VecZeroEntries(M_inv));
256a175e481SJames Wright     PetscCall(VecSet(ones, 1));
257a175e481SJames Wright     PetscCall(MatMult(mat_mass, ones, M_inv));
258a175e481SJames Wright     PetscCall(VecReciprocal(M_inv));
259a175e481SJames Wright     user->spanstats.M_inv = M_inv;
260a175e481SJames Wright     PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &ones));
261a175e481SJames Wright 
262*f967ad79SJames Wright     PetscCall(KSPCreate(comm, &ksp));
263b7d66439SJames Wright     PetscCall(KSPSetOptionsPrefix(ksp, "turbulence_spanstats_"));
264a175e481SJames Wright     {
265a175e481SJames Wright       PC pc;
266a175e481SJames Wright       PetscCall(KSPGetPC(ksp, &pc));
267a175e481SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
268a175e481SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
269a175e481SJames Wright       PetscCall(KSPSetType(ksp, KSPCG));
270a175e481SJames Wright       PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
271a175e481SJames Wright       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
272a175e481SJames Wright     }
273a175e481SJames Wright     PetscCall(KSPSetOperators(ksp, mat_mass, mat_mass));
274a175e481SJames Wright     PetscCall(KSPSetFromOptions(ksp));
275a175e481SJames Wright     user->spanstats.ksp = ksp;
276a175e481SJames Wright   }
277a175e481SJames Wright 
278a175e481SJames Wright   // Cleanup
279a175e481SJames Wright   CeedQFunctionDestroy(&qf_mass);
280a175e481SJames Wright   PetscFunctionReturn(0);
281a175e481SJames Wright }
282a175e481SJames Wright 
283a175e481SJames Wright // Create CeedOperators and KSP for the statistics collection and processing
284a175e481SJames Wright PetscErrorCode CreateStatisticsOperators(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
285a175e481SJames Wright   CeedInt                     num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q;
286a175e481SJames Wright   CeedOperator                op_setup_sur;
287495b9769SJames Wright   Turbulence_SpanStatsContext collect_ctx;
288495b9769SJames Wright   NewtonianIdealGasContext    newtonian_ig_ctx;
289495b9769SJames Wright   CeedQFunctionContext        collect_context;
290a175e481SJames Wright   PetscFunctionBeginUser;
291a175e481SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q);
292a175e481SJames Wright 
293a175e481SJames Wright   // Create Operator for statistics collection
294a175e481SJames Wright   switch (user->phys->state_var) {
295a175e481SJames Wright     case STATEVAR_PRIMITIVE:
296a175e481SJames Wright       CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &ceed_data->spanstats.qf_stats_collect);
297a175e481SJames Wright       break;
298a175e481SJames Wright     case STATEVAR_CONSERVATIVE:
299a175e481SJames Wright       CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &ceed_data->spanstats.qf_stats_collect);
300a175e481SJames Wright       break;
301a175e481SJames Wright   }
302a175e481SJames Wright 
303b7d66439SJames Wright   if (user->spanstats.do_mms_test) {
304a175e481SJames Wright     CeedQFunctionDestroy(&ceed_data->spanstats.qf_stats_collect);
305b7d66439SJames Wright     CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest, ChildStatsCollectionMMSTest_loc, &ceed_data->spanstats.qf_stats_collect);
306a175e481SJames Wright   }
307a175e481SJames Wright 
308495b9769SJames Wright   // Setup Collection Context
309495b9769SJames Wright   {
310495b9769SJames Wright     PetscCall(PetscNew(&collect_ctx));
311495b9769SJames Wright     CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
312495b9769SJames Wright     collect_ctx->gas = *newtonian_ig_ctx;
313495b9769SJames Wright 
314495b9769SJames Wright     CeedQFunctionContextCreate(user->ceed, &collect_context);
315495b9769SJames Wright     CeedQFunctionContextSetData(collect_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx);
316495b9769SJames Wright     CeedQFunctionContextSetDataDestroy(collect_context, CEED_MEM_HOST, FreeContextPetsc);
317495b9769SJames Wright 
318495b9769SJames Wright     CeedQFunctionContextRegisterDouble(collect_context, "solution time", offsetof(struct Turbulence_SpanStatsContext_, solution_time), 1,
319495b9769SJames Wright                                        "Current solution time");
320495b9769SJames Wright     CeedQFunctionContextRegisterDouble(collect_context, "previous time", offsetof(struct Turbulence_SpanStatsContext_, previous_time), 1,
321495b9769SJames Wright                                        "Previous time statistics collection was done");
322495b9769SJames Wright 
323495b9769SJames Wright     CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
324495b9769SJames Wright   }
325495b9769SJames Wright 
326495b9769SJames Wright   CeedQFunctionSetContext(ceed_data->spanstats.qf_stats_collect, collect_context);
327495b9769SJames Wright   CeedQFunctionContextDestroy(&collect_context);
328a175e481SJames Wright   CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP);
329a175e481SJames Wright   CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE);
330a175e481SJames Wright   CeedQFunctionAddInput(ceed_data->spanstats.qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP);
331a175e481SJames Wright   CeedQFunctionAddOutput(ceed_data->spanstats.qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE);
332a175e481SJames Wright 
333a175e481SJames Wright   CeedOperatorCreate(ceed, ceed_data->spanstats.qf_stats_collect, NULL, NULL, &user->spanstats.op_stats_collect);
334a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
335a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
336a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
337a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_collect, "v", ceed_data->spanstats.elem_restr_child_colloc, CEED_BASIS_COLLOCATED,
338a175e481SJames Wright                        CEED_VECTOR_ACTIVE);
339a175e481SJames Wright 
340495b9769SJames Wright   CeedOperatorContextGetFieldLabel(user->spanstats.op_stats_collect, "solution time", &user->spanstats.solution_time_label);
341495b9769SJames Wright   CeedOperatorContextGetFieldLabel(user->spanstats.op_stats_collect, "previous time", &user->spanstats.previous_time_label);
342495b9769SJames Wright 
343823a1283SJames Wright   // Create Operator for L^2 projection of statistics
344a175e481SJames Wright   // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function.
345a175e481SJames Wright   // Therefore, an Identity QF is sufficient
346a175e481SJames Wright   CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &ceed_data->spanstats.qf_stats_proj);
347a175e481SJames Wright 
348a175e481SJames Wright   CeedOperatorCreate(ceed, ceed_data->spanstats.qf_stats_proj, NULL, NULL, &user->spanstats.op_stats_proj);
349a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_proj, "input", ceed_data->spanstats.elem_restr_parent_colloc, CEED_BASIS_COLLOCATED,
350a175e481SJames Wright                        CEED_VECTOR_ACTIVE);
351a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_proj, "output", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats,
352a175e481SJames Wright                        CEED_VECTOR_ACTIVE);
353a175e481SJames Wright 
354a175e481SJames Wright   // Get q_data for lumped mass matrix formation
355a175e481SJames Wright   CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur);
356a175e481SJames Wright   CeedOperatorSetField(op_setup_sur, "dx", ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, CEED_VECTOR_ACTIVE);
357a175e481SJames Wright   CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->spanstats.basis_x, CEED_VECTOR_NONE);
358a175e481SJames Wright   CeedOperatorSetField(op_setup_sur, "surface qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
359a175e481SJames Wright   CeedOperatorApply(op_setup_sur, ceed_data->spanstats.x_coord, ceed_data->spanstats.q_data, CEED_REQUEST_IMMEDIATE);
360a175e481SJames Wright 
361a175e481SJames Wright   CeedOperatorDestroy(&op_setup_sur);
362a175e481SJames Wright   PetscFunctionReturn(0);
363a175e481SJames Wright }
364a175e481SJames Wright 
365b7d66439SJames Wright // Creates operator for calculating error of method of manufactured solution (MMS) test
366b7d66439SJames Wright PetscErrorCode SetupMMSErrorChecking(Ceed ceed, User user, CeedData ceed_data) {
367823a1283SJames Wright   CeedInt       num_comp_stats = user->spanstats.num_comp_stats, num_comp_x;
368823a1283SJames Wright   CeedQFunction qf_error;
369823a1283SJames Wright   CeedOperator  op_error;
370823a1283SJames Wright   CeedInt       q_data_size;
371823a1283SJames Wright   CeedVector    x_ceed, y_ceed;
372823a1283SJames Wright   PetscFunctionBeginUser;
373823a1283SJames Wright 
374823a1283SJames Wright   CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_qd, &q_data_size);
375823a1283SJames Wright   CeedBasisGetNumComponents(ceed_data->spanstats.basis_x, &num_comp_x);
376823a1283SJames Wright 
377b7d66439SJames Wright   CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest_Error, ChildStatsCollectionMMSTest_Error_loc, &qf_error);
378823a1283SJames Wright   CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP);
379823a1283SJames Wright   CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE);
380823a1283SJames Wright   CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP);
381823a1283SJames Wright   CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP);
382823a1283SJames Wright 
383823a1283SJames Wright   CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
384823a1283SJames Wright   CeedOperatorSetField(op_error, "q", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE);
385823a1283SJames Wright   CeedOperatorSetField(op_error, "qdata", ceed_data->spanstats.elem_restr_parent_qd, CEED_BASIS_COLLOCATED, ceed_data->spanstats.q_data);
386823a1283SJames Wright   CeedOperatorSetField(op_error, "x", ceed_data->spanstats.elem_restr_parent_x, ceed_data->spanstats.basis_x, ceed_data->spanstats.x_coord);
387823a1283SJames Wright   CeedOperatorSetField(op_error, "v", ceed_data->spanstats.elem_restr_parent_stats, ceed_data->spanstats.basis_stats, CEED_VECTOR_ACTIVE);
388823a1283SJames Wright 
389823a1283SJames Wright   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &x_ceed, NULL);
390823a1283SJames Wright   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &y_ceed, NULL);
391823a1283SJames Wright 
392b7d66439SJames Wright   PetscCall(PetscCalloc1(1, &user->spanstats.mms_error_ctx));
393*f967ad79SJames Wright   PetscCall(SetupMatopApplyCtx(PetscObjectComm((PetscObject)user->spanstats.dm), user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL,
394*f967ad79SJames Wright                                user->spanstats.mms_error_ctx));
395823a1283SJames Wright 
396823a1283SJames Wright   PetscFunctionReturn(0);
397823a1283SJames Wright }
398823a1283SJames Wright 
399a175e481SJames Wright // Setup for statistics collection
40019706a06SJames Wright PetscErrorCode SetupStatsCollection(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
40119706a06SJames Wright   DM                 dm   = user->spanstats.dm;
40219706a06SJames Wright   MPI_Comm           comm = PetscObjectComm((PetscObject)dm);
40319706a06SJames Wright   CeedInt            dim, P, Q, num_comp_x;
40419706a06SJames Wright   Vec                X_loc;
40519706a06SJames Wright   PetscMemType       X_loc_memtype;
40619706a06SJames Wright   const PetscScalar *X_loc_array;
4074eed8630SJames Wright   PetscLogStage      stage_stats_setup;
40819706a06SJames Wright   PetscFunctionBeginUser;
40919706a06SJames Wright 
4104eed8630SJames Wright   PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup));
4114eed8630SJames Wright   if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup));
4124eed8630SJames Wright   PetscCall(PetscLogStagePush(stage_stats_setup));
4134eed8630SJames Wright 
41419706a06SJames Wright   PetscCall(DMGetDimension(dm, &dim));
415a175e481SJames Wright   CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &Q);
416a175e481SJames Wright   CeedBasisGetNumNodes1D(ceed_data->basis_q, &P);
41719706a06SJames Wright 
4181737222fSJames Wright   PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, problem->q_data_size_sur, &ceed_data->spanstats.elem_restr_parent_stats,
41919706a06SJames Wright                                     &ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.elem_restr_parent_qd));
42019706a06SJames Wright   CeedElemRestrictionGetNumComponents(ceed_data->spanstats.elem_restr_parent_x, &num_comp_x);
42119706a06SJames Wright   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_x, &ceed_data->spanstats.x_coord, NULL);
422a175e481SJames Wright   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_stats, &user->spanstats.rhs_ceed, NULL);
4231737222fSJames Wright   CeedElemRestrictionCreateVector(ceed_data->spanstats.elem_restr_parent_qd, &ceed_data->spanstats.q_data, NULL);
42419706a06SJames Wright 
425a175e481SJames Wright   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS, &ceed_data->spanstats.basis_x);
42619706a06SJames Wright   CeedBasisCreateTensorH1Lagrange(ceed, dim, user->spanstats.num_comp_stats, P, Q, CEED_GAUSS, &ceed_data->spanstats.basis_stats);
42719706a06SJames Wright 
4281737222fSJames Wright   PetscCall(CreateElemRestrColloc(ceed, user->spanstats.num_comp_stats, ceed_data->spanstats.basis_stats,
4291737222fSJames Wright                                   ceed_data->spanstats.elem_restr_parent_stats, &ceed_data->spanstats.elem_restr_parent_colloc,
4301737222fSJames Wright                                   &user->spanstats.parent_stats, NULL));
4311737222fSJames Wright   PetscCall(CreateElemRestrColloc(ceed, user->spanstats.num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q,
4321737222fSJames Wright                                   &ceed_data->spanstats.elem_restr_child_colloc, &user->spanstats.child_stats, NULL));
433a175e481SJames Wright   CeedVectorSetValue(user->spanstats.child_stats, 0);
4341737222fSJames Wright 
435*f967ad79SJames Wright   {  // -- Copy DM coordinates into CeedVector
43619706a06SJames Wright     DM cdm;
43719706a06SJames Wright     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
43819706a06SJames Wright     if (cdm) {
43919706a06SJames Wright       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
44019706a06SJames Wright     } else {
44119706a06SJames Wright       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
44219706a06SJames Wright     }
44319706a06SJames Wright   }
44419706a06SJames Wright   PetscCall(VecScale(X_loc, problem->dm_scale));
44519706a06SJames Wright   PetscCall(VecGetArrayReadAndMemType(X_loc, &X_loc_array, &X_loc_memtype));
44619706a06SJames Wright   CeedVectorSetArray(ceed_data->spanstats.x_coord, MemTypeP2C(X_loc_memtype), CEED_COPY_VALUES, (PetscScalar *)X_loc_array);
44719706a06SJames Wright   PetscCall(VecRestoreArrayRead(X_loc, &X_loc_array));
44819706a06SJames Wright 
44919706a06SJames Wright   // Create SF for communicating child data back their respective parents
45019706a06SJames Wright   PetscCall(PetscSFCreate(comm, &user->spanstats.sf));
45119706a06SJames Wright   PetscCall(CreateStatsSF(ceed, ceed_data, user->dm, user->spanstats.dm, user->spanstats.sf));
45251ee423eSJames Wright 
453a175e481SJames Wright   // Create CeedOperators for statistics collection
454a175e481SJames Wright   PetscCall(CreateStatisticsOperators(ceed, user, ceed_data, problem));
455a175e481SJames Wright 
456a175e481SJames Wright   // Setup KSP and Mat for L^2 projection of statistics
457a175e481SJames Wright   PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data));
458a175e481SJames Wright 
459b7d66439SJames Wright   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_turbulence_spanstats_mms", &user->spanstats.do_mms_test, NULL));
460b7d66439SJames Wright   if (user->spanstats.do_mms_test) {
461b7d66439SJames Wright     PetscCall(SetupMMSErrorChecking(ceed, user, ceed_data));
462823a1283SJames Wright   }
463823a1283SJames Wright 
464*f967ad79SJames Wright   {  // Setup stats viewer with prefix
4658ed52730SJames Wright     PetscViewerType viewer_type;
466b7d66439SJames Wright     PetscCall(PetscViewerGetType(user->app_ctx->turb_spanstats_viewer, &viewer_type));
467b7d66439SJames Wright     PetscCall(PetscOptionsSetValue(NULL, "-ts_monitor_turbulence_spanstats_viewer_type", viewer_type));
4688ed52730SJames Wright 
469b7d66439SJames Wright     PetscCall(PetscViewerSetOptionsPrefix(user->app_ctx->turb_spanstats_viewer, "ts_monitor_turbulence_spanstats_"));
470b7d66439SJames Wright     PetscCall(PetscViewerSetFromOptions(user->app_ctx->turb_spanstats_viewer));
4718ed52730SJames Wright   }
4724eed8630SJames Wright 
4734eed8630SJames Wright   PetscCall(PetscLogStagePop());
474a175e481SJames Wright   PetscFunctionReturn(0);
475a175e481SJames Wright }
476a175e481SJames Wright 
477a175e481SJames Wright // Collect statistics based on the solution Q
478a175e481SJames Wright PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) {
479a175e481SJames Wright   PetscMemType       q_mem_type;
480a175e481SJames Wright   const PetscScalar *q_arr;
481a175e481SJames Wright   PetscFunctionBeginUser;
482a175e481SJames Wright 
4836dcea3beSJames Wright   PetscLogStage stage_stats_collect;
4846dcea3beSJames Wright   PetscCall(PetscLogStageGetId("Stats Collect", &stage_stats_collect));
4856dcea3beSJames Wright   if (stage_stats_collect == -1) PetscCall(PetscLogStageRegister("Stats Collect", &stage_stats_collect));
4866dcea3beSJames Wright   PetscCall(PetscLogStagePush(stage_stats_collect));
4876dcea3beSJames Wright 
488a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, user->Q_loc, solution_time));
489495b9769SJames Wright   CeedOperatorContextSetDouble(user->spanstats.op_stats_collect, user->spanstats.solution_time_label, &solution_time);
490a0b9a424SJames Wright   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, user->Q_loc));
491a0b9a424SJames Wright   PetscCall(VecGetArrayReadAndMemType(user->Q_loc, &q_arr, &q_mem_type));
492a175e481SJames Wright   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q_arr);
493a175e481SJames Wright 
494495b9769SJames Wright   CeedOperatorApplyAdd(user->spanstats.op_stats_collect, user->q_ceed, user->spanstats.child_stats, CEED_REQUEST_IMMEDIATE);
495a175e481SJames Wright 
496a175e481SJames Wright   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
497a0b9a424SJames Wright   PetscCall(VecRestoreArrayReadAndMemType(user->Q_loc, &q_arr));
498a175e481SJames Wright 
499495b9769SJames Wright   CeedOperatorContextSetDouble(user->spanstats.op_stats_collect, user->spanstats.previous_time_label, &solution_time);
500a175e481SJames Wright 
5016dcea3beSJames Wright   PetscCall(PetscLogStagePop());
502a175e481SJames Wright   PetscFunctionReturn(0);
503a175e481SJames Wright }
504a175e481SJames Wright 
505a175e481SJames Wright // Process the child statistics into parent statistics and project them onto stats
50678bbfb6fSJed Brown PetscErrorCode ProcessStatistics(User user, Vec stats) {
507a175e481SJames Wright   Span_Stats         user_stats = user->spanstats;
508a175e481SJames Wright   const PetscScalar *child_stats;
509a175e481SJames Wright   PetscScalar       *parent_stats;
510a175e481SJames Wright   MPI_Datatype       unit;
511a175e481SJames Wright   Vec                rhs_loc, rhs;
512a175e481SJames Wright   PetscMemType       rhs_mem_type;
513a175e481SJames Wright   CeedScalar        *rhs_arr;
514a175e481SJames Wright   CeedMemType        ceed_mem_type;
515a175e481SJames Wright   PetscFunctionBeginUser;
516a175e481SJames Wright 
5176dcea3beSJames Wright   PetscLogStage stage_stats_process;
5186dcea3beSJames Wright   PetscCall(PetscLogStageGetId("Stats Process", &stage_stats_process));
5196dcea3beSJames Wright   if (stage_stats_process == -1) PetscCall(PetscLogStageRegister("Stats Process", &stage_stats_process));
5206dcea3beSJames Wright   PetscCall(PetscLogStagePush(stage_stats_process));
5216dcea3beSJames Wright 
522a175e481SJames Wright   CeedGetPreferredMemType(user->ceed, &ceed_mem_type);
523a175e481SJames Wright   CeedVectorSetValue(user_stats.parent_stats, 0);
524a175e481SJames Wright 
525a175e481SJames Wright   CeedVectorGetArrayRead(user_stats.child_stats, ceed_mem_type, &child_stats);
526a175e481SJames Wright   CeedVectorGetArray(user_stats.parent_stats, ceed_mem_type, &parent_stats);
527a175e481SJames Wright 
528a175e481SJames Wright   if (user_stats.num_comp_stats == 1) unit = MPIU_REAL;
529a175e481SJames Wright   else {
530a175e481SJames Wright     PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit));
531a175e481SJames Wright     PetscCallMPI(MPI_Type_commit(&unit));
532a175e481SJames Wright   }
533a175e481SJames Wright 
534a175e481SJames Wright   PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
535a175e481SJames Wright   PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
536a175e481SJames Wright 
537a175e481SJames Wright   CeedVectorRestoreArrayRead(user_stats.child_stats, &child_stats);
538a175e481SJames Wright   CeedVectorRestoreArray(user_stats.parent_stats, &parent_stats);
539a175e481SJames Wright   PetscCallMPI(MPI_Type_free(&unit));
540a175e481SJames Wright 
54178bbfb6fSJed Brown   PetscReal solution_time;
54278bbfb6fSJed Brown   PetscCall(DMGetOutputSequenceNumber(user_stats.dm, NULL, &solution_time));
54378bbfb6fSJed Brown   PetscReal summing_duration = solution_time - user->app_ctx->cont_time;
54478bbfb6fSJed Brown   CeedVectorScale(user_stats.parent_stats, 1 / (summing_duration * user_stats.span_width));
545a175e481SJames Wright 
546a175e481SJames Wright   // L^2 projection with the parent_data
547a175e481SJames Wright   PetscCall(DMGetLocalVector(user_stats.dm, &rhs_loc));
548a175e481SJames Wright   PetscCall(VecZeroEntries(rhs_loc));
549a175e481SJames Wright   PetscCall(VecGetArrayWriteAndMemType(rhs_loc, &rhs_arr, &rhs_mem_type));
550a175e481SJames Wright   CeedVectorSetArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), CEED_USE_POINTER, (PetscScalar *)rhs_arr);
551a175e481SJames Wright 
552a175e481SJames Wright   CeedOperatorApply(user_stats.op_stats_proj, user_stats.parent_stats, user_stats.rhs_ceed, CEED_REQUEST_IMMEDIATE);
553a175e481SJames Wright 
554a175e481SJames Wright   CeedVectorTakeArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), &rhs_arr);
555a175e481SJames Wright   PetscCall(VecRestoreArrayAndMemType(rhs_loc, &rhs_arr));
55678bbfb6fSJed Brown 
55778bbfb6fSJed Brown   PetscCall(DMGetGlobalVector(user_stats.dm, &rhs));
55878bbfb6fSJed Brown   PetscCall(VecZeroEntries(rhs));
559a175e481SJames Wright   PetscCall(DMLocalToGlobal(user_stats.dm, rhs_loc, ADD_VALUES, rhs));
56078bbfb6fSJed Brown   PetscCall(DMRestoreLocalVector(user_stats.dm, &rhs_loc));
561a175e481SJames Wright 
56278bbfb6fSJed Brown   PetscCall(KSPSolve(user_stats.ksp, rhs, stats));
563a175e481SJames Wright 
56478bbfb6fSJed Brown   PetscCall(DMRestoreGlobalVector(user_stats.dm, &rhs));
5656dcea3beSJames Wright   PetscCall(PetscLogStagePop());
566a175e481SJames Wright   PetscFunctionReturn(0);
567a175e481SJames Wright }
568a175e481SJames Wright 
569a175e481SJames Wright // TSMonitor for the statistics collection and processing
570a175e481SJames Wright PetscErrorCode TSMonitor_Statistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) {
571a175e481SJames Wright   User              user = (User)ctx;
572a175e481SJames Wright   Vec               stats;
57378bbfb6fSJed Brown   TSConvergedReason reason;
574b7d66439SJames Wright   PetscInt collect_interval = user->app_ctx->turb_spanstats_collect_interval, viewer_interval = user->app_ctx->turb_spanstats_viewer_interval;
575a175e481SJames Wright   PetscFunctionBeginUser;
57678bbfb6fSJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
577a175e481SJames Wright   // Do not collect or process on the first step of the run (ie. on the initial condition)
57878bbfb6fSJed Brown   if (steps == user->app_ctx->cont_steps && reason == TS_CONVERGED_ITERATING) PetscFunctionReturn(0);
579a175e481SJames Wright 
58027240365SJames Wright   PetscBool run_processing_and_viewer = (steps % viewer_interval == 0 && viewer_interval != -1) || reason != TS_CONVERGED_ITERATING;
581a175e481SJames Wright 
58227240365SJames Wright   if (steps % collect_interval == 0 || run_processing_and_viewer) {
58327240365SJames Wright     PetscCall(CollectStatistics(user, solution_time, Q));
58427240365SJames Wright 
58527240365SJames Wright     if (run_processing_and_viewer) {
5868ed52730SJames Wright       PetscCall(DMSetOutputSequenceNumber(user->spanstats.dm, steps, solution_time));
58778bbfb6fSJed Brown       PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats));
58878bbfb6fSJed Brown       PetscCall(ProcessStatistics(user, stats));
5898fb33541SJames Wright       if (user->app_ctx->test_type == TESTTYPE_NONE) {
590b7d66439SJames Wright         PetscCall(PetscViewerPushFormat(user->app_ctx->turb_spanstats_viewer, user->app_ctx->turb_spanstats_viewer_format));
591b7d66439SJames Wright         PetscCall(VecView(stats, user->app_ctx->turb_spanstats_viewer));
592b7d66439SJames Wright         PetscCall(PetscViewerPopFormat(user->app_ctx->turb_spanstats_viewer));
5938fb33541SJames Wright       }
5948fb33541SJames Wright       if (user->app_ctx->test_type == TESTTYPE_TURB_SPANSTATS && reason != TS_CONVERGED_ITERATING) {
5958fb33541SJames Wright         PetscCall(RegressionTests_NS(user->app_ctx, stats));
5968fb33541SJames Wright       }
597b7d66439SJames Wright       if (user->spanstats.do_mms_test && reason != TS_CONVERGED_ITERATING) {
598823a1283SJames Wright         Vec error;
599823a1283SJames Wright         PetscCall(VecDuplicate(stats, &error));
600b7d66439SJames Wright         PetscCall(ApplyLocal_Ceed(stats, error, user->spanstats.mms_error_ctx));
601823a1283SJames Wright         PetscScalar error_sq = 0;
602823a1283SJames Wright         PetscCall(VecSum(error, &error_sq));
603823a1283SJames Wright         PetscScalar l2_error = sqrt(error_sq);
604823a1283SJames Wright         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error));
605823a1283SJames Wright       }
606a175e481SJames Wright       PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats));
607522ee345SJames Wright     }
60827240365SJames Wright   }
60951ee423eSJames Wright   PetscFunctionReturn(0);
61051ee423eSJames Wright }
611fd170fd0SJames Wright 
612fd170fd0SJames Wright PetscErrorCode CleanupStats(User user, CeedData ceed_data) {
613fd170fd0SJames Wright   PetscFunctionBeginUser;
614fd170fd0SJames Wright 
615fd170fd0SJames Wright   // -- CeedVectors
616fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.child_stats);
617fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.parent_stats);
618fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.rhs_ceed);
619fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.x_ceed);
620fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.y_ceed);
621fd170fd0SJames Wright   CeedVectorDestroy(&ceed_data->spanstats.x_coord);
622fd170fd0SJames Wright   CeedVectorDestroy(&ceed_data->spanstats.q_data);
623fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.M_ctx->x_ceed);
624fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.M_ctx->y_ceed);
625b7d66439SJames Wright   if (user->spanstats.do_mms_test) {
626b7d66439SJames Wright     CeedVectorDestroy(&user->spanstats.mms_error_ctx->x_ceed);
627b7d66439SJames Wright     CeedVectorDestroy(&user->spanstats.mms_error_ctx->y_ceed);
628fd170fd0SJames Wright   }
629fd170fd0SJames Wright 
630fd170fd0SJames Wright   // -- QFunctions
631fd170fd0SJames Wright   CeedQFunctionDestroy(&ceed_data->spanstats.qf_stats_collect);
632fd170fd0SJames Wright   CeedQFunctionDestroy(&ceed_data->spanstats.qf_stats_proj);
633fd170fd0SJames Wright 
634fd170fd0SJames Wright   // -- CeedBasis
635fd170fd0SJames Wright   CeedBasisDestroy(&ceed_data->spanstats.basis_stats);
636fd170fd0SJames Wright   CeedBasisDestroy(&ceed_data->spanstats.basis_x);
637fd170fd0SJames Wright 
638fd170fd0SJames Wright   // -- CeedElemRestriction
639fd170fd0SJames Wright   CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_parent_stats);
640fd170fd0SJames Wright   CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_parent_qd);
641fd170fd0SJames Wright   CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_parent_x);
642fd170fd0SJames Wright   CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_parent_colloc);
643fd170fd0SJames Wright   CeedElemRestrictionDestroy(&ceed_data->spanstats.elem_restr_child_colloc);
644fd170fd0SJames Wright 
645fd170fd0SJames Wright   // -- CeedOperators
646fd170fd0SJames Wright   CeedOperatorDestroy(&user->spanstats.op_stats_collect);
647fd170fd0SJames Wright   CeedOperatorDestroy(&user->spanstats.op_stats_proj);
648fd170fd0SJames Wright   CeedOperatorDestroy(&user->spanstats.M_ctx->op);
649b7d66439SJames Wright   if (user->spanstats.do_mms_test) CeedOperatorDestroy(&user->spanstats.mms_error_ctx->op);
650fd170fd0SJames Wright 
651fd170fd0SJames Wright   // -- Vec
652fd170fd0SJames Wright   PetscCall(VecDestroy(&user->spanstats.M_inv));
653fd170fd0SJames Wright 
654fd170fd0SJames Wright   // -- KSP
655fd170fd0SJames Wright   PetscCall(KSPDestroy(&user->spanstats.ksp));
656fd170fd0SJames Wright 
657fd170fd0SJames Wright   // -- SF
658fd170fd0SJames Wright   PetscCall(PetscSFDestroy(&user->spanstats.sf));
659fd170fd0SJames Wright 
660fd170fd0SJames Wright   // -- DM
661fd170fd0SJames Wright   PetscCall(DMDestroy(&user->spanstats.dm));
662fd170fd0SJames Wright 
663fd170fd0SJames Wright   PetscFunctionReturn(0);
664fd170fd0SJames Wright }
665