1*5aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 251ee423eSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 351ee423eSJames Wright // 451ee423eSJames Wright // SPDX-License-Identifier: BSD-2-Clause 551ee423eSJames Wright // 651ee423eSJames Wright // This file is part of CEED: http://github.com/ceed 751ee423eSJames Wright /// @file 8701bc832SJames Wright /// Functions for setting up and performing spanwise-statistics collection 9701bc832SJames Wright /// 10701bc832SJames Wright /// "Parent" refers to the 2D plane on which statistics are collected *onto*. 11701bc832SJames Wright /// "Child" refers to the 3D domain where statistics are gathered *from*. 12701bc832SJames Wright /// Each quadrature point on the parent plane has several children in the child domain that it performs spanwise averaging with. 13a175e481SJames Wright 14a175e481SJames Wright #include "../qfunctions/turb_spanstats.h" 1551ee423eSJames Wright 1649aac155SJeremy L Thompson #include <ceed.h> 1749aac155SJeremy L Thompson #include <petscdmplex.h> 1849aac155SJeremy L Thompson #include <petscsf.h> 1949aac155SJeremy L Thompson 204021610dSJames Wright #include "../include/petsc_ops.h" 2151ee423eSJames Wright #include "../navierstokes.h" 2251ee423eSJames Wright 23269a910fSJames Wright typedef struct { 24269a910fSJames Wright CeedElemRestriction elem_restr_parent_x, elem_restr_parent_stats, elem_restr_parent_qd, elem_restr_parent_colloc, elem_restr_child_colloc; 25269a910fSJames Wright CeedBasis basis_x, basis_stats; 26269a910fSJames Wright CeedVector x_coord, q_data; 27269a910fSJames Wright } *SpanStatsSetupData; 28269a910fSJames Wright 29f5452247SJames Wright PetscErrorCode CreateStatsDM(User user, ProblemData *problem, PetscInt degree) { 303a4208e6SJames Wright user->spanstats.num_comp_stats = TURB_NUM_COMPONENTS; 31a175e481SJames Wright PetscReal domain_min[3], domain_max[3]; 3278837792SJames Wright PetscSection section; 334eed8630SJames Wright PetscLogStage stage_stats_setup; 341f595ac1SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)user->dm); 3551ee423eSJames Wright 36f17d818dSJames Wright PetscFunctionBeginUser; 374eed8630SJames Wright PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup)); 384eed8630SJames Wright if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup)); 394eed8630SJames Wright PetscCall(PetscLogStagePush(stage_stats_setup)); 404eed8630SJames Wright 41a175e481SJames Wright // Get spanwise length 42a175e481SJames Wright PetscCall(DMGetBoundingBox(user->dm, domain_min, domain_max)); 43a175e481SJames Wright user->spanstats.span_width = domain_max[2] - domain_min[1]; 44a175e481SJames Wright 45f967ad79SJames Wright { // Get DM from surface 461f595ac1SJames Wright DM parent_distributed_dm; 47c9198418SJames Wright PetscSF isoperiodicface; 4851ee423eSJames Wright DMLabel label; 491f595ac1SJames Wright PetscMPIInt size; 50c9198418SJames Wright 51c9198418SJames Wright PetscCall(DMPlexGetIsoperiodicFaceSF(user->dm, &isoperiodicface)); 52c9198418SJames Wright 53c9198418SJames Wright if (isoperiodicface) { 54c9198418SJames Wright PetscSF inv_isoperiodicface; 55c9198418SJames Wright PetscInt nleaves; 56c9198418SJames Wright const PetscInt *ilocal; 57c9198418SJames Wright 58c9198418SJames Wright PetscCall(PetscSFCreateInverseSF(isoperiodicface, &inv_isoperiodicface)); 59c9198418SJames Wright PetscCall(PetscSFGetGraph(inv_isoperiodicface, NULL, &nleaves, &ilocal, NULL)); 60c9198418SJames Wright PetscCall(DMCreateLabel(user->dm, "Periodic Face")); 61c9198418SJames Wright PetscCall(DMGetLabel(user->dm, "Periodic Face", &label)); 62c9198418SJames Wright for (PetscInt i = 0; i < nleaves; i++) { 63c9198418SJames Wright PetscCall(DMLabelSetValue(label, ilocal[i], 1)); 64c9198418SJames Wright } 65c9198418SJames Wright } else { 6651ee423eSJames Wright PetscCall(DMGetLabel(user->dm, "Face Sets", &label)); 67c9198418SJames Wright } 68c9198418SJames Wright 6951ee423eSJames Wright PetscCall(DMPlexLabelComplete(user->dm, label)); 70342228c7SAlex Lindsay PetscCall(DMPlexFilter(user->dm, label, 1, PETSC_FALSE, PETSC_FALSE, NULL, &user->spanstats.dm)); 7149a40c8aSKenneth E. Jansen PetscCall(DMSetCoordinateDisc(user->spanstats.dm, NULL, PETSC_TRUE)); // Ensure that a coordinate FE exists 721f595ac1SJames Wright 731f595ac1SJames Wright PetscCall(DMPlexDistribute(user->spanstats.dm, 0, NULL, &parent_distributed_dm)); 741f595ac1SJames Wright PetscCallMPI(MPI_Comm_size(comm, &size)); 751f595ac1SJames Wright if (parent_distributed_dm) { 761f595ac1SJames Wright PetscCall(DMDestroy(&user->spanstats.dm)); 771f595ac1SJames Wright user->spanstats.dm = parent_distributed_dm; 781f595ac1SJames Wright } else if (size > 1) { 791f595ac1SJames Wright PetscCall(PetscPrintf(comm, "WARNING: Turbulent spanwise statistics: parent DM could not be distributed accross %d ranks.\n", size)); 801f595ac1SJames Wright } 8151ee423eSJames Wright } 82d68a66c4SJames Wright { 83d68a66c4SJames Wright PetscBool is_simplex = PETSC_FALSE; 84d68a66c4SJames Wright PetscCall(DMPlexIsSimplex(user->spanstats.dm, &is_simplex)); 85d68a66c4SJames Wright PetscCheck(is_simplex != PETSC_TRUE, comm, PETSC_ERR_ARG_WRONGSTATE, "Spanwise Statistics is not implemented for non-tensor product grids"); 86d68a66c4SJames Wright } 8751ee423eSJames Wright 8851ee423eSJames Wright PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats")); 89b7d66439SJames Wright PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "turbulence_spanstats_")); 9051ee423eSJames Wright PetscCall(DMSetFromOptions(user->spanstats.dm)); 91b57b8e72SJames Wright PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view")); 92c9198418SJames Wright 9378837792SJames Wright // Create FE space for parent DM 94d68a66c4SJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, user->app_ctx->degree, 1, user->app_ctx->q_extra, 1, &user->spanstats.num_comp_stats, 95d68a66c4SJames Wright user->spanstats.dm)); 9651ee423eSJames Wright 9778837792SJames Wright // Create Section for data 9851ee423eSJames Wright PetscCall(DMGetLocalSection(user->spanstats.dm, §ion)); 9951ee423eSJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 1003a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY, "MeanDensity")); 1013a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE, "MeanPressure")); 1023a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_SQUARED, "MeanPressureSquared")); 1033a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_X, "MeanPressureVelocityX")); 1043a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Y, "MeanPressureVelocityY")); 1053a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Z, "MeanPressureVelocityZ")); 1063a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE, "MeanDensityTemperature")); 1073a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_X, "MeanDensityTemperatureFluxX")); 1083a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Y, "MeanDensityTemperatureFluxY")); 1093a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Z, "MeanDensityTemperatureFluxZ")); 1103a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_X, "MeanMomentumX")); 1113a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Y, "MeanMomentumY")); 1123a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Z, "MeanMomentumZ")); 1133a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XX, "MeanMomentumFluxXX")); 1143a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YY, "MeanMomentumFluxYY")); 1153a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_ZZ, "MeanMomentumFluxZZ")); 1163a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YZ, "MeanMomentumFluxYZ")); 1173a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XZ, "MeanMomentumFluxXZ")); 1183a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XY, "MeanMomentumFluxXY")); 1193a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_X, "MeanVelocityX")); 1203a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Y, "MeanVelocityY")); 1213a4208e6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Z, "MeanVelocityZ")); 12219706a06SJames Wright 1234eed8630SJames Wright PetscCall(PetscLogStagePop()); 124ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 12519706a06SJames Wright } 12619706a06SJames Wright 127701bc832SJames Wright /** @brief Create CeedElemRestriction for collocated data in component-major order. 128701bc832SJames Wright a. Sets the strides of the restriction to component-major order 129701bc832SJames Wright Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction. 130701bc832SJames Wright */ 131701bc832SJames Wright static PetscErrorCode CreateElemRestrColloc_CompMajor(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base, 132ed331efdSJames Wright CeedElemRestriction *elem_restr_collocated) { 1331737222fSJames Wright CeedInt num_elem_qpts, loc_num_elem; 1341737222fSJames Wright 135f17d818dSJames Wright PetscFunctionBeginUser; 136a424bcd0SJames Wright PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts)); 137a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem)); 1381737222fSJames Wright 1391737222fSJames Wright const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp}; 140a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides, 141a424bcd0SJames Wright elem_restr_collocated)); 142ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 1431737222fSJames Wright } 1441737222fSJames Wright 145a175e481SJames Wright // Get coordinates of quadrature points 146ed331efdSJames Wright PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, Vec *Qx_coords, 14719706a06SJames Wright PetscInt *total_nqpnts) { 148269a910fSJames Wright CeedElemRestriction elem_restr_qx; 14919706a06SJames Wright CeedQFunction qf_quad_coords; 15019706a06SJames Wright CeedOperator op_quad_coords; 151701bc832SJames Wright CeedInt num_comp_x; 152701bc832SJames Wright CeedSize l_vec_size; 153ed331efdSJames Wright OperatorApplyContext op_quad_coords_ctx; 15419706a06SJames Wright 155ed331efdSJames Wright PetscFunctionBeginUser; 156a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x)); 157701bc832SJames Wright PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx)); 158701bc832SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetLVectorSize(elem_restr_qx, &l_vec_size)); 159701bc832SJames Wright *total_nqpnts = l_vec_size / num_comp_x; 16019706a06SJames Wright 16119706a06SJames Wright // Create QFunction 162a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords)); 16319706a06SJames Wright 16419706a06SJames Wright // Create Operator 165a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords)); 166a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, x_coords)); 167356036faSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 16819706a06SJames Wright 169ed331efdSJames Wright PetscCall(CeedOperatorCreateLocalVecs(op_quad_coords, DMReturnVecType(dm), PetscObjectComm((PetscObject)dm), NULL, Qx_coords)); 170ed331efdSJames Wright PetscCall(OperatorApplyContextCreate(NULL, NULL, ceed, op_quad_coords, CEED_VECTOR_NONE, NULL, NULL, NULL, &op_quad_coords_ctx)); 17119706a06SJames Wright 172ed331efdSJames Wright PetscCall(ApplyCeedOperatorLocalToLocal(NULL, *Qx_coords, op_quad_coords_ctx)); 173ed331efdSJames Wright 174ed331efdSJames Wright PetscCall(OperatorApplyContextDestroy(op_quad_coords_ctx)); 175a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qx)); 176a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_quad_coords)); 177a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_quad_coords)); 178ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 17919706a06SJames Wright } 18019706a06SJames Wright 181269a910fSJames Wright PetscErrorCode SpanStatsSetupDataCreate(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem, SpanStatsSetupData *stats_data) { 182269a910fSJames Wright DM dm = user->spanstats.dm; 183457e73b2SJames Wright PetscInt dim; 184bb85d312SJames Wright CeedInt num_comp_x, num_comp_stats = user->spanstats.num_comp_stats; 185269a910fSJames Wright Vec X_loc; 186bb85d312SJames Wright DMLabel domain_label = NULL; 187bb85d312SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 188269a910fSJames Wright 189f17d818dSJames Wright PetscFunctionBeginUser; 190269a910fSJames Wright PetscCall(PetscNew(stats_data)); 191269a910fSJames Wright 192269a910fSJames Wright PetscCall(DMGetDimension(dm, &dim)); 193bb85d312SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &(*stats_data)->elem_restr_parent_stats)); 194bb85d312SJames Wright PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &(*stats_data)->elem_restr_parent_x)); 195bb85d312SJames Wright PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, problem->q_data_size_sur, 196bb85d312SJames Wright &(*stats_data)->elem_restr_parent_qd)); 197a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents((*stats_data)->elem_restr_parent_x, &num_comp_x)); 198a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_x, &(*stats_data)->x_coord, NULL)); 199a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_qd, &(*stats_data)->q_data, NULL)); 200269a910fSJames Wright 2010814d5a7SKenneth E. Jansen { 2020814d5a7SKenneth E. Jansen DM dm_coord; 2030814d5a7SKenneth E. Jansen PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 204bb85d312SJames Wright PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &(*stats_data)->basis_x)); 205bb85d312SJames Wright PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &(*stats_data)->basis_stats)); 2060814d5a7SKenneth E. Jansen } 207269a910fSJames Wright 208701bc832SJames Wright PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_stats, (*stats_data)->basis_stats, (*stats_data)->elem_restr_parent_stats, 209ed331efdSJames Wright &(*stats_data)->elem_restr_parent_colloc)); 210701bc832SJames Wright PetscCall( 211701bc832SJames Wright CreateElemRestrColloc_CompMajor(ceed, num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q, &(*stats_data)->elem_restr_child_colloc)); 212269a910fSJames Wright 213269a910fSJames Wright { // -- Copy DM coordinates into CeedVector 214269a910fSJames Wright DM cdm; 215269a910fSJames Wright PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 216269a910fSJames Wright if (cdm) { 217269a910fSJames Wright PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 218269a910fSJames Wright } else { 219269a910fSJames Wright PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 220269a910fSJames Wright } 221269a910fSJames Wright } 222269a910fSJames Wright PetscCall(VecScale(X_loc, problem->dm_scale)); 223d0593705SJames Wright PetscCall(VecCopyPetscToCeed(X_loc, (*stats_data)->x_coord)); 224ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 225269a910fSJames Wright } 226269a910fSJames Wright 227269a910fSJames Wright PetscErrorCode SpanStatsSetupDataDestroy(SpanStatsSetupData data) { 228a424bcd0SJames Wright Ceed ceed; 229a424bcd0SJames Wright 230269a910fSJames Wright PetscFunctionBeginUser; 231a424bcd0SJames Wright PetscCall(CeedElemRestrictionGetCeed(data->elem_restr_parent_x, &ceed)); 232a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_x)); 233a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_stats)); 234a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_qd)); 235a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_colloc)); 236a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_child_colloc)); 237269a910fSJames Wright 238a424bcd0SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_x)); 239a424bcd0SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_stats)); 240269a910fSJames Wright 241a424bcd0SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&data->x_coord)); 242a424bcd0SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&data->q_data)); 243269a910fSJames Wright 244269a910fSJames Wright PetscCall(PetscFree(data)); 245ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 246269a910fSJames Wright } 247269a910fSJames Wright 248a175e481SJames Wright // Create PetscSF for child-to-parent communication 249269a910fSJames Wright PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, SpanStatsSetupData stats_data, DM parentdm, DM childdm, PetscSF *statssf) { 250457e73b2SJames Wright PetscInt child_num_qpnts, parent_num_qpnts; 251457e73b2SJames Wright CeedInt num_comp_x; 252ed331efdSJames Wright Vec Child_qx_coords, Parent_qx_coords; 25319706a06SJames Wright 254ed331efdSJames Wright PetscFunctionBeginUser; 255269a910fSJames Wright PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)childdm), statssf)); 256269a910fSJames Wright 25719706a06SJames Wright // Assume that child and parent have the same number of components 258a424bcd0SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x)); 25919706a06SJames Wright const PetscInt num_comp_sf = num_comp_x - 1; // Number of coord components used in the creation of the SF 26019706a06SJames Wright 261ed331efdSJames Wright // Get quad_coords for child and parent DM 262ed331efdSJames Wright PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &Child_qx_coords, &child_num_qpnts)); 263ed331efdSJames Wright PetscCall(GetQuadratureCoords(ceed, parentdm, stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord, &Parent_qx_coords, 264269a910fSJames Wright &parent_num_qpnts)); 26519706a06SJames Wright 266ed331efdSJames Wright { // Remove z component of coordinates for matching 26719706a06SJames Wright const PetscReal *child_quad_coords, *parent_quad_coords; 268ed331efdSJames Wright PetscReal *child_coords, *parent_coords; 26919706a06SJames Wright 270ed331efdSJames Wright PetscCall(VecGetArrayRead(Child_qx_coords, &child_quad_coords)); 271ed331efdSJames Wright PetscCall(VecGetArrayRead(Parent_qx_coords, &parent_quad_coords)); 27219706a06SJames Wright 27319706a06SJames Wright PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords)); 27419706a06SJames Wright for (int i = 0; i < child_num_qpnts; i++) { 27519706a06SJames Wright child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x]; 27619706a06SJames Wright child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x]; 27719706a06SJames Wright } 27819706a06SJames Wright for (int i = 0; i < parent_num_qpnts; i++) { 27919706a06SJames Wright parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x]; 28019706a06SJames Wright parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x]; 28119706a06SJames Wright } 282ed331efdSJames Wright PetscCall(VecRestoreArrayRead(Child_qx_coords, &child_quad_coords)); 283ed331efdSJames Wright PetscCall(VecRestoreArrayRead(Parent_qx_coords, &parent_quad_coords)); 28419706a06SJames Wright 285269a910fSJames Wright PetscCall(PetscSFSetGraphFromCoordinates(*statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords)); 286ed331efdSJames Wright PetscCall(PetscFree2(child_coords, parent_coords)); 287ed331efdSJames Wright } 28819706a06SJames Wright 289269a910fSJames Wright PetscCall(PetscSFViewFromOptions(*statssf, NULL, "-spanstats_sf_view")); 290a175e481SJames Wright 291ed331efdSJames Wright PetscCall(VecDestroy(&Child_qx_coords)); 292ed331efdSJames Wright PetscCall(VecDestroy(&Parent_qx_coords)); 293ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 29419706a06SJames Wright } 29519706a06SJames Wright 296ed331efdSJames Wright // @brief Setup RHS and LHS for L^2 projection of statistics 297269a910fSJames Wright PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) { 298ed331efdSJames Wright CeedOperator op_mass, op_setup_sur, op_proj_rhs; 299269a910fSJames Wright CeedQFunction qf_mass, qf_stats_proj; 300269a910fSJames Wright CeedInt q_data_size, num_comp_stats = user->spanstats.num_comp_stats; 301f967ad79SJames Wright MPI_Comm comm = PetscObjectComm((PetscObject)user->spanstats.dm); 302a175e481SJames Wright 303ed331efdSJames Wright PetscFunctionBeginUser; 304ed331efdSJames Wright // -- Create Operator for RHS of L^2 projection of statistics 305269a910fSJames Wright // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function. 306269a910fSJames Wright // Therefore, an Identity QF is sufficient 307a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &qf_stats_proj)); 308269a910fSJames Wright 309a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_proj, NULL, NULL, &op_proj_rhs)); 310356036faSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "input", stats_data->elem_restr_parent_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 311a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "output", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 312269a910fSJames Wright 313ed331efdSJames Wright PetscCall(OperatorApplyContextCreate(NULL, user->spanstats.dm, ceed, op_proj_rhs, NULL, NULL, NULL, NULL, &user->spanstats.op_proj_rhs_ctx)); 314ed331efdSJames Wright PetscCall(CeedOperatorCreateLocalVecs(op_proj_rhs, DMReturnVecType(user->spanstats.dm), comm, &user->spanstats.Parent_Stats_loc, NULL)); 315ed331efdSJames Wright 316ed331efdSJames Wright // -- Setup LHS of L^2 projection 317269a910fSJames Wright // Get q_data for mass matrix operator 318a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur)); 319a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", stats_data->elem_restr_parent_x, stats_data->basis_x, CEED_VECTOR_ACTIVE)); 320a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, stats_data->basis_x, CEED_VECTOR_NONE)); 321356036faSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 322a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, stats_data->x_coord, stats_data->q_data, CEED_REQUEST_IMMEDIATE)); 323269a910fSJames Wright 324a175e481SJames Wright // CEED Restriction 325a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size)); 326a175e481SJames Wright 327a175e481SJames Wright // Create Mass CeedOperator 328269a910fSJames Wright PetscCall(CreateMassQFunction(ceed, num_comp_stats, q_data_size, &qf_mass)); 329a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 330a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 331356036faSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_NONE, stats_data->q_data)); 332a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 333a175e481SJames Wright 334f967ad79SJames Wright { // Setup KSP for L^2 projection 3357f2a9303SJames Wright Mat mat_mass; 336a175e481SJames Wright KSP ksp; 337a175e481SJames Wright 3387f2a9303SJames Wright PetscCall(MatCeedCreate(user->spanstats.dm, user->spanstats.dm, op_mass, NULL, &mat_mass)); 339a175e481SJames Wright 340f967ad79SJames Wright PetscCall(KSPCreate(comm, &ksp)); 341b7d66439SJames Wright PetscCall(KSPSetOptionsPrefix(ksp, "turbulence_spanstats_")); 342a175e481SJames Wright { 343a175e481SJames Wright PC pc; 344a175e481SJames Wright PetscCall(KSPGetPC(ksp, &pc)); 345a175e481SJames Wright PetscCall(PCSetType(pc, PCJACOBI)); 346a175e481SJames Wright PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL)); 347a175e481SJames Wright PetscCall(KSPSetType(ksp, KSPCG)); 348a175e481SJames Wright PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL)); 349a175e481SJames Wright PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 350a175e481SJames Wright } 3517f2a9303SJames Wright PetscCall(KSPSetFromOptions_WithMatCeed(ksp, mat_mass)); 352a175e481SJames Wright user->spanstats.ksp = ksp; 353cde30410SJames Wright PetscCall(MatDestroy(&mat_mass)); 354a175e481SJames Wright } 355a175e481SJames Wright 356a175e481SJames Wright // Cleanup 357a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 358a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_proj)); 359a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 360a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur)); 361a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_proj_rhs)); 362ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 363a175e481SJames Wright } 364a175e481SJames Wright 365269a910fSJames Wright // Create CeedOperator for statistics collection 366269a910fSJames Wright PetscErrorCode CreateStatisticCollectionOperator(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data, ProblemData *problem) { 367a175e481SJames Wright CeedInt num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q; 368495b9769SJames Wright Turbulence_SpanStatsContext collect_ctx; 369495b9769SJames Wright NewtonianIdealGasContext newtonian_ig_ctx; 370495b9769SJames Wright CeedQFunctionContext collect_context; 371269a910fSJames Wright CeedQFunction qf_stats_collect; 372ed331efdSJames Wright CeedOperator op_stats_collect; 373ed331efdSJames Wright 374a175e481SJames Wright PetscFunctionBeginUser; 375a424bcd0SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q)); 376a175e481SJames Wright 377a175e481SJames Wright // Create Operator for statistics collection 378a175e481SJames Wright switch (user->phys->state_var) { 379a175e481SJames Wright case STATEVAR_PRIMITIVE: 380a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &qf_stats_collect)); 381a175e481SJames Wright break; 382a175e481SJames Wright case STATEVAR_CONSERVATIVE: 383a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &qf_stats_collect)); 384a175e481SJames Wright break; 385269a910fSJames Wright default: 386269a910fSJames Wright SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No statisics collection available for chosen state variable"); 387a175e481SJames Wright } 388a175e481SJames Wright 389b7d66439SJames Wright if (user->spanstats.do_mms_test) { 390a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect)); 391a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest, ChildStatsCollectionMMSTest_loc, &qf_stats_collect)); 392a175e481SJames Wright } 393a175e481SJames Wright 394269a910fSJames Wright { // Setup Collection Context 395495b9769SJames Wright PetscCall(PetscNew(&collect_ctx)); 396a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx)); 397495b9769SJames Wright collect_ctx->gas = *newtonian_ig_ctx; 398495b9769SJames Wright 399a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &collect_context)); 400a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(collect_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx)); 401a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(collect_context, CEED_MEM_HOST, FreeContextPetsc)); 402495b9769SJames Wright 403a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(collect_context, "solution time", 404a424bcd0SJames Wright offsetof(struct Turbulence_SpanStatsContext_, solution_time), 1, "Current solution time")); 405a424bcd0SJames Wright PetscCallCeed( 406a424bcd0SJames Wright ceed, CeedQFunctionContextRegisterDouble(collect_context, "previous time", offsetof(struct Turbulence_SpanStatsContext_, previous_time), 1, 407a424bcd0SJames Wright "Previous time statistics collection was done")); 408495b9769SJames Wright 409a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx)); 410495b9769SJames Wright } 411495b9769SJames Wright 412a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_stats_collect, collect_context)); 413a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&collect_context)); 414a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP)); 415a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE)); 416a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP)); 417a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE)); 418a175e481SJames Wright 419a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_collect, NULL, NULL, &op_stats_collect)); 420a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 421356036faSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 422a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord)); 423356036faSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "v", stats_data->elem_restr_child_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 424a175e481SJames Wright 425a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "solution time", &user->spanstats.solution_time_label)); 426a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "previous time", &user->spanstats.previous_time_label)); 427ed331efdSJames Wright 428ed331efdSJames Wright PetscCall(OperatorApplyContextCreate(user->dm, user->spanstats.dm, user->ceed, op_stats_collect, user->q_ceed, NULL, NULL, NULL, 429ed331efdSJames Wright &user->spanstats.op_stats_collect_ctx)); 430ed331efdSJames Wright 431ed331efdSJames Wright PetscCall(CeedOperatorCreateLocalVecs(op_stats_collect, DMReturnVecType(user->spanstats.dm), PetscObjectComm((PetscObject)user->spanstats.dm), NULL, 432ed331efdSJames Wright &user->spanstats.Child_Stats_loc)); 433ed331efdSJames Wright PetscCall(VecZeroEntries(user->spanstats.Child_Stats_loc)); 434495b9769SJames Wright 435a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect)); 436a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_stats_collect)); 437ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 438a175e481SJames Wright } 439a175e481SJames Wright 440b7d66439SJames Wright // Creates operator for calculating error of method of manufactured solution (MMS) test 441269a910fSJames Wright PetscErrorCode SetupMMSErrorChecking(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) { 442269a910fSJames Wright CeedInt num_comp_stats = user->spanstats.num_comp_stats, num_comp_x, q_data_size; 443823a1283SJames Wright CeedQFunction qf_error; 444823a1283SJames Wright CeedOperator op_error; 445823a1283SJames Wright CeedVector x_ceed, y_ceed; 446823a1283SJames Wright 447f17d818dSJames Wright PetscFunctionBeginUser; 448a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size)); 449a424bcd0SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(stats_data->basis_x, &num_comp_x)); 450823a1283SJames Wright 451a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest_Error, ChildStatsCollectionMMSTest_Error_loc, &qf_error)); 452a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP)); 453a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE)); 454a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP)); 455a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP)); 456823a1283SJames Wright 457a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error)); 458a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_error, "q", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 459356036faSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_error, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_NONE, stats_data->q_data)); 460a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_error, "x", stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord)); 461a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_error, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 462823a1283SJames Wright 463a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &x_ceed, NULL)); 464a424bcd0SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &y_ceed, NULL)); 4654021610dSJames Wright PetscCall(OperatorApplyContextCreate(user->spanstats.dm, user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, NULL, 466ea168fcfSJames Wright &user->spanstats.mms_error_ctx)); 467823a1283SJames Wright 468a424bcd0SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_error)); 469a424bcd0SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_error)); 470a424bcd0SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&x_ceed)); 471a424bcd0SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&y_ceed)); 472ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 473823a1283SJames Wright } 474823a1283SJames Wright 475a175e481SJames Wright // Setup for statistics collection 476f5452247SJames Wright PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) { 477269a910fSJames Wright SpanStatsSetupData stats_data; 4784eed8630SJames Wright PetscLogStage stage_stats_setup; 479f5452247SJames Wright 48019706a06SJames Wright PetscFunctionBeginUser; 4814eed8630SJames Wright PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup)); 4824eed8630SJames Wright if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup)); 4834eed8630SJames Wright PetscCall(PetscLogStagePush(stage_stats_setup)); 4844eed8630SJames Wright 485f5452247SJames Wright // Create parent DM 486f5452247SJames Wright PetscCall(CreateStatsDM(user, problem, user->app_ctx->degree)); 487f5452247SJames Wright 488269a910fSJames Wright // Create necessary CeedObjects for setting up statistics 489269a910fSJames Wright PetscCall(SpanStatsSetupDataCreate(ceed, user, ceed_data, problem, &stats_data)); 4901737222fSJames Wright 49119706a06SJames Wright // Create SF for communicating child data back their respective parents 492269a910fSJames Wright PetscCall(CreateStatsSF(ceed, ceed_data, stats_data, user->dm, user->spanstats.dm, &user->spanstats.sf)); 49351ee423eSJames Wright 494a175e481SJames Wright // Create CeedOperators for statistics collection 495269a910fSJames Wright PetscCall(CreateStatisticCollectionOperator(ceed, user, ceed_data, stats_data, problem)); 496a175e481SJames Wright 497a175e481SJames Wright // Setup KSP and Mat for L^2 projection of statistics 498269a910fSJames Wright PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data, stats_data)); 499a175e481SJames Wright 500b7d66439SJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_turbulence_spanstats_mms", &user->spanstats.do_mms_test, NULL)); 501b7d66439SJames Wright if (user->spanstats.do_mms_test) { 502269a910fSJames Wright PetscCall(SetupMMSErrorChecking(ceed, user, ceed_data, stats_data)); 503823a1283SJames Wright } 504823a1283SJames Wright 505f967ad79SJames Wright { // Setup stats viewer with prefix 5068ed52730SJames Wright PetscViewerType viewer_type; 507b7d66439SJames Wright PetscCall(PetscViewerGetType(user->app_ctx->turb_spanstats_viewer, &viewer_type)); 508b7d66439SJames Wright PetscCall(PetscOptionsSetValue(NULL, "-ts_monitor_turbulence_spanstats_viewer_type", viewer_type)); 5098ed52730SJames Wright 510b7d66439SJames Wright PetscCall(PetscViewerSetOptionsPrefix(user->app_ctx->turb_spanstats_viewer, "ts_monitor_turbulence_spanstats_")); 511b7d66439SJames Wright PetscCall(PetscViewerSetFromOptions(user->app_ctx->turb_spanstats_viewer)); 5128ed52730SJames Wright } 5134eed8630SJames Wright 514269a910fSJames Wright PetscCall(SpanStatsSetupDataDestroy(stats_data)); 5154eed8630SJames Wright PetscCall(PetscLogStagePop()); 516ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 517a175e481SJames Wright } 518a175e481SJames Wright 519a175e481SJames Wright // Collect statistics based on the solution Q 520a175e481SJames Wright PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) { 521cbef7084SJames Wright SpanStatsData user_stats = user->spanstats; 522a175e481SJames Wright 523ed331efdSJames Wright PetscFunctionBeginUser; 5246dcea3beSJames Wright PetscLogStage stage_stats_collect; 5256dcea3beSJames Wright PetscCall(PetscLogStageGetId("Stats Collect", &stage_stats_collect)); 5266dcea3beSJames Wright if (stage_stats_collect == -1) PetscCall(PetscLogStageRegister("Stats Collect", &stage_stats_collect)); 5276dcea3beSJames Wright PetscCall(PetscLogStagePush(stage_stats_collect)); 5286dcea3beSJames Wright 529a0b9a424SJames Wright PetscCall(UpdateBoundaryValues(user, user->Q_loc, solution_time)); 530a424bcd0SJames Wright PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.solution_time_label, &solution_time)); 531a0b9a424SJames Wright PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, user->Q_loc)); 532ed331efdSJames Wright PetscCall(ApplyAddCeedOperatorLocalToLocal(user->Q_loc, user_stats.Child_Stats_loc, user_stats.op_stats_collect_ctx)); 533a175e481SJames Wright 534a424bcd0SJames Wright PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.previous_time_label, &solution_time)); 535a175e481SJames Wright 5366dcea3beSJames Wright PetscCall(PetscLogStagePop()); 537ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 538a175e481SJames Wright } 539a175e481SJames Wright 540a175e481SJames Wright // Process the child statistics into parent statistics and project them onto stats 54178bbfb6fSJed Brown PetscErrorCode ProcessStatistics(User user, Vec stats) { 542cbef7084SJames Wright SpanStatsData user_stats = user->spanstats; 543a175e481SJames Wright const PetscScalar *child_stats; 544a175e481SJames Wright PetscScalar *parent_stats; 545a175e481SJames Wright MPI_Datatype unit; 546ed331efdSJames Wright Vec RHS; 547a175e481SJames Wright 548ed331efdSJames Wright PetscFunctionBeginUser; 5496dcea3beSJames Wright PetscLogStage stage_stats_process; 5506dcea3beSJames Wright PetscCall(PetscLogStageGetId("Stats Process", &stage_stats_process)); 5516dcea3beSJames Wright if (stage_stats_process == -1) PetscCall(PetscLogStageRegister("Stats Process", &stage_stats_process)); 5526dcea3beSJames Wright PetscCall(PetscLogStagePush(stage_stats_process)); 5536dcea3beSJames Wright 554ed331efdSJames Wright PetscCall(VecZeroEntries(user_stats.Parent_Stats_loc)); 555a175e481SJames Wright 556ed331efdSJames Wright PetscCall(VecGetArrayRead(user_stats.Child_Stats_loc, &child_stats)); 557ed331efdSJames Wright PetscCall(VecGetArray(user_stats.Parent_Stats_loc, &parent_stats)); 558a175e481SJames Wright 559a175e481SJames Wright if (user_stats.num_comp_stats == 1) unit = MPIU_REAL; 560a175e481SJames Wright else { 561a175e481SJames Wright PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit)); 562a175e481SJames Wright PetscCallMPI(MPI_Type_commit(&unit)); 563a175e481SJames Wright } 564a175e481SJames Wright 565a175e481SJames Wright PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM)); 566a175e481SJames Wright PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM)); 567a175e481SJames Wright 568ed331efdSJames Wright PetscCall(VecRestoreArrayRead(user_stats.Child_Stats_loc, &child_stats)); 569ed331efdSJames Wright PetscCall(VecRestoreArray(user_stats.Parent_Stats_loc, &parent_stats)); 570a175e481SJames Wright PetscCallMPI(MPI_Type_free(&unit)); 571a175e481SJames Wright 57278bbfb6fSJed Brown PetscReal solution_time; 57378bbfb6fSJed Brown PetscCall(DMGetOutputSequenceNumber(user_stats.dm, NULL, &solution_time)); 57478bbfb6fSJed Brown PetscReal summing_duration = solution_time - user->app_ctx->cont_time; 575ed331efdSJames Wright PetscCall(VecScale(user_stats.Parent_Stats_loc, 1 / (summing_duration * user_stats.span_width))); 576a175e481SJames Wright 577a175e481SJames Wright // L^2 projection with the parent_data 578ed331efdSJames Wright PetscCall(DMGetGlobalVector(user_stats.dm, &RHS)); 579ed331efdSJames Wright PetscCall(ApplyCeedOperatorLocalToGlobal(user_stats.Parent_Stats_loc, RHS, user_stats.op_proj_rhs_ctx)); 580a175e481SJames Wright 581ed331efdSJames Wright PetscCall(KSPSolve(user_stats.ksp, RHS, stats)); 582a175e481SJames Wright 583ed331efdSJames Wright PetscCall(DMRestoreGlobalVector(user_stats.dm, &RHS)); 5846dcea3beSJames Wright PetscCall(PetscLogStagePop()); 585ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 586a175e481SJames Wright } 587a175e481SJames Wright 588a175e481SJames Wright // TSMonitor for the statistics collection and processing 589f5452247SJames Wright PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) { 590a175e481SJames Wright User user = (User)ctx; 591a175e481SJames Wright Vec stats; 59278bbfb6fSJed Brown TSConvergedReason reason; 593b7d66439SJames Wright PetscInt collect_interval = user->app_ctx->turb_spanstats_collect_interval, viewer_interval = user->app_ctx->turb_spanstats_viewer_interval; 594f17d818dSJames Wright 595a175e481SJames Wright PetscFunctionBeginUser; 59678bbfb6fSJed Brown PetscCall(TSGetConvergedReason(ts, &reason)); 597a175e481SJames Wright // Do not collect or process on the first step of the run (ie. on the initial condition) 5986852f6f6SJames Wright if (steps == user->app_ctx->cont_steps) PetscFunctionReturn(PETSC_SUCCESS); 599a175e481SJames Wright 60027240365SJames Wright PetscBool run_processing_and_viewer = (steps % viewer_interval == 0 && viewer_interval != -1) || reason != TS_CONVERGED_ITERATING; 601a175e481SJames Wright 60227240365SJames Wright if (steps % collect_interval == 0 || run_processing_and_viewer) { 60327240365SJames Wright PetscCall(CollectStatistics(user, solution_time, Q)); 60427240365SJames Wright 60527240365SJames Wright if (run_processing_and_viewer) { 6068ed52730SJames Wright PetscCall(DMSetOutputSequenceNumber(user->spanstats.dm, steps, solution_time)); 60778bbfb6fSJed Brown PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats)); 60878bbfb6fSJed Brown PetscCall(ProcessStatistics(user, stats)); 6098fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_NONE) { 610b7d66439SJames Wright PetscCall(PetscViewerPushFormat(user->app_ctx->turb_spanstats_viewer, user->app_ctx->turb_spanstats_viewer_format)); 611b7d66439SJames Wright PetscCall(VecView(stats, user->app_ctx->turb_spanstats_viewer)); 612b7d66439SJames Wright PetscCall(PetscViewerPopFormat(user->app_ctx->turb_spanstats_viewer)); 6138fb33541SJames Wright } 6148fb33541SJames Wright if (user->app_ctx->test_type == TESTTYPE_TURB_SPANSTATS && reason != TS_CONVERGED_ITERATING) { 6153c9e7ad1SJames Wright PetscCall(RegressionTest(user->app_ctx, stats)); 6168fb33541SJames Wright } 617b7d66439SJames Wright if (user->spanstats.do_mms_test && reason != TS_CONVERGED_ITERATING) { 618823a1283SJames Wright Vec error; 619823a1283SJames Wright PetscCall(VecDuplicate(stats, &error)); 6204021610dSJames Wright PetscCall(ApplyCeedOperatorGlobalToGlobal(stats, error, user->spanstats.mms_error_ctx)); 621823a1283SJames Wright PetscScalar error_sq = 0; 622823a1283SJames Wright PetscCall(VecSum(error, &error_sq)); 623823a1283SJames Wright PetscScalar l2_error = sqrt(error_sq); 624823a1283SJames Wright PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error)); 625823a1283SJames Wright } 626a175e481SJames Wright PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats)); 627522ee345SJames Wright } 62827240365SJames Wright } 629ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 63051ee423eSJames Wright } 631fd170fd0SJames Wright 632f5452247SJames Wright PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data) { 633fd170fd0SJames Wright PetscFunctionBeginUser; 634ed331efdSJames Wright PetscCall(VecDestroy(&user->spanstats.Child_Stats_loc)); 635ed331efdSJames Wright PetscCall(VecDestroy(&user->spanstats.Parent_Stats_loc)); 636fd170fd0SJames Wright 637ed331efdSJames Wright PetscCall(OperatorApplyContextDestroy(user->spanstats.op_stats_collect_ctx)); 638ed331efdSJames Wright PetscCall(OperatorApplyContextDestroy(user->spanstats.op_proj_rhs_ctx)); 6394021610dSJames Wright PetscCall(OperatorApplyContextDestroy(user->spanstats.mms_error_ctx)); 640fd170fd0SJames Wright 641fd170fd0SJames Wright PetscCall(KSPDestroy(&user->spanstats.ksp)); 642fd170fd0SJames Wright PetscCall(PetscSFDestroy(&user->spanstats.sf)); 643fd170fd0SJames Wright PetscCall(DMDestroy(&user->spanstats.dm)); 644ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 645fd170fd0SJames Wright } 646