xref: /libCEED/examples/fluids/src/turb_spanstats.c (revision d4cc18453651bd0f94c1a2e078b2646a92dafdcc)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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 
CreateStatsDM(User user,ProblemData problem,PetscInt degree)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));
434aa5fbb4SJames 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       }
6838690fecSJames Wright       PetscCall(PetscSFDestroy(&inv_isoperiodicface));
69c9198418SJames Wright     } else {
7051ee423eSJames Wright       PetscCall(DMGetLabel(user->dm, "Face Sets", &label));
71c9198418SJames Wright     }
72c9198418SJames Wright 
7351ee423eSJames Wright     PetscCall(DMPlexLabelComplete(user->dm, label));
74342228c7SAlex Lindsay     PetscCall(DMPlexFilter(user->dm, label, 1, PETSC_FALSE, PETSC_FALSE, NULL, &user->spanstats.dm));
7549a40c8aSKenneth E. Jansen     PetscCall(DMSetCoordinateDisc(user->spanstats.dm, NULL, PETSC_TRUE));  // Ensure that a coordinate FE exists
761f595ac1SJames Wright 
771f595ac1SJames Wright     PetscCall(DMPlexDistribute(user->spanstats.dm, 0, NULL, &parent_distributed_dm));
781f595ac1SJames Wright     PetscCallMPI(MPI_Comm_size(comm, &size));
791f595ac1SJames Wright     if (parent_distributed_dm) {
801f595ac1SJames Wright       PetscCall(DMDestroy(&user->spanstats.dm));
811f595ac1SJames Wright       user->spanstats.dm = parent_distributed_dm;
821f595ac1SJames Wright     } else if (size > 1) {
831f595ac1SJames Wright       PetscCall(PetscPrintf(comm, "WARNING: Turbulent spanwise statistics: parent DM could not be distributed accross %d ranks.\n", size));
841f595ac1SJames Wright     }
8551ee423eSJames Wright   }
86d68a66c4SJames Wright   {
87d68a66c4SJames Wright     PetscBool is_simplex = PETSC_FALSE;
88d68a66c4SJames Wright     PetscCall(DMPlexIsSimplex(user->spanstats.dm, &is_simplex));
89d68a66c4SJames Wright     PetscCheck(is_simplex != PETSC_TRUE, comm, PETSC_ERR_ARG_WRONGSTATE, "Spanwise Statistics is not implemented for non-tensor product grids");
90d68a66c4SJames Wright   }
9151ee423eSJames Wright 
9251ee423eSJames Wright   PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats"));
93b7d66439SJames Wright   PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "turbulence_spanstats_"));
9451ee423eSJames Wright   PetscCall(DMSetFromOptions(user->spanstats.dm));
95b57b8e72SJames Wright   PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view"));
96c9198418SJames Wright 
9778837792SJames Wright   // Create FE space for parent DM
98d68a66c4SJames Wright   PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, user->app_ctx->degree, 1, user->app_ctx->q_extra, 1, &user->spanstats.num_comp_stats,
99d68a66c4SJames Wright                                user->spanstats.dm));
10051ee423eSJames Wright 
10178837792SJames Wright   // Create Section for data
10251ee423eSJames Wright   PetscCall(DMGetLocalSection(user->spanstats.dm, &section));
10351ee423eSJames Wright   PetscCall(PetscSectionSetFieldName(section, 0, ""));
1043a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY, "MeanDensity"));
1053a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE, "MeanPressure"));
1063a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_SQUARED, "MeanPressureSquared"));
1073a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_X, "MeanPressureVelocityX"));
1083a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Y, "MeanPressureVelocityY"));
1093a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Z, "MeanPressureVelocityZ"));
1103a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE, "MeanDensityTemperature"));
1113a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_X, "MeanDensityTemperatureFluxX"));
1123a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Y, "MeanDensityTemperatureFluxY"));
1133a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Z, "MeanDensityTemperatureFluxZ"));
1143a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_X, "MeanMomentumX"));
1153a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Y, "MeanMomentumY"));
1163a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Z, "MeanMomentumZ"));
1173a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XX, "MeanMomentumFluxXX"));
1183a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YY, "MeanMomentumFluxYY"));
1193a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_ZZ, "MeanMomentumFluxZZ"));
1203a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YZ, "MeanMomentumFluxYZ"));
1213a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XZ, "MeanMomentumFluxXZ"));
1223a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XY, "MeanMomentumFluxXY"));
1233a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_X, "MeanVelocityX"));
1243a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Y, "MeanVelocityY"));
1253a4208e6SJames Wright   PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Z, "MeanVelocityZ"));
12619706a06SJames Wright 
1274eed8630SJames Wright   PetscCall(PetscLogStagePop());
128ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
12919706a06SJames Wright }
13019706a06SJames Wright 
131701bc832SJames Wright /** @brief Create CeedElemRestriction for collocated data in component-major order.
132701bc832SJames Wright a. Sets the strides of the restriction to component-major order
133701bc832SJames Wright  Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction.
134701bc832SJames Wright */
CreateElemRestrColloc_CompMajor(Ceed ceed,CeedInt num_comp,CeedBasis basis,CeedElemRestriction elem_restr_base,CeedElemRestriction * elem_restr_collocated)135701bc832SJames Wright static PetscErrorCode CreateElemRestrColloc_CompMajor(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base,
136ed331efdSJames Wright                                                       CeedElemRestriction *elem_restr_collocated) {
1371737222fSJames Wright   CeedInt num_elem_qpts, loc_num_elem;
1381737222fSJames Wright 
139f17d818dSJames Wright   PetscFunctionBeginUser;
140a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts));
141a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem));
1421737222fSJames Wright 
1431737222fSJames Wright   const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp};
144a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides,
145a424bcd0SJames Wright                                                        elem_restr_collocated));
146ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1471737222fSJames Wright }
1481737222fSJames Wright 
149a175e481SJames Wright // Get coordinates of quadrature points
GetQuadratureCoords(Ceed ceed,DM dm,CeedElemRestriction elem_restr_x,CeedBasis basis_x,CeedVector x_coords,Vec * Qx_coords,PetscInt * total_nqpnts)150ed331efdSJames Wright PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, Vec *Qx_coords,
15119706a06SJames Wright                                    PetscInt *total_nqpnts) {
152269a910fSJames Wright   CeedElemRestriction  elem_restr_qx;
15319706a06SJames Wright   CeedQFunction        qf_quad_coords;
15419706a06SJames Wright   CeedOperator         op_quad_coords;
155701bc832SJames Wright   CeedInt              num_comp_x;
156701bc832SJames Wright   CeedSize             l_vec_size;
157ed331efdSJames Wright   OperatorApplyContext op_quad_coords_ctx;
15819706a06SJames Wright 
159ed331efdSJames Wright   PetscFunctionBeginUser;
160a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
161701bc832SJames Wright   PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx));
162701bc832SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetLVectorSize(elem_restr_qx, &l_vec_size));
163701bc832SJames Wright   *total_nqpnts = l_vec_size / num_comp_x;
16419706a06SJames Wright 
16519706a06SJames Wright   // Create QFunction
166a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords));
16719706a06SJames Wright 
16819706a06SJames Wright   // Create Operator
169a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords));
170a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, x_coords));
171356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
17219706a06SJames Wright 
1732788647cSJames Wright   PetscCall(CeedOperatorCreateLocalVecs(op_quad_coords, DMReturnVecType(dm), PETSC_COMM_SELF, NULL, Qx_coords));
174ed331efdSJames Wright   PetscCall(OperatorApplyContextCreate(NULL, NULL, ceed, op_quad_coords, CEED_VECTOR_NONE, NULL, NULL, NULL, &op_quad_coords_ctx));
17519706a06SJames Wright 
176ed331efdSJames Wright   PetscCall(ApplyCeedOperatorLocalToLocal(NULL, *Qx_coords, op_quad_coords_ctx));
177ed331efdSJames Wright 
178ed331efdSJames Wright   PetscCall(OperatorApplyContextDestroy(op_quad_coords_ctx));
179a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qx));
180a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_quad_coords));
181a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_quad_coords));
182ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
18319706a06SJames Wright }
18419706a06SJames Wright 
SpanStatsSetupDataCreate(Ceed ceed,User user,CeedData ceed_data,ProblemData problem,SpanStatsSetupData * stats_data)185731c13d7SJames Wright PetscErrorCode SpanStatsSetupDataCreate(Ceed ceed, User user, CeedData ceed_data, ProblemData problem, SpanStatsSetupData *stats_data) {
186269a910fSJames Wright   DM       dm = user->spanstats.dm;
187bb85d312SJames Wright   CeedInt  num_comp_x, num_comp_stats = user->spanstats.num_comp_stats;
188269a910fSJames Wright   Vec      X_loc;
189bb85d312SJames Wright   DMLabel  domain_label = NULL;
190bb85d312SJames Wright   PetscInt label_value = 0, height = 0, dm_field = 0;
191269a910fSJames Wright 
192f17d818dSJames Wright   PetscFunctionBeginUser;
193269a910fSJames Wright   PetscCall(PetscNew(stats_data));
194269a910fSJames Wright 
195bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &(*stats_data)->elem_restr_parent_stats));
196bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &(*stats_data)->elem_restr_parent_x));
197a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents((*stats_data)->elem_restr_parent_x, &num_comp_x));
198a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_x, &(*stats_data)->x_coord, NULL));
199269a910fSJames Wright 
2000814d5a7SKenneth E. Jansen   {
2010814d5a7SKenneth E. Jansen     DM dm_coord;
2020814d5a7SKenneth E. Jansen     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
203bb85d312SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &(*stats_data)->basis_x));
204bb85d312SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &(*stats_data)->basis_stats));
2050814d5a7SKenneth E. Jansen   }
206269a910fSJames Wright 
207701bc832SJames Wright   PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_stats, (*stats_data)->basis_stats, (*stats_data)->elem_restr_parent_stats,
208ed331efdSJames Wright                                             &(*stats_data)->elem_restr_parent_colloc));
2091a8516d0SJames Wright   PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q,
2101a8516d0SJames Wright                                             &(*stats_data)->elem_restr_child_colloc));
211269a910fSJames Wright 
212269a910fSJames Wright   {  // -- Copy DM coordinates into CeedVector
213269a910fSJames Wright     DM cdm;
214269a910fSJames Wright     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
215269a910fSJames Wright     if (cdm) {
216269a910fSJames Wright       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
217269a910fSJames Wright     } else {
218269a910fSJames Wright       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
219269a910fSJames Wright     }
220269a910fSJames Wright   }
221269a910fSJames Wright   PetscCall(VecScale(X_loc, problem->dm_scale));
222d0593705SJames Wright   PetscCall(VecCopyPetscToCeed(X_loc, (*stats_data)->x_coord));
223ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
224269a910fSJames Wright }
225269a910fSJames Wright 
SpanStatsSetupDataDestroy(SpanStatsSetupData data)226269a910fSJames Wright PetscErrorCode SpanStatsSetupDataDestroy(SpanStatsSetupData data) {
227a424bcd0SJames Wright   Ceed ceed;
228a424bcd0SJames Wright 
229269a910fSJames Wright   PetscFunctionBeginUser;
230a424bcd0SJames Wright   PetscCall(CeedElemRestrictionGetCeed(data->elem_restr_parent_x, &ceed));
231a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_x));
232a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_stats));
233a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_colloc));
234a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_child_colloc));
235269a910fSJames Wright 
236a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_x));
237a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_stats));
238269a910fSJames Wright 
239a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&data->x_coord));
240269a910fSJames Wright 
2419bc66399SJeremy L Thompson   PetscCheck(CeedDestroy(&ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_WORLD, PETSC_ERR_LIB, "Destroying Ceed object failed");
2429bc66399SJeremy L Thompson 
243269a910fSJames Wright   PetscCall(PetscFree(data));
244ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
245269a910fSJames Wright }
246269a910fSJames Wright 
247a175e481SJames Wright // Create PetscSF for child-to-parent communication
CreateStatsSF(Ceed ceed,CeedData ceed_data,SpanStatsSetupData stats_data,DM parentdm,DM childdm,PetscSF * statssf)248269a910fSJames Wright PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, SpanStatsSetupData stats_data, DM parentdm, DM childdm, PetscSF *statssf) {
249457e73b2SJames Wright   PetscInt child_num_qpnts, parent_num_qpnts;
250457e73b2SJames Wright   CeedInt  num_comp_x;
251ed331efdSJames Wright   Vec      Child_qx_coords, Parent_qx_coords;
25219706a06SJames Wright 
253ed331efdSJames Wright   PetscFunctionBeginUser;
254269a910fSJames Wright   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)childdm), statssf));
255269a910fSJames Wright 
25619706a06SJames Wright   // Assume that child and parent have the same number of components
257a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x));
25819706a06SJames Wright   const PetscInt num_comp_sf = num_comp_x - 1;  // Number of coord components used in the creation of the SF
25919706a06SJames Wright 
260ed331efdSJames Wright   // Get quad_coords for child and parent DM
261ed331efdSJames Wright   PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &Child_qx_coords, &child_num_qpnts));
262ed331efdSJames Wright   PetscCall(GetQuadratureCoords(ceed, parentdm, stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord, &Parent_qx_coords,
263269a910fSJames Wright                                 &parent_num_qpnts));
26419706a06SJames Wright 
265ed331efdSJames Wright   {  // Remove z component of coordinates for matching
26619706a06SJames Wright     const PetscReal *child_quad_coords, *parent_quad_coords;
267ed331efdSJames Wright     PetscReal       *child_coords, *parent_coords;
26819706a06SJames Wright 
269ed331efdSJames Wright     PetscCall(VecGetArrayRead(Child_qx_coords, &child_quad_coords));
270ed331efdSJames Wright     PetscCall(VecGetArrayRead(Parent_qx_coords, &parent_quad_coords));
27119706a06SJames Wright 
27219706a06SJames Wright     PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords));
27319706a06SJames Wright     for (int i = 0; i < child_num_qpnts; i++) {
27419706a06SJames Wright       child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x];
27519706a06SJames Wright       child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x];
27619706a06SJames Wright     }
27719706a06SJames Wright     for (int i = 0; i < parent_num_qpnts; i++) {
27819706a06SJames Wright       parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x];
27919706a06SJames Wright       parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x];
28019706a06SJames Wright     }
281ed331efdSJames Wright     PetscCall(VecRestoreArrayRead(Child_qx_coords, &child_quad_coords));
282ed331efdSJames Wright     PetscCall(VecRestoreArrayRead(Parent_qx_coords, &parent_quad_coords));
28319706a06SJames Wright 
284269a910fSJames Wright     PetscCall(PetscSFSetGraphFromCoordinates(*statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords));
285ed331efdSJames Wright     PetscCall(PetscFree2(child_coords, parent_coords));
286ed331efdSJames Wright   }
28719706a06SJames Wright 
288269a910fSJames Wright   PetscCall(PetscSFViewFromOptions(*statssf, NULL, "-spanstats_sf_view"));
289a175e481SJames Wright 
290ed331efdSJames Wright   PetscCall(VecDestroy(&Child_qx_coords));
291ed331efdSJames Wright   PetscCall(VecDestroy(&Parent_qx_coords));
292ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
29319706a06SJames Wright }
29419706a06SJames Wright 
295ed331efdSJames Wright // @brief Setup RHS and LHS for L^2 projection of statistics
SetupL2ProjectionStats(Ceed ceed,User user,CeedData ceed_data,SpanStatsSetupData stats_data)296269a910fSJames Wright PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) {
297cfb075a4SJames Wright   CeedOperator        op_mass, op_proj_rhs;
298269a910fSJames Wright   CeedQFunction       qf_mass, qf_stats_proj;
299269a910fSJames Wright   CeedInt             q_data_size, num_comp_stats = user->spanstats.num_comp_stats;
300cfb075a4SJames Wright   CeedElemRestriction elem_restr_qd;
301cfb075a4SJames Wright   CeedVector          q_data;
302cfb075a4SJames Wright   DMLabel             domain_label = NULL;
303cfb075a4SJames Wright   PetscInt            label_value  = 0;
304a175e481SJames Wright 
305ed331efdSJames Wright   PetscFunctionBeginUser;
306ed331efdSJames Wright   // -- Create Operator for RHS of L^2 projection of statistics
307269a910fSJames Wright   // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function.
308269a910fSJames Wright   // Therefore, an Identity QF is sufficient
309a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &qf_stats_proj));
310269a910fSJames Wright 
311a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_proj, NULL, NULL, &op_proj_rhs));
312356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "input", stats_data->elem_restr_parent_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
313a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "output", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
314269a910fSJames Wright 
315ed331efdSJames Wright   PetscCall(OperatorApplyContextCreate(NULL, user->spanstats.dm, ceed, op_proj_rhs, NULL, NULL, NULL, NULL, &user->spanstats.op_proj_rhs_ctx));
3162788647cSJames Wright   PetscCall(CeedOperatorCreateLocalVecs(op_proj_rhs, DMReturnVecType(user->spanstats.dm), PETSC_COMM_SELF, &user->spanstats.Parent_Stats_loc, NULL));
317cfb075a4SJames 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,
318cfb075a4SJames Wright                      &elem_restr_qd, &q_data, &q_data_size));
319a175e481SJames Wright 
320a175e481SJames Wright   // Create Mass CeedOperator
321269a910fSJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_stats, q_data_size, &qf_mass));
322a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass));
323a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
324cfb075a4SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
325a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
326a175e481SJames Wright 
327f967ad79SJames Wright   {  // Setup KSP for L^2 projection
3287f2a9303SJames Wright     Mat mat_mass;
329a175e481SJames Wright     KSP ksp;
330a175e481SJames Wright 
331243afec9SJames Wright     PetscCall(MatCreateCeed(user->spanstats.dm, user->spanstats.dm, op_mass, NULL, &mat_mass));
332a175e481SJames Wright 
3332788647cSJames Wright     PetscCall(KSPCreate(PetscObjectComm((PetscObject)user->spanstats.dm), &ksp));
334b7d66439SJames Wright     PetscCall(KSPSetOptionsPrefix(ksp, "turbulence_spanstats_"));
335a175e481SJames Wright     {
336a175e481SJames Wright       PC pc;
337a175e481SJames Wright       PetscCall(KSPGetPC(ksp, &pc));
338a175e481SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
339a175e481SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
340a175e481SJames Wright       PetscCall(KSPSetType(ksp, KSPCG));
341a175e481SJames Wright       PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
342a175e481SJames Wright       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
343a175e481SJames Wright     }
3447f2a9303SJames Wright     PetscCall(KSPSetFromOptions_WithMatCeed(ksp, mat_mass));
345a175e481SJames Wright     user->spanstats.ksp = ksp;
346cde30410SJames Wright     PetscCall(MatDestroy(&mat_mass));
347a175e481SJames Wright   }
348a175e481SJames Wright 
349a175e481SJames Wright   // Cleanup
350cfb075a4SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
351cfb075a4SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
352a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
353a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_proj));
354a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
355a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_proj_rhs));
356ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
357a175e481SJames Wright }
358a175e481SJames Wright 
359269a910fSJames Wright // Create CeedOperator for statistics collection
CreateStatisticCollectionOperator(Ceed ceed,User user,CeedData ceed_data,SpanStatsSetupData stats_data,ProblemData problem)360731c13d7SJames Wright PetscErrorCode CreateStatisticCollectionOperator(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data, ProblemData problem) {
361a175e481SJames Wright   CeedInt                     num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q;
362495b9769SJames Wright   Turbulence_SpanStatsContext collect_ctx;
363495b9769SJames Wright   NewtonianIdealGasContext    newtonian_ig_ctx;
364495b9769SJames Wright   CeedQFunctionContext        collect_context;
365269a910fSJames Wright   CeedQFunction               qf_stats_collect;
366ed331efdSJames Wright   CeedOperator                op_stats_collect;
367ed331efdSJames Wright 
368a175e481SJames Wright   PetscFunctionBeginUser;
369a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q));
370a175e481SJames Wright 
371a175e481SJames Wright   // Create Operator for statistics collection
372a175e481SJames Wright   switch (user->phys->state_var) {
373a175e481SJames Wright     case STATEVAR_PRIMITIVE:
374a424bcd0SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &qf_stats_collect));
375a175e481SJames Wright       break;
376a175e481SJames Wright     case STATEVAR_CONSERVATIVE:
377a424bcd0SJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &qf_stats_collect));
378a175e481SJames Wright       break;
379a2d72b6fSJames Wright     case STATEVAR_ENTROPY:
380a2d72b6fSJames Wright       PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Entropy, ChildStatsCollection_Entropy_loc, &qf_stats_collect));
381a2d72b6fSJames Wright       break;
382a175e481SJames Wright   }
383a175e481SJames Wright 
384b7d66439SJames Wright   if (user->spanstats.do_mms_test) {
385a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect));
386a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest, ChildStatsCollectionMMSTest_loc, &qf_stats_collect));
387a175e481SJames Wright   }
388a175e481SJames Wright 
389269a910fSJames Wright   {  // Setup Collection Context
390495b9769SJames Wright     PetscCall(PetscNew(&collect_ctx));
391a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx));
392495b9769SJames Wright     collect_ctx->gas = *newtonian_ig_ctx;
393495b9769SJames Wright 
394a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &collect_context));
395a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextSetData(collect_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx));
396a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(collect_context, CEED_MEM_HOST, FreeContextPetsc));
397495b9769SJames Wright 
398a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(collect_context, "solution time",
399a424bcd0SJames Wright                                                            offsetof(struct Turbulence_SpanStatsContext_, solution_time), 1, "Current solution time"));
4001a8516d0SJames Wright     PetscCallCeed(ceed,
4011a8516d0SJames Wright                   CeedQFunctionContextRegisterDouble(collect_context, "previous time", offsetof(struct Turbulence_SpanStatsContext_, previous_time),
4021a8516d0SJames Wright                                                      1, "Previous time statistics collection was done"));
403495b9769SJames Wright 
404a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx));
405495b9769SJames Wright   }
406495b9769SJames Wright 
407a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_stats_collect, collect_context));
408a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&collect_context));
409a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP));
410a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE));
411a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP));
412a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE));
413a175e481SJames Wright 
414a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_collect, NULL, NULL, &op_stats_collect));
415a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
416356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
417a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
418356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "v", stats_data->elem_restr_child_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
419a175e481SJames Wright 
420a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "solution time", &user->spanstats.solution_time_label));
421a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "previous time", &user->spanstats.previous_time_label));
422ed331efdSJames Wright 
423ed331efdSJames Wright   PetscCall(OperatorApplyContextCreate(user->dm, user->spanstats.dm, user->ceed, op_stats_collect, user->q_ceed, NULL, NULL, NULL,
424ed331efdSJames Wright                                        &user->spanstats.op_stats_collect_ctx));
425ed331efdSJames Wright 
4261a8516d0SJames Wright   PetscCall(CeedOperatorCreateLocalVecs(op_stats_collect, DMReturnVecType(user->spanstats.dm), PETSC_COMM_SELF, NULL,
4271a8516d0SJames Wright                                         &user->spanstats.Child_Stats_loc));
428ed331efdSJames Wright   PetscCall(VecZeroEntries(user->spanstats.Child_Stats_loc));
429495b9769SJames Wright 
430a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect));
431a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_stats_collect));
432ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
433a175e481SJames Wright }
434a175e481SJames Wright 
435b7d66439SJames Wright // Creates operator for calculating error of method of manufactured solution (MMS) test
SetupMMSErrorChecking(Ceed ceed,User user,CeedData ceed_data,SpanStatsSetupData stats_data)436269a910fSJames Wright PetscErrorCode SetupMMSErrorChecking(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) {
437269a910fSJames Wright   CeedInt             num_comp_stats = user->spanstats.num_comp_stats, num_comp_x, q_data_size;
438823a1283SJames Wright   CeedQFunction       qf_error;
439823a1283SJames Wright   CeedOperator        op_error;
440823a1283SJames Wright   CeedVector          x_ceed, y_ceed;
441cfb075a4SJames Wright   DMLabel             domain_label = NULL;
442cfb075a4SJames Wright   PetscInt            label_value  = 0;
443cfb075a4SJames Wright   CeedVector          q_data;
444cfb075a4SJames Wright   CeedElemRestriction elem_restr_parent_qd;
445823a1283SJames Wright 
446f17d818dSJames Wright   PetscFunctionBeginUser;
447cfb075a4SJames 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,
448cfb075a4SJames Wright                      &elem_restr_parent_qd, &q_data, &q_data_size));
449a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumComponents(stats_data->basis_x, &num_comp_x));
450823a1283SJames Wright 
451a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest_Error, ChildStatsCollectionMMSTest_Error_loc, &qf_error));
452a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP));
453a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE));
454a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP));
455a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP));
456823a1283SJames Wright 
457a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error));
458a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "q", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
459cfb075a4SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "qdata", elem_restr_parent_qd, CEED_BASIS_NONE, q_data));
460a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "x", stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord));
461a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_error, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE));
462823a1283SJames Wright 
463a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &x_ceed, NULL));
464a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &y_ceed, NULL));
4654021610dSJames Wright   PetscCall(OperatorApplyContextCreate(user->spanstats.dm, user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, NULL,
466ea168fcfSJames Wright                                        &user->spanstats.mms_error_ctx));
467823a1283SJames Wright 
468cfb075a4SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
469a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&x_ceed));
470a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&y_ceed));
471cfb075a4SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_parent_qd));
472cfb075a4SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_error));
473cfb075a4SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_error));
474ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
475823a1283SJames Wright }
476823a1283SJames Wright 
477a175e481SJames Wright // Setup for statistics collection
TurbulenceStatisticsSetup(Ceed ceed,User user,CeedData ceed_data,ProblemData problem)478731c13d7SJames Wright PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) {
479269a910fSJames Wright   SpanStatsSetupData stats_data;
4804eed8630SJames Wright   PetscLogStage      stage_stats_setup;
481f5452247SJames Wright 
48219706a06SJames Wright   PetscFunctionBeginUser;
4834eed8630SJames Wright   PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup));
4844eed8630SJames Wright   if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup));
4854eed8630SJames Wright   PetscCall(PetscLogStagePush(stage_stats_setup));
4864eed8630SJames Wright 
487f5452247SJames Wright   // Create parent DM
488f5452247SJames Wright   PetscCall(CreateStatsDM(user, problem, user->app_ctx->degree));
489f5452247SJames Wright 
490269a910fSJames Wright   // Create necessary CeedObjects for setting up statistics
491269a910fSJames Wright   PetscCall(SpanStatsSetupDataCreate(ceed, user, ceed_data, problem, &stats_data));
4921737222fSJames Wright 
49319706a06SJames Wright   // Create SF for communicating child data back their respective parents
494269a910fSJames Wright   PetscCall(CreateStatsSF(ceed, ceed_data, stats_data, user->dm, user->spanstats.dm, &user->spanstats.sf));
49551ee423eSJames Wright 
496a175e481SJames Wright   // Create CeedOperators for statistics collection
497269a910fSJames Wright   PetscCall(CreateStatisticCollectionOperator(ceed, user, ceed_data, stats_data, problem));
498a175e481SJames Wright 
499a175e481SJames Wright   // Setup KSP and Mat for L^2 projection of statistics
500269a910fSJames Wright   PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data, stats_data));
501a175e481SJames Wright 
502b7d66439SJames Wright   PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_turbulence_spanstats_mms", &user->spanstats.do_mms_test, NULL));
503b7d66439SJames Wright   if (user->spanstats.do_mms_test) {
504269a910fSJames Wright     PetscCall(SetupMMSErrorChecking(ceed, user, ceed_data, stats_data));
505823a1283SJames Wright   }
506823a1283SJames Wright 
507f967ad79SJames Wright   {  // Setup stats viewer with prefix
5088ed52730SJames Wright     PetscViewerType viewer_type;
509b7d66439SJames Wright     PetscCall(PetscViewerGetType(user->app_ctx->turb_spanstats_viewer, &viewer_type));
510b7d66439SJames Wright     PetscCall(PetscOptionsSetValue(NULL, "-ts_monitor_turbulence_spanstats_viewer_type", viewer_type));
5118ed52730SJames Wright 
512b7d66439SJames Wright     PetscCall(PetscViewerSetOptionsPrefix(user->app_ctx->turb_spanstats_viewer, "ts_monitor_turbulence_spanstats_"));
513b7d66439SJames Wright     PetscCall(PetscViewerSetFromOptions(user->app_ctx->turb_spanstats_viewer));
5148ed52730SJames Wright   }
5154eed8630SJames Wright 
516269a910fSJames Wright   PetscCall(SpanStatsSetupDataDestroy(stats_data));
5174eed8630SJames Wright   PetscCall(PetscLogStagePop());
518ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
519a175e481SJames Wright }
520a175e481SJames Wright 
521a175e481SJames Wright // Collect statistics based on the solution Q
CollectStatistics(User user,PetscScalar solution_time,Vec Q)522a175e481SJames Wright PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) {
523cbef7084SJames Wright   SpanStatsData user_stats = user->spanstats;
524a175e481SJames Wright 
525ed331efdSJames Wright   PetscFunctionBeginUser;
5266dcea3beSJames Wright   PetscLogStage stage_stats_collect;
5276dcea3beSJames Wright   PetscCall(PetscLogStageGetId("Stats Collect", &stage_stats_collect));
5286dcea3beSJames Wright   if (stage_stats_collect == -1) PetscCall(PetscLogStageRegister("Stats Collect", &stage_stats_collect));
5296dcea3beSJames Wright   PetscCall(PetscLogStagePush(stage_stats_collect));
5306dcea3beSJames Wright 
531a0b9a424SJames Wright   PetscCall(UpdateBoundaryValues(user, user->Q_loc, solution_time));
532a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.solution_time_label, &solution_time));
533a0b9a424SJames Wright   PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, user->Q_loc));
534ed331efdSJames Wright   PetscCall(ApplyAddCeedOperatorLocalToLocal(user->Q_loc, user_stats.Child_Stats_loc, user_stats.op_stats_collect_ctx));
535a175e481SJames Wright 
536a424bcd0SJames Wright   PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.previous_time_label, &solution_time));
537a175e481SJames Wright 
5386dcea3beSJames Wright   PetscCall(PetscLogStagePop());
539ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
540a175e481SJames Wright }
541a175e481SJames Wright 
542a175e481SJames Wright // Process the child statistics into parent statistics and project them onto stats
ProcessStatistics(User user,Vec stats)54378bbfb6fSJed Brown PetscErrorCode ProcessStatistics(User user, Vec stats) {
544cbef7084SJames Wright   SpanStatsData      user_stats = user->spanstats;
545a175e481SJames Wright   const PetscScalar *child_stats;
546a175e481SJames Wright   PetscScalar       *parent_stats;
547a175e481SJames Wright   MPI_Datatype       unit;
548ed331efdSJames Wright   Vec                RHS;
549a175e481SJames Wright 
550ed331efdSJames Wright   PetscFunctionBeginUser;
5516dcea3beSJames Wright   PetscLogStage stage_stats_process;
5526dcea3beSJames Wright   PetscCall(PetscLogStageGetId("Stats Process", &stage_stats_process));
5536dcea3beSJames Wright   if (stage_stats_process == -1) PetscCall(PetscLogStageRegister("Stats Process", &stage_stats_process));
5546dcea3beSJames Wright   PetscCall(PetscLogStagePush(stage_stats_process));
5556dcea3beSJames Wright 
556ed331efdSJames Wright   PetscCall(VecZeroEntries(user_stats.Parent_Stats_loc));
557a175e481SJames Wright 
558ed331efdSJames Wright   PetscCall(VecGetArrayRead(user_stats.Child_Stats_loc, &child_stats));
559ed331efdSJames Wright   PetscCall(VecGetArray(user_stats.Parent_Stats_loc, &parent_stats));
560a175e481SJames Wright 
561a175e481SJames Wright   if (user_stats.num_comp_stats == 1) unit = MPIU_REAL;
562a175e481SJames Wright   else {
563a175e481SJames Wright     PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit));
564a175e481SJames Wright     PetscCallMPI(MPI_Type_commit(&unit));
565a175e481SJames Wright   }
566a175e481SJames Wright 
567a175e481SJames Wright   PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
568a175e481SJames Wright   PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM));
569a175e481SJames Wright 
570ed331efdSJames Wright   PetscCall(VecRestoreArrayRead(user_stats.Child_Stats_loc, &child_stats));
571ed331efdSJames Wright   PetscCall(VecRestoreArray(user_stats.Parent_Stats_loc, &parent_stats));
572a175e481SJames Wright   PetscCallMPI(MPI_Type_free(&unit));
573a175e481SJames Wright 
57478bbfb6fSJed Brown   PetscReal solution_time;
57578bbfb6fSJed Brown   PetscCall(DMGetOutputSequenceNumber(user_stats.dm, NULL, &solution_time));
57678bbfb6fSJed Brown   PetscReal summing_duration = solution_time - user->app_ctx->cont_time;
577ed331efdSJames Wright   PetscCall(VecScale(user_stats.Parent_Stats_loc, 1 / (summing_duration * user_stats.span_width)));
578a175e481SJames Wright 
579a175e481SJames Wright   // L^2 projection with the parent_data
580ed331efdSJames Wright   PetscCall(DMGetGlobalVector(user_stats.dm, &RHS));
581ed331efdSJames Wright   PetscCall(ApplyCeedOperatorLocalToGlobal(user_stats.Parent_Stats_loc, RHS, user_stats.op_proj_rhs_ctx));
582a175e481SJames Wright 
583ed331efdSJames Wright   PetscCall(KSPSolve(user_stats.ksp, RHS, stats));
584a175e481SJames Wright 
585ed331efdSJames Wright   PetscCall(DMRestoreGlobalVector(user_stats.dm, &RHS));
5866dcea3beSJames Wright   PetscCall(PetscLogStagePop());
587ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
588a175e481SJames Wright }
589a175e481SJames Wright 
590a175e481SJames Wright // TSMonitor for the statistics collection and processing
TSMonitor_TurbulenceStatistics(TS ts,PetscInt steps,PetscReal solution_time,Vec Q,void * ctx)591f5452247SJames Wright PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) {
592a175e481SJames Wright   User              user = (User)ctx;
593a175e481SJames Wright   Vec               stats;
59478bbfb6fSJed Brown   TSConvergedReason reason;
595b7d66439SJames Wright   PetscInt collect_interval = user->app_ctx->turb_spanstats_collect_interval, viewer_interval = user->app_ctx->turb_spanstats_viewer_interval;
596f17d818dSJames Wright 
597a175e481SJames Wright   PetscFunctionBeginUser;
59878bbfb6fSJed Brown   PetscCall(TSGetConvergedReason(ts, &reason));
599a175e481SJames Wright   // Do not collect or process on the first step of the run (ie. on the initial condition)
6006852f6f6SJames Wright   if (steps == user->app_ctx->cont_steps) PetscFunctionReturn(PETSC_SUCCESS);
601a175e481SJames Wright 
60227240365SJames Wright   PetscBool run_processing_and_viewer = (steps % viewer_interval == 0 && viewer_interval != -1) || reason != TS_CONVERGED_ITERATING;
603a175e481SJames Wright 
60427240365SJames Wright   if (steps % collect_interval == 0 || run_processing_and_viewer) {
60527240365SJames Wright     PetscCall(CollectStatistics(user, solution_time, Q));
60627240365SJames Wright 
60727240365SJames Wright     if (run_processing_and_viewer) {
6088ed52730SJames Wright       PetscCall(DMSetOutputSequenceNumber(user->spanstats.dm, steps, solution_time));
60978bbfb6fSJed Brown       PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats));
61078bbfb6fSJed Brown       PetscCall(ProcessStatistics(user, stats));
6118fb33541SJames Wright       if (user->app_ctx->test_type == TESTTYPE_NONE) {
612b7d66439SJames Wright         PetscCall(PetscViewerPushFormat(user->app_ctx->turb_spanstats_viewer, user->app_ctx->turb_spanstats_viewer_format));
613b7d66439SJames Wright         PetscCall(VecView(stats, user->app_ctx->turb_spanstats_viewer));
614b7d66439SJames Wright         PetscCall(PetscViewerPopFormat(user->app_ctx->turb_spanstats_viewer));
6158fb33541SJames Wright       }
6168fb33541SJames Wright       if (user->app_ctx->test_type == TESTTYPE_TURB_SPANSTATS && reason != TS_CONVERGED_ITERATING) {
6173c9e7ad1SJames Wright         PetscCall(RegressionTest(user->app_ctx, stats));
6188fb33541SJames Wright       }
619b7d66439SJames Wright       if (user->spanstats.do_mms_test && reason != TS_CONVERGED_ITERATING) {
620823a1283SJames Wright         Vec error;
621823a1283SJames Wright         PetscCall(VecDuplicate(stats, &error));
6224021610dSJames Wright         PetscCall(ApplyCeedOperatorGlobalToGlobal(stats, error, user->spanstats.mms_error_ctx));
623823a1283SJames Wright         PetscScalar error_sq = 0;
624823a1283SJames Wright         PetscCall(VecSum(error, &error_sq));
625823a1283SJames Wright         PetscScalar l2_error = sqrt(error_sq);
626823a1283SJames Wright         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error));
627823a1283SJames Wright       }
628a175e481SJames Wright       PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats));
629522ee345SJames Wright     }
63027240365SJames Wright   }
631ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
63251ee423eSJames Wright }
633fd170fd0SJames Wright 
TurbulenceStatisticsDestroy(User user,CeedData ceed_data)634f5452247SJames Wright PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data) {
635fd170fd0SJames Wright   PetscFunctionBeginUser;
636ed331efdSJames Wright   PetscCall(VecDestroy(&user->spanstats.Child_Stats_loc));
637ed331efdSJames Wright   PetscCall(VecDestroy(&user->spanstats.Parent_Stats_loc));
638fd170fd0SJames Wright 
639ed331efdSJames Wright   PetscCall(OperatorApplyContextDestroy(user->spanstats.op_stats_collect_ctx));
640ed331efdSJames Wright   PetscCall(OperatorApplyContextDestroy(user->spanstats.op_proj_rhs_ctx));
6414021610dSJames Wright   PetscCall(OperatorApplyContextDestroy(user->spanstats.mms_error_ctx));
642fd170fd0SJames Wright 
643fd170fd0SJames Wright   PetscCall(KSPDestroy(&user->spanstats.ksp));
644fd170fd0SJames Wright   PetscCall(PetscSFDestroy(&user->spanstats.sf));
645fd170fd0SJames Wright   PetscCall(DMDestroy(&user->spanstats.dm));
646ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
647fd170fd0SJames Wright }
648