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 2551ee423eSJames Wright PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree, SimpleBC bc) { 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, §ion)); 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()); 12219706a06SJames Wright PetscFunctionReturn(0); 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, 128*ed331efdSJames 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); 1381737222fSJames Wright PetscFunctionReturn(0); 1391737222fSJames Wright } 1401737222fSJames Wright 141a175e481SJames Wright // Get coordinates of quadrature points 142*ed331efdSJames 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; 14719706a06SJames Wright PetscInt num_comp_x, loc_num_elem, num_elem_qpts; 148*ed331efdSJames Wright OperatorApplyContext op_quad_coords_ctx; 14919706a06SJames Wright 150*ed331efdSJames 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; 156*ed331efdSJames 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); 163*ed331efdSJames 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 166*ed331efdSJames Wright PetscCall(CeedOperatorCreateLocalVecs(op_quad_coords, DMReturnVecType(dm), PetscObjectComm((PetscObject)dm), NULL, Qx_coords)); 167*ed331efdSJames Wright PetscCall(OperatorApplyContextCreate(NULL, NULL, ceed, op_quad_coords, CEED_VECTOR_NONE, NULL, NULL, NULL, &op_quad_coords_ctx)); 16819706a06SJames Wright 169*ed331efdSJames Wright PetscCall(ApplyCeedOperatorLocalToLocal(NULL, *Qx_coords, op_quad_coords_ctx)); 170*ed331efdSJames Wright 171*ed331efdSJames 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); 17519706a06SJames Wright PetscFunctionReturn(0); 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; 180269a910fSJames Wright CeedInt dim, P, Q, num_comp_x, num_comp_stats = user->spanstats.num_comp_stats; 181269a910fSJames Wright Vec X_loc; 182269a910fSJames Wright PetscFunctionBeginUser; 183269a910fSJames Wright 184269a910fSJames Wright PetscCall(PetscNew(stats_data)); 185269a910fSJames Wright 186269a910fSJames Wright PetscCall(DMGetDimension(dm, &dim)); 187269a910fSJames Wright CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &Q); 188269a910fSJames Wright CeedBasisGetNumNodes1D(ceed_data->basis_q, &P); 189269a910fSJames Wright 190269a910fSJames Wright PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, problem->q_data_size_sur, &(*stats_data)->elem_restr_parent_stats, 191269a910fSJames Wright &(*stats_data)->elem_restr_parent_x, &(*stats_data)->elem_restr_parent_qd)); 192269a910fSJames Wright CeedElemRestrictionGetNumComponents((*stats_data)->elem_restr_parent_x, &num_comp_x); 193269a910fSJames Wright CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_x, &(*stats_data)->x_coord, NULL); 194269a910fSJames Wright CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_qd, &(*stats_data)->q_data, NULL); 195269a910fSJames Wright 196269a910fSJames Wright CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS, &(*stats_data)->basis_x); 197269a910fSJames Wright CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_stats, P, Q, CEED_GAUSS, &(*stats_data)->basis_stats); 198269a910fSJames Wright 199269a910fSJames Wright PetscCall(CreateElemRestrColloc(ceed, num_comp_stats, (*stats_data)->basis_stats, (*stats_data)->elem_restr_parent_stats, 200*ed331efdSJames Wright &(*stats_data)->elem_restr_parent_colloc)); 201*ed331efdSJames Wright PetscCall(CreateElemRestrColloc(ceed, num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q, &(*stats_data)->elem_restr_child_colloc)); 202269a910fSJames Wright 203269a910fSJames Wright { // -- Copy DM coordinates into CeedVector 204269a910fSJames Wright DM cdm; 205269a910fSJames Wright PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 206269a910fSJames Wright if (cdm) { 207269a910fSJames Wright PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 208269a910fSJames Wright } else { 209269a910fSJames Wright PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 210269a910fSJames Wright } 211269a910fSJames Wright } 212269a910fSJames Wright PetscCall(VecScale(X_loc, problem->dm_scale)); 213fe69b334SJames Wright PetscCall(VecCopyP2C(X_loc, (*stats_data)->x_coord)); 214269a910fSJames Wright 215269a910fSJames Wright PetscFunctionReturn(0); 216269a910fSJames Wright } 217269a910fSJames Wright 218269a910fSJames Wright PetscErrorCode SpanStatsSetupDataDestroy(SpanStatsSetupData data) { 219269a910fSJames Wright PetscFunctionBeginUser; 220269a910fSJames Wright 221269a910fSJames Wright CeedElemRestrictionDestroy(&data->elem_restr_parent_x); 222269a910fSJames Wright CeedElemRestrictionDestroy(&data->elem_restr_parent_stats); 223269a910fSJames Wright CeedElemRestrictionDestroy(&data->elem_restr_parent_qd); 224269a910fSJames Wright CeedElemRestrictionDestroy(&data->elem_restr_parent_colloc); 225269a910fSJames Wright CeedElemRestrictionDestroy(&data->elem_restr_child_colloc); 226269a910fSJames Wright 227269a910fSJames Wright CeedBasisDestroy(&data->basis_x); 228269a910fSJames Wright CeedBasisDestroy(&data->basis_stats); 229269a910fSJames Wright 230269a910fSJames Wright CeedVectorDestroy(&data->x_coord); 231269a910fSJames Wright CeedVectorDestroy(&data->q_data); 232269a910fSJames Wright 233269a910fSJames Wright PetscCall(PetscFree(data)); 234269a910fSJames Wright PetscFunctionReturn(0); 235269a910fSJames Wright } 236269a910fSJames Wright 237a175e481SJames Wright // Create PetscSF for child-to-parent communication 238269a910fSJames Wright PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, SpanStatsSetupData stats_data, DM parentdm, DM childdm, PetscSF *statssf) { 23919706a06SJames Wright PetscInt child_num_qpnts, parent_num_qpnts, num_comp_x; 240*ed331efdSJames Wright Vec Child_qx_coords, Parent_qx_coords; 24119706a06SJames Wright 242*ed331efdSJames Wright PetscFunctionBeginUser; 243269a910fSJames Wright PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)childdm), statssf)); 244269a910fSJames Wright 24519706a06SJames Wright // Assume that child and parent have the same number of components 24619706a06SJames Wright CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x); 24719706a06SJames Wright const PetscInt num_comp_sf = num_comp_x - 1; // Number of coord components used in the creation of the SF 24819706a06SJames Wright 249*ed331efdSJames Wright // Get quad_coords for child and parent DM 250*ed331efdSJames Wright PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &Child_qx_coords, &child_num_qpnts)); 251*ed331efdSJames Wright PetscCall(GetQuadratureCoords(ceed, parentdm, stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord, &Parent_qx_coords, 252269a910fSJames Wright &parent_num_qpnts)); 25319706a06SJames Wright 254*ed331efdSJames Wright { // Remove z component of coordinates for matching 25519706a06SJames Wright const PetscReal *child_quad_coords, *parent_quad_coords; 256*ed331efdSJames Wright PetscReal *child_coords, *parent_coords; 25719706a06SJames Wright 258*ed331efdSJames Wright PetscCall(VecGetArrayRead(Child_qx_coords, &child_quad_coords)); 259*ed331efdSJames Wright PetscCall(VecGetArrayRead(Parent_qx_coords, &parent_quad_coords)); 26019706a06SJames Wright 26119706a06SJames Wright PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords)); 26219706a06SJames Wright for (int i = 0; i < child_num_qpnts; i++) { 26319706a06SJames Wright child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x]; 26419706a06SJames Wright child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x]; 26519706a06SJames Wright } 26619706a06SJames Wright for (int i = 0; i < parent_num_qpnts; i++) { 26719706a06SJames Wright parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x]; 26819706a06SJames Wright parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x]; 26919706a06SJames Wright } 270*ed331efdSJames Wright PetscCall(VecRestoreArrayRead(Child_qx_coords, &child_quad_coords)); 271*ed331efdSJames Wright PetscCall(VecRestoreArrayRead(Parent_qx_coords, &parent_quad_coords)); 27219706a06SJames Wright 273269a910fSJames Wright PetscCall(PetscSFSetGraphFromCoordinates(*statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords)); 274*ed331efdSJames Wright PetscCall(PetscFree2(child_coords, parent_coords)); 275*ed331efdSJames Wright } 27619706a06SJames Wright 277269a910fSJames Wright PetscCall(PetscSFViewFromOptions(*statssf, NULL, "-spanstats_sf_view")); 278a175e481SJames Wright 279*ed331efdSJames Wright PetscCall(VecDestroy(&Child_qx_coords)); 280*ed331efdSJames Wright PetscCall(VecDestroy(&Parent_qx_coords)); 28119706a06SJames Wright PetscFunctionReturn(0); 28219706a06SJames Wright } 28319706a06SJames Wright 284*ed331efdSJames Wright // @brief Setup RHS and LHS for L^2 projection of statistics 285269a910fSJames Wright PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) { 286*ed331efdSJames Wright CeedOperator op_mass, op_setup_sur, op_proj_rhs; 287269a910fSJames Wright CeedQFunction qf_mass, qf_stats_proj; 288269a910fSJames Wright CeedInt q_data_size, num_comp_stats = user->spanstats.num_comp_stats; 289f967ad79SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)user->spanstats.dm); 290a175e481SJames Wright 291*ed331efdSJames Wright PetscFunctionBeginUser; 292*ed331efdSJames Wright // -- Create Operator for RHS of L^2 projection of statistics 293269a910fSJames Wright // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function. 294269a910fSJames Wright // Therefore, an Identity QF is sufficient 295269a910fSJames Wright CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &qf_stats_proj); 296269a910fSJames Wright 297*ed331efdSJames Wright CeedOperatorCreate(ceed, qf_stats_proj, NULL, NULL, &op_proj_rhs); 298*ed331efdSJames Wright CeedOperatorSetField(op_proj_rhs, "input", stats_data->elem_restr_parent_colloc, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 299*ed331efdSJames Wright CeedOperatorSetField(op_proj_rhs, "output", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE); 300269a910fSJames Wright 301*ed331efdSJames Wright PetscCall(OperatorApplyContextCreate(NULL, user->spanstats.dm, ceed, op_proj_rhs, NULL, NULL, NULL, NULL, &user->spanstats.op_proj_rhs_ctx)); 302*ed331efdSJames Wright PetscCall(CeedOperatorCreateLocalVecs(op_proj_rhs, DMReturnVecType(user->spanstats.dm), comm, &user->spanstats.Parent_Stats_loc, NULL)); 303*ed331efdSJames Wright 304*ed331efdSJames Wright // -- Setup LHS of L^2 projection 305269a910fSJames Wright // Get q_data for mass matrix operator 306269a910fSJames Wright CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur); 307269a910fSJames Wright CeedOperatorSetField(op_setup_sur, "dx", stats_data->elem_restr_parent_x, stats_data->basis_x, CEED_VECTOR_ACTIVE); 308269a910fSJames Wright CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, stats_data->basis_x, CEED_VECTOR_NONE); 309269a910fSJames Wright CeedOperatorSetField(op_setup_sur, "surface qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 310269a910fSJames Wright CeedOperatorApply(op_setup_sur, stats_data->x_coord, stats_data->q_data, CEED_REQUEST_IMMEDIATE); 311269a910fSJames Wright 312a175e481SJames Wright // CEED Restriction 313269a910fSJames Wright CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size); 314a175e481SJames Wright 315a175e481SJames Wright // Create Mass CeedOperator 316269a910fSJames Wright PetscCall(CreateMassQFunction(ceed, num_comp_stats, q_data_size, &qf_mass)); 317a175e481SJames Wright CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass); 318269a910fSJames Wright CeedOperatorSetField(op_mass, "u", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE); 319269a910fSJames Wright CeedOperatorSetField(op_mass, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_COLLOCATED, stats_data->q_data); 320269a910fSJames Wright CeedOperatorSetField(op_mass, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE); 321a175e481SJames Wright 322f967ad79SJames Wright { // Setup KSP for L^2 projection 3234021610dSJames Wright OperatorApplyContext M_ctx; 324a175e481SJames Wright Mat mat_mass; 325a175e481SJames Wright KSP ksp; 326a175e481SJames Wright 327*ed331efdSJames Wright PetscCall(OperatorApplyContextCreate(user->spanstats.dm, user->spanstats.dm, user->ceed, op_mass, NULL, NULL, NULL, NULL, &M_ctx)); 328*ed331efdSJames Wright PetscCall(CreateMatShell_Ceed(M_ctx, &mat_mass)); 329a175e481SJames Wright 330f967ad79SJames Wright PetscCall(KSPCreate(comm, &ksp)); 331b7d66439SJames Wright PetscCall(KSPSetOptionsPrefix(ksp, "turbulence_spanstats_")); 332a175e481SJames Wright { 333a175e481SJames Wright PC pc; 334a175e481SJames Wright PetscCall(KSPGetPC(ksp, &pc)); 335a175e481SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 336a175e481SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL)); 337a175e481SJames Wright PetscCall(KSPSetType(ksp, KSPCG)); 338a175e481SJames Wright PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL)); 339a175e481SJames Wright PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 340a175e481SJames Wright } 341a175e481SJames Wright PetscCall(KSPSetOperators(ksp, mat_mass, mat_mass)); 342a175e481SJames Wright PetscCall(KSPSetFromOptions(ksp)); 343a175e481SJames Wright user->spanstats.ksp = ksp; 344a175e481SJames Wright } 345a175e481SJames Wright 346a175e481SJames Wright // Cleanup 347a175e481SJames Wright CeedQFunctionDestroy(&qf_mass); 348269a910fSJames Wright CeedQFunctionDestroy(&qf_stats_proj); 3496665b873SJames Wright CeedOperatorDestroy(&op_mass); 350269a910fSJames Wright CeedOperatorDestroy(&op_setup_sur); 351*ed331efdSJames Wright CeedOperatorDestroy(&op_proj_rhs); 352a175e481SJames Wright PetscFunctionReturn(0); 353a175e481SJames Wright } 354a175e481SJames Wright 355269a910fSJames Wright // Create CeedOperator for statistics collection 356269a910fSJames Wright PetscErrorCode CreateStatisticCollectionOperator(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data, ProblemData *problem) { 357a175e481SJames Wright CeedInt num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q; 358495b9769SJames Wright Turbulence_SpanStatsContext collect_ctx; 359495b9769SJames Wright NewtonianIdealGasContext newtonian_ig_ctx; 360495b9769SJames Wright CeedQFunctionContext collect_context; 361269a910fSJames Wright CeedQFunction qf_stats_collect; 362*ed331efdSJames Wright CeedOperator op_stats_collect; 363*ed331efdSJames Wright 364a175e481SJames Wright PetscFunctionBeginUser; 365a175e481SJames Wright CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q); 366a175e481SJames Wright 367a175e481SJames Wright // Create Operator for statistics collection 368a175e481SJames Wright switch (user->phys->state_var) { 369a175e481SJames Wright case STATEVAR_PRIMITIVE: 370269a910fSJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &qf_stats_collect); 371a175e481SJames Wright break; 372a175e481SJames Wright case STATEVAR_CONSERVATIVE: 373269a910fSJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &qf_stats_collect); 374a175e481SJames Wright break; 375269a910fSJames Wright default: 376269a910fSJames Wright SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No statisics collection available for chosen state variable"); 377a175e481SJames Wright } 378a175e481SJames Wright 379b7d66439SJames Wright if (user->spanstats.do_mms_test) { 380269a910fSJames Wright CeedQFunctionDestroy(&qf_stats_collect); 381269a910fSJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest, ChildStatsCollectionMMSTest_loc, &qf_stats_collect); 382a175e481SJames Wright } 383a175e481SJames Wright 384269a910fSJames Wright { // Setup Collection Context 385495b9769SJames Wright PetscCall(PetscNew(&collect_ctx)); 386495b9769SJames Wright CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx); 387495b9769SJames Wright collect_ctx->gas = *newtonian_ig_ctx; 388495b9769SJames Wright 389495b9769SJames Wright CeedQFunctionContextCreate(user->ceed, &collect_context); 390495b9769SJames Wright CeedQFunctionContextSetData(collect_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx); 391495b9769SJames Wright CeedQFunctionContextSetDataDestroy(collect_context, CEED_MEM_HOST, FreeContextPetsc); 392495b9769SJames Wright 393495b9769SJames Wright CeedQFunctionContextRegisterDouble(collect_context, "solution time", offsetof(struct Turbulence_SpanStatsContext_, solution_time), 1, 394495b9769SJames Wright "Current solution time"); 395495b9769SJames Wright CeedQFunctionContextRegisterDouble(collect_context, "previous time", offsetof(struct Turbulence_SpanStatsContext_, previous_time), 1, 396495b9769SJames Wright "Previous time statistics collection was done"); 397495b9769SJames Wright 398495b9769SJames Wright CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx); 399495b9769SJames Wright } 400495b9769SJames Wright 401269a910fSJames Wright CeedQFunctionSetContext(qf_stats_collect, collect_context); 402495b9769SJames Wright CeedQFunctionContextDestroy(&collect_context); 403269a910fSJames Wright CeedQFunctionAddInput(qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP); 404269a910fSJames Wright CeedQFunctionAddInput(qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE); 405269a910fSJames Wright CeedQFunctionAddInput(qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP); 406269a910fSJames Wright CeedQFunctionAddOutput(qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE); 407a175e481SJames Wright 408*ed331efdSJames Wright CeedOperatorCreate(ceed, qf_stats_collect, NULL, NULL, &op_stats_collect); 409*ed331efdSJames Wright CeedOperatorSetField(op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE); 410*ed331efdSJames Wright CeedOperatorSetField(op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 411*ed331efdSJames Wright CeedOperatorSetField(op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord); 412*ed331efdSJames Wright CeedOperatorSetField(op_stats_collect, "v", stats_data->elem_restr_child_colloc, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 413a175e481SJames Wright 414*ed331efdSJames Wright CeedOperatorGetContextFieldLabel(op_stats_collect, "solution time", &user->spanstats.solution_time_label); 415*ed331efdSJames Wright CeedOperatorGetContextFieldLabel(op_stats_collect, "previous time", &user->spanstats.previous_time_label); 416*ed331efdSJames Wright 417*ed331efdSJames Wright PetscCall(OperatorApplyContextCreate(user->dm, user->spanstats.dm, user->ceed, op_stats_collect, user->q_ceed, NULL, NULL, NULL, 418*ed331efdSJames Wright &user->spanstats.op_stats_collect_ctx)); 419*ed331efdSJames Wright 420*ed331efdSJames Wright PetscCall(CeedOperatorCreateLocalVecs(op_stats_collect, DMReturnVecType(user->spanstats.dm), PetscObjectComm((PetscObject)user->spanstats.dm), NULL, 421*ed331efdSJames Wright &user->spanstats.Child_Stats_loc)); 422*ed331efdSJames Wright PetscCall(VecZeroEntries(user->spanstats.Child_Stats_loc)); 423495b9769SJames Wright 424269a910fSJames Wright CeedQFunctionDestroy(&qf_stats_collect); 425*ed331efdSJames Wright CeedOperatorDestroy(&op_stats_collect); 426a175e481SJames Wright PetscFunctionReturn(0); 427a175e481SJames Wright } 428a175e481SJames Wright 429b7d66439SJames Wright // Creates operator for calculating error of method of manufactured solution (MMS) test 430269a910fSJames Wright PetscErrorCode SetupMMSErrorChecking(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) { 431269a910fSJames Wright CeedInt num_comp_stats = user->spanstats.num_comp_stats, num_comp_x, q_data_size; 432823a1283SJames Wright CeedQFunction qf_error; 433823a1283SJames Wright CeedOperator op_error; 434823a1283SJames Wright CeedVector x_ceed, y_ceed; 435823a1283SJames Wright PetscFunctionBeginUser; 436823a1283SJames Wright 437269a910fSJames Wright CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size); 438269a910fSJames Wright CeedBasisGetNumComponents(stats_data->basis_x, &num_comp_x); 439823a1283SJames Wright 440b7d66439SJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest_Error, ChildStatsCollectionMMSTest_Error_loc, &qf_error); 441823a1283SJames Wright CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP); 442823a1283SJames Wright CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE); 443823a1283SJames Wright CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP); 444823a1283SJames Wright CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP); 445823a1283SJames Wright 446823a1283SJames Wright CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error); 447269a910fSJames Wright CeedOperatorSetField(op_error, "q", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE); 448269a910fSJames Wright CeedOperatorSetField(op_error, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_COLLOCATED, stats_data->q_data); 449269a910fSJames Wright CeedOperatorSetField(op_error, "x", stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord); 450269a910fSJames Wright CeedOperatorSetField(op_error, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE); 451823a1283SJames Wright 452269a910fSJames Wright CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &x_ceed, NULL); 453269a910fSJames Wright CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &y_ceed, NULL); 4544021610dSJames Wright PetscCall(OperatorApplyContextCreate(user->spanstats.dm, user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, NULL, 455ea168fcfSJames Wright &user->spanstats.mms_error_ctx)); 456823a1283SJames Wright 4576665b873SJames Wright CeedOperatorDestroy(&op_error); 458269a910fSJames Wright CeedQFunctionDestroy(&qf_error); 4596665b873SJames Wright CeedVectorDestroy(&x_ceed); 4606665b873SJames Wright CeedVectorDestroy(&y_ceed); 461823a1283SJames Wright PetscFunctionReturn(0); 462823a1283SJames Wright } 463823a1283SJames Wright 464a175e481SJames Wright // Setup for statistics collection 46519706a06SJames Wright PetscErrorCode SetupStatsCollection(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) { 466269a910fSJames Wright SpanStatsSetupData stats_data; 4674eed8630SJames Wright PetscLogStage stage_stats_setup; 46819706a06SJames Wright PetscFunctionBeginUser; 4694eed8630SJames Wright PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup)); 4704eed8630SJames Wright if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup)); 4714eed8630SJames Wright PetscCall(PetscLogStagePush(stage_stats_setup)); 4724eed8630SJames Wright 473269a910fSJames Wright // Create necessary CeedObjects for setting up statistics 474269a910fSJames Wright PetscCall(SpanStatsSetupDataCreate(ceed, user, ceed_data, problem, &stats_data)); 4751737222fSJames Wright 47619706a06SJames Wright // Create SF for communicating child data back their respective parents 477269a910fSJames Wright PetscCall(CreateStatsSF(ceed, ceed_data, stats_data, user->dm, user->spanstats.dm, &user->spanstats.sf)); 47851ee423eSJames Wright 479a175e481SJames Wright // Create CeedOperators for statistics collection 480269a910fSJames Wright PetscCall(CreateStatisticCollectionOperator(ceed, user, ceed_data, stats_data, problem)); 481a175e481SJames Wright 482a175e481SJames Wright // Setup KSP and Mat for L^2 projection of statistics 483269a910fSJames Wright PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data, stats_data)); 484a175e481SJames Wright 485b7d66439SJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_turbulence_spanstats_mms", &user->spanstats.do_mms_test, NULL)); 486b7d66439SJames Wright if (user->spanstats.do_mms_test) { 487269a910fSJames Wright PetscCall(SetupMMSErrorChecking(ceed, user, ceed_data, stats_data)); 488823a1283SJames Wright } 489823a1283SJames Wright 490f967ad79SJames Wright { // Setup stats viewer with prefix 4918ed52730SJames Wright PetscViewerType viewer_type; 492b7d66439SJames Wright PetscCall(PetscViewerGetType(user->app_ctx->turb_spanstats_viewer, &viewer_type)); 493b7d66439SJames Wright PetscCall(PetscOptionsSetValue(NULL, "-ts_monitor_turbulence_spanstats_viewer_type", viewer_type)); 4948ed52730SJames Wright 495b7d66439SJames Wright PetscCall(PetscViewerSetOptionsPrefix(user->app_ctx->turb_spanstats_viewer, "ts_monitor_turbulence_spanstats_")); 496b7d66439SJames Wright PetscCall(PetscViewerSetFromOptions(user->app_ctx->turb_spanstats_viewer)); 4978ed52730SJames Wright } 4984eed8630SJames Wright 499269a910fSJames Wright PetscCall(SpanStatsSetupDataDestroy(stats_data)); 5004eed8630SJames Wright PetscCall(PetscLogStagePop()); 501a175e481SJames Wright PetscFunctionReturn(0); 502a175e481SJames Wright } 503a175e481SJames Wright 504a175e481SJames Wright // Collect statistics based on the solution Q 505a175e481SJames Wright PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) { 506*ed331efdSJames Wright Span_Stats user_stats = user->spanstats; 507a175e481SJames Wright 508*ed331efdSJames Wright PetscFunctionBeginUser; 5096dcea3beSJames Wright PetscLogStage stage_stats_collect; 5106dcea3beSJames Wright PetscCall(PetscLogStageGetId("Stats Collect", &stage_stats_collect)); 5116dcea3beSJames Wright if (stage_stats_collect == -1) PetscCall(PetscLogStageRegister("Stats Collect", &stage_stats_collect)); 5126dcea3beSJames Wright PetscCall(PetscLogStagePush(stage_stats_collect)); 5136dcea3beSJames Wright 514a0b9a424SJames Wright PetscCall(UpdateBoundaryValues(user, user->Q_loc, solution_time)); 515*ed331efdSJames Wright CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.solution_time_label, &solution_time); 516a0b9a424SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, user->Q_loc)); 517*ed331efdSJames Wright PetscCall(ApplyAddCeedOperatorLocalToLocal(user->Q_loc, user_stats.Child_Stats_loc, user_stats.op_stats_collect_ctx)); 518a175e481SJames Wright 519*ed331efdSJames Wright CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.previous_time_label, &solution_time); 520a175e481SJames Wright 5216dcea3beSJames Wright PetscCall(PetscLogStagePop()); 522a175e481SJames Wright PetscFunctionReturn(0); 523a175e481SJames Wright } 524a175e481SJames Wright 525a175e481SJames Wright // Process the child statistics into parent statistics and project them onto stats 52678bbfb6fSJed Brown PetscErrorCode ProcessStatistics(User user, Vec stats) { 527a175e481SJames Wright Span_Stats user_stats = user->spanstats; 528a175e481SJames Wright const PetscScalar *child_stats; 529a175e481SJames Wright PetscScalar *parent_stats; 530a175e481SJames Wright MPI_Datatype unit; 531*ed331efdSJames Wright Vec RHS; 532a175e481SJames Wright 533*ed331efdSJames Wright PetscFunctionBeginUser; 5346dcea3beSJames Wright PetscLogStage stage_stats_process; 5356dcea3beSJames Wright PetscCall(PetscLogStageGetId("Stats Process", &stage_stats_process)); 5366dcea3beSJames Wright if (stage_stats_process == -1) PetscCall(PetscLogStageRegister("Stats Process", &stage_stats_process)); 5376dcea3beSJames Wright PetscCall(PetscLogStagePush(stage_stats_process)); 5386dcea3beSJames Wright 539*ed331efdSJames Wright PetscCall(VecZeroEntries(user_stats.Parent_Stats_loc)); 540a175e481SJames Wright 541*ed331efdSJames Wright PetscCall(VecGetArrayRead(user_stats.Child_Stats_loc, &child_stats)); 542*ed331efdSJames Wright PetscCall(VecGetArray(user_stats.Parent_Stats_loc, &parent_stats)); 543a175e481SJames Wright 544a175e481SJames Wright if (user_stats.num_comp_stats == 1) unit = MPIU_REAL; 545a175e481SJames Wright else { 546a175e481SJames Wright PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit)); 547a175e481SJames Wright PetscCallMPI(MPI_Type_commit(&unit)); 548a175e481SJames Wright } 549a175e481SJames Wright 550a175e481SJames Wright PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM)); 551a175e481SJames Wright PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM)); 552a175e481SJames Wright 553*ed331efdSJames Wright PetscCall(VecRestoreArrayRead(user_stats.Child_Stats_loc, &child_stats)); 554*ed331efdSJames Wright PetscCall(VecRestoreArray(user_stats.Parent_Stats_loc, &parent_stats)); 555a175e481SJames Wright PetscCallMPI(MPI_Type_free(&unit)); 556a175e481SJames Wright 55778bbfb6fSJed Brown PetscReal solution_time; 55878bbfb6fSJed Brown PetscCall(DMGetOutputSequenceNumber(user_stats.dm, NULL, &solution_time)); 55978bbfb6fSJed Brown PetscReal summing_duration = solution_time - user->app_ctx->cont_time; 560*ed331efdSJames Wright PetscCall(VecScale(user_stats.Parent_Stats_loc, 1 / (summing_duration * user_stats.span_width))); 561a175e481SJames Wright 562a175e481SJames Wright // L^2 projection with the parent_data 563*ed331efdSJames Wright PetscCall(DMGetGlobalVector(user_stats.dm, &RHS)); 564*ed331efdSJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(user_stats.Parent_Stats_loc, RHS, user_stats.op_proj_rhs_ctx)); 565a175e481SJames Wright 566*ed331efdSJames Wright PetscCall(KSPSolve(user_stats.ksp, RHS, stats)); 567a175e481SJames Wright 568*ed331efdSJames Wright PetscCall(DMRestoreGlobalVector(user_stats.dm, &RHS)); 5696dcea3beSJames Wright PetscCall(PetscLogStagePop()); 570a175e481SJames Wright PetscFunctionReturn(0); 571a175e481SJames Wright } 572a175e481SJames Wright 573a175e481SJames Wright // TSMonitor for the statistics collection and processing 574a175e481SJames Wright PetscErrorCode TSMonitor_Statistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) { 575a175e481SJames Wright User user = (User)ctx; 576a175e481SJames Wright Vec stats; 57778bbfb6fSJed Brown TSConvergedReason reason; 578b7d66439SJames Wright PetscInt collect_interval = user->app_ctx->turb_spanstats_collect_interval, viewer_interval = user->app_ctx->turb_spanstats_viewer_interval; 579a175e481SJames Wright PetscFunctionBeginUser; 58078bbfb6fSJed Brown PetscCall(TSGetConvergedReason(ts, &reason)); 581a175e481SJames Wright // Do not collect or process on the first step of the run (ie. on the initial condition) 58278bbfb6fSJed Brown if (steps == user->app_ctx->cont_steps && reason == TS_CONVERGED_ITERATING) PetscFunctionReturn(0); 583a175e481SJames Wright 58427240365SJames Wright PetscBool run_processing_and_viewer = (steps % viewer_interval == 0 && viewer_interval != -1) || reason != TS_CONVERGED_ITERATING; 585a175e481SJames Wright 58627240365SJames Wright if (steps % collect_interval == 0 || run_processing_and_viewer) { 58727240365SJames Wright PetscCall(CollectStatistics(user, solution_time, Q)); 58827240365SJames Wright 58927240365SJames Wright if (run_processing_and_viewer) { 5908ed52730SJames Wright PetscCall(DMSetOutputSequenceNumber(user->spanstats.dm, steps, solution_time)); 59178bbfb6fSJed Brown PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats)); 59278bbfb6fSJed Brown PetscCall(ProcessStatistics(user, stats)); 5938fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_NONE) { 594b7d66439SJames Wright PetscCall(PetscViewerPushFormat(user->app_ctx->turb_spanstats_viewer, user->app_ctx->turb_spanstats_viewer_format)); 595b7d66439SJames Wright PetscCall(VecView(stats, user->app_ctx->turb_spanstats_viewer)); 596b7d66439SJames Wright PetscCall(PetscViewerPopFormat(user->app_ctx->turb_spanstats_viewer)); 5978fb33541SJames Wright } 5988fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_TURB_SPANSTATS && reason != TS_CONVERGED_ITERATING) { 5998fb33541SJames Wright PetscCall(RegressionTests_NS(user->app_ctx, stats)); 6008fb33541SJames Wright } 601b7d66439SJames Wright if (user->spanstats.do_mms_test && reason != TS_CONVERGED_ITERATING) { 602823a1283SJames Wright Vec error; 603823a1283SJames Wright PetscCall(VecDuplicate(stats, &error)); 6044021610dSJames Wright PetscCall(ApplyCeedOperatorGlobalToGlobal(stats, error, user->spanstats.mms_error_ctx)); 605823a1283SJames Wright PetscScalar error_sq = 0; 606823a1283SJames Wright PetscCall(VecSum(error, &error_sq)); 607823a1283SJames Wright PetscScalar l2_error = sqrt(error_sq); 608823a1283SJames Wright PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error)); 609823a1283SJames Wright } 610a175e481SJames Wright PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats)); 611522ee345SJames Wright } 61227240365SJames Wright } 61351ee423eSJames Wright PetscFunctionReturn(0); 61451ee423eSJames Wright } 615fd170fd0SJames Wright 6166665b873SJames Wright PetscErrorCode DestroyStats(User user, CeedData ceed_data) { 617fd170fd0SJames Wright PetscFunctionBeginUser; 618fd170fd0SJames Wright 619fd170fd0SJames Wright // -- CeedVectors 620*ed331efdSJames Wright PetscCall(VecDestroy(&user->spanstats.Child_Stats_loc)); 621*ed331efdSJames Wright PetscCall(VecDestroy(&user->spanstats.Parent_Stats_loc)); 622fd170fd0SJames Wright 623fd170fd0SJames Wright // -- CeedOperators 624*ed331efdSJames Wright PetscCall(OperatorApplyContextDestroy(user->spanstats.op_stats_collect_ctx)); 625*ed331efdSJames Wright PetscCall(OperatorApplyContextDestroy(user->spanstats.op_proj_rhs_ctx)); 6264021610dSJames Wright PetscCall(OperatorApplyContextDestroy(user->spanstats.mms_error_ctx)); 627fd170fd0SJames Wright 628fd170fd0SJames Wright // -- KSP 629fd170fd0SJames Wright PetscCall(KSPDestroy(&user->spanstats.ksp)); 630fd170fd0SJames Wright 631fd170fd0SJames Wright // -- SF 632fd170fd0SJames Wright PetscCall(PetscSFDestroy(&user->spanstats.sf)); 633fd170fd0SJames Wright 634fd170fd0SJames Wright // -- DM 635fd170fd0SJames Wright PetscCall(DMDestroy(&user->spanstats.dm)); 636fd170fd0SJames Wright 637fd170fd0SJames Wright PetscFunctionReturn(0); 638fd170fd0SJames Wright } 639