xref: /libCEED/examples/fluids/src/turb_spanstats.c (revision ea168fcf12290090da9ec39ea6d0fe45bfe59a86)
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 
12a175e481SJames Wright #include "../include/matops.h"
1351ee423eSJames Wright #include "../navierstokes.h"
1451ee423eSJames Wright 
15269a910fSJames Wright typedef struct {
16269a910fSJames Wright   CeedElemRestriction elem_restr_parent_x, elem_restr_parent_stats, elem_restr_parent_qd, elem_restr_parent_colloc, elem_restr_child_colloc;
17269a910fSJames Wright   CeedBasis           basis_x, basis_stats;
18269a910fSJames Wright   CeedVector          x_coord, q_data;
19269a910fSJames Wright } *SpanStatsSetupData;
20269a910fSJames Wright 
2151ee423eSJames Wright PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree, SimpleBC bc) {
223a4208e6SJames Wright   user->spanstats.num_comp_stats = TURB_NUM_COMPONENTS;
23a175e481SJames Wright   PetscReal     domain_min[3], domain_max[3];
2478837792SJames Wright   PetscFE       fe;
2578837792SJames Wright   PetscSection  section;
264eed8630SJames Wright   PetscLogStage stage_stats_setup;
271f595ac1SJames Wright   MPI_Comm      comm = PetscObjectComm((PetscObject)user->dm);
2851ee423eSJames Wright   PetscFunctionBeginUser;
2951ee423eSJames Wright 
304eed8630SJames Wright   PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup));
314eed8630SJames Wright   if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup));
324eed8630SJames Wright   PetscCall(PetscLogStagePush(stage_stats_setup));
334eed8630SJames Wright 
34a175e481SJames Wright   // Get spanwise length
35a175e481SJames Wright   PetscCall(DMGetBoundingBox(user->dm, domain_min, domain_max));
36a175e481SJames Wright   user->spanstats.span_width = domain_max[2] - domain_min[1];
37a175e481SJames Wright 
38f967ad79SJames Wright   {  // Get DM from surface
391f595ac1SJames Wright     DM          parent_distributed_dm;
40c9198418SJames Wright     PetscSF     isoperiodicface;
4151ee423eSJames Wright     DMLabel     label;
421f595ac1SJames Wright     PetscMPIInt size;
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
651f595ac1SJames Wright 
661f595ac1SJames Wright     PetscCall(DMPlexDistribute(user->spanstats.dm, 0, NULL, &parent_distributed_dm));
671f595ac1SJames Wright     PetscCallMPI(MPI_Comm_size(comm, &size));
681f595ac1SJames Wright     if (parent_distributed_dm) {
691f595ac1SJames Wright       PetscCall(DMDestroy(&user->spanstats.dm));
701f595ac1SJames Wright       user->spanstats.dm = parent_distributed_dm;
711f595ac1SJames Wright     } else if (size > 1) {
721f595ac1SJames Wright       PetscCall(PetscPrintf(comm, "WARNING: Turbulent spanwise statistics: parent DM could not be distributed accross %d ranks.\n", size));
731f595ac1SJames Wright     }
7451ee423eSJames Wright   }
7551ee423eSJames Wright 
7651ee423eSJames Wright   PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats"));
77b7d66439SJames Wright   PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "turbulence_spanstats_"));
7851ee423eSJames Wright   PetscCall(DMSetFromOptions(user->spanstats.dm));
79b57b8e72SJames Wright   PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view"));
80c9198418SJames Wright 
8178837792SJames Wright   // Create FE space for parent DM
8251ee423eSJames Wright   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim - 1, user->spanstats.num_comp_stats, PETSC_FALSE, degree, PETSC_DECIDE, &fe));
8351ee423eSJames Wright   PetscCall(PetscObjectSetName((PetscObject)fe, "stats"));
8451ee423eSJames Wright   PetscCall(DMAddField(user->spanstats.dm, NULL, (PetscObject)fe));
8551ee423eSJames Wright   PetscCall(DMCreateDS(user->spanstats.dm));
8651ee423eSJames Wright   PetscCall(DMPlexSetClosurePermutationTensor(user->spanstats.dm, PETSC_DETERMINE, NULL));
8751ee423eSJames Wright 
8878837792SJames Wright   // Create Section for data
8951ee423eSJames Wright   PetscCall(DMGetLocalSection(user->spanstats.dm, &section));
9051ee423eSJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
913a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY, "MeanDensity"));
923a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE, "MeanPressure"));
933a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_SQUARED, "MeanPressureSquared"));
943a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_X, "MeanPressureVelocityX"));
953a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Y, "MeanPressureVelocityY"));
963a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Z, "MeanPressureVelocityZ"));
973a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE, "MeanDensityTemperature"));
983a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_X, "MeanDensityTemperatureFluxX"));
993a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Y, "MeanDensityTemperatureFluxY"));
1003a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Z, "MeanDensityTemperatureFluxZ"));
1013a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_X, "MeanMomentumX"));
1023a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Y, "MeanMomentumY"));
1033a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Z, "MeanMomentumZ"));
1043a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XX, "MeanMomentumFluxXX"));
1053a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YY, "MeanMomentumFluxYY"));
1063a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_ZZ, "MeanMomentumFluxZZ"));
1073a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YZ, "MeanMomentumFluxYZ"));
1083a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XZ, "MeanMomentumFluxXZ"));
1093a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XY, "MeanMomentumFluxXY"));
1103a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_X, "MeanVelocityX"));
1113a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Y, "MeanVelocityY"));
1123a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Z, "MeanVelocityZ"));
11319706a06SJames Wright 
11478837792SJames Wright   // Cleanup
11578837792SJames Wright   PetscCall(PetscFEDestroy(&fe));
11678837792SJames Wright 
1174eed8630SJames Wright   PetscCall(PetscLogStagePop());
11819706a06SJames Wright   PetscFunctionReturn(0);
11919706a06SJames Wright }
12019706a06SJames Wright 
1211737222fSJames Wright // Create CeedElemRestriction for collocated data based on associated CeedBasis and CeedElemRestriction
1221737222fSJames Wright // Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction
1231737222fSJames Wright PetscErrorCode CreateElemRestrColloc(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base,
1241737222fSJames Wright                                      CeedElemRestriction *elem_restr_collocated, CeedVector *l_vec, CeedVector *e_vec) {
1251737222fSJames Wright   CeedInt num_elem_qpts, loc_num_elem;
1261737222fSJames Wright   PetscFunctionBeginUser;
1271737222fSJames Wright 
1281737222fSJames Wright   CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts);
1291737222fSJames Wright   CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem);
1301737222fSJames Wright 
1311737222fSJames Wright   const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp};
1321737222fSJames Wright   CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides,
1331737222fSJames Wright                                    elem_restr_collocated);
1341737222fSJames Wright   CeedElemRestrictionCreateVector(*elem_restr_collocated, l_vec, e_vec);
1351737222fSJames Wright   PetscFunctionReturn(0);
1361737222fSJames Wright }
1371737222fSJames Wright 
138a175e481SJames Wright // Get coordinates of quadrature points
13919706a06SJames Wright PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, CeedVector *qx_coords,
14019706a06SJames Wright                                    PetscInt *total_nqpnts) {
141269a910fSJames Wright   CeedElemRestriction elem_restr_qx;
14219706a06SJames Wright   CeedQFunction       qf_quad_coords;
14319706a06SJames Wright   CeedOperator        op_quad_coords;
14419706a06SJames Wright   PetscInt            num_comp_x, loc_num_elem, num_elem_qpts;
14519706a06SJames Wright   PetscFunctionBeginUser;
14619706a06SJames Wright 
14719706a06SJames Wright   // Create Element Restriction and CeedVector for quadrature coordinates
14819706a06SJames Wright   CeedBasisGetNumQuadraturePoints(basis_x, &num_elem_qpts);
14919706a06SJames Wright   CeedElemRestrictionGetNumElements(elem_restr_x, &loc_num_elem);
15019706a06SJames Wright   CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x);
15119706a06SJames Wright   *total_nqpnts = num_elem_qpts * loc_num_elem;
1521737222fSJames Wright   PetscCall(CreateElemRestrColloc(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx, qx_coords, NULL));
15319706a06SJames Wright 
15419706a06SJames Wright   // Create QFunction
15519706a06SJames Wright   CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords);
15619706a06SJames Wright 
15719706a06SJames Wright   // Create Operator
15819706a06SJames Wright   CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords);
15919706a06SJames Wright   CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
16019706a06SJames Wright   CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
16119706a06SJames Wright 
16219706a06SJames Wright   CeedOperatorApply(op_quad_coords, x_coords, *qx_coords, CEED_REQUEST_IMMEDIATE);
16319706a06SJames Wright 
164269a910fSJames Wright   CeedElemRestrictionDestroy(&elem_restr_qx);
16519706a06SJames Wright   CeedQFunctionDestroy(&qf_quad_coords);
16619706a06SJames Wright   CeedOperatorDestroy(&op_quad_coords);
16719706a06SJames Wright   PetscFunctionReturn(0);
16819706a06SJames Wright }
16919706a06SJames Wright 
170269a910fSJames Wright PetscErrorCode SpanStatsSetupDataCreate(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem, SpanStatsSetupData *stats_data) {
171269a910fSJames Wright   DM                 dm = user->spanstats.dm;
172269a910fSJames Wright   CeedInt            dim, P, Q, num_comp_x, num_comp_stats = user->spanstats.num_comp_stats;
173269a910fSJames Wright   Vec                X_loc;
174269a910fSJames Wright   PetscMemType       X_loc_memtype;
175269a910fSJames Wright   const PetscScalar *X_loc_array;
176269a910fSJames Wright   PetscFunctionBeginUser;
177269a910fSJames Wright 
178269a910fSJames Wright   PetscCall(PetscNew(stats_data));
179269a910fSJames Wright 
180269a910fSJames Wright   PetscCall(DMGetDimension(dm, &dim));
181269a910fSJames Wright   CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &Q);
182269a910fSJames Wright   CeedBasisGetNumNodes1D(ceed_data->basis_q, &P);
183269a910fSJames Wright 
184269a910fSJames Wright   PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, problem->q_data_size_sur, &(*stats_data)->elem_restr_parent_stats,
185269a910fSJames Wright                                     &(*stats_data)->elem_restr_parent_x, &(*stats_data)->elem_restr_parent_qd));
186269a910fSJames Wright   CeedElemRestrictionGetNumComponents((*stats_data)->elem_restr_parent_x, &num_comp_x);
187269a910fSJames Wright   CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_x, &(*stats_data)->x_coord, NULL);
188269a910fSJames Wright   CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_qd, &(*stats_data)->q_data, NULL);
189269a910fSJames Wright 
190269a910fSJames Wright   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS, &(*stats_data)->basis_x);
191269a910fSJames Wright   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_stats, P, Q, CEED_GAUSS, &(*stats_data)->basis_stats);
192269a910fSJames Wright 
193269a910fSJames Wright   PetscCall(CreateElemRestrColloc(ceed, num_comp_stats, (*stats_data)->basis_stats, (*stats_data)->elem_restr_parent_stats,
194269a910fSJames Wright                                   &(*stats_data)->elem_restr_parent_colloc, NULL, NULL));
195269a910fSJames Wright   PetscCall(
196269a910fSJames Wright       CreateElemRestrColloc(ceed, num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q, &(*stats_data)->elem_restr_child_colloc, NULL, NULL));
197269a910fSJames Wright 
198269a910fSJames Wright   {  // -- Copy DM coordinates into CeedVector
199269a910fSJames Wright     DM cdm;
200269a910fSJames Wright     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
201269a910fSJames Wright     if (cdm) {
202269a910fSJames Wright       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
203269a910fSJames Wright     } else {
204269a910fSJames Wright       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
205269a910fSJames Wright     }
206269a910fSJames Wright   }
207269a910fSJames Wright   PetscCall(VecScale(X_loc, problem->dm_scale));
208269a910fSJames Wright   PetscCall(VecGetArrayReadAndMemType(X_loc, &X_loc_array, &X_loc_memtype));
209269a910fSJames Wright   CeedVectorSetArray((*stats_data)->x_coord, MemTypeP2C(X_loc_memtype), CEED_COPY_VALUES, (PetscScalar *)X_loc_array);
210269a910fSJames Wright   PetscCall(VecRestoreArrayRead(X_loc, &X_loc_array));
211269a910fSJames Wright 
212269a910fSJames Wright   PetscFunctionReturn(0);
213269a910fSJames Wright }
214269a910fSJames Wright 
215269a910fSJames Wright PetscErrorCode SpanStatsSetupDataDestroy(SpanStatsSetupData data) {
216269a910fSJames Wright   PetscFunctionBeginUser;
217269a910fSJames Wright 
218269a910fSJames Wright   CeedElemRestrictionDestroy(&data->elem_restr_parent_x);
219269a910fSJames Wright   CeedElemRestrictionDestroy(&data->elem_restr_parent_stats);
220269a910fSJames Wright   CeedElemRestrictionDestroy(&data->elem_restr_parent_qd);
221269a910fSJames Wright   CeedElemRestrictionDestroy(&data->elem_restr_parent_colloc);
222269a910fSJames Wright   CeedElemRestrictionDestroy(&data->elem_restr_child_colloc);
223269a910fSJames Wright 
224269a910fSJames Wright   CeedBasisDestroy(&data->basis_x);
225269a910fSJames Wright   CeedBasisDestroy(&data->basis_stats);
226269a910fSJames Wright 
227269a910fSJames Wright   CeedVectorDestroy(&data->x_coord);
228269a910fSJames Wright   CeedVectorDestroy(&data->q_data);
229269a910fSJames Wright 
230269a910fSJames Wright   PetscCall(PetscFree(data));
231269a910fSJames Wright   PetscFunctionReturn(0);
232269a910fSJames Wright }
233269a910fSJames Wright 
234a175e481SJames Wright // Create PetscSF for child-to-parent communication
235269a910fSJames Wright PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, SpanStatsSetupData stats_data, DM parentdm, DM childdm, PetscSF *statssf) {
23619706a06SJames Wright   PetscInt   child_num_qpnts, parent_num_qpnts, num_comp_x;
23719706a06SJames Wright   CeedVector child_qx_coords, parent_qx_coords;
23819706a06SJames Wright   PetscReal *child_coords, *parent_coords;
23919706a06SJames Wright   PetscFunctionBeginUser;
24019706a06SJames Wright 
241269a910fSJames Wright   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)childdm), statssf));
242269a910fSJames Wright 
24319706a06SJames Wright   // Assume that child and parent have the same number of components
24419706a06SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x);
24519706a06SJames Wright   const PetscInt num_comp_sf = num_comp_x - 1;  // Number of coord components used in the creation of the SF
24619706a06SJames Wright 
24719706a06SJames Wright   // Get quad_coords for child DM
248a175e481SJames Wright   PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &child_qx_coords, &child_num_qpnts));
24919706a06SJames Wright 
25019706a06SJames Wright   // Get quad_coords for parent DM
251269a910fSJames Wright   PetscCall(GetQuadratureCoords(ceed, parentdm, stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord, &parent_qx_coords,
252269a910fSJames Wright                                 &parent_num_qpnts));
25319706a06SJames Wright 
25419706a06SJames Wright   // Remove z component of coordinates for matching
25519706a06SJames Wright   {
25619706a06SJames Wright     const PetscReal *child_quad_coords, *parent_quad_coords;
25719706a06SJames Wright 
25819706a06SJames Wright     CeedVectorGetArrayRead(child_qx_coords, CEED_MEM_HOST, &child_quad_coords);
25919706a06SJames Wright     CeedVectorGetArrayRead(parent_qx_coords, CEED_MEM_HOST, &parent_quad_coords);
26019706a06SJames Wright 
26119706a06SJames Wright     PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords));
26219706a06SJames Wright     for (int i = 0; i < child_num_qpnts; i++) {
26319706a06SJames Wright       child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x];
26419706a06SJames Wright       child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x];
26519706a06SJames Wright     }
26619706a06SJames Wright     for (int i = 0; i < parent_num_qpnts; i++) {
26719706a06SJames Wright       parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x];
26819706a06SJames Wright       parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x];
26919706a06SJames Wright     }
27019706a06SJames Wright     CeedVectorRestoreArrayRead(child_qx_coords, &child_quad_coords);
27119706a06SJames Wright     CeedVectorRestoreArrayRead(parent_qx_coords, &parent_quad_coords);
27219706a06SJames Wright   }
27319706a06SJames Wright 
274269a910fSJames Wright   PetscCall(PetscSFSetGraphFromCoordinates(*statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords));
27519706a06SJames Wright 
276269a910fSJames Wright   PetscCall(PetscSFViewFromOptions(*statssf, NULL, "-spanstats_sf_view"));
277a175e481SJames Wright 
27819706a06SJames Wright   PetscCall(PetscFree2(child_coords, parent_coords));
27919706a06SJames Wright   CeedVectorDestroy(&child_qx_coords);
28019706a06SJames Wright   CeedVectorDestroy(&parent_qx_coords);
28119706a06SJames Wright   PetscFunctionReturn(0);
28219706a06SJames Wright }
28319706a06SJames Wright 
284a175e481SJames Wright // Compute mass matrix for statistics projection
285269a910fSJames Wright PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) {
286269a910fSJames Wright   CeedQFunction qf_mass, qf_stats_proj;
287269a910fSJames Wright   CeedOperator  op_mass, op_setup_sur;
288269a910fSJames Wright   CeedInt       q_data_size, num_comp_stats = user->spanstats.num_comp_stats;
289f967ad79SJames Wright   MPI_Comm      comm = PetscObjectComm((PetscObject)user->spanstats.dm);
290a175e481SJames Wright   PetscFunctionBeginUser;
291a175e481SJames Wright 
292269a910fSJames Wright   // Create Operator for L^2 projection of statistics
293269a910fSJames Wright   // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function.
294269a910fSJames Wright   // Therefore, an Identity QF is sufficient
295269a910fSJames Wright   CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &qf_stats_proj);
296269a910fSJames Wright 
297269a910fSJames Wright   CeedOperatorCreate(ceed, qf_stats_proj, NULL, NULL, &user->spanstats.op_stats_proj);
298269a910fSJames Wright   CeedOperatorSetField(user->spanstats.op_stats_proj, "input", stats_data->elem_restr_parent_colloc, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
299269a910fSJames Wright   CeedOperatorSetField(user->spanstats.op_stats_proj, "output", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE);
300269a910fSJames Wright 
301269a910fSJames Wright   // Get q_data for mass matrix operator
302269a910fSJames Wright   CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur);
303269a910fSJames Wright   CeedOperatorSetField(op_setup_sur, "dx", stats_data->elem_restr_parent_x, stats_data->basis_x, CEED_VECTOR_ACTIVE);
304269a910fSJames Wright   CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, stats_data->basis_x, CEED_VECTOR_NONE);
305269a910fSJames Wright   CeedOperatorSetField(op_setup_sur, "surface qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
306269a910fSJames Wright   CeedOperatorApply(op_setup_sur, stats_data->x_coord, stats_data->q_data, CEED_REQUEST_IMMEDIATE);
307269a910fSJames Wright 
308a175e481SJames Wright   // CEED Restriction
309269a910fSJames Wright   CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size);
310a175e481SJames Wright 
311a175e481SJames Wright   // Create Mass CeedOperator
312269a910fSJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_stats, q_data_size, &qf_mass));
313a175e481SJames Wright   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
314269a910fSJames Wright   CeedOperatorSetField(op_mass, "u", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE);
315269a910fSJames Wright   CeedOperatorSetField(op_mass, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_COLLOCATED, stats_data->q_data);
316269a910fSJames Wright   CeedOperatorSetField(op_mass, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE);
317a175e481SJames Wright 
318f967ad79SJames Wright   {  // Setup KSP for L^2 projection
319a175e481SJames Wright     MatopApplyContext M_ctx;
320a175e481SJames Wright     PetscInt          l_size, g_size;
321a175e481SJames Wright     Mat               mat_mass;
322a175e481SJames Wright     VecType           vec_type;
323a175e481SJames Wright     KSP               ksp;
324a175e481SJames Wright     Vec               ones, M_inv;
325a175e481SJames Wright     CeedVector        x_ceed, y_ceed;
326a175e481SJames Wright 
327a175e481SJames Wright     PetscCall(DMCreateGlobalVector(user->spanstats.dm, &M_inv));
328a175e481SJames Wright     PetscCall(VecGetLocalSize(M_inv, &l_size));
329a175e481SJames Wright     PetscCall(VecGetSize(M_inv, &g_size));
330a175e481SJames Wright     PetscCall(VecGetType(M_inv, &vec_type));
331a175e481SJames Wright 
332269a910fSJames Wright     CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &x_ceed, NULL);
333269a910fSJames Wright     CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &y_ceed, NULL);
334*ea168fcfSJames Wright     PetscCall(MatopApplyContextCreate(user->spanstats.dm, user->spanstats.dm, user->ceed, op_mass, x_ceed, y_ceed, NULL, NULL, &M_ctx));
3356665b873SJames Wright     CeedVectorDestroy(&x_ceed);
3366665b873SJames Wright     CeedVectorDestroy(&y_ceed);
3376665b873SJames Wright 
338f967ad79SJames Wright     PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, M_ctx, &mat_mass));
3396665b873SJames Wright     PetscCall(MatShellSetContextDestroy(mat_mass, (PetscErrorCode(*)(void *))MatopApplyContextDestroy));
340a175e481SJames Wright     PetscCall(MatShellSetOperation(mat_mass, MATOP_MULT, (void (*)(void))MatMult_Ceed));
341a175e481SJames Wright     PetscCall(MatShellSetOperation(mat_mass, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag_Ceed));
342a175e481SJames Wright     PetscCall(MatShellSetVecType(mat_mass, vec_type));
343a175e481SJames Wright 
344a175e481SJames Wright     // Create lumped mass matrix inverse
345a175e481SJames Wright     PetscCall(DMGetGlobalVector(user->spanstats.dm, &ones));
346a175e481SJames Wright     PetscCall(VecZeroEntries(M_inv));
347a175e481SJames Wright     PetscCall(VecSet(ones, 1));
348a175e481SJames Wright     PetscCall(MatMult(mat_mass, ones, M_inv));
349a175e481SJames Wright     PetscCall(VecReciprocal(M_inv));
350a175e481SJames Wright     user->spanstats.M_inv = M_inv;
351a175e481SJames Wright     PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &ones));
352a175e481SJames Wright 
353f967ad79SJames Wright     PetscCall(KSPCreate(comm, &ksp));
354b7d66439SJames Wright     PetscCall(KSPSetOptionsPrefix(ksp, "turbulence_spanstats_"));
355a175e481SJames Wright     {
356a175e481SJames Wright       PC pc;
357a175e481SJames Wright       PetscCall(KSPGetPC(ksp, &pc));
358a175e481SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
359a175e481SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
360a175e481SJames Wright       PetscCall(KSPSetType(ksp, KSPCG));
361a175e481SJames Wright       PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
362a175e481SJames Wright       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
363a175e481SJames Wright     }
364a175e481SJames Wright     PetscCall(KSPSetOperators(ksp, mat_mass, mat_mass));
365a175e481SJames Wright     PetscCall(KSPSetFromOptions(ksp));
366a175e481SJames Wright     user->spanstats.ksp = ksp;
367a175e481SJames Wright   }
368a175e481SJames Wright 
369a175e481SJames Wright   // Cleanup
370a175e481SJames Wright   CeedQFunctionDestroy(&qf_mass);
371269a910fSJames Wright   CeedQFunctionDestroy(&qf_stats_proj);
3726665b873SJames Wright   CeedOperatorDestroy(&op_mass);
373269a910fSJames Wright   CeedOperatorDestroy(&op_setup_sur);
374a175e481SJames Wright   PetscFunctionReturn(0);
375a175e481SJames Wright }
376a175e481SJames Wright 
377269a910fSJames Wright // Create CeedOperator for statistics collection
378269a910fSJames Wright PetscErrorCode CreateStatisticCollectionOperator(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data, ProblemData *problem) {
379a175e481SJames Wright   CeedInt                     num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q;
380495b9769SJames Wright   Turbulence_SpanStatsContext collect_ctx;
381495b9769SJames Wright   NewtonianIdealGasContext    newtonian_ig_ctx;
382495b9769SJames Wright   CeedQFunctionContext        collect_context;
383269a910fSJames Wright   CeedQFunction               qf_stats_collect;
384a175e481SJames Wright   PetscFunctionBeginUser;
385a175e481SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q);
386a175e481SJames Wright 
387a175e481SJames Wright   // Create Operator for statistics collection
388a175e481SJames Wright   switch (user->phys->state_var) {
389a175e481SJames Wright     case STATEVAR_PRIMITIVE:
390269a910fSJames Wright       CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &qf_stats_collect);
391a175e481SJames Wright       break;
392a175e481SJames Wright     case STATEVAR_CONSERVATIVE:
393269a910fSJames Wright       CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &qf_stats_collect);
394a175e481SJames Wright       break;
395269a910fSJames Wright     default:
396269a910fSJames Wright       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No statisics collection available for chosen state variable");
397a175e481SJames Wright   }
398a175e481SJames Wright 
399b7d66439SJames Wright   if (user->spanstats.do_mms_test) {
400269a910fSJames Wright     CeedQFunctionDestroy(&qf_stats_collect);
401269a910fSJames Wright     CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest, ChildStatsCollectionMMSTest_loc, &qf_stats_collect);
402a175e481SJames Wright   }
403a175e481SJames Wright 
404269a910fSJames Wright   {  // Setup Collection Context
405495b9769SJames Wright     PetscCall(PetscNew(&collect_ctx));
406495b9769SJames Wright     CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
407495b9769SJames Wright     collect_ctx->gas = *newtonian_ig_ctx;
408495b9769SJames Wright 
409495b9769SJames Wright     CeedQFunctionContextCreate(user->ceed, &collect_context);
410495b9769SJames Wright     CeedQFunctionContextSetData(collect_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx);
411495b9769SJames Wright     CeedQFunctionContextSetDataDestroy(collect_context, CEED_MEM_HOST, FreeContextPetsc);
412495b9769SJames Wright 
413495b9769SJames Wright     CeedQFunctionContextRegisterDouble(collect_context, "solution time", offsetof(struct Turbulence_SpanStatsContext_, solution_time), 1,
414495b9769SJames Wright                                        "Current solution time");
415495b9769SJames Wright     CeedQFunctionContextRegisterDouble(collect_context, "previous time", offsetof(struct Turbulence_SpanStatsContext_, previous_time), 1,
416495b9769SJames Wright                                        "Previous time statistics collection was done");
417495b9769SJames Wright 
418495b9769SJames Wright     CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
419495b9769SJames Wright   }
420495b9769SJames Wright 
421269a910fSJames Wright   CeedQFunctionSetContext(qf_stats_collect, collect_context);
422495b9769SJames Wright   CeedQFunctionContextDestroy(&collect_context);
423269a910fSJames Wright   CeedQFunctionAddInput(qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP);
424269a910fSJames Wright   CeedQFunctionAddInput(qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE);
425269a910fSJames Wright   CeedQFunctionAddInput(qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP);
426269a910fSJames Wright   CeedQFunctionAddOutput(qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE);
427a175e481SJames Wright 
428269a910fSJames Wright   CeedOperatorCreate(ceed, qf_stats_collect, NULL, NULL, &user->spanstats.op_stats_collect);
429a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
430a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
431a175e481SJames Wright   CeedOperatorSetField(user->spanstats.op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
432269a910fSJames Wright   CeedOperatorSetField(user->spanstats.op_stats_collect, "v", stats_data->elem_restr_child_colloc, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
433a175e481SJames Wright 
434495b9769SJames Wright   CeedOperatorContextGetFieldLabel(user->spanstats.op_stats_collect, "solution time", &user->spanstats.solution_time_label);
435495b9769SJames Wright   CeedOperatorContextGetFieldLabel(user->spanstats.op_stats_collect, "previous time", &user->spanstats.previous_time_label);
436495b9769SJames Wright 
437269a910fSJames Wright   CeedQFunctionDestroy(&qf_stats_collect);
438a175e481SJames Wright   PetscFunctionReturn(0);
439a175e481SJames Wright }
440a175e481SJames Wright 
441b7d66439SJames Wright // Creates operator for calculating error of method of manufactured solution (MMS) test
442269a910fSJames Wright PetscErrorCode SetupMMSErrorChecking(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) {
443269a910fSJames Wright   CeedInt       num_comp_stats = user->spanstats.num_comp_stats, num_comp_x, q_data_size;
444823a1283SJames Wright   CeedQFunction qf_error;
445823a1283SJames Wright   CeedOperator  op_error;
446823a1283SJames Wright   CeedVector    x_ceed, y_ceed;
447823a1283SJames Wright   PetscFunctionBeginUser;
448823a1283SJames Wright 
449269a910fSJames Wright   CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size);
450269a910fSJames Wright   CeedBasisGetNumComponents(stats_data->basis_x, &num_comp_x);
451823a1283SJames Wright 
452b7d66439SJames Wright   CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest_Error, ChildStatsCollectionMMSTest_Error_loc, &qf_error);
453823a1283SJames Wright   CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP);
454823a1283SJames Wright   CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE);
455823a1283SJames Wright   CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP);
456823a1283SJames Wright   CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP);
457823a1283SJames Wright 
458823a1283SJames Wright   CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
459269a910fSJames Wright   CeedOperatorSetField(op_error, "q", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE);
460269a910fSJames Wright   CeedOperatorSetField(op_error, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_COLLOCATED, stats_data->q_data);
461269a910fSJames Wright   CeedOperatorSetField(op_error, "x", stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord);
462269a910fSJames Wright   CeedOperatorSetField(op_error, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE);
463823a1283SJames Wright 
464269a910fSJames Wright   CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &x_ceed, NULL);
465269a910fSJames Wright   CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &y_ceed, NULL);
466*ea168fcfSJames Wright   PetscCall(MatopApplyContextCreate(user->spanstats.dm, user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, NULL,
467*ea168fcfSJames Wright                                     &user->spanstats.mms_error_ctx));
468823a1283SJames Wright 
4696665b873SJames Wright   CeedOperatorDestroy(&op_error);
470269a910fSJames Wright   CeedQFunctionDestroy(&qf_error);
4716665b873SJames Wright   CeedVectorDestroy(&x_ceed);
4726665b873SJames Wright   CeedVectorDestroy(&y_ceed);
473823a1283SJames Wright   PetscFunctionReturn(0);
474823a1283SJames Wright }
475823a1283SJames Wright 
476a175e481SJames Wright // Setup for statistics collection
47719706a06SJames Wright PetscErrorCode SetupStatsCollection(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
478269a910fSJames Wright   SpanStatsSetupData stats_data;
4794eed8630SJames Wright   PetscLogStage      stage_stats_setup;
48019706a06SJames Wright   PetscFunctionBeginUser;
4814eed8630SJames Wright   PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup));
4824eed8630SJames Wright   if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup));
4834eed8630SJames Wright   PetscCall(PetscLogStagePush(stage_stats_setup));
4844eed8630SJames Wright 
485269a910fSJames Wright   // Create necessary CeedObjects for setting up statistics
486269a910fSJames Wright   PetscCall(SpanStatsSetupDataCreate(ceed, user, ceed_data, problem, &stats_data));
487269a910fSJames Wright   CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &user->spanstats.rhs_ceed, NULL);
488269a910fSJames Wright   CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_colloc, &user->spanstats.parent_stats, NULL);
489269a910fSJames Wright   CeedElemRestrictionCreateVector(stats_data->elem_restr_child_colloc, &user->spanstats.child_stats, NULL);
490a175e481SJames Wright   CeedVectorSetValue(user->spanstats.child_stats, 0);
4911737222fSJames Wright 
49219706a06SJames Wright   // Create SF for communicating child data back their respective parents
493269a910fSJames Wright   PetscCall(CreateStatsSF(ceed, ceed_data, stats_data, user->dm, user->spanstats.dm, &user->spanstats.sf));
49451ee423eSJames Wright 
495a175e481SJames Wright   // Create CeedOperators for statistics collection
496269a910fSJames Wright   PetscCall(CreateStatisticCollectionOperator(ceed, user, ceed_data, stats_data, problem));
497a175e481SJames Wright 
498a175e481SJames Wright   // Setup KSP and Mat for L^2 projection of statistics
499269a910fSJames Wright   PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data, stats_data));
500a175e481SJames Wright 
501b7d66439SJames Wright   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_turbulence_spanstats_mms", &user->spanstats.do_mms_test, NULL));
502b7d66439SJames Wright   if (user->spanstats.do_mms_test) {
503269a910fSJames Wright     PetscCall(SetupMMSErrorChecking(ceed, user, ceed_data, stats_data));
504823a1283SJames Wright   }
505823a1283SJames Wright 
506f967ad79SJames Wright   {  // Setup stats viewer with prefix
5078ed52730SJames Wright     PetscViewerType viewer_type;
508b7d66439SJames Wright     PetscCall(PetscViewerGetType(user->app_ctx->turb_spanstats_viewer, &viewer_type));
509b7d66439SJames Wright     PetscCall(PetscOptionsSetValue(NULL, "-ts_monitor_turbulence_spanstats_viewer_type", viewer_type));
5108ed52730SJames Wright 
511b7d66439SJames Wright     PetscCall(PetscViewerSetOptionsPrefix(user->app_ctx->turb_spanstats_viewer, "ts_monitor_turbulence_spanstats_"));
512b7d66439SJames Wright     PetscCall(PetscViewerSetFromOptions(user->app_ctx->turb_spanstats_viewer));
5138ed52730SJames Wright   }
5144eed8630SJames Wright 
515269a910fSJames Wright   PetscCall(SpanStatsSetupDataDestroy(stats_data));
5164eed8630SJames Wright   PetscCall(PetscLogStagePop());
517a175e481SJames Wright   PetscFunctionReturn(0);
518a175e481SJames Wright }
519a175e481SJames Wright 
520a175e481SJames Wright // Collect statistics based on the solution Q
521a175e481SJames Wright PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) {
522a175e481SJames Wright   PetscMemType       q_mem_type;
523a175e481SJames Wright   const PetscScalar *q_arr;
524a175e481SJames Wright   PetscFunctionBeginUser;
525a175e481SJames Wright 
5266dcea3beSJames Wright   PetscLogStage stage_stats_collect;
5276dcea3beSJames Wright   PetscCall(PetscLogStageGetId("Stats Collect", &stage_stats_collect));
5286dcea3beSJames Wright   if (stage_stats_collect == -1) PetscCall(PetscLogStageRegister("Stats Collect", &stage_stats_collect));
5296dcea3beSJames Wright   PetscCall(PetscLogStagePush(stage_stats_collect));
5306dcea3beSJames Wright 
531a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, user->Q_loc, solution_time));
532495b9769SJames Wright   CeedOperatorContextSetDouble(user->spanstats.op_stats_collect, user->spanstats.solution_time_label, &solution_time);
533a0b9a424SJames Wright   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, user->Q_loc));
534a0b9a424SJames Wright   PetscCall(VecGetArrayReadAndMemType(user->Q_loc, &q_arr, &q_mem_type));
535a175e481SJames Wright   CeedVectorSetArray(user->q_ceed, MemTypeP2C(q_mem_type), CEED_USE_POINTER, (PetscScalar *)q_arr);
536a175e481SJames Wright 
537495b9769SJames Wright   CeedOperatorApplyAdd(user->spanstats.op_stats_collect, user->q_ceed, user->spanstats.child_stats, CEED_REQUEST_IMMEDIATE);
538a175e481SJames Wright 
539a175e481SJames Wright   CeedVectorTakeArray(user->q_ceed, MemTypeP2C(q_mem_type), NULL);
540a0b9a424SJames Wright   PetscCall(VecRestoreArrayReadAndMemType(user->Q_loc, &q_arr));
541a175e481SJames Wright 
542495b9769SJames Wright   CeedOperatorContextSetDouble(user->spanstats.op_stats_collect, user->spanstats.previous_time_label, &solution_time);
543a175e481SJames Wright 
5446dcea3beSJames Wright   PetscCall(PetscLogStagePop());
545a175e481SJames Wright   PetscFunctionReturn(0);
546a175e481SJames Wright }
547a175e481SJames Wright 
548a175e481SJames Wright // Process the child statistics into parent statistics and project them onto stats
54978bbfb6fSJed Brown PetscErrorCode ProcessStatistics(User user, Vec stats) {
550a175e481SJames Wright   Span_Stats         user_stats = user->spanstats;
551a175e481SJames Wright   const PetscScalar *child_stats;
552a175e481SJames Wright   PetscScalar       *parent_stats;
553a175e481SJames Wright   MPI_Datatype       unit;
554a175e481SJames Wright   Vec                rhs_loc, rhs;
555a175e481SJames Wright   PetscMemType       rhs_mem_type;
556a175e481SJames Wright   CeedScalar        *rhs_arr;
557a175e481SJames Wright   CeedMemType        ceed_mem_type;
558a175e481SJames Wright   PetscFunctionBeginUser;
559a175e481SJames Wright 
5606dcea3beSJames Wright   PetscLogStage stage_stats_process;
5616dcea3beSJames Wright   PetscCall(PetscLogStageGetId("Stats Process", &stage_stats_process));
5626dcea3beSJames Wright   if (stage_stats_process == -1) PetscCall(PetscLogStageRegister("Stats Process", &stage_stats_process));
5636dcea3beSJames Wright   PetscCall(PetscLogStagePush(stage_stats_process));
5646dcea3beSJames Wright 
565a175e481SJames Wright   CeedGetPreferredMemType(user->ceed, &ceed_mem_type);
566a175e481SJames Wright   CeedVectorSetValue(user_stats.parent_stats, 0);
567a175e481SJames Wright 
568a175e481SJames Wright   CeedVectorGetArrayRead(user_stats.child_stats, ceed_mem_type, &child_stats);
569a175e481SJames Wright   CeedVectorGetArray(user_stats.parent_stats, ceed_mem_type, &parent_stats);
570a175e481SJames Wright 
571a175e481SJames Wright   if (user_stats.num_comp_stats == 1) unit = MPIU_REAL;
572a175e481SJames Wright   else {
573a175e481SJames Wright     PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit));
574a175e481SJames Wright     PetscCallMPI(MPI_Type_commit(&unit));
575a175e481SJames Wright   }
576a175e481SJames Wright 
577a175e481SJames Wright   PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
578a175e481SJames Wright   PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
579a175e481SJames Wright 
580a175e481SJames Wright   CeedVectorRestoreArrayRead(user_stats.child_stats, &child_stats);
581a175e481SJames Wright   CeedVectorRestoreArray(user_stats.parent_stats, &parent_stats);
582a175e481SJames Wright   PetscCallMPI(MPI_Type_free(&unit));
583a175e481SJames Wright 
58478bbfb6fSJed Brown   PetscReal solution_time;
58578bbfb6fSJed Brown   PetscCall(DMGetOutputSequenceNumber(user_stats.dm, NULL, &solution_time));
58678bbfb6fSJed Brown   PetscReal summing_duration = solution_time - user->app_ctx->cont_time;
58778bbfb6fSJed Brown   CeedVectorScale(user_stats.parent_stats, 1 / (summing_duration * user_stats.span_width));
588a175e481SJames Wright 
589a175e481SJames Wright   // L^2 projection with the parent_data
590a175e481SJames Wright   PetscCall(DMGetLocalVector(user_stats.dm, &rhs_loc));
591a175e481SJames Wright   PetscCall(VecZeroEntries(rhs_loc));
592a175e481SJames Wright   PetscCall(VecGetArrayWriteAndMemType(rhs_loc, &rhs_arr, &rhs_mem_type));
593a175e481SJames Wright   CeedVectorSetArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), CEED_USE_POINTER, (PetscScalar *)rhs_arr);
594a175e481SJames Wright 
595a175e481SJames Wright   CeedOperatorApply(user_stats.op_stats_proj, user_stats.parent_stats, user_stats.rhs_ceed, CEED_REQUEST_IMMEDIATE);
596a175e481SJames Wright 
597a175e481SJames Wright   CeedVectorTakeArray(user_stats.rhs_ceed, MemTypeP2C(rhs_mem_type), &rhs_arr);
598a175e481SJames Wright   PetscCall(VecRestoreArrayAndMemType(rhs_loc, &rhs_arr));
59978bbfb6fSJed Brown 
60078bbfb6fSJed Brown   PetscCall(DMGetGlobalVector(user_stats.dm, &rhs));
60178bbfb6fSJed Brown   PetscCall(VecZeroEntries(rhs));
602a175e481SJames Wright   PetscCall(DMLocalToGlobal(user_stats.dm, rhs_loc, ADD_VALUES, rhs));
60378bbfb6fSJed Brown   PetscCall(DMRestoreLocalVector(user_stats.dm, &rhs_loc));
604a175e481SJames Wright 
60578bbfb6fSJed Brown   PetscCall(KSPSolve(user_stats.ksp, rhs, stats));
606a175e481SJames Wright 
60778bbfb6fSJed Brown   PetscCall(DMRestoreGlobalVector(user_stats.dm, &rhs));
6086dcea3beSJames Wright   PetscCall(PetscLogStagePop());
609a175e481SJames Wright   PetscFunctionReturn(0);
610a175e481SJames Wright }
611a175e481SJames Wright 
612a175e481SJames Wright // TSMonitor for the statistics collection and processing
613a175e481SJames Wright PetscErrorCode TSMonitor_Statistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) {
614a175e481SJames Wright   User              user = (User)ctx;
615a175e481SJames Wright   Vec               stats;
61678bbfb6fSJed Brown   TSConvergedReason reason;
617b7d66439SJames Wright   PetscInt collect_interval = user->app_ctx->turb_spanstats_collect_interval, viewer_interval = user->app_ctx->turb_spanstats_viewer_interval;
618a175e481SJames Wright   PetscFunctionBeginUser;
61978bbfb6fSJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
620a175e481SJames Wright   // Do not collect or process on the first step of the run (ie. on the initial condition)
62178bbfb6fSJed Brown   if (steps == user->app_ctx->cont_steps && reason == TS_CONVERGED_ITERATING) PetscFunctionReturn(0);
622a175e481SJames Wright 
62327240365SJames Wright   PetscBool run_processing_and_viewer = (steps % viewer_interval == 0 && viewer_interval != -1) || reason != TS_CONVERGED_ITERATING;
624a175e481SJames Wright 
62527240365SJames Wright   if (steps % collect_interval == 0 || run_processing_and_viewer) {
62627240365SJames Wright     PetscCall(CollectStatistics(user, solution_time, Q));
62727240365SJames Wright 
62827240365SJames Wright     if (run_processing_and_viewer) {
6298ed52730SJames Wright       PetscCall(DMSetOutputSequenceNumber(user->spanstats.dm, steps, solution_time));
63078bbfb6fSJed Brown       PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats));
63178bbfb6fSJed Brown       PetscCall(ProcessStatistics(user, stats));
6328fb33541SJames Wright       if (user->app_ctx->test_type == TESTTYPE_NONE) {
633b7d66439SJames Wright         PetscCall(PetscViewerPushFormat(user->app_ctx->turb_spanstats_viewer, user->app_ctx->turb_spanstats_viewer_format));
634b7d66439SJames Wright         PetscCall(VecView(stats, user->app_ctx->turb_spanstats_viewer));
635b7d66439SJames Wright         PetscCall(PetscViewerPopFormat(user->app_ctx->turb_spanstats_viewer));
6368fb33541SJames Wright       }
6378fb33541SJames Wright       if (user->app_ctx->test_type == TESTTYPE_TURB_SPANSTATS && reason != TS_CONVERGED_ITERATING) {
6388fb33541SJames Wright         PetscCall(RegressionTests_NS(user->app_ctx, stats));
6398fb33541SJames Wright       }
640b7d66439SJames Wright       if (user->spanstats.do_mms_test && reason != TS_CONVERGED_ITERATING) {
641823a1283SJames Wright         Vec error;
642823a1283SJames Wright         PetscCall(VecDuplicate(stats, &error));
643b7d66439SJames Wright         PetscCall(ApplyLocal_Ceed(stats, error, user->spanstats.mms_error_ctx));
644823a1283SJames Wright         PetscScalar error_sq = 0;
645823a1283SJames Wright         PetscCall(VecSum(error, &error_sq));
646823a1283SJames Wright         PetscScalar l2_error = sqrt(error_sq);
647823a1283SJames Wright         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error));
648823a1283SJames Wright       }
649a175e481SJames Wright       PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats));
650522ee345SJames Wright     }
65127240365SJames Wright   }
65251ee423eSJames Wright   PetscFunctionReturn(0);
65351ee423eSJames Wright }
654fd170fd0SJames Wright 
6556665b873SJames Wright PetscErrorCode DestroyStats(User user, CeedData ceed_data) {
656fd170fd0SJames Wright   PetscFunctionBeginUser;
657fd170fd0SJames Wright 
658fd170fd0SJames Wright   // -- CeedVectors
659fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.child_stats);
660fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.parent_stats);
661fd170fd0SJames Wright   CeedVectorDestroy(&user->spanstats.rhs_ceed);
662fd170fd0SJames Wright 
663fd170fd0SJames Wright   // -- CeedOperators
664fd170fd0SJames Wright   CeedOperatorDestroy(&user->spanstats.op_stats_collect);
665fd170fd0SJames Wright   CeedOperatorDestroy(&user->spanstats.op_stats_proj);
6666665b873SJames Wright   PetscCall(MatopApplyContextDestroy(user->spanstats.mms_error_ctx));
667fd170fd0SJames Wright 
668fd170fd0SJames Wright   // -- Vec
669fd170fd0SJames Wright   PetscCall(VecDestroy(&user->spanstats.M_inv));
670fd170fd0SJames Wright 
671fd170fd0SJames Wright   // -- KSP
672fd170fd0SJames Wright   PetscCall(KSPDestroy(&user->spanstats.ksp));
673fd170fd0SJames Wright 
674fd170fd0SJames Wright   // -- SF
675fd170fd0SJames Wright   PetscCall(PetscSFDestroy(&user->spanstats.sf));
676fd170fd0SJames Wright 
677fd170fd0SJames Wright   // -- DM
678fd170fd0SJames Wright   PetscCall(DMDestroy(&user->spanstats.dm));
679fd170fd0SJames Wright 
680fd170fd0SJames Wright   PetscFunctionReturn(0);
681fd170fd0SJames Wright }
682