xref: /libCEED/examples/fluids/src/turb_spanstats.c (revision 457e73b2156816fd50927125e96139f6cbacb95c)
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 
1249aac155SJeremy L Thompson #include <ceed.h>
1349aac155SJeremy L Thompson #include <petscdmplex.h>
1449aac155SJeremy L Thompson #include <petscsf.h>
1549aac155SJeremy L Thompson 
164021610dSJames Wright #include "../include/petsc_ops.h"
1751ee423eSJames Wright #include "../navierstokes.h"
1851ee423eSJames Wright 
19269a910fSJames Wright typedef struct {
20269a910fSJames Wright   CeedElemRestriction elem_restr_parent_x, elem_restr_parent_stats, elem_restr_parent_qd, elem_restr_parent_colloc, elem_restr_child_colloc;
21269a910fSJames Wright   CeedBasis           basis_x, basis_stats;
22269a910fSJames Wright   CeedVector          x_coord, q_data;
23269a910fSJames Wright } *SpanStatsSetupData;
24269a910fSJames Wright 
25f5452247SJames Wright PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree) {
263a4208e6SJames Wright   user->spanstats.num_comp_stats = TURB_NUM_COMPONENTS;
27a175e481SJames Wright   PetscReal     domain_min[3], domain_max[3];
2878837792SJames Wright   PetscFE       fe;
2978837792SJames Wright   PetscSection  section;
304eed8630SJames Wright   PetscLogStage stage_stats_setup;
311f595ac1SJames Wright   MPI_Comm      comm = PetscObjectComm((PetscObject)user->dm);
3251ee423eSJames Wright   PetscFunctionBeginUser;
3351ee423eSJames Wright 
344eed8630SJames Wright   PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup));
354eed8630SJames Wright   if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup));
364eed8630SJames Wright   PetscCall(PetscLogStagePush(stage_stats_setup));
374eed8630SJames Wright 
38a175e481SJames Wright   // Get spanwise length
39a175e481SJames Wright   PetscCall(DMGetBoundingBox(user->dm, domain_min, domain_max));
40a175e481SJames Wright   user->spanstats.span_width = domain_max[2] - domain_min[1];
41a175e481SJames Wright 
42f967ad79SJames Wright   {  // Get DM from surface
431f595ac1SJames Wright     DM          parent_distributed_dm;
44c9198418SJames Wright     PetscSF     isoperiodicface;
4551ee423eSJames Wright     DMLabel     label;
461f595ac1SJames Wright     PetscMPIInt size;
47c9198418SJames Wright 
48c9198418SJames Wright     PetscCall(DMPlexGetIsoperiodicFaceSF(user->dm, &isoperiodicface));
49c9198418SJames Wright 
50c9198418SJames Wright     if (isoperiodicface) {
51c9198418SJames Wright       PetscSF         inv_isoperiodicface;
52c9198418SJames Wright       PetscInt        nleaves;
53c9198418SJames Wright       const PetscInt *ilocal;
54c9198418SJames Wright 
55c9198418SJames Wright       PetscCall(PetscSFCreateInverseSF(isoperiodicface, &inv_isoperiodicface));
56c9198418SJames Wright       PetscCall(PetscSFGetGraph(inv_isoperiodicface, NULL, &nleaves, &ilocal, NULL));
57c9198418SJames Wright       PetscCall(DMCreateLabel(user->dm, "Periodic Face"));
58c9198418SJames Wright       PetscCall(DMGetLabel(user->dm, "Periodic Face", &label));
59c9198418SJames Wright       for (PetscInt i = 0; i < nleaves; i++) {
60c9198418SJames Wright         PetscCall(DMLabelSetValue(label, ilocal[i], 1));
61c9198418SJames Wright       }
62c9198418SJames Wright     } else {
6351ee423eSJames Wright       PetscCall(DMGetLabel(user->dm, "Face Sets", &label));
64c9198418SJames Wright     }
65c9198418SJames Wright 
6651ee423eSJames Wright     PetscCall(DMPlexLabelComplete(user->dm, label));
6751ee423eSJames Wright     PetscCall(DMPlexFilter(user->dm, label, 1, &user->spanstats.dm));
6851ee423eSJames Wright     PetscCall(DMProjectCoordinates(user->spanstats.dm, NULL));  // Ensure that a coordinate FE exists
691f595ac1SJames Wright 
701f595ac1SJames Wright     PetscCall(DMPlexDistribute(user->spanstats.dm, 0, NULL, &parent_distributed_dm));
711f595ac1SJames Wright     PetscCallMPI(MPI_Comm_size(comm, &size));
721f595ac1SJames Wright     if (parent_distributed_dm) {
731f595ac1SJames Wright       PetscCall(DMDestroy(&user->spanstats.dm));
741f595ac1SJames Wright       user->spanstats.dm = parent_distributed_dm;
751f595ac1SJames Wright     } else if (size > 1) {
761f595ac1SJames Wright       PetscCall(PetscPrintf(comm, "WARNING: Turbulent spanwise statistics: parent DM could not be distributed accross %d ranks.\n", size));
771f595ac1SJames Wright     }
7851ee423eSJames Wright   }
7951ee423eSJames Wright 
8051ee423eSJames Wright   PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats"));
81b7d66439SJames Wright   PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "turbulence_spanstats_"));
8251ee423eSJames Wright   PetscCall(DMSetFromOptions(user->spanstats.dm));
83b57b8e72SJames Wright   PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view"));
84c9198418SJames Wright 
8578837792SJames Wright   // Create FE space for parent DM
8651ee423eSJames Wright   PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, problem->dim - 1, user->spanstats.num_comp_stats, PETSC_FALSE, degree, PETSC_DECIDE, &fe));
8751ee423eSJames Wright   PetscCall(PetscObjectSetName((PetscObject)fe, "stats"));
8851ee423eSJames Wright   PetscCall(DMAddField(user->spanstats.dm, NULL, (PetscObject)fe));
8951ee423eSJames Wright   PetscCall(DMCreateDS(user->spanstats.dm));
9051ee423eSJames Wright   PetscCall(DMPlexSetClosurePermutationTensor(user->spanstats.dm, PETSC_DETERMINE, NULL));
9151ee423eSJames Wright 
9278837792SJames Wright   // Create Section for data
9351ee423eSJames Wright   PetscCall(DMGetLocalSection(user->spanstats.dm, &section));
9451ee423eSJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
953a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY, "MeanDensity"));
963a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE, "MeanPressure"));
973a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_SQUARED, "MeanPressureSquared"));
983a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_X, "MeanPressureVelocityX"));
993a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Y, "MeanPressureVelocityY"));
1003a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Z, "MeanPressureVelocityZ"));
1013a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE, "MeanDensityTemperature"));
1023a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_X, "MeanDensityTemperatureFluxX"));
1033a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Y, "MeanDensityTemperatureFluxY"));
1043a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Z, "MeanDensityTemperatureFluxZ"));
1053a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_X, "MeanMomentumX"));
1063a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Y, "MeanMomentumY"));
1073a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Z, "MeanMomentumZ"));
1083a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XX, "MeanMomentumFluxXX"));
1093a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YY, "MeanMomentumFluxYY"));
1103a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_ZZ, "MeanMomentumFluxZZ"));
1113a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YZ, "MeanMomentumFluxYZ"));
1123a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XZ, "MeanMomentumFluxXZ"));
1133a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XY, "MeanMomentumFluxXY"));
1143a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_X, "MeanVelocityX"));
1153a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Y, "MeanVelocityY"));
1163a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Z, "MeanVelocityZ"));
11719706a06SJames Wright 
11878837792SJames Wright   // Cleanup
11978837792SJames Wright   PetscCall(PetscFEDestroy(&fe));
12078837792SJames Wright 
1214eed8630SJames Wright   PetscCall(PetscLogStagePop());
122ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
12319706a06SJames Wright }
12419706a06SJames Wright 
1251737222fSJames Wright // Create CeedElemRestriction for collocated data based on associated CeedBasis and CeedElemRestriction
1261737222fSJames Wright // Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction
1271737222fSJames Wright PetscErrorCode CreateElemRestrColloc(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base,
128ed331efdSJames Wright                                      CeedElemRestriction *elem_restr_collocated) {
1291737222fSJames Wright   CeedInt num_elem_qpts, loc_num_elem;
1301737222fSJames Wright   PetscFunctionBeginUser;
1311737222fSJames Wright 
1321737222fSJames Wright   CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts);
1331737222fSJames Wright   CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem);
1341737222fSJames Wright 
1351737222fSJames Wright   const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp};
1361737222fSJames Wright   CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides,
1371737222fSJames Wright                                    elem_restr_collocated);
138ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1391737222fSJames Wright }
1401737222fSJames Wright 
141a175e481SJames Wright // Get coordinates of quadrature points
142ed331efdSJames Wright PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, Vec *Qx_coords,
14319706a06SJames Wright                                    PetscInt *total_nqpnts) {
144269a910fSJames Wright   CeedElemRestriction  elem_restr_qx;
14519706a06SJames Wright   CeedQFunction        qf_quad_coords;
14619706a06SJames Wright   CeedOperator         op_quad_coords;
147*457e73b2SJames Wright   CeedInt              num_comp_x, loc_num_elem, num_elem_qpts;
148ed331efdSJames Wright   OperatorApplyContext op_quad_coords_ctx;
14919706a06SJames Wright 
150ed331efdSJames Wright   PetscFunctionBeginUser;
15119706a06SJames Wright   // Create Element Restriction and CeedVector for quadrature coordinates
15219706a06SJames Wright   CeedBasisGetNumQuadraturePoints(basis_x, &num_elem_qpts);
15319706a06SJames Wright   CeedElemRestrictionGetNumElements(elem_restr_x, &loc_num_elem);
15419706a06SJames Wright   CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x);
15519706a06SJames Wright   *total_nqpnts = num_elem_qpts * loc_num_elem;
156ed331efdSJames Wright   PetscCall(CreateElemRestrColloc(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx));
15719706a06SJames Wright 
15819706a06SJames Wright   // Create QFunction
15919706a06SJames Wright   CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords);
16019706a06SJames Wright 
16119706a06SJames Wright   // Create Operator
16219706a06SJames Wright   CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords);
163ed331efdSJames Wright   CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, x_coords);
16419706a06SJames Wright   CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
16519706a06SJames Wright 
166ed331efdSJames Wright   PetscCall(CeedOperatorCreateLocalVecs(op_quad_coords, DMReturnVecType(dm), PetscObjectComm((PetscObject)dm), NULL, Qx_coords));
167ed331efdSJames Wright   PetscCall(OperatorApplyContextCreate(NULL, NULL, ceed, op_quad_coords, CEED_VECTOR_NONE, NULL, NULL, NULL, &op_quad_coords_ctx));
16819706a06SJames Wright 
169ed331efdSJames Wright   PetscCall(ApplyCeedOperatorLocalToLocal(NULL, *Qx_coords, op_quad_coords_ctx));
170ed331efdSJames Wright 
171ed331efdSJames Wright   PetscCall(OperatorApplyContextDestroy(op_quad_coords_ctx));
172269a910fSJames Wright   CeedElemRestrictionDestroy(&elem_restr_qx);
17319706a06SJames Wright   CeedQFunctionDestroy(&qf_quad_coords);
17419706a06SJames Wright   CeedOperatorDestroy(&op_quad_coords);
175ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
17619706a06SJames Wright }
17719706a06SJames Wright 
178269a910fSJames Wright PetscErrorCode SpanStatsSetupDataCreate(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem, SpanStatsSetupData *stats_data) {
179269a910fSJames Wright   DM       dm = user->spanstats.dm;
180*457e73b2SJames Wright   PetscInt dim;
181*457e73b2SJames Wright   CeedInt  P, Q, num_comp_x, num_comp_stats = user->spanstats.num_comp_stats;
182269a910fSJames Wright   Vec      X_loc;
183269a910fSJames Wright   PetscFunctionBeginUser;
184269a910fSJames Wright 
185269a910fSJames Wright   PetscCall(PetscNew(stats_data));
186269a910fSJames Wright 
187269a910fSJames Wright   PetscCall(DMGetDimension(dm, &dim));
188269a910fSJames Wright   CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &Q);
189269a910fSJames Wright   CeedBasisGetNumNodes1D(ceed_data->basis_q, &P);
190269a910fSJames Wright 
191431cd09aSJames Wright   PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, 0, Q, problem->q_data_size_sur, &(*stats_data)->elem_restr_parent_stats,
192269a910fSJames Wright                                     &(*stats_data)->elem_restr_parent_x, &(*stats_data)->elem_restr_parent_qd));
193269a910fSJames Wright   CeedElemRestrictionGetNumComponents((*stats_data)->elem_restr_parent_x, &num_comp_x);
194269a910fSJames Wright   CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_x, &(*stats_data)->x_coord, NULL);
195269a910fSJames Wright   CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_qd, &(*stats_data)->q_data, NULL);
196269a910fSJames Wright 
197269a910fSJames Wright   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS, &(*stats_data)->basis_x);
198269a910fSJames Wright   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_stats, P, Q, CEED_GAUSS, &(*stats_data)->basis_stats);
199269a910fSJames Wright 
200269a910fSJames Wright   PetscCall(CreateElemRestrColloc(ceed, num_comp_stats, (*stats_data)->basis_stats, (*stats_data)->elem_restr_parent_stats,
201ed331efdSJames Wright                                   &(*stats_data)->elem_restr_parent_colloc));
202ed331efdSJames Wright   PetscCall(CreateElemRestrColloc(ceed, num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q, &(*stats_data)->elem_restr_child_colloc));
203269a910fSJames Wright 
204269a910fSJames Wright   {  // -- Copy DM coordinates into CeedVector
205269a910fSJames Wright     DM cdm;
206269a910fSJames Wright     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
207269a910fSJames Wright     if (cdm) {
208269a910fSJames Wright       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
209269a910fSJames Wright     } else {
210269a910fSJames Wright       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
211269a910fSJames Wright     }
212269a910fSJames Wright   }
213269a910fSJames Wright   PetscCall(VecScale(X_loc, problem->dm_scale));
214fe69b334SJames Wright   PetscCall(VecCopyP2C(X_loc, (*stats_data)->x_coord));
215269a910fSJames Wright 
216ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
217269a910fSJames Wright }
218269a910fSJames Wright 
219269a910fSJames Wright PetscErrorCode SpanStatsSetupDataDestroy(SpanStatsSetupData data) {
220269a910fSJames Wright   PetscFunctionBeginUser;
221269a910fSJames Wright 
222269a910fSJames Wright   CeedElemRestrictionDestroy(&data->elem_restr_parent_x);
223269a910fSJames Wright   CeedElemRestrictionDestroy(&data->elem_restr_parent_stats);
224269a910fSJames Wright   CeedElemRestrictionDestroy(&data->elem_restr_parent_qd);
225269a910fSJames Wright   CeedElemRestrictionDestroy(&data->elem_restr_parent_colloc);
226269a910fSJames Wright   CeedElemRestrictionDestroy(&data->elem_restr_child_colloc);
227269a910fSJames Wright 
228269a910fSJames Wright   CeedBasisDestroy(&data->basis_x);
229269a910fSJames Wright   CeedBasisDestroy(&data->basis_stats);
230269a910fSJames Wright 
231269a910fSJames Wright   CeedVectorDestroy(&data->x_coord);
232269a910fSJames Wright   CeedVectorDestroy(&data->q_data);
233269a910fSJames Wright 
234269a910fSJames Wright   PetscCall(PetscFree(data));
235ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
236269a910fSJames Wright }
237269a910fSJames Wright 
238a175e481SJames Wright // Create PetscSF for child-to-parent communication
239269a910fSJames Wright PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, SpanStatsSetupData stats_data, DM parentdm, DM childdm, PetscSF *statssf) {
240*457e73b2SJames Wright   PetscInt child_num_qpnts, parent_num_qpnts;
241*457e73b2SJames Wright   CeedInt  num_comp_x;
242ed331efdSJames Wright   Vec      Child_qx_coords, Parent_qx_coords;
24319706a06SJames Wright 
244ed331efdSJames Wright   PetscFunctionBeginUser;
245269a910fSJames Wright   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)childdm), statssf));
246269a910fSJames Wright 
24719706a06SJames Wright   // Assume that child and parent have the same number of components
24819706a06SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x);
24919706a06SJames Wright   const PetscInt num_comp_sf = num_comp_x - 1;  // Number of coord components used in the creation of the SF
25019706a06SJames Wright 
251ed331efdSJames Wright   // Get quad_coords for child and parent DM
252ed331efdSJames Wright   PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &Child_qx_coords, &child_num_qpnts));
253ed331efdSJames Wright   PetscCall(GetQuadratureCoords(ceed, parentdm, stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord, &Parent_qx_coords,
254269a910fSJames Wright                                 &parent_num_qpnts));
25519706a06SJames Wright 
256ed331efdSJames Wright   {  // Remove z component of coordinates for matching
25719706a06SJames Wright     const PetscReal *child_quad_coords, *parent_quad_coords;
258ed331efdSJames Wright     PetscReal       *child_coords, *parent_coords;
25919706a06SJames Wright 
260ed331efdSJames Wright     PetscCall(VecGetArrayRead(Child_qx_coords, &child_quad_coords));
261ed331efdSJames Wright     PetscCall(VecGetArrayRead(Parent_qx_coords, &parent_quad_coords));
26219706a06SJames Wright 
26319706a06SJames Wright     PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords));
26419706a06SJames Wright     for (int i = 0; i < child_num_qpnts; i++) {
26519706a06SJames Wright       child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x];
26619706a06SJames Wright       child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x];
26719706a06SJames Wright     }
26819706a06SJames Wright     for (int i = 0; i < parent_num_qpnts; i++) {
26919706a06SJames Wright       parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x];
27019706a06SJames Wright       parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x];
27119706a06SJames Wright     }
272ed331efdSJames Wright     PetscCall(VecRestoreArrayRead(Child_qx_coords, &child_quad_coords));
273ed331efdSJames Wright     PetscCall(VecRestoreArrayRead(Parent_qx_coords, &parent_quad_coords));
27419706a06SJames Wright 
275269a910fSJames Wright     PetscCall(PetscSFSetGraphFromCoordinates(*statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords));
276ed331efdSJames Wright     PetscCall(PetscFree2(child_coords, parent_coords));
277ed331efdSJames Wright   }
27819706a06SJames Wright 
279269a910fSJames Wright   PetscCall(PetscSFViewFromOptions(*statssf, NULL, "-spanstats_sf_view"));
280a175e481SJames Wright 
281ed331efdSJames Wright   PetscCall(VecDestroy(&Child_qx_coords));
282ed331efdSJames Wright   PetscCall(VecDestroy(&Parent_qx_coords));
283ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
28419706a06SJames Wright }
28519706a06SJames Wright 
286ed331efdSJames Wright // @brief Setup RHS and LHS for L^2 projection of statistics
287269a910fSJames Wright PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) {
288ed331efdSJames Wright   CeedOperator  op_mass, op_setup_sur, op_proj_rhs;
289269a910fSJames Wright   CeedQFunction qf_mass, qf_stats_proj;
290269a910fSJames Wright   CeedInt       q_data_size, num_comp_stats = user->spanstats.num_comp_stats;
291f967ad79SJames Wright   MPI_Comm      comm = PetscObjectComm((PetscObject)user->spanstats.dm);
292a175e481SJames Wright 
293ed331efdSJames Wright   PetscFunctionBeginUser;
294ed331efdSJames Wright   // -- Create Operator for RHS of L^2 projection of statistics
295269a910fSJames Wright   // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function.
296269a910fSJames Wright   // Therefore, an Identity QF is sufficient
297269a910fSJames Wright   CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &qf_stats_proj);
298269a910fSJames Wright 
299ed331efdSJames Wright   CeedOperatorCreate(ceed, qf_stats_proj, NULL, NULL, &op_proj_rhs);
300ed331efdSJames Wright   CeedOperatorSetField(op_proj_rhs, "input", stats_data->elem_restr_parent_colloc, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
301ed331efdSJames Wright   CeedOperatorSetField(op_proj_rhs, "output", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE);
302269a910fSJames Wright 
303ed331efdSJames Wright   PetscCall(OperatorApplyContextCreate(NULL, user->spanstats.dm, ceed, op_proj_rhs, NULL, NULL, NULL, NULL, &user->spanstats.op_proj_rhs_ctx));
304ed331efdSJames Wright   PetscCall(CeedOperatorCreateLocalVecs(op_proj_rhs, DMReturnVecType(user->spanstats.dm), comm, &user->spanstats.Parent_Stats_loc, NULL));
305ed331efdSJames Wright 
306ed331efdSJames Wright   // -- Setup LHS of L^2 projection
307269a910fSJames Wright   // Get q_data for mass matrix operator
308269a910fSJames Wright   CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur);
309269a910fSJames Wright   CeedOperatorSetField(op_setup_sur, "dx", stats_data->elem_restr_parent_x, stats_data->basis_x, CEED_VECTOR_ACTIVE);
310269a910fSJames Wright   CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, stats_data->basis_x, CEED_VECTOR_NONE);
311269a910fSJames Wright   CeedOperatorSetField(op_setup_sur, "surface qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
312269a910fSJames Wright   CeedOperatorApply(op_setup_sur, stats_data->x_coord, stats_data->q_data, CEED_REQUEST_IMMEDIATE);
313269a910fSJames Wright 
314a175e481SJames Wright   // CEED Restriction
315269a910fSJames Wright   CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size);
316a175e481SJames Wright 
317a175e481SJames Wright   // Create Mass CeedOperator
318269a910fSJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_stats, q_data_size, &qf_mass));
319a175e481SJames Wright   CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass);
320269a910fSJames Wright   CeedOperatorSetField(op_mass, "u", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE);
321269a910fSJames Wright   CeedOperatorSetField(op_mass, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_COLLOCATED, stats_data->q_data);
322269a910fSJames Wright   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
3254021610dSJames Wright     OperatorApplyContext M_ctx;
326a175e481SJames Wright     Mat                  mat_mass;
327a175e481SJames Wright     KSP                  ksp;
328a175e481SJames Wright 
329ed331efdSJames Wright     PetscCall(OperatorApplyContextCreate(user->spanstats.dm, user->spanstats.dm, user->ceed, op_mass, NULL, NULL, NULL, NULL, &M_ctx));
330ed331efdSJames Wright     PetscCall(CreateMatShell_Ceed(M_ctx, &mat_mass));
331a175e481SJames Wright 
332f967ad79SJames Wright     PetscCall(KSPCreate(comm, &ksp));
333b7d66439SJames Wright     PetscCall(KSPSetOptionsPrefix(ksp, "turbulence_spanstats_"));
334a175e481SJames Wright     {
335a175e481SJames Wright       PC pc;
336a175e481SJames Wright       PetscCall(KSPGetPC(ksp, &pc));
337a175e481SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
338a175e481SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
339a175e481SJames Wright       PetscCall(KSPSetType(ksp, KSPCG));
340a175e481SJames Wright       PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
341a175e481SJames Wright       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
342a175e481SJames Wright     }
343a175e481SJames Wright     PetscCall(KSPSetOperators(ksp, mat_mass, mat_mass));
344a175e481SJames Wright     PetscCall(KSPSetFromOptions(ksp));
345a175e481SJames Wright     user->spanstats.ksp = ksp;
346a175e481SJames Wright   }
347a175e481SJames Wright 
348a175e481SJames Wright   // Cleanup
349a175e481SJames Wright   CeedQFunctionDestroy(&qf_mass);
350269a910fSJames Wright   CeedQFunctionDestroy(&qf_stats_proj);
3516665b873SJames Wright   CeedOperatorDestroy(&op_mass);
352269a910fSJames Wright   CeedOperatorDestroy(&op_setup_sur);
353ed331efdSJames Wright   CeedOperatorDestroy(&op_proj_rhs);
354ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
355a175e481SJames Wright }
356a175e481SJames Wright 
357269a910fSJames Wright // Create CeedOperator for statistics collection
358269a910fSJames Wright PetscErrorCode CreateStatisticCollectionOperator(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data, ProblemData *problem) {
359a175e481SJames Wright   CeedInt                     num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q;
360495b9769SJames Wright   Turbulence_SpanStatsContext collect_ctx;
361495b9769SJames Wright   NewtonianIdealGasContext    newtonian_ig_ctx;
362495b9769SJames Wright   CeedQFunctionContext        collect_context;
363269a910fSJames Wright   CeedQFunction               qf_stats_collect;
364ed331efdSJames Wright   CeedOperator                op_stats_collect;
365ed331efdSJames Wright 
366a175e481SJames Wright   PetscFunctionBeginUser;
367a175e481SJames Wright   CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q);
368a175e481SJames Wright 
369a175e481SJames Wright   // Create Operator for statistics collection
370a175e481SJames Wright   switch (user->phys->state_var) {
371a175e481SJames Wright     case STATEVAR_PRIMITIVE:
372269a910fSJames Wright       CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &qf_stats_collect);
373a175e481SJames Wright       break;
374a175e481SJames Wright     case STATEVAR_CONSERVATIVE:
375269a910fSJames Wright       CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &qf_stats_collect);
376a175e481SJames Wright       break;
377269a910fSJames Wright     default:
378269a910fSJames Wright       SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No statisics collection available for chosen state variable");
379a175e481SJames Wright   }
380a175e481SJames Wright 
381b7d66439SJames Wright   if (user->spanstats.do_mms_test) {
382269a910fSJames Wright     CeedQFunctionDestroy(&qf_stats_collect);
383269a910fSJames Wright     CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest, ChildStatsCollectionMMSTest_loc, &qf_stats_collect);
384a175e481SJames Wright   }
385a175e481SJames Wright 
386269a910fSJames Wright   {  // Setup Collection Context
387495b9769SJames Wright     PetscCall(PetscNew(&collect_ctx));
388495b9769SJames Wright     CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx);
389495b9769SJames Wright     collect_ctx->gas = *newtonian_ig_ctx;
390495b9769SJames Wright 
391495b9769SJames Wright     CeedQFunctionContextCreate(user->ceed, &collect_context);
392495b9769SJames Wright     CeedQFunctionContextSetData(collect_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx);
393495b9769SJames Wright     CeedQFunctionContextSetDataDestroy(collect_context, CEED_MEM_HOST, FreeContextPetsc);
394495b9769SJames Wright 
395495b9769SJames Wright     CeedQFunctionContextRegisterDouble(collect_context, "solution time", offsetof(struct Turbulence_SpanStatsContext_, solution_time), 1,
396495b9769SJames Wright                                        "Current solution time");
397495b9769SJames Wright     CeedQFunctionContextRegisterDouble(collect_context, "previous time", offsetof(struct Turbulence_SpanStatsContext_, previous_time), 1,
398495b9769SJames Wright                                        "Previous time statistics collection was done");
399495b9769SJames Wright 
400495b9769SJames Wright     CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx);
401495b9769SJames Wright   }
402495b9769SJames Wright 
403269a910fSJames Wright   CeedQFunctionSetContext(qf_stats_collect, collect_context);
404495b9769SJames Wright   CeedQFunctionContextDestroy(&collect_context);
405269a910fSJames Wright   CeedQFunctionAddInput(qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP);
406269a910fSJames Wright   CeedQFunctionAddInput(qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE);
407269a910fSJames Wright   CeedQFunctionAddInput(qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP);
408269a910fSJames Wright   CeedQFunctionAddOutput(qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE);
409a175e481SJames Wright 
410ed331efdSJames Wright   CeedOperatorCreate(ceed, qf_stats_collect, NULL, NULL, &op_stats_collect);
411ed331efdSJames Wright   CeedOperatorSetField(op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
412ed331efdSJames Wright   CeedOperatorSetField(op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
413ed331efdSJames Wright   CeedOperatorSetField(op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
414ed331efdSJames Wright   CeedOperatorSetField(op_stats_collect, "v", stats_data->elem_restr_child_colloc, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
415a175e481SJames Wright 
416ed331efdSJames Wright   CeedOperatorGetContextFieldLabel(op_stats_collect, "solution time", &user->spanstats.solution_time_label);
417ed331efdSJames Wright   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 
422ed331efdSJames Wright   PetscCall(CeedOperatorCreateLocalVecs(op_stats_collect, DMReturnVecType(user->spanstats.dm), PetscObjectComm((PetscObject)user->spanstats.dm), NULL,
423ed331efdSJames Wright                                         &user->spanstats.Child_Stats_loc));
424ed331efdSJames Wright   PetscCall(VecZeroEntries(user->spanstats.Child_Stats_loc));
425495b9769SJames Wright 
426269a910fSJames Wright   CeedQFunctionDestroy(&qf_stats_collect);
427ed331efdSJames Wright   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;
437823a1283SJames Wright   PetscFunctionBeginUser;
438823a1283SJames Wright 
439269a910fSJames Wright   CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size);
440269a910fSJames Wright   CeedBasisGetNumComponents(stats_data->basis_x, &num_comp_x);
441823a1283SJames Wright 
442b7d66439SJames Wright   CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest_Error, ChildStatsCollectionMMSTest_Error_loc, &qf_error);
443823a1283SJames Wright   CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP);
444823a1283SJames Wright   CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE);
445823a1283SJames Wright   CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP);
446823a1283SJames Wright   CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP);
447823a1283SJames Wright 
448823a1283SJames Wright   CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
449269a910fSJames Wright   CeedOperatorSetField(op_error, "q", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE);
450269a910fSJames Wright   CeedOperatorSetField(op_error, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_COLLOCATED, stats_data->q_data);
451269a910fSJames Wright   CeedOperatorSetField(op_error, "x", stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord);
452269a910fSJames Wright   CeedOperatorSetField(op_error, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE);
453823a1283SJames Wright 
454269a910fSJames Wright   CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &x_ceed, NULL);
455269a910fSJames Wright   CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &y_ceed, NULL);
4564021610dSJames Wright   PetscCall(OperatorApplyContextCreate(user->spanstats.dm, user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, NULL,
457ea168fcfSJames Wright                                        &user->spanstats.mms_error_ctx));
458823a1283SJames Wright 
4596665b873SJames Wright   CeedOperatorDestroy(&op_error);
460269a910fSJames Wright   CeedQFunctionDestroy(&qf_error);
4616665b873SJames Wright   CeedVectorDestroy(&x_ceed);
4626665b873SJames Wright   CeedVectorDestroy(&y_ceed);
463ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
464823a1283SJames Wright }
465823a1283SJames Wright 
466a175e481SJames Wright // Setup for statistics collection
467f5452247SJames Wright PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) {
468269a910fSJames Wright   SpanStatsSetupData stats_data;
4694eed8630SJames Wright   PetscLogStage      stage_stats_setup;
470f5452247SJames Wright 
47119706a06SJames Wright   PetscFunctionBeginUser;
4724eed8630SJames Wright   PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup));
4734eed8630SJames Wright   if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup));
4744eed8630SJames Wright   PetscCall(PetscLogStagePush(stage_stats_setup));
4754eed8630SJames Wright 
476f5452247SJames Wright   // Create parent DM
477f5452247SJames Wright   PetscCall(CreateStatsDM(user, problem, user->app_ctx->degree));
478f5452247SJames Wright 
479269a910fSJames Wright   // Create necessary CeedObjects for setting up statistics
480269a910fSJames Wright   PetscCall(SpanStatsSetupDataCreate(ceed, user, ceed_data, problem, &stats_data));
4811737222fSJames Wright 
48219706a06SJames Wright   // Create SF for communicating child data back their respective parents
483269a910fSJames Wright   PetscCall(CreateStatsSF(ceed, ceed_data, stats_data, user->dm, user->spanstats.dm, &user->spanstats.sf));
48451ee423eSJames Wright 
485a175e481SJames Wright   // Create CeedOperators for statistics collection
486269a910fSJames Wright   PetscCall(CreateStatisticCollectionOperator(ceed, user, ceed_data, stats_data, problem));
487a175e481SJames Wright 
488a175e481SJames Wright   // Setup KSP and Mat for L^2 projection of statistics
489269a910fSJames Wright   PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data, stats_data));
490a175e481SJames Wright 
491b7d66439SJames Wright   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_turbulence_spanstats_mms", &user->spanstats.do_mms_test, NULL));
492b7d66439SJames Wright   if (user->spanstats.do_mms_test) {
493269a910fSJames Wright     PetscCall(SetupMMSErrorChecking(ceed, user, ceed_data, stats_data));
494823a1283SJames Wright   }
495823a1283SJames Wright 
496f967ad79SJames Wright   {  // Setup stats viewer with prefix
4978ed52730SJames Wright     PetscViewerType viewer_type;
498b7d66439SJames Wright     PetscCall(PetscViewerGetType(user->app_ctx->turb_spanstats_viewer, &viewer_type));
499b7d66439SJames Wright     PetscCall(PetscOptionsSetValue(NULL, "-ts_monitor_turbulence_spanstats_viewer_type", viewer_type));
5008ed52730SJames Wright 
501b7d66439SJames Wright     PetscCall(PetscViewerSetOptionsPrefix(user->app_ctx->turb_spanstats_viewer, "ts_monitor_turbulence_spanstats_"));
502b7d66439SJames Wright     PetscCall(PetscViewerSetFromOptions(user->app_ctx->turb_spanstats_viewer));
5038ed52730SJames Wright   }
5044eed8630SJames Wright 
505269a910fSJames Wright   PetscCall(SpanStatsSetupDataDestroy(stats_data));
5064eed8630SJames Wright   PetscCall(PetscLogStagePop());
507ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
508a175e481SJames Wright }
509a175e481SJames Wright 
510a175e481SJames Wright // Collect statistics based on the solution Q
511a175e481SJames Wright PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) {
512ed331efdSJames Wright   Span_Stats user_stats = user->spanstats;
513a175e481SJames Wright 
514ed331efdSJames Wright   PetscFunctionBeginUser;
5156dcea3beSJames Wright   PetscLogStage stage_stats_collect;
5166dcea3beSJames Wright   PetscCall(PetscLogStageGetId("Stats Collect", &stage_stats_collect));
5176dcea3beSJames Wright   if (stage_stats_collect == -1) PetscCall(PetscLogStageRegister("Stats Collect", &stage_stats_collect));
5186dcea3beSJames Wright   PetscCall(PetscLogStagePush(stage_stats_collect));
5196dcea3beSJames Wright 
520a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, user->Q_loc, solution_time));
521ed331efdSJames Wright   CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.solution_time_label, &solution_time);
522a0b9a424SJames Wright   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, user->Q_loc));
523ed331efdSJames Wright   PetscCall(ApplyAddCeedOperatorLocalToLocal(user->Q_loc, user_stats.Child_Stats_loc, user_stats.op_stats_collect_ctx));
524a175e481SJames Wright 
525ed331efdSJames Wright   CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.previous_time_label, &solution_time);
526a175e481SJames Wright 
5276dcea3beSJames Wright   PetscCall(PetscLogStagePop());
528ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
529a175e481SJames Wright }
530a175e481SJames Wright 
531a175e481SJames Wright // Process the child statistics into parent statistics and project them onto stats
53278bbfb6fSJed Brown PetscErrorCode ProcessStatistics(User user, Vec stats) {
533a175e481SJames Wright   Span_Stats         user_stats = user->spanstats;
534a175e481SJames Wright   const PetscScalar *child_stats;
535a175e481SJames Wright   PetscScalar       *parent_stats;
536a175e481SJames Wright   MPI_Datatype       unit;
537ed331efdSJames Wright   Vec                RHS;
538a175e481SJames Wright 
539ed331efdSJames Wright   PetscFunctionBeginUser;
5406dcea3beSJames Wright   PetscLogStage stage_stats_process;
5416dcea3beSJames Wright   PetscCall(PetscLogStageGetId("Stats Process", &stage_stats_process));
5426dcea3beSJames Wright   if (stage_stats_process == -1) PetscCall(PetscLogStageRegister("Stats Process", &stage_stats_process));
5436dcea3beSJames Wright   PetscCall(PetscLogStagePush(stage_stats_process));
5446dcea3beSJames Wright 
545ed331efdSJames Wright   PetscCall(VecZeroEntries(user_stats.Parent_Stats_loc));
546a175e481SJames Wright 
547ed331efdSJames Wright   PetscCall(VecGetArrayRead(user_stats.Child_Stats_loc, &child_stats));
548ed331efdSJames Wright   PetscCall(VecGetArray(user_stats.Parent_Stats_loc, &parent_stats));
549a175e481SJames Wright 
550a175e481SJames Wright   if (user_stats.num_comp_stats == 1) unit = MPIU_REAL;
551a175e481SJames Wright   else {
552a175e481SJames Wright     PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit));
553a175e481SJames Wright     PetscCallMPI(MPI_Type_commit(&unit));
554a175e481SJames Wright   }
555a175e481SJames Wright 
556a175e481SJames Wright   PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
557a175e481SJames Wright   PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
558a175e481SJames Wright 
559ed331efdSJames Wright   PetscCall(VecRestoreArrayRead(user_stats.Child_Stats_loc, &child_stats));
560ed331efdSJames Wright   PetscCall(VecRestoreArray(user_stats.Parent_Stats_loc, &parent_stats));
561a175e481SJames Wright   PetscCallMPI(MPI_Type_free(&unit));
562a175e481SJames Wright 
56378bbfb6fSJed Brown   PetscReal solution_time;
56478bbfb6fSJed Brown   PetscCall(DMGetOutputSequenceNumber(user_stats.dm, NULL, &solution_time));
56578bbfb6fSJed Brown   PetscReal summing_duration = solution_time - user->app_ctx->cont_time;
566ed331efdSJames Wright   PetscCall(VecScale(user_stats.Parent_Stats_loc, 1 / (summing_duration * user_stats.span_width)));
567a175e481SJames Wright 
568a175e481SJames Wright   // L^2 projection with the parent_data
569ed331efdSJames Wright   PetscCall(DMGetGlobalVector(user_stats.dm, &RHS));
570ed331efdSJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(user_stats.Parent_Stats_loc, RHS, user_stats.op_proj_rhs_ctx));
571a175e481SJames Wright 
572ed331efdSJames Wright   PetscCall(KSPSolve(user_stats.ksp, RHS, stats));
573a175e481SJames Wright 
574ed331efdSJames Wright   PetscCall(DMRestoreGlobalVector(user_stats.dm, &RHS));
5756dcea3beSJames Wright   PetscCall(PetscLogStagePop());
576ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
577a175e481SJames Wright }
578a175e481SJames Wright 
579a175e481SJames Wright // TSMonitor for the statistics collection and processing
580f5452247SJames Wright PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) {
581a175e481SJames Wright   User              user = (User)ctx;
582a175e481SJames Wright   Vec               stats;
58378bbfb6fSJed Brown   TSConvergedReason reason;
584b7d66439SJames Wright   PetscInt collect_interval = user->app_ctx->turb_spanstats_collect_interval, viewer_interval = user->app_ctx->turb_spanstats_viewer_interval;
585a175e481SJames Wright   PetscFunctionBeginUser;
58678bbfb6fSJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
587a175e481SJames Wright   // Do not collect or process on the first step of the run (ie. on the initial condition)
588ee4ca9cbSJames Wright   if (steps == user->app_ctx->cont_steps && reason == TS_CONVERGED_ITERATING) PetscFunctionReturn(PETSC_SUCCESS);
589a175e481SJames Wright 
59027240365SJames Wright   PetscBool run_processing_and_viewer = (steps % viewer_interval == 0 && viewer_interval != -1) || reason != TS_CONVERGED_ITERATING;
591a175e481SJames Wright 
59227240365SJames Wright   if (steps % collect_interval == 0 || run_processing_and_viewer) {
59327240365SJames Wright     PetscCall(CollectStatistics(user, solution_time, Q));
59427240365SJames Wright 
59527240365SJames Wright     if (run_processing_and_viewer) {
5968ed52730SJames Wright       PetscCall(DMSetOutputSequenceNumber(user->spanstats.dm, steps, solution_time));
59778bbfb6fSJed Brown       PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats));
59878bbfb6fSJed Brown       PetscCall(ProcessStatistics(user, stats));
5998fb33541SJames Wright       if (user->app_ctx->test_type == TESTTYPE_NONE) {
600b7d66439SJames Wright         PetscCall(PetscViewerPushFormat(user->app_ctx->turb_spanstats_viewer, user->app_ctx->turb_spanstats_viewer_format));
601b7d66439SJames Wright         PetscCall(VecView(stats, user->app_ctx->turb_spanstats_viewer));
602b7d66439SJames Wright         PetscCall(PetscViewerPopFormat(user->app_ctx->turb_spanstats_viewer));
6038fb33541SJames Wright       }
6048fb33541SJames Wright       if (user->app_ctx->test_type == TESTTYPE_TURB_SPANSTATS && reason != TS_CONVERGED_ITERATING) {
6058fb33541SJames Wright         PetscCall(RegressionTests_NS(user->app_ctx, stats));
6068fb33541SJames Wright       }
607b7d66439SJames Wright       if (user->spanstats.do_mms_test && reason != TS_CONVERGED_ITERATING) {
608823a1283SJames Wright         Vec error;
609823a1283SJames Wright         PetscCall(VecDuplicate(stats, &error));
6104021610dSJames Wright         PetscCall(ApplyCeedOperatorGlobalToGlobal(stats, error, user->spanstats.mms_error_ctx));
611823a1283SJames Wright         PetscScalar error_sq = 0;
612823a1283SJames Wright         PetscCall(VecSum(error, &error_sq));
613823a1283SJames Wright         PetscScalar l2_error = sqrt(error_sq);
614823a1283SJames Wright         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error));
615823a1283SJames Wright       }
616a175e481SJames Wright       PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats));
617522ee345SJames Wright     }
61827240365SJames Wright   }
619ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
62051ee423eSJames Wright }
621fd170fd0SJames Wright 
622f5452247SJames Wright PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data) {
623fd170fd0SJames Wright   PetscFunctionBeginUser;
624fd170fd0SJames Wright 
625fd170fd0SJames Wright   // -- CeedVectors
626ed331efdSJames Wright   PetscCall(VecDestroy(&user->spanstats.Child_Stats_loc));
627ed331efdSJames Wright   PetscCall(VecDestroy(&user->spanstats.Parent_Stats_loc));
628fd170fd0SJames Wright 
629fd170fd0SJames Wright   // -- CeedOperators
630ed331efdSJames Wright   PetscCall(OperatorApplyContextDestroy(user->spanstats.op_stats_collect_ctx));
631ed331efdSJames Wright   PetscCall(OperatorApplyContextDestroy(user->spanstats.op_proj_rhs_ctx));
6324021610dSJames Wright   PetscCall(OperatorApplyContextDestroy(user->spanstats.mms_error_ctx));
633fd170fd0SJames Wright 
634fd170fd0SJames Wright   // -- KSP
635fd170fd0SJames Wright   PetscCall(KSPDestroy(&user->spanstats.ksp));
636fd170fd0SJames Wright 
637fd170fd0SJames Wright   // -- SF
638fd170fd0SJames Wright   PetscCall(PetscSFDestroy(&user->spanstats.sf));
639fd170fd0SJames Wright 
640fd170fd0SJames Wright   // -- DM
641fd170fd0SJames Wright   PetscCall(DMDestroy(&user->spanstats.dm));
642fd170fd0SJames Wright 
643ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
644fd170fd0SJames Wright }
645