151ee423eSJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 251ee423eSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 351ee423eSJames Wright // 451ee423eSJames Wright // SPDX-License-Identifier: BSD-2-Clause 551ee423eSJames Wright // 651ee423eSJames Wright // This file is part of CEED: http://github.com/ceed 751ee423eSJames Wright /// @file 8a175e481SJames Wright /// Functions for setting up and performing statistics collection 9a175e481SJames Wright 10a175e481SJames Wright #include "../qfunctions/turb_spanstats.h" 1151ee423eSJames Wright 1249aac155SJeremy L Thompson #include <ceed.h> 1349aac155SJeremy L Thompson #include <petscdmplex.h> 1449aac155SJeremy L Thompson #include <petscsf.h> 1549aac155SJeremy L Thompson 164021610dSJames Wright #include "../include/petsc_ops.h" 1751ee423eSJames Wright #include "../navierstokes.h" 1851ee423eSJames Wright 19269a910fSJames Wright typedef struct { 20269a910fSJames Wright CeedElemRestriction elem_restr_parent_x, elem_restr_parent_stats, elem_restr_parent_qd, elem_restr_parent_colloc, elem_restr_child_colloc; 21269a910fSJames Wright CeedBasis basis_x, basis_stats; 22269a910fSJames Wright CeedVector x_coord, q_data; 23269a910fSJames Wright } *SpanStatsSetupData; 24269a910fSJames Wright 25f5452247SJames Wright PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree) { 263a4208e6SJames Wright user->spanstats.num_comp_stats = TURB_NUM_COMPONENTS; 27a175e481SJames Wright PetscReal domain_min[3], domain_max[3]; 2878837792SJames Wright PetscSection section; 294eed8630SJames Wright PetscLogStage stage_stats_setup; 301f595ac1SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)user->dm); 3151ee423eSJames Wright PetscFunctionBeginUser; 3251ee423eSJames Wright 334eed8630SJames Wright PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup)); 344eed8630SJames Wright if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup)); 354eed8630SJames Wright PetscCall(PetscLogStagePush(stage_stats_setup)); 364eed8630SJames Wright 37a175e481SJames Wright // Get spanwise length 38a175e481SJames Wright PetscCall(DMGetBoundingBox(user->dm, domain_min, domain_max)); 39a175e481SJames Wright user->spanstats.span_width = domain_max[2] - domain_min[1]; 40a175e481SJames Wright 41f967ad79SJames Wright { // Get DM from surface 421f595ac1SJames Wright DM parent_distributed_dm; 43c9198418SJames Wright PetscSF isoperiodicface; 4451ee423eSJames Wright DMLabel label; 451f595ac1SJames Wright PetscMPIInt size; 46c9198418SJames Wright 47c9198418SJames Wright PetscCall(DMPlexGetIsoperiodicFaceSF(user->dm, &isoperiodicface)); 48c9198418SJames Wright 49c9198418SJames Wright if (isoperiodicface) { 50c9198418SJames Wright PetscSF inv_isoperiodicface; 51c9198418SJames Wright PetscInt nleaves; 52c9198418SJames Wright const PetscInt *ilocal; 53c9198418SJames Wright 54c9198418SJames Wright PetscCall(PetscSFCreateInverseSF(isoperiodicface, &inv_isoperiodicface)); 55c9198418SJames Wright PetscCall(PetscSFGetGraph(inv_isoperiodicface, NULL, &nleaves, &ilocal, NULL)); 56c9198418SJames Wright PetscCall(DMCreateLabel(user->dm, "Periodic Face")); 57c9198418SJames Wright PetscCall(DMGetLabel(user->dm, "Periodic Face", &label)); 58c9198418SJames Wright for (PetscInt i = 0; i < nleaves; i++) { 59c9198418SJames Wright PetscCall(DMLabelSetValue(label, ilocal[i], 1)); 60c9198418SJames Wright } 61c9198418SJames Wright } else { 6251ee423eSJames Wright PetscCall(DMGetLabel(user->dm, "Face Sets", &label)); 63c9198418SJames Wright } 64c9198418SJames Wright 6551ee423eSJames Wright PetscCall(DMPlexLabelComplete(user->dm, label)); 6651ee423eSJames Wright PetscCall(DMPlexFilter(user->dm, label, 1, &user->spanstats.dm)); 6751ee423eSJames Wright PetscCall(DMProjectCoordinates(user->spanstats.dm, NULL)); // Ensure that a coordinate FE exists 681f595ac1SJames Wright 691f595ac1SJames Wright PetscCall(DMPlexDistribute(user->spanstats.dm, 0, NULL, &parent_distributed_dm)); 701f595ac1SJames Wright PetscCallMPI(MPI_Comm_size(comm, &size)); 711f595ac1SJames Wright if (parent_distributed_dm) { 721f595ac1SJames Wright PetscCall(DMDestroy(&user->spanstats.dm)); 731f595ac1SJames Wright user->spanstats.dm = parent_distributed_dm; 741f595ac1SJames Wright } else if (size > 1) { 751f595ac1SJames Wright PetscCall(PetscPrintf(comm, "WARNING: Turbulent spanwise statistics: parent DM could not be distributed accross %d ranks.\n", size)); 761f595ac1SJames Wright } 7751ee423eSJames Wright } 78*d68a66c4SJames Wright { 79*d68a66c4SJames Wright PetscBool is_simplex = PETSC_FALSE; 80*d68a66c4SJames Wright PetscCall(DMPlexIsSimplex(user->spanstats.dm, &is_simplex)); 81*d68a66c4SJames Wright PetscCheck(is_simplex != PETSC_TRUE, comm, PETSC_ERR_ARG_WRONGSTATE, "Spanwise Statistics is not implemented for non-tensor product grids"); 82*d68a66c4SJames Wright } 8351ee423eSJames Wright 8451ee423eSJames Wright PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats")); 85b7d66439SJames Wright PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "turbulence_spanstats_")); 8651ee423eSJames Wright PetscCall(DMSetFromOptions(user->spanstats.dm)); 87b57b8e72SJames Wright PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view")); 88c9198418SJames Wright 8978837792SJames Wright // Create FE space for parent DM 90*d68a66c4SJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, user->app_ctx->degree, 1, user->app_ctx->q_extra, 1, &user->spanstats.num_comp_stats, 91*d68a66c4SJames Wright user->spanstats.dm)); 9251ee423eSJames Wright 9378837792SJames Wright // Create Section for data 9451ee423eSJames Wright PetscCall(DMGetLocalSection(user->spanstats.dm, §ion)); 9551ee423eSJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 963a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY, "MeanDensity")); 973a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE, "MeanPressure")); 983a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_SQUARED, "MeanPressureSquared")); 993a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_X, "MeanPressureVelocityX")); 1003a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Y, "MeanPressureVelocityY")); 1013a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Z, "MeanPressureVelocityZ")); 1023a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE, "MeanDensityTemperature")); 1033a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_X, "MeanDensityTemperatureFluxX")); 1043a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Y, "MeanDensityTemperatureFluxY")); 1053a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Z, "MeanDensityTemperatureFluxZ")); 1063a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_X, "MeanMomentumX")); 1073a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Y, "MeanMomentumY")); 1083a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Z, "MeanMomentumZ")); 1093a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XX, "MeanMomentumFluxXX")); 1103a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YY, "MeanMomentumFluxYY")); 1113a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_ZZ, "MeanMomentumFluxZZ")); 1123a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YZ, "MeanMomentumFluxYZ")); 1133a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XZ, "MeanMomentumFluxXZ")); 1143a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XY, "MeanMomentumFluxXY")); 1153a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_X, "MeanVelocityX")); 1163a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Y, "MeanVelocityY")); 1173a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Z, "MeanVelocityZ")); 11819706a06SJames Wright 1194eed8630SJames Wright PetscCall(PetscLogStagePop()); 120ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 12119706a06SJames Wright } 12219706a06SJames Wright 1231737222fSJames Wright // Create CeedElemRestriction for collocated data based on associated CeedBasis and CeedElemRestriction 1241737222fSJames Wright // Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction 1251737222fSJames Wright PetscErrorCode CreateElemRestrColloc(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base, 126ed331efdSJames Wright CeedElemRestriction *elem_restr_collocated) { 1271737222fSJames Wright CeedInt num_elem_qpts, loc_num_elem; 1281737222fSJames Wright PetscFunctionBeginUser; 1291737222fSJames Wright 130a424bcd0SJames Wright PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts)); 131a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem)); 1321737222fSJames Wright 1331737222fSJames Wright const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp}; 134a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides, 135a424bcd0SJames Wright elem_restr_collocated)); 136ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1371737222fSJames Wright } 1381737222fSJames Wright 139a175e481SJames Wright // Get coordinates of quadrature points 140ed331efdSJames Wright PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, Vec *Qx_coords, 14119706a06SJames Wright PetscInt *total_nqpnts) { 142269a910fSJames Wright CeedElemRestriction elem_restr_qx; 14319706a06SJames Wright CeedQFunction qf_quad_coords; 14419706a06SJames Wright CeedOperator op_quad_coords; 145457e73b2SJames Wright CeedInt num_comp_x, loc_num_elem, num_elem_qpts; 146ed331efdSJames Wright OperatorApplyContext op_quad_coords_ctx; 14719706a06SJames Wright 148ed331efdSJames Wright PetscFunctionBeginUser; 14919706a06SJames Wright // Create Element Restriction and CeedVector for quadrature coordinates 150a424bcd0SJames Wright PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(basis_x, &num_elem_qpts)); 151a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_x, &loc_num_elem)); 152a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x)); 15319706a06SJames Wright *total_nqpnts = num_elem_qpts * loc_num_elem; 154ed331efdSJames Wright PetscCall(CreateElemRestrColloc(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx)); 15519706a06SJames Wright 15619706a06SJames Wright // Create QFunction 157a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords)); 15819706a06SJames Wright 15919706a06SJames Wright // Create Operator 160a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords)); 161a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, x_coords)); 162a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE)); 16319706a06SJames Wright 164ed331efdSJames Wright PetscCall(CeedOperatorCreateLocalVecs(op_quad_coords, DMReturnVecType(dm), PetscObjectComm((PetscObject)dm), NULL, Qx_coords)); 165ed331efdSJames Wright PetscCall(OperatorApplyContextCreate(NULL, NULL, ceed, op_quad_coords, CEED_VECTOR_NONE, NULL, NULL, NULL, &op_quad_coords_ctx)); 16619706a06SJames Wright 167ed331efdSJames Wright PetscCall(ApplyCeedOperatorLocalToLocal(NULL, *Qx_coords, op_quad_coords_ctx)); 168ed331efdSJames Wright 169ed331efdSJames Wright PetscCall(OperatorApplyContextDestroy(op_quad_coords_ctx)); 170a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qx)); 171a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_quad_coords)); 172a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_quad_coords)); 173ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 17419706a06SJames Wright } 17519706a06SJames Wright 176269a910fSJames Wright PetscErrorCode SpanStatsSetupDataCreate(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem, SpanStatsSetupData *stats_data) { 177269a910fSJames Wright DM dm = user->spanstats.dm; 178457e73b2SJames Wright PetscInt dim; 1790814d5a7SKenneth E. Jansen CeedInt num_qpts_child1d, num_comp_x, num_comp_stats = user->spanstats.num_comp_stats; 180269a910fSJames Wright Vec X_loc; 181269a910fSJames Wright PetscFunctionBeginUser; 182269a910fSJames Wright 183269a910fSJames Wright PetscCall(PetscNew(stats_data)); 184269a910fSJames Wright 185269a910fSJames Wright PetscCall(DMGetDimension(dm, &dim)); 186a424bcd0SJames Wright PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints1D(ceed_data->basis_q, &num_qpts_child1d)); 1870814d5a7SKenneth E. Jansen CeedInt num_qpts_parent = CeedIntPow(num_qpts_child1d, dim); 1880814d5a7SKenneth E. Jansen PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, 0, num_qpts_parent, problem->q_data_size_sur, &(*stats_data)->elem_restr_parent_stats, 189269a910fSJames Wright &(*stats_data)->elem_restr_parent_x, &(*stats_data)->elem_restr_parent_qd)); 190a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents((*stats_data)->elem_restr_parent_x, &num_comp_x)); 191a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_x, &(*stats_data)->x_coord, NULL)); 192a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_qd, &(*stats_data)->q_data, NULL)); 193269a910fSJames Wright 1940814d5a7SKenneth E. Jansen { 1950814d5a7SKenneth E. Jansen DM dm_coord; 1960814d5a7SKenneth E. Jansen PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 1970814d5a7SKenneth E. Jansen PetscCall(CreateBasisFromPlex(ceed, dm_coord, 0, 0, 0, 0, &(*stats_data)->basis_x)); 1980814d5a7SKenneth E. Jansen PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, &(*stats_data)->basis_stats)); 1990814d5a7SKenneth E. Jansen } 200269a910fSJames Wright 201269a910fSJames Wright PetscCall(CreateElemRestrColloc(ceed, num_comp_stats, (*stats_data)->basis_stats, (*stats_data)->elem_restr_parent_stats, 202ed331efdSJames Wright &(*stats_data)->elem_restr_parent_colloc)); 203ed331efdSJames Wright PetscCall(CreateElemRestrColloc(ceed, num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q, &(*stats_data)->elem_restr_child_colloc)); 204269a910fSJames Wright 205269a910fSJames Wright { // -- Copy DM coordinates into CeedVector 206269a910fSJames Wright DM cdm; 207269a910fSJames Wright PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 208269a910fSJames Wright if (cdm) { 209269a910fSJames Wright PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 210269a910fSJames Wright } else { 211269a910fSJames Wright PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 212269a910fSJames Wright } 213269a910fSJames Wright } 214269a910fSJames Wright PetscCall(VecScale(X_loc, problem->dm_scale)); 215fe69b334SJames Wright PetscCall(VecCopyP2C(X_loc, (*stats_data)->x_coord)); 216269a910fSJames Wright 217ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 218269a910fSJames Wright } 219269a910fSJames Wright 220269a910fSJames Wright PetscErrorCode SpanStatsSetupDataDestroy(SpanStatsSetupData data) { 221a424bcd0SJames Wright Ceed ceed; 222a424bcd0SJames Wright 223269a910fSJames Wright PetscFunctionBeginUser; 224a424bcd0SJames Wright PetscCall(CeedElemRestrictionGetCeed(data->elem_restr_parent_x, &ceed)); 225a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_x)); 226a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_stats)); 227a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_qd)); 228a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_colloc)); 229a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_child_colloc)); 230269a910fSJames Wright 231a424bcd0SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_x)); 232a424bcd0SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_stats)); 233269a910fSJames Wright 234a424bcd0SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&data->x_coord)); 235a424bcd0SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&data->q_data)); 236269a910fSJames Wright 237269a910fSJames Wright PetscCall(PetscFree(data)); 238ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 239269a910fSJames Wright } 240269a910fSJames Wright 241a175e481SJames Wright // Create PetscSF for child-to-parent communication 242269a910fSJames Wright PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, SpanStatsSetupData stats_data, DM parentdm, DM childdm, PetscSF *statssf) { 243457e73b2SJames Wright PetscInt child_num_qpnts, parent_num_qpnts; 244457e73b2SJames Wright CeedInt num_comp_x; 245ed331efdSJames Wright Vec Child_qx_coords, Parent_qx_coords; 24619706a06SJames Wright 247ed331efdSJames Wright PetscFunctionBeginUser; 248269a910fSJames Wright PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)childdm), statssf)); 249269a910fSJames Wright 25019706a06SJames Wright // Assume that child and parent have the same number of components 251a424bcd0SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x)); 25219706a06SJames Wright const PetscInt num_comp_sf = num_comp_x - 1; // Number of coord components used in the creation of the SF 25319706a06SJames Wright 254ed331efdSJames Wright // Get quad_coords for child and parent DM 255ed331efdSJames Wright PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &Child_qx_coords, &child_num_qpnts)); 256ed331efdSJames Wright PetscCall(GetQuadratureCoords(ceed, parentdm, stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord, &Parent_qx_coords, 257269a910fSJames Wright &parent_num_qpnts)); 25819706a06SJames Wright 259ed331efdSJames Wright { // Remove z component of coordinates for matching 26019706a06SJames Wright const PetscReal *child_quad_coords, *parent_quad_coords; 261ed331efdSJames Wright PetscReal *child_coords, *parent_coords; 26219706a06SJames Wright 263ed331efdSJames Wright PetscCall(VecGetArrayRead(Child_qx_coords, &child_quad_coords)); 264ed331efdSJames Wright PetscCall(VecGetArrayRead(Parent_qx_coords, &parent_quad_coords)); 26519706a06SJames Wright 26619706a06SJames Wright PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords)); 26719706a06SJames Wright for (int i = 0; i < child_num_qpnts; i++) { 26819706a06SJames Wright child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x]; 26919706a06SJames Wright child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x]; 27019706a06SJames Wright } 27119706a06SJames Wright for (int i = 0; i < parent_num_qpnts; i++) { 27219706a06SJames Wright parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x]; 27319706a06SJames Wright parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x]; 27419706a06SJames Wright } 275ed331efdSJames Wright PetscCall(VecRestoreArrayRead(Child_qx_coords, &child_quad_coords)); 276ed331efdSJames Wright PetscCall(VecRestoreArrayRead(Parent_qx_coords, &parent_quad_coords)); 27719706a06SJames Wright 278269a910fSJames Wright PetscCall(PetscSFSetGraphFromCoordinates(*statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords)); 279ed331efdSJames Wright PetscCall(PetscFree2(child_coords, parent_coords)); 280ed331efdSJames Wright } 28119706a06SJames Wright 282269a910fSJames Wright PetscCall(PetscSFViewFromOptions(*statssf, NULL, "-spanstats_sf_view")); 283a175e481SJames Wright 284ed331efdSJames Wright PetscCall(VecDestroy(&Child_qx_coords)); 285ed331efdSJames Wright PetscCall(VecDestroy(&Parent_qx_coords)); 286ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 28719706a06SJames Wright } 28819706a06SJames Wright 289ed331efdSJames Wright // @brief Setup RHS and LHS for L^2 projection of statistics 290269a910fSJames Wright PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) { 291ed331efdSJames Wright CeedOperator op_mass, op_setup_sur, op_proj_rhs; 292269a910fSJames Wright CeedQFunction qf_mass, qf_stats_proj; 293269a910fSJames Wright CeedInt q_data_size, num_comp_stats = user->spanstats.num_comp_stats; 294f967ad79SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)user->spanstats.dm); 295a175e481SJames Wright 296ed331efdSJames Wright PetscFunctionBeginUser; 297ed331efdSJames Wright // -- Create Operator for RHS of L^2 projection of statistics 298269a910fSJames Wright // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function. 299269a910fSJames Wright // Therefore, an Identity QF is sufficient 300a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &qf_stats_proj)); 301269a910fSJames Wright 302a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_proj, NULL, NULL, &op_proj_rhs)); 303a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "input", stats_data->elem_restr_parent_colloc, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE)); 304a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "output", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 305269a910fSJames Wright 306ed331efdSJames Wright PetscCall(OperatorApplyContextCreate(NULL, user->spanstats.dm, ceed, op_proj_rhs, NULL, NULL, NULL, NULL, &user->spanstats.op_proj_rhs_ctx)); 307ed331efdSJames Wright PetscCall(CeedOperatorCreateLocalVecs(op_proj_rhs, DMReturnVecType(user->spanstats.dm), comm, &user->spanstats.Parent_Stats_loc, NULL)); 308ed331efdSJames Wright 309ed331efdSJames Wright // -- Setup LHS of L^2 projection 310269a910fSJames Wright // Get q_data for mass matrix operator 311a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur)); 312a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", stats_data->elem_restr_parent_x, stats_data->basis_x, CEED_VECTOR_ACTIVE)); 313a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, stats_data->basis_x, CEED_VECTOR_NONE)); 314a424bcd0SJames Wright PetscCallCeed(ceed, 315a424bcd0SJames Wright CeedOperatorSetField(op_setup_sur, "surface qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE)); 316a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, stats_data->x_coord, stats_data->q_data, CEED_REQUEST_IMMEDIATE)); 317269a910fSJames Wright 318a175e481SJames Wright // CEED Restriction 319a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size)); 320a175e481SJames Wright 321a175e481SJames Wright // Create Mass CeedOperator 322269a910fSJames Wright PetscCall(CreateMassQFunction(ceed, num_comp_stats, q_data_size, &qf_mass)); 323a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 324a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 325a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_COLLOCATED, stats_data->q_data)); 326a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 327a175e481SJames Wright 328f967ad79SJames Wright { // Setup KSP for L^2 projection 3294021610dSJames Wright OperatorApplyContext M_ctx; 330a175e481SJames Wright Mat mat_mass; 331a175e481SJames Wright KSP ksp; 332a175e481SJames Wright 333ed331efdSJames Wright PetscCall(OperatorApplyContextCreate(user->spanstats.dm, user->spanstats.dm, user->ceed, op_mass, NULL, NULL, NULL, NULL, &M_ctx)); 334ed331efdSJames Wright PetscCall(CreateMatShell_Ceed(M_ctx, &mat_mass)); 335a175e481SJames Wright 336f967ad79SJames Wright PetscCall(KSPCreate(comm, &ksp)); 337b7d66439SJames Wright PetscCall(KSPSetOptionsPrefix(ksp, "turbulence_spanstats_")); 338a175e481SJames Wright { 339a175e481SJames Wright PC pc; 340a175e481SJames Wright PetscCall(KSPGetPC(ksp, &pc)); 341a175e481SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 342a175e481SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL)); 343a175e481SJames Wright PetscCall(KSPSetType(ksp, KSPCG)); 344a175e481SJames Wright PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL)); 345a175e481SJames Wright PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 346a175e481SJames Wright } 347a175e481SJames Wright PetscCall(KSPSetOperators(ksp, mat_mass, mat_mass)); 348a175e481SJames Wright PetscCall(KSPSetFromOptions(ksp)); 349a175e481SJames Wright user->spanstats.ksp = ksp; 350a175e481SJames Wright } 351a175e481SJames Wright 352a175e481SJames Wright // Cleanup 353a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 354a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_proj)); 355a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 356a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur)); 357a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_proj_rhs)); 358ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 359a175e481SJames Wright } 360a175e481SJames Wright 361269a910fSJames Wright // Create CeedOperator for statistics collection 362269a910fSJames Wright PetscErrorCode CreateStatisticCollectionOperator(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data, ProblemData *problem) { 363a175e481SJames Wright CeedInt num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q; 364495b9769SJames Wright Turbulence_SpanStatsContext collect_ctx; 365495b9769SJames Wright NewtonianIdealGasContext newtonian_ig_ctx; 366495b9769SJames Wright CeedQFunctionContext collect_context; 367269a910fSJames Wright CeedQFunction qf_stats_collect; 368ed331efdSJames Wright CeedOperator op_stats_collect; 369ed331efdSJames Wright 370a175e481SJames Wright PetscFunctionBeginUser; 371a424bcd0SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q)); 372a175e481SJames Wright 373a175e481SJames Wright // Create Operator for statistics collection 374a175e481SJames Wright switch (user->phys->state_var) { 375a175e481SJames Wright case STATEVAR_PRIMITIVE: 376a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &qf_stats_collect)); 377a175e481SJames Wright break; 378a175e481SJames Wright case STATEVAR_CONSERVATIVE: 379a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &qf_stats_collect)); 380a175e481SJames Wright break; 381269a910fSJames Wright default: 382269a910fSJames Wright SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No statisics collection available for chosen state variable"); 383a175e481SJames Wright } 384a175e481SJames Wright 385b7d66439SJames Wright if (user->spanstats.do_mms_test) { 386a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect)); 387a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest, ChildStatsCollectionMMSTest_loc, &qf_stats_collect)); 388a175e481SJames Wright } 389a175e481SJames Wright 390269a910fSJames Wright { // Setup Collection Context 391495b9769SJames Wright PetscCall(PetscNew(&collect_ctx)); 392a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx)); 393495b9769SJames Wright collect_ctx->gas = *newtonian_ig_ctx; 394495b9769SJames Wright 395a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &collect_context)); 396a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(collect_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx)); 397a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(collect_context, CEED_MEM_HOST, FreeContextPetsc)); 398495b9769SJames Wright 399a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(collect_context, "solution time", 400a424bcd0SJames Wright offsetof(struct Turbulence_SpanStatsContext_, solution_time), 1, "Current solution time")); 401a424bcd0SJames Wright PetscCallCeed( 402a424bcd0SJames Wright ceed, CeedQFunctionContextRegisterDouble(collect_context, "previous time", offsetof(struct Turbulence_SpanStatsContext_, previous_time), 1, 403a424bcd0SJames Wright "Previous time statistics collection was done")); 404495b9769SJames Wright 405a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx)); 406495b9769SJames Wright } 407495b9769SJames Wright 408a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_stats_collect, collect_context)); 409a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&collect_context)); 410a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP)); 411a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE)); 412a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP)); 413a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE)); 414a175e481SJames Wright 415a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_collect, NULL, NULL, &op_stats_collect)); 416a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 417a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data)); 418a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord)); 419a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "v", stats_data->elem_restr_child_colloc, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE)); 420a175e481SJames Wright 421a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "solution time", &user->spanstats.solution_time_label)); 422a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "previous time", &user->spanstats.previous_time_label)); 423ed331efdSJames Wright 424ed331efdSJames Wright PetscCall(OperatorApplyContextCreate(user->dm, user->spanstats.dm, user->ceed, op_stats_collect, user->q_ceed, NULL, NULL, NULL, 425ed331efdSJames Wright &user->spanstats.op_stats_collect_ctx)); 426ed331efdSJames Wright 427ed331efdSJames Wright PetscCall(CeedOperatorCreateLocalVecs(op_stats_collect, DMReturnVecType(user->spanstats.dm), PetscObjectComm((PetscObject)user->spanstats.dm), NULL, 428ed331efdSJames Wright &user->spanstats.Child_Stats_loc)); 429ed331efdSJames Wright PetscCall(VecZeroEntries(user->spanstats.Child_Stats_loc)); 430495b9769SJames Wright 431a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect)); 432a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_stats_collect)); 433ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 434a175e481SJames Wright } 435a175e481SJames Wright 436b7d66439SJames Wright // Creates operator for calculating error of method of manufactured solution (MMS) test 437269a910fSJames Wright PetscErrorCode SetupMMSErrorChecking(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) { 438269a910fSJames Wright CeedInt num_comp_stats = user->spanstats.num_comp_stats, num_comp_x, q_data_size; 439823a1283SJames Wright CeedQFunction qf_error; 440823a1283SJames Wright CeedOperator op_error; 441823a1283SJames Wright CeedVector x_ceed, y_ceed; 442823a1283SJames Wright PetscFunctionBeginUser; 443823a1283SJames Wright 444a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size)); 445a424bcd0SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(stats_data->basis_x, &num_comp_x)); 446823a1283SJames Wright 447a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest_Error, ChildStatsCollectionMMSTest_Error_loc, &qf_error)); 448a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP)); 449a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE)); 450a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP)); 451a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP)); 452823a1283SJames Wright 453a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error)); 454a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_error, "q", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 455a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_error, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_COLLOCATED, stats_data->q_data)); 456a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_error, "x", stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord)); 457a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_error, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 458823a1283SJames Wright 459a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &x_ceed, NULL)); 460a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &y_ceed, NULL)); 4614021610dSJames Wright PetscCall(OperatorApplyContextCreate(user->spanstats.dm, user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, NULL, 462ea168fcfSJames Wright &user->spanstats.mms_error_ctx)); 463823a1283SJames Wright 464a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_error)); 465a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_error)); 466a424bcd0SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&x_ceed)); 467a424bcd0SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&y_ceed)); 468ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 469823a1283SJames Wright } 470823a1283SJames Wright 471a175e481SJames Wright // Setup for statistics collection 472f5452247SJames Wright PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) { 473269a910fSJames Wright SpanStatsSetupData stats_data; 4744eed8630SJames Wright PetscLogStage stage_stats_setup; 475f5452247SJames Wright 47619706a06SJames Wright PetscFunctionBeginUser; 4774eed8630SJames Wright PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup)); 4784eed8630SJames Wright if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup)); 4794eed8630SJames Wright PetscCall(PetscLogStagePush(stage_stats_setup)); 4804eed8630SJames Wright 481f5452247SJames Wright // Create parent DM 482f5452247SJames Wright PetscCall(CreateStatsDM(user, problem, user->app_ctx->degree)); 483f5452247SJames Wright 484269a910fSJames Wright // Create necessary CeedObjects for setting up statistics 485269a910fSJames Wright PetscCall(SpanStatsSetupDataCreate(ceed, user, ceed_data, problem, &stats_data)); 4861737222fSJames Wright 48719706a06SJames Wright // Create SF for communicating child data back their respective parents 488269a910fSJames Wright PetscCall(CreateStatsSF(ceed, ceed_data, stats_data, user->dm, user->spanstats.dm, &user->spanstats.sf)); 48951ee423eSJames Wright 490a175e481SJames Wright // Create CeedOperators for statistics collection 491269a910fSJames Wright PetscCall(CreateStatisticCollectionOperator(ceed, user, ceed_data, stats_data, problem)); 492a175e481SJames Wright 493a175e481SJames Wright // Setup KSP and Mat for L^2 projection of statistics 494269a910fSJames Wright PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data, stats_data)); 495a175e481SJames Wright 496b7d66439SJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_turbulence_spanstats_mms", &user->spanstats.do_mms_test, NULL)); 497b7d66439SJames Wright if (user->spanstats.do_mms_test) { 498269a910fSJames Wright PetscCall(SetupMMSErrorChecking(ceed, user, ceed_data, stats_data)); 499823a1283SJames Wright } 500823a1283SJames Wright 501f967ad79SJames Wright { // Setup stats viewer with prefix 5028ed52730SJames Wright PetscViewerType viewer_type; 503b7d66439SJames Wright PetscCall(PetscViewerGetType(user->app_ctx->turb_spanstats_viewer, &viewer_type)); 504b7d66439SJames Wright PetscCall(PetscOptionsSetValue(NULL, "-ts_monitor_turbulence_spanstats_viewer_type", viewer_type)); 5058ed52730SJames Wright 506b7d66439SJames Wright PetscCall(PetscViewerSetOptionsPrefix(user->app_ctx->turb_spanstats_viewer, "ts_monitor_turbulence_spanstats_")); 507b7d66439SJames Wright PetscCall(PetscViewerSetFromOptions(user->app_ctx->turb_spanstats_viewer)); 5088ed52730SJames Wright } 5094eed8630SJames Wright 510269a910fSJames Wright PetscCall(SpanStatsSetupDataDestroy(stats_data)); 5114eed8630SJames Wright PetscCall(PetscLogStagePop()); 512ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 513a175e481SJames Wright } 514a175e481SJames Wright 515a175e481SJames Wright // Collect statistics based on the solution Q 516a175e481SJames Wright PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) { 517ed331efdSJames Wright Span_Stats user_stats = user->spanstats; 518a175e481SJames Wright 519ed331efdSJames Wright PetscFunctionBeginUser; 5206dcea3beSJames Wright PetscLogStage stage_stats_collect; 5216dcea3beSJames Wright PetscCall(PetscLogStageGetId("Stats Collect", &stage_stats_collect)); 5226dcea3beSJames Wright if (stage_stats_collect == -1) PetscCall(PetscLogStageRegister("Stats Collect", &stage_stats_collect)); 5236dcea3beSJames Wright PetscCall(PetscLogStagePush(stage_stats_collect)); 5246dcea3beSJames Wright 525a0b9a424SJames Wright PetscCall(UpdateBoundaryValues(user, user->Q_loc, solution_time)); 526a424bcd0SJames Wright PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.solution_time_label, &solution_time)); 527a0b9a424SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, user->Q_loc)); 528ed331efdSJames Wright PetscCall(ApplyAddCeedOperatorLocalToLocal(user->Q_loc, user_stats.Child_Stats_loc, user_stats.op_stats_collect_ctx)); 529a175e481SJames Wright 530a424bcd0SJames Wright PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.previous_time_label, &solution_time)); 531a175e481SJames Wright 5326dcea3beSJames Wright PetscCall(PetscLogStagePop()); 533ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 534a175e481SJames Wright } 535a175e481SJames Wright 536a175e481SJames Wright // Process the child statistics into parent statistics and project them onto stats 53778bbfb6fSJed Brown PetscErrorCode ProcessStatistics(User user, Vec stats) { 538a175e481SJames Wright Span_Stats user_stats = user->spanstats; 539a175e481SJames Wright const PetscScalar *child_stats; 540a175e481SJames Wright PetscScalar *parent_stats; 541a175e481SJames Wright MPI_Datatype unit; 542ed331efdSJames Wright Vec RHS; 543a175e481SJames Wright 544ed331efdSJames Wright PetscFunctionBeginUser; 5456dcea3beSJames Wright PetscLogStage stage_stats_process; 5466dcea3beSJames Wright PetscCall(PetscLogStageGetId("Stats Process", &stage_stats_process)); 5476dcea3beSJames Wright if (stage_stats_process == -1) PetscCall(PetscLogStageRegister("Stats Process", &stage_stats_process)); 5486dcea3beSJames Wright PetscCall(PetscLogStagePush(stage_stats_process)); 5496dcea3beSJames Wright 550ed331efdSJames Wright PetscCall(VecZeroEntries(user_stats.Parent_Stats_loc)); 551a175e481SJames Wright 552ed331efdSJames Wright PetscCall(VecGetArrayRead(user_stats.Child_Stats_loc, &child_stats)); 553ed331efdSJames Wright PetscCall(VecGetArray(user_stats.Parent_Stats_loc, &parent_stats)); 554a175e481SJames Wright 555a175e481SJames Wright if (user_stats.num_comp_stats == 1) unit = MPIU_REAL; 556a175e481SJames Wright else { 557a175e481SJames Wright PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit)); 558a175e481SJames Wright PetscCallMPI(MPI_Type_commit(&unit)); 559a175e481SJames Wright } 560a175e481SJames Wright 561a175e481SJames Wright PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM)); 562a175e481SJames Wright PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM)); 563a175e481SJames Wright 564ed331efdSJames Wright PetscCall(VecRestoreArrayRead(user_stats.Child_Stats_loc, &child_stats)); 565ed331efdSJames Wright PetscCall(VecRestoreArray(user_stats.Parent_Stats_loc, &parent_stats)); 566a175e481SJames Wright PetscCallMPI(MPI_Type_free(&unit)); 567a175e481SJames Wright 56878bbfb6fSJed Brown PetscReal solution_time; 56978bbfb6fSJed Brown PetscCall(DMGetOutputSequenceNumber(user_stats.dm, NULL, &solution_time)); 57078bbfb6fSJed Brown PetscReal summing_duration = solution_time - user->app_ctx->cont_time; 571ed331efdSJames Wright PetscCall(VecScale(user_stats.Parent_Stats_loc, 1 / (summing_duration * user_stats.span_width))); 572a175e481SJames Wright 573a175e481SJames Wright // L^2 projection with the parent_data 574ed331efdSJames Wright PetscCall(DMGetGlobalVector(user_stats.dm, &RHS)); 575ed331efdSJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(user_stats.Parent_Stats_loc, RHS, user_stats.op_proj_rhs_ctx)); 576a175e481SJames Wright 577ed331efdSJames Wright PetscCall(KSPSolve(user_stats.ksp, RHS, stats)); 578a175e481SJames Wright 579ed331efdSJames Wright PetscCall(DMRestoreGlobalVector(user_stats.dm, &RHS)); 5806dcea3beSJames Wright PetscCall(PetscLogStagePop()); 581ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 582a175e481SJames Wright } 583a175e481SJames Wright 584a175e481SJames Wright // TSMonitor for the statistics collection and processing 585f5452247SJames Wright PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) { 586a175e481SJames Wright User user = (User)ctx; 587a175e481SJames Wright Vec stats; 58878bbfb6fSJed Brown TSConvergedReason reason; 589b7d66439SJames Wright PetscInt collect_interval = user->app_ctx->turb_spanstats_collect_interval, viewer_interval = user->app_ctx->turb_spanstats_viewer_interval; 590a175e481SJames Wright PetscFunctionBeginUser; 59178bbfb6fSJed Brown PetscCall(TSGetConvergedReason(ts, &reason)); 592a175e481SJames Wright // Do not collect or process on the first step of the run (ie. on the initial condition) 593ee4ca9cbSJames Wright if (steps == user->app_ctx->cont_steps && reason == TS_CONVERGED_ITERATING) PetscFunctionReturn(PETSC_SUCCESS); 594a175e481SJames Wright 59527240365SJames Wright PetscBool run_processing_and_viewer = (steps % viewer_interval == 0 && viewer_interval != -1) || reason != TS_CONVERGED_ITERATING; 596a175e481SJames Wright 59727240365SJames Wright if (steps % collect_interval == 0 || run_processing_and_viewer) { 59827240365SJames Wright PetscCall(CollectStatistics(user, solution_time, Q)); 59927240365SJames Wright 60027240365SJames Wright if (run_processing_and_viewer) { 6018ed52730SJames Wright PetscCall(DMSetOutputSequenceNumber(user->spanstats.dm, steps, solution_time)); 60278bbfb6fSJed Brown PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats)); 60378bbfb6fSJed Brown PetscCall(ProcessStatistics(user, stats)); 6048fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_NONE) { 605b7d66439SJames Wright PetscCall(PetscViewerPushFormat(user->app_ctx->turb_spanstats_viewer, user->app_ctx->turb_spanstats_viewer_format)); 606b7d66439SJames Wright PetscCall(VecView(stats, user->app_ctx->turb_spanstats_viewer)); 607b7d66439SJames Wright PetscCall(PetscViewerPopFormat(user->app_ctx->turb_spanstats_viewer)); 6088fb33541SJames Wright } 6098fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_TURB_SPANSTATS && reason != TS_CONVERGED_ITERATING) { 6108fb33541SJames Wright PetscCall(RegressionTests_NS(user->app_ctx, stats)); 6118fb33541SJames Wright } 612b7d66439SJames Wright if (user->spanstats.do_mms_test && reason != TS_CONVERGED_ITERATING) { 613823a1283SJames Wright Vec error; 614823a1283SJames Wright PetscCall(VecDuplicate(stats, &error)); 6154021610dSJames Wright PetscCall(ApplyCeedOperatorGlobalToGlobal(stats, error, user->spanstats.mms_error_ctx)); 616823a1283SJames Wright PetscScalar error_sq = 0; 617823a1283SJames Wright PetscCall(VecSum(error, &error_sq)); 618823a1283SJames Wright PetscScalar l2_error = sqrt(error_sq); 619823a1283SJames Wright PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error)); 620823a1283SJames Wright } 621a175e481SJames Wright PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats)); 622522ee345SJames Wright } 62327240365SJames Wright } 624ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 62551ee423eSJames Wright } 626fd170fd0SJames Wright 627f5452247SJames Wright PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data) { 628fd170fd0SJames Wright PetscFunctionBeginUser; 629fd170fd0SJames Wright 630fd170fd0SJames Wright // -- CeedVectors 631ed331efdSJames Wright PetscCall(VecDestroy(&user->spanstats.Child_Stats_loc)); 632ed331efdSJames Wright PetscCall(VecDestroy(&user->spanstats.Parent_Stats_loc)); 633fd170fd0SJames Wright 634fd170fd0SJames Wright // -- CeedOperators 635ed331efdSJames Wright PetscCall(OperatorApplyContextDestroy(user->spanstats.op_stats_collect_ctx)); 636ed331efdSJames Wright PetscCall(OperatorApplyContextDestroy(user->spanstats.op_proj_rhs_ctx)); 6374021610dSJames Wright PetscCall(OperatorApplyContextDestroy(user->spanstats.mms_error_ctx)); 638fd170fd0SJames Wright 639fd170fd0SJames Wright // -- KSP 640fd170fd0SJames Wright PetscCall(KSPDestroy(&user->spanstats.ksp)); 641fd170fd0SJames Wright 642fd170fd0SJames Wright // -- SF 643fd170fd0SJames Wright PetscCall(PetscSFDestroy(&user->spanstats.sf)); 644fd170fd0SJames Wright 645fd170fd0SJames Wright // -- DM 646fd170fd0SJames Wright PetscCall(DMDestroy(&user->spanstats.dm)); 647fd170fd0SJames Wright 648ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 649fd170fd0SJames Wright } 650