xref: /libCEED/examples/fluids/src/turb_spanstats.c (revision 4aa5fbb4c55e771431ba815a9dddd623e4e3362a)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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
8701bc832SJames Wright /// Functions for setting up and performing spanwise-statistics collection
9701bc832SJames Wright ///
10701bc832SJames Wright /// "Parent" refers to the 2D plane on which statistics are collected *onto*.
11701bc832SJames Wright /// "Child" refers to the 3D domain where statistics are gathered *from*.
12701bc832SJames Wright /// Each quadrature point on the parent plane has several children in the child domain that it performs spanwise averaging with.
13a175e481SJames Wright 
14a175e481SJames Wright #include "../qfunctions/turb_spanstats.h"
1551ee423eSJames Wright 
1649aac155SJeremy L Thompson #include <ceed.h>
1749aac155SJeremy L Thompson #include <petscdmplex.h>
1849aac155SJeremy L Thompson #include <petscsf.h>
1949aac155SJeremy L Thompson 
204021610dSJames Wright #include "../include/petsc_ops.h"
2151ee423eSJames Wright #include "../navierstokes.h"
2251ee423eSJames Wright 
23269a910fSJames Wright typedef struct {
24cfb075a4SJames Wright   CeedElemRestriction elem_restr_parent_x, elem_restr_parent_stats, elem_restr_parent_colloc, elem_restr_child_colloc;
25269a910fSJames Wright   CeedBasis           basis_x, basis_stats;
26cfb075a4SJames Wright   CeedVector          x_coord;
27269a910fSJames Wright } *SpanStatsSetupData;
28269a910fSJames Wright 
29731c13d7SJames Wright PetscErrorCode CreateStatsDM(User user, ProblemData problem, PetscInt degree) {
303a4208e6SJames Wright   user->spanstats.num_comp_stats = TURB_NUM_COMPONENTS;
31a175e481SJames Wright   PetscReal     domain_min[3], domain_max[3];
3278837792SJames Wright   PetscSection  section;
334eed8630SJames Wright   PetscLogStage stage_stats_setup;
341f595ac1SJames Wright   MPI_Comm      comm = PetscObjectComm((PetscObject)user->dm);
3551ee423eSJames Wright 
36f17d818dSJames Wright   PetscFunctionBeginUser;
374eed8630SJames Wright   PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup));
384eed8630SJames Wright   if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup));
394eed8630SJames Wright   PetscCall(PetscLogStagePush(stage_stats_setup));
404eed8630SJames Wright 
41a175e481SJames Wright   // Get spanwise length
42a175e481SJames Wright   PetscCall(DMGetBoundingBox(user->dm, domain_min, domain_max));
43*4aa5fbb4SJames Wright   user->spanstats.span_width = domain_max[2] - domain_min[2];
44a175e481SJames Wright 
45f967ad79SJames Wright   {  // Get DM from surface
461f595ac1SJames Wright     DM             parent_distributed_dm;
47fd9d666bSJames Wright     const PetscSF *isoperiodicface;
48fd9d666bSJames Wright     PetscInt       num_isoperiodicface;
4951ee423eSJames Wright     DMLabel        label;
501f595ac1SJames Wright     PetscMPIInt    size;
51c9198418SJames Wright 
52fd9d666bSJames Wright     PetscCall(DMPlexGetIsoperiodicFaceSF(user->dm, &num_isoperiodicface, &isoperiodicface));
53c9198418SJames Wright 
54c9198418SJames Wright     if (isoperiodicface) {
55c9198418SJames Wright       PetscSF         inv_isoperiodicface;
56fd9d666bSJames Wright       PetscInt        nleaves, isoperiodicface_index = -1;
57c9198418SJames Wright       const PetscInt *ilocal;
58c9198418SJames Wright 
59fd9d666bSJames Wright       PetscCall(PetscOptionsGetInt(NULL, NULL, "-ts_monitor_turbulence_spanstats_isoperiodic_facesf", &isoperiodicface_index, NULL));
60fd9d666bSJames Wright       isoperiodicface_index = isoperiodicface_index == -1 ? num_isoperiodicface - 1 : isoperiodicface_index;
61fd9d666bSJames Wright       PetscCall(PetscSFCreateInverseSF(isoperiodicface[isoperiodicface_index], &inv_isoperiodicface));
62c9198418SJames Wright       PetscCall(PetscSFGetGraph(inv_isoperiodicface, NULL, &nleaves, &ilocal, NULL));
63c9198418SJames Wright       PetscCall(DMCreateLabel(user->dm, "Periodic Face"));
64c9198418SJames Wright       PetscCall(DMGetLabel(user->dm, "Periodic Face", &label));
65c9198418SJames Wright       for (PetscInt i = 0; i < nleaves; i++) {
66c9198418SJames Wright         PetscCall(DMLabelSetValue(label, ilocal[i], 1));
67c9198418SJames Wright       }
68c9198418SJames Wright     } else {
6951ee423eSJames Wright       PetscCall(DMGetLabel(user->dm, "Face Sets", &label));
70c9198418SJames Wright     }
71c9198418SJames Wright 
7251ee423eSJames Wright     PetscCall(DMPlexLabelComplete(user->dm, label));
73342228c7SAlex Lindsay     PetscCall(DMPlexFilter(user->dm, label, 1, PETSC_FALSE, PETSC_FALSE, NULL, &user->spanstats.dm));
7449a40c8aSKenneth E. Jansen     PetscCall(DMSetCoordinateDisc(user->spanstats.dm, NULL, PETSC_TRUE));  // Ensure that a coordinate FE exists
751f595ac1SJames Wright 
761f595ac1SJames Wright     PetscCall(DMPlexDistribute(user->spanstats.dm, 0, NULL, &parent_distributed_dm));
771f595ac1SJames Wright     PetscCallMPI(MPI_Comm_size(comm, &size));
781f595ac1SJames Wright     if (parent_distributed_dm) {
791f595ac1SJames Wright       PetscCall(DMDestroy(&user->spanstats.dm));
801f595ac1SJames Wright       user->spanstats.dm = parent_distributed_dm;
811f595ac1SJames Wright     } else if (size > 1) {
821f595ac1SJames Wright       PetscCall(PetscPrintf(comm, "WARNING: Turbulent spanwise statistics: parent DM could not be distributed accross %d ranks.\n", size));
831f595ac1SJames Wright     }
8451ee423eSJames Wright   }
85d68a66c4SJames Wright   {
86d68a66c4SJames Wright     PetscBool is_simplex = PETSC_FALSE;
87d68a66c4SJames Wright     PetscCall(DMPlexIsSimplex(user->spanstats.dm, &is_simplex));
88d68a66c4SJames Wright     PetscCheck(is_simplex != PETSC_TRUE, comm, PETSC_ERR_ARG_WRONGSTATE, "Spanwise Statistics is not implemented for non-tensor product grids");
89d68a66c4SJames Wright   }
9051ee423eSJames Wright 
9151ee423eSJames Wright   PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats"));
92b7d66439SJames Wright   PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "turbulence_spanstats_"));
9351ee423eSJames Wright   PetscCall(DMSetFromOptions(user->spanstats.dm));
94b57b8e72SJames Wright   PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view"));
95c9198418SJames Wright 
9678837792SJames Wright   // Create FE space for parent DM
97d68a66c4SJames Wright   PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, user->app_ctx->degree, 1, user->app_ctx->q_extra, 1, &user->spanstats.num_comp_stats,
98d68a66c4SJames Wright                                user->spanstats.dm));
9951ee423eSJames Wright 
10078837792SJames Wright   // Create Section for data
10151ee423eSJames Wright   PetscCall(DMGetLocalSection(user->spanstats.dm, &section));
10251ee423eSJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
1033a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY, "MeanDensity"));
1043a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE, "MeanPressure"));
1053a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_SQUARED, "MeanPressureSquared"));
1063a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_X, "MeanPressureVelocityX"));
1073a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Y, "MeanPressureVelocityY"));
1083a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Z, "MeanPressureVelocityZ"));
1093a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE, "MeanDensityTemperature"));
1103a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_X, "MeanDensityTemperatureFluxX"));
1113a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Y, "MeanDensityTemperatureFluxY"));
1123a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Z, "MeanDensityTemperatureFluxZ"));
1133a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_X, "MeanMomentumX"));
1143a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Y, "MeanMomentumY"));
1153a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Z, "MeanMomentumZ"));
1163a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XX, "MeanMomentumFluxXX"));
1173a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YY, "MeanMomentumFluxYY"));
1183a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_ZZ, "MeanMomentumFluxZZ"));
1193a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YZ, "MeanMomentumFluxYZ"));
1203a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XZ, "MeanMomentumFluxXZ"));
1213a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XY, "MeanMomentumFluxXY"));
1223a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_X, "MeanVelocityX"));
1233a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Y, "MeanVelocityY"));
1243a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Z, "MeanVelocityZ"));
12519706a06SJames Wright 
1264eed8630SJames Wright   PetscCall(PetscLogStagePop());
127ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
12819706a06SJames Wright }
12919706a06SJames Wright 
130701bc832SJames Wright /** @brief Create CeedElemRestriction for collocated data in component-major order.
131701bc832SJames Wright a. Sets the strides of the restriction to component-major order
132701bc832SJames Wright  Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction.
133701bc832SJames Wright */
134701bc832SJames Wright static PetscErrorCode CreateElemRestrColloc_CompMajor(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base,
135ed331efdSJames Wright                                                       CeedElemRestriction *elem_restr_collocated) {
1361737222fSJames Wright   CeedInt num_elem_qpts, loc_num_elem;
1371737222fSJames Wright 
138f17d818dSJames Wright   PetscFunctionBeginUser;
139a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts));
140a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem));
1411737222fSJames Wright 
1421737222fSJames Wright   const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp};
143a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides,
144a424bcd0SJames Wright                                                        elem_restr_collocated));
145ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1461737222fSJames Wright }
1471737222fSJames Wright 
148a175e481SJames Wright // Get coordinates of quadrature points
149ed331efdSJames Wright PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, Vec *Qx_coords,
15019706a06SJames Wright                                    PetscInt *total_nqpnts) {
151269a910fSJames Wright   CeedElemRestriction  elem_restr_qx;
15219706a06SJames Wright   CeedQFunction        qf_quad_coords;
15319706a06SJames Wright   CeedOperator         op_quad_coords;
154701bc832SJames Wright   CeedInt              num_comp_x;
155701bc832SJames Wright   CeedSize             l_vec_size;
156ed331efdSJames Wright   OperatorApplyContext op_quad_coords_ctx;
15719706a06SJames Wright 
158ed331efdSJames Wright   PetscFunctionBeginUser;
159a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
160701bc832SJames Wright   PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx));
161701bc832SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetLVectorSize(elem_restr_qx, &l_vec_size));
162701bc832SJames Wright   *total_nqpnts = l_vec_size / num_comp_x;
16319706a06SJames Wright 
16419706a06SJames Wright   // Create QFunction
165a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords));
16619706a06SJames Wright 
16719706a06SJames Wright   // Create Operator
168a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords));
169a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, x_coords));
170356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
17119706a06SJames Wright 
1722788647cSJames Wright   PetscCall(CeedOperatorCreateLocalVecs(op_quad_coords, DMReturnVecType(dm), PETSC_COMM_SELF, NULL, Qx_coords));
173ed331efdSJames Wright   PetscCall(OperatorApplyContextCreate(NULL, NULL, ceed, op_quad_coords, CEED_VECTOR_NONE, NULL, NULL, NULL, &op_quad_coords_ctx));
17419706a06SJames Wright 
175ed331efdSJames Wright   PetscCall(ApplyCeedOperatorLocalToLocal(NULL, *Qx_coords, op_quad_coords_ctx));
176ed331efdSJames Wright 
177ed331efdSJames Wright   PetscCall(OperatorApplyContextDestroy(op_quad_coords_ctx));
178a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qx));
179a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_quad_coords));
180a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_quad_coords));
181ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
18219706a06SJames Wright }
18319706a06SJames Wright 
184731c13d7SJames Wright PetscErrorCode SpanStatsSetupDataCreate(Ceed ceed, User user, CeedData ceed_data, ProblemData problem, SpanStatsSetupData *stats_data) {
185269a910fSJames Wright   DM       dm = user->spanstats.dm;
186bb85d312SJames Wright   CeedInt  num_comp_x, num_comp_stats = user->spanstats.num_comp_stats;
187269a910fSJames Wright   Vec      X_loc;
188bb85d312SJames Wright   DMLabel  domain_label = NULL;
189bb85d312SJames Wright   PetscInt label_value = 0, height = 0, dm_field = 0;
190269a910fSJames Wright 
191f17d818dSJames Wright   PetscFunctionBeginUser;
192269a910fSJames Wright   PetscCall(PetscNew(stats_data));
193269a910fSJames Wright 
194bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &(*stats_data)->elem_restr_parent_stats));
195bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &(*stats_data)->elem_restr_parent_x));
196a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents((*stats_data)->elem_restr_parent_x, &num_comp_x));
197a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_x, &(*stats_data)->x_coord, NULL));
198269a910fSJames Wright 
1990814d5a7SKenneth E. Jansen   {
2000814d5a7SKenneth E. Jansen     DM dm_coord;
2010814d5a7SKenneth E. Jansen     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
202bb85d312SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &(*stats_data)->basis_x));
203bb85d312SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &(*stats_data)->basis_stats));
2040814d5a7SKenneth E. Jansen   }
205269a910fSJames Wright 
206701bc832SJames Wright   PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_stats, (*stats_data)->basis_stats, (*stats_data)->elem_restr_parent_stats,
207ed331efdSJames Wright                                             &(*stats_data)->elem_restr_parent_colloc));
208701bc832SJames Wright   PetscCall(
209701bc832SJames Wright       CreateElemRestrColloc_CompMajor(ceed, num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q, &(*stats_data)->elem_restr_child_colloc));
210269a910fSJames Wright 
211269a910fSJames Wright   {  // -- Copy DM coordinates into CeedVector
212269a910fSJames Wright     DM cdm;
213269a910fSJames Wright     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
214269a910fSJames Wright     if (cdm) {
215269a910fSJames Wright       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
216269a910fSJames Wright     } else {
217269a910fSJames Wright       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
218269a910fSJames Wright     }
219269a910fSJames Wright   }
220269a910fSJames Wright   PetscCall(VecScale(X_loc, problem->dm_scale));
221d0593705SJames Wright   PetscCall(VecCopyPetscToCeed(X_loc, (*stats_data)->x_coord));
222ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
223269a910fSJames Wright }
224269a910fSJames Wright 
225269a910fSJames Wright PetscErrorCode SpanStatsSetupDataDestroy(SpanStatsSetupData data) {
226a424bcd0SJames Wright   Ceed ceed;
227a424bcd0SJames Wright 
228269a910fSJames Wright   PetscFunctionBeginUser;
229a424bcd0SJames Wright   PetscCall(CeedElemRestrictionGetCeed(data->elem_restr_parent_x, &ceed));
230a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_x));
231a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_stats));
232a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_colloc));
233a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_child_colloc));
234269a910fSJames Wright 
235a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_x));
236a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_stats));
237269a910fSJames Wright 
238a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&data->x_coord));
239269a910fSJames Wright 
240269a910fSJames Wright   PetscCall(PetscFree(data));
241ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
242269a910fSJames Wright }
243269a910fSJames Wright 
244a175e481SJames Wright // Create PetscSF for child-to-parent communication
245269a910fSJames Wright PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, SpanStatsSetupData stats_data, DM parentdm, DM childdm, PetscSF *statssf) {
246457e73b2SJames Wright   PetscInt child_num_qpnts, parent_num_qpnts;
247457e73b2SJames Wright   CeedInt  num_comp_x;
248ed331efdSJames Wright   Vec      Child_qx_coords, Parent_qx_coords;
24919706a06SJames Wright 
250ed331efdSJames Wright   PetscFunctionBeginUser;
251269a910fSJames Wright   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)childdm), statssf));
252269a910fSJames Wright 
25319706a06SJames Wright   // Assume that child and parent have the same number of components
254a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x));
25519706a06SJames Wright   const PetscInt num_comp_sf = num_comp_x - 1;  // Number of coord components used in the creation of the SF
25619706a06SJames Wright 
257ed331efdSJames Wright   // Get quad_coords for child and parent DM
258ed331efdSJames Wright   PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &Child_qx_coords, &child_num_qpnts));
259ed331efdSJames Wright   PetscCall(GetQuadratureCoords(ceed, parentdm, stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord, &Parent_qx_coords,
260269a910fSJames Wright                                 &parent_num_qpnts));
26119706a06SJames Wright 
262ed331efdSJames Wright   {  // Remove z component of coordinates for matching
26319706a06SJames Wright     const PetscReal *child_quad_coords, *parent_quad_coords;
264ed331efdSJames Wright     PetscReal       *child_coords, *parent_coords;
26519706a06SJames Wright 
266ed331efdSJames Wright     PetscCall(VecGetArrayRead(Child_qx_coords, &child_quad_coords));
267ed331efdSJames Wright     PetscCall(VecGetArrayRead(Parent_qx_coords, &parent_quad_coords));
26819706a06SJames Wright 
26919706a06SJames Wright     PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords));
27019706a06SJames Wright     for (int i = 0; i < child_num_qpnts; i++) {
27119706a06SJames Wright       child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x];
27219706a06SJames Wright       child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x];
27319706a06SJames Wright     }
27419706a06SJames Wright     for (int i = 0; i < parent_num_qpnts; i++) {
27519706a06SJames Wright       parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x];
27619706a06SJames Wright       parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x];
27719706a06SJames Wright     }
278ed331efdSJames Wright     PetscCall(VecRestoreArrayRead(Child_qx_coords, &child_quad_coords));
279ed331efdSJames Wright     PetscCall(VecRestoreArrayRead(Parent_qx_coords, &parent_quad_coords));
28019706a06SJames Wright 
281269a910fSJames Wright     PetscCall(PetscSFSetGraphFromCoordinates(*statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords));
282ed331efdSJames Wright     PetscCall(PetscFree2(child_coords, parent_coords));
283ed331efdSJames Wright   }
28419706a06SJames Wright 
285269a910fSJames Wright   PetscCall(PetscSFViewFromOptions(*statssf, NULL, "-spanstats_sf_view"));
286a175e481SJames Wright 
287ed331efdSJames Wright   PetscCall(VecDestroy(&Child_qx_coords));
288ed331efdSJames Wright   PetscCall(VecDestroy(&Parent_qx_coords));
289ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
29019706a06SJames Wright }
29119706a06SJames Wright 
292ed331efdSJames Wright // @brief Setup RHS and LHS for L^2 projection of statistics
293269a910fSJames Wright PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) {
294cfb075a4SJames Wright   CeedOperator        op_mass, op_proj_rhs;
295269a910fSJames Wright   CeedQFunction       qf_mass, qf_stats_proj;
296269a910fSJames Wright   CeedInt             q_data_size, num_comp_stats = user->spanstats.num_comp_stats;
297cfb075a4SJames Wright   CeedElemRestriction elem_restr_qd;
298cfb075a4SJames Wright   CeedVector          q_data;
299cfb075a4SJames Wright   DMLabel             domain_label = NULL;
300cfb075a4SJames Wright   PetscInt            label_value  = 0;
301a175e481SJames Wright 
302ed331efdSJames Wright   PetscFunctionBeginUser;
303ed331efdSJames Wright   // -- Create Operator for RHS of L^2 projection of statistics
304269a910fSJames Wright   // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function.
305269a910fSJames Wright   // Therefore, an Identity QF is sufficient
306a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &qf_stats_proj));
307269a910fSJames Wright 
308a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_proj, NULL, NULL, &op_proj_rhs));
309356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "input", stats_data->elem_restr_parent_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
310a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "output", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
311269a910fSJames Wright 
312ed331efdSJames Wright   PetscCall(OperatorApplyContextCreate(NULL, user->spanstats.dm, ceed, op_proj_rhs, NULL, NULL, NULL, NULL, &user->spanstats.op_proj_rhs_ctx));
3132788647cSJames Wright   PetscCall(CeedOperatorCreateLocalVecs(op_proj_rhs, DMReturnVecType(user->spanstats.dm), PETSC_COMM_SELF, &user->spanstats.Parent_Stats_loc, NULL));
314cfb075a4SJames Wright   PetscCall(QDataGet(ceed, user->spanstats.dm, domain_label, label_value, stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord,
315cfb075a4SJames Wright                      &elem_restr_qd, &q_data, &q_data_size));
316a175e481SJames Wright 
317a175e481SJames Wright   // Create Mass CeedOperator
318269a910fSJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_stats, q_data_size, &qf_mass));
319a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
320a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
321cfb075a4SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
322a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
323a175e481SJames Wright 
324f967ad79SJames Wright   {  // Setup KSP for L^2 projection
3257f2a9303SJames Wright     Mat mat_mass;
326a175e481SJames Wright     KSP ksp;
327a175e481SJames Wright 
3287f2a9303SJames Wright     PetscCall(MatCeedCreate(user->spanstats.dm, user->spanstats.dm, op_mass, NULL, &mat_mass));
329a175e481SJames Wright 
3302788647cSJames Wright     PetscCall(KSPCreate(PetscObjectComm((PetscObject)user->spanstats.dm), &ksp));
331b7d66439SJames Wright     PetscCall(KSPSetOptionsPrefix(ksp, "turbulence_spanstats_"));
332a175e481SJames Wright     {
333a175e481SJames Wright       PC pc;
334a175e481SJames Wright       PetscCall(KSPGetPC(ksp, &pc));
335a175e481SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
336a175e481SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
337a175e481SJames Wright       PetscCall(KSPSetType(ksp, KSPCG));
338a175e481SJames Wright       PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
339a175e481SJames Wright       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
340a175e481SJames Wright     }
3417f2a9303SJames Wright     PetscCall(KSPSetFromOptions_WithMatCeed(ksp, mat_mass));
342a175e481SJames Wright     user->spanstats.ksp = ksp;
343cde30410SJames Wright     PetscCall(MatDestroy(&mat_mass));
344a175e481SJames Wright   }
345a175e481SJames Wright 
346a175e481SJames Wright   // Cleanup
347cfb075a4SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
348cfb075a4SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
349a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
350a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_proj));
351a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
352a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_proj_rhs));
353ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
354a175e481SJames Wright }
355a175e481SJames Wright 
356269a910fSJames Wright // Create CeedOperator for statistics collection
357731c13d7SJames Wright PetscErrorCode CreateStatisticCollectionOperator(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data, ProblemData problem) {
358a175e481SJames Wright   CeedInt                     num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q;
359495b9769SJames Wright   Turbulence_SpanStatsContext collect_ctx;
360495b9769SJames Wright   NewtonianIdealGasContext    newtonian_ig_ctx;
361495b9769SJames Wright   CeedQFunctionContext        collect_context;
362269a910fSJames Wright   CeedQFunction               qf_stats_collect;
363ed331efdSJames Wright   CeedOperator                op_stats_collect;
364ed331efdSJames Wright 
365a175e481SJames Wright   PetscFunctionBeginUser;
366a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q));
367a175e481SJames Wright 
368a175e481SJames Wright   // Create Operator for statistics collection
369a175e481SJames Wright   switch (user->phys->state_var) {
370a175e481SJames Wright     case STATEVAR_PRIMITIVE:
371a424bcd0SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &qf_stats_collect));
372a175e481SJames Wright       break;
373a175e481SJames Wright     case STATEVAR_CONSERVATIVE:
374a424bcd0SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &qf_stats_collect));
375a175e481SJames Wright       break;
376269a910fSJames Wright     default:
377269a910fSJames Wright       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No statisics collection available for chosen state variable");
378a175e481SJames Wright   }
379a175e481SJames Wright 
380b7d66439SJames Wright   if (user->spanstats.do_mms_test) {
381a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect));
382a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest, ChildStatsCollectionMMSTest_loc, &qf_stats_collect));
383a175e481SJames Wright   }
384a175e481SJames Wright 
385269a910fSJames Wright   {  // Setup Collection Context
386495b9769SJames Wright     PetscCall(PetscNew(&collect_ctx));
387a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx));
388495b9769SJames Wright     collect_ctx->gas = *newtonian_ig_ctx;
389495b9769SJames Wright 
390a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &collect_context));
391a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextSetData(collect_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx));
392a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(collect_context, CEED_MEM_HOST, FreeContextPetsc));
393495b9769SJames Wright 
394a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(collect_context, "solution time",
395a424bcd0SJames Wright                                                            offsetof(struct Turbulence_SpanStatsContext_, solution_time), 1, "Current solution time"));
396a424bcd0SJames Wright     PetscCallCeed(
397a424bcd0SJames Wright         ceed, CeedQFunctionContextRegisterDouble(collect_context, "previous time", offsetof(struct Turbulence_SpanStatsContext_, previous_time), 1,
398a424bcd0SJames Wright                                                  "Previous time statistics collection was done"));
399495b9769SJames Wright 
400a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx));
401495b9769SJames Wright   }
402495b9769SJames Wright 
403a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_stats_collect, collect_context));
404a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&collect_context));
405a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP));
406a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE));
407a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP));
408a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE));
409a175e481SJames Wright 
410a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_collect, NULL, NULL, &op_stats_collect));
411a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
412356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
413a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
414356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "v", stats_data->elem_restr_child_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
415a175e481SJames Wright 
416a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "solution time", &user->spanstats.solution_time_label));
417a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "previous time", &user->spanstats.previous_time_label));
418ed331efdSJames Wright 
419ed331efdSJames Wright   PetscCall(OperatorApplyContextCreate(user->dm, user->spanstats.dm, user->ceed, op_stats_collect, user->q_ceed, NULL, NULL, NULL,
420ed331efdSJames Wright                                        &user->spanstats.op_stats_collect_ctx));
421ed331efdSJames Wright 
4222788647cSJames Wright   PetscCall(
4232788647cSJames Wright       CeedOperatorCreateLocalVecs(op_stats_collect, DMReturnVecType(user->spanstats.dm), PETSC_COMM_SELF, NULL, &user->spanstats.Child_Stats_loc));
424ed331efdSJames Wright   PetscCall(VecZeroEntries(user->spanstats.Child_Stats_loc));
425495b9769SJames Wright 
426a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect));
427a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_stats_collect));
428ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
429a175e481SJames Wright }
430a175e481SJames Wright 
431b7d66439SJames Wright // Creates operator for calculating error of method of manufactured solution (MMS) test
432269a910fSJames Wright PetscErrorCode SetupMMSErrorChecking(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) {
433269a910fSJames Wright   CeedInt             num_comp_stats = user->spanstats.num_comp_stats, num_comp_x, q_data_size;
434823a1283SJames Wright   CeedQFunction       qf_error;
435823a1283SJames Wright   CeedOperator        op_error;
436823a1283SJames Wright   CeedVector          x_ceed, y_ceed;
437cfb075a4SJames Wright   DMLabel             domain_label = NULL;
438cfb075a4SJames Wright   PetscInt            label_value  = 0;
439cfb075a4SJames Wright   CeedVector          q_data;
440cfb075a4SJames Wright   CeedElemRestriction elem_restr_parent_qd;
441823a1283SJames Wright 
442f17d818dSJames Wright   PetscFunctionBeginUser;
443cfb075a4SJames Wright   PetscCall(QDataGet(ceed, user->spanstats.dm, domain_label, label_value, stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord,
444cfb075a4SJames Wright                      &elem_restr_parent_qd, &q_data, &q_data_size));
445a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(stats_data->basis_x, &num_comp_x));
446823a1283SJames Wright 
447a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest_Error, ChildStatsCollectionMMSTest_Error_loc, &qf_error));
448a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP));
449a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE));
450a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP));
451a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP));
452823a1283SJames Wright 
453a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error));
454a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "q", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
455cfb075a4SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "qdata", elem_restr_parent_qd, CEED_BASIS_NONE, q_data));
456a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "x", stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord));
457a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
458823a1283SJames Wright 
459a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &x_ceed, NULL));
460a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &y_ceed, NULL));
4614021610dSJames Wright   PetscCall(OperatorApplyContextCreate(user->spanstats.dm, user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, NULL,
462ea168fcfSJames Wright                                        &user->spanstats.mms_error_ctx));
463823a1283SJames Wright 
464cfb075a4SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
465a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&x_ceed));
466a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&y_ceed));
467cfb075a4SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_parent_qd));
468cfb075a4SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_error));
469cfb075a4SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_error));
470ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
471823a1283SJames Wright }
472823a1283SJames Wright 
473a175e481SJames Wright // Setup for statistics collection
474731c13d7SJames Wright PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) {
475269a910fSJames Wright   SpanStatsSetupData stats_data;
4764eed8630SJames Wright   PetscLogStage      stage_stats_setup;
477f5452247SJames Wright 
47819706a06SJames Wright   PetscFunctionBeginUser;
4794eed8630SJames Wright   PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup));
4804eed8630SJames Wright   if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup));
4814eed8630SJames Wright   PetscCall(PetscLogStagePush(stage_stats_setup));
4824eed8630SJames Wright 
483f5452247SJames Wright   // Create parent DM
484f5452247SJames Wright   PetscCall(CreateStatsDM(user, problem, user->app_ctx->degree));
485f5452247SJames Wright 
486269a910fSJames Wright   // Create necessary CeedObjects for setting up statistics
487269a910fSJames Wright   PetscCall(SpanStatsSetupDataCreate(ceed, user, ceed_data, problem, &stats_data));
4881737222fSJames Wright 
48919706a06SJames Wright   // Create SF for communicating child data back their respective parents
490269a910fSJames Wright   PetscCall(CreateStatsSF(ceed, ceed_data, stats_data, user->dm, user->spanstats.dm, &user->spanstats.sf));
49151ee423eSJames Wright 
492a175e481SJames Wright   // Create CeedOperators for statistics collection
493269a910fSJames Wright   PetscCall(CreateStatisticCollectionOperator(ceed, user, ceed_data, stats_data, problem));
494a175e481SJames Wright 
495a175e481SJames Wright   // Setup KSP and Mat for L^2 projection of statistics
496269a910fSJames Wright   PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data, stats_data));
497a175e481SJames Wright 
498b7d66439SJames Wright   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_turbulence_spanstats_mms", &user->spanstats.do_mms_test, NULL));
499b7d66439SJames Wright   if (user->spanstats.do_mms_test) {
500269a910fSJames Wright     PetscCall(SetupMMSErrorChecking(ceed, user, ceed_data, stats_data));
501823a1283SJames Wright   }
502823a1283SJames Wright 
503f967ad79SJames Wright   {  // Setup stats viewer with prefix
5048ed52730SJames Wright     PetscViewerType viewer_type;
505b7d66439SJames Wright     PetscCall(PetscViewerGetType(user->app_ctx->turb_spanstats_viewer, &viewer_type));
506b7d66439SJames Wright     PetscCall(PetscOptionsSetValue(NULL, "-ts_monitor_turbulence_spanstats_viewer_type", viewer_type));
5078ed52730SJames Wright 
508b7d66439SJames Wright     PetscCall(PetscViewerSetOptionsPrefix(user->app_ctx->turb_spanstats_viewer, "ts_monitor_turbulence_spanstats_"));
509b7d66439SJames Wright     PetscCall(PetscViewerSetFromOptions(user->app_ctx->turb_spanstats_viewer));
5108ed52730SJames Wright   }
5114eed8630SJames Wright 
512269a910fSJames Wright   PetscCall(SpanStatsSetupDataDestroy(stats_data));
5134eed8630SJames Wright   PetscCall(PetscLogStagePop());
514ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
515a175e481SJames Wright }
516a175e481SJames Wright 
517a175e481SJames Wright // Collect statistics based on the solution Q
518a175e481SJames Wright PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) {
519cbef7084SJames Wright   SpanStatsData user_stats = user->spanstats;
520a175e481SJames Wright 
521ed331efdSJames Wright   PetscFunctionBeginUser;
5226dcea3beSJames Wright   PetscLogStage stage_stats_collect;
5236dcea3beSJames Wright   PetscCall(PetscLogStageGetId("Stats Collect", &stage_stats_collect));
5246dcea3beSJames Wright   if (stage_stats_collect == -1) PetscCall(PetscLogStageRegister("Stats Collect", &stage_stats_collect));
5256dcea3beSJames Wright   PetscCall(PetscLogStagePush(stage_stats_collect));
5266dcea3beSJames Wright 
527a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, user->Q_loc, solution_time));
528a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.solution_time_label, &solution_time));
529a0b9a424SJames Wright   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, user->Q_loc));
530ed331efdSJames Wright   PetscCall(ApplyAddCeedOperatorLocalToLocal(user->Q_loc, user_stats.Child_Stats_loc, user_stats.op_stats_collect_ctx));
531a175e481SJames Wright 
532a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.previous_time_label, &solution_time));
533a175e481SJames Wright 
5346dcea3beSJames Wright   PetscCall(PetscLogStagePop());
535ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
536a175e481SJames Wright }
537a175e481SJames Wright 
538a175e481SJames Wright // Process the child statistics into parent statistics and project them onto stats
53978bbfb6fSJed Brown PetscErrorCode ProcessStatistics(User user, Vec stats) {
540cbef7084SJames Wright   SpanStatsData      user_stats = user->spanstats;
541a175e481SJames Wright   const PetscScalar *child_stats;
542a175e481SJames Wright   PetscScalar       *parent_stats;
543a175e481SJames Wright   MPI_Datatype       unit;
544ed331efdSJames Wright   Vec                RHS;
545a175e481SJames Wright 
546ed331efdSJames Wright   PetscFunctionBeginUser;
5476dcea3beSJames Wright   PetscLogStage stage_stats_process;
5486dcea3beSJames Wright   PetscCall(PetscLogStageGetId("Stats Process", &stage_stats_process));
5496dcea3beSJames Wright   if (stage_stats_process == -1) PetscCall(PetscLogStageRegister("Stats Process", &stage_stats_process));
5506dcea3beSJames Wright   PetscCall(PetscLogStagePush(stage_stats_process));
5516dcea3beSJames Wright 
552ed331efdSJames Wright   PetscCall(VecZeroEntries(user_stats.Parent_Stats_loc));
553a175e481SJames Wright 
554ed331efdSJames Wright   PetscCall(VecGetArrayRead(user_stats.Child_Stats_loc, &child_stats));
555ed331efdSJames Wright   PetscCall(VecGetArray(user_stats.Parent_Stats_loc, &parent_stats));
556a175e481SJames Wright 
557a175e481SJames Wright   if (user_stats.num_comp_stats == 1) unit = MPIU_REAL;
558a175e481SJames Wright   else {
559a175e481SJames Wright     PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit));
560a175e481SJames Wright     PetscCallMPI(MPI_Type_commit(&unit));
561a175e481SJames Wright   }
562a175e481SJames Wright 
563a175e481SJames Wright   PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
564a175e481SJames Wright   PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
565a175e481SJames Wright 
566ed331efdSJames Wright   PetscCall(VecRestoreArrayRead(user_stats.Child_Stats_loc, &child_stats));
567ed331efdSJames Wright   PetscCall(VecRestoreArray(user_stats.Parent_Stats_loc, &parent_stats));
568a175e481SJames Wright   PetscCallMPI(MPI_Type_free(&unit));
569a175e481SJames Wright 
57078bbfb6fSJed Brown   PetscReal solution_time;
57178bbfb6fSJed Brown   PetscCall(DMGetOutputSequenceNumber(user_stats.dm, NULL, &solution_time));
57278bbfb6fSJed Brown   PetscReal summing_duration = solution_time - user->app_ctx->cont_time;
573ed331efdSJames Wright   PetscCall(VecScale(user_stats.Parent_Stats_loc, 1 / (summing_duration * user_stats.span_width)));
574a175e481SJames Wright 
575a175e481SJames Wright   // L^2 projection with the parent_data
576ed331efdSJames Wright   PetscCall(DMGetGlobalVector(user_stats.dm, &RHS));
577ed331efdSJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(user_stats.Parent_Stats_loc, RHS, user_stats.op_proj_rhs_ctx));
578a175e481SJames Wright 
579ed331efdSJames Wright   PetscCall(KSPSolve(user_stats.ksp, RHS, stats));
580a175e481SJames Wright 
581ed331efdSJames Wright   PetscCall(DMRestoreGlobalVector(user_stats.dm, &RHS));
5826dcea3beSJames Wright   PetscCall(PetscLogStagePop());
583ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
584a175e481SJames Wright }
585a175e481SJames Wright 
586a175e481SJames Wright // TSMonitor for the statistics collection and processing
587f5452247SJames Wright PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) {
588a175e481SJames Wright   User              user = (User)ctx;
589a175e481SJames Wright   Vec               stats;
59078bbfb6fSJed Brown   TSConvergedReason reason;
591b7d66439SJames Wright   PetscInt collect_interval = user->app_ctx->turb_spanstats_collect_interval, viewer_interval = user->app_ctx->turb_spanstats_viewer_interval;
592f17d818dSJames Wright 
593a175e481SJames Wright   PetscFunctionBeginUser;
59478bbfb6fSJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
595a175e481SJames Wright   // Do not collect or process on the first step of the run (ie. on the initial condition)
5966852f6f6SJames Wright   if (steps == user->app_ctx->cont_steps) PetscFunctionReturn(PETSC_SUCCESS);
597a175e481SJames Wright 
59827240365SJames Wright   PetscBool run_processing_and_viewer = (steps % viewer_interval == 0 && viewer_interval != -1) || reason != TS_CONVERGED_ITERATING;
599a175e481SJames Wright 
60027240365SJames Wright   if (steps % collect_interval == 0 || run_processing_and_viewer) {
60127240365SJames Wright     PetscCall(CollectStatistics(user, solution_time, Q));
60227240365SJames Wright 
60327240365SJames Wright     if (run_processing_and_viewer) {
6048ed52730SJames Wright       PetscCall(DMSetOutputSequenceNumber(user->spanstats.dm, steps, solution_time));
60578bbfb6fSJed Brown       PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats));
60678bbfb6fSJed Brown       PetscCall(ProcessStatistics(user, stats));
6078fb33541SJames Wright       if (user->app_ctx->test_type == TESTTYPE_NONE) {
608b7d66439SJames Wright         PetscCall(PetscViewerPushFormat(user->app_ctx->turb_spanstats_viewer, user->app_ctx->turb_spanstats_viewer_format));
609b7d66439SJames Wright         PetscCall(VecView(stats, user->app_ctx->turb_spanstats_viewer));
610b7d66439SJames Wright         PetscCall(PetscViewerPopFormat(user->app_ctx->turb_spanstats_viewer));
6118fb33541SJames Wright       }
6128fb33541SJames Wright       if (user->app_ctx->test_type == TESTTYPE_TURB_SPANSTATS && reason != TS_CONVERGED_ITERATING) {
6133c9e7ad1SJames Wright         PetscCall(RegressionTest(user->app_ctx, stats));
6148fb33541SJames Wright       }
615b7d66439SJames Wright       if (user->spanstats.do_mms_test && reason != TS_CONVERGED_ITERATING) {
616823a1283SJames Wright         Vec error;
617823a1283SJames Wright         PetscCall(VecDuplicate(stats, &error));
6184021610dSJames Wright         PetscCall(ApplyCeedOperatorGlobalToGlobal(stats, error, user->spanstats.mms_error_ctx));
619823a1283SJames Wright         PetscScalar error_sq = 0;
620823a1283SJames Wright         PetscCall(VecSum(error, &error_sq));
621823a1283SJames Wright         PetscScalar l2_error = sqrt(error_sq);
622823a1283SJames Wright         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error));
623823a1283SJames Wright       }
624a175e481SJames Wright       PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats));
625522ee345SJames Wright     }
62627240365SJames Wright   }
627ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
62851ee423eSJames Wright }
629fd170fd0SJames Wright 
630f5452247SJames Wright PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data) {
631fd170fd0SJames Wright   PetscFunctionBeginUser;
632ed331efdSJames Wright   PetscCall(VecDestroy(&user->spanstats.Child_Stats_loc));
633ed331efdSJames Wright   PetscCall(VecDestroy(&user->spanstats.Parent_Stats_loc));
634fd170fd0SJames Wright 
635ed331efdSJames Wright   PetscCall(OperatorApplyContextDestroy(user->spanstats.op_stats_collect_ctx));
636ed331efdSJames Wright   PetscCall(OperatorApplyContextDestroy(user->spanstats.op_proj_rhs_ctx));
6374021610dSJames Wright   PetscCall(OperatorApplyContextDestroy(user->spanstats.mms_error_ctx));
638fd170fd0SJames Wright 
639fd170fd0SJames Wright   PetscCall(KSPDestroy(&user->spanstats.ksp));
640fd170fd0SJames Wright   PetscCall(PetscSFDestroy(&user->spanstats.sf));
641fd170fd0SJames Wright   PetscCall(DMDestroy(&user->spanstats.dm));
642ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
643fd170fd0SJames Wright }
644