1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 /// @file 8 /// Functions for setting up and performing spanwise-statistics collection 9 /// 10 /// "Parent" refers to the 2D plane on which statistics are collected *onto*. 11 /// "Child" refers to the 3D domain where statistics are gathered *from*. 12 /// Each quadrature point on the parent plane has several children in the child domain that it performs spanwise averaging with. 13 14 #include "../qfunctions/turb_spanstats.h" 15 16 #include <ceed.h> 17 #include <petscdmplex.h> 18 #include <petscsf.h> 19 20 #include "../include/petsc_ops.h" 21 #include "../navierstokes.h" 22 23 typedef struct { 24 CeedElemRestriction elem_restr_parent_x, elem_restr_parent_stats, elem_restr_parent_qd, elem_restr_parent_colloc, elem_restr_child_colloc; 25 CeedBasis basis_x, basis_stats; 26 CeedVector x_coord, q_data; 27 } *SpanStatsSetupData; 28 29 PetscErrorCode CreateStatsDM(User user, ProblemData problem, PetscInt degree) { 30 user->spanstats.num_comp_stats = TURB_NUM_COMPONENTS; 31 PetscReal domain_min[3], domain_max[3]; 32 PetscSection section; 33 PetscLogStage stage_stats_setup; 34 MPI_Comm comm = PetscObjectComm((PetscObject)user->dm); 35 36 PetscFunctionBeginUser; 37 PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup)); 38 if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup)); 39 PetscCall(PetscLogStagePush(stage_stats_setup)); 40 41 // Get spanwise length 42 PetscCall(DMGetBoundingBox(user->dm, domain_min, domain_max)); 43 user->spanstats.span_width = domain_max[2] - domain_min[1]; 44 45 { // Get DM from surface 46 DM parent_distributed_dm; 47 const PetscSF *isoperiodicface; 48 PetscInt num_isoperiodicface; 49 DMLabel label; 50 PetscMPIInt size; 51 52 PetscCall(DMPlexGetIsoperiodicFaceSF(user->dm, &num_isoperiodicface, &isoperiodicface)); 53 54 if (isoperiodicface) { 55 PetscSF inv_isoperiodicface; 56 PetscInt nleaves, isoperiodicface_index = -1; 57 const PetscInt *ilocal; 58 59 PetscCall(PetscOptionsGetInt(NULL, NULL, "-ts_monitor_turbulence_spanstats_isoperiodic_facesf", &isoperiodicface_index, NULL)); 60 isoperiodicface_index = isoperiodicface_index == -1 ? num_isoperiodicface - 1 : isoperiodicface_index; 61 PetscCall(PetscSFCreateInverseSF(isoperiodicface[isoperiodicface_index], &inv_isoperiodicface)); 62 PetscCall(PetscSFGetGraph(inv_isoperiodicface, NULL, &nleaves, &ilocal, NULL)); 63 PetscCall(DMCreateLabel(user->dm, "Periodic Face")); 64 PetscCall(DMGetLabel(user->dm, "Periodic Face", &label)); 65 for (PetscInt i = 0; i < nleaves; i++) { 66 PetscCall(DMLabelSetValue(label, ilocal[i], 1)); 67 } 68 } else { 69 PetscCall(DMGetLabel(user->dm, "Face Sets", &label)); 70 } 71 72 PetscCall(DMPlexLabelComplete(user->dm, label)); 73 PetscCall(DMPlexFilter(user->dm, label, 1, PETSC_FALSE, PETSC_FALSE, NULL, &user->spanstats.dm)); 74 PetscCall(DMSetCoordinateDisc(user->spanstats.dm, NULL, PETSC_TRUE)); // Ensure that a coordinate FE exists 75 76 PetscCall(DMPlexDistribute(user->spanstats.dm, 0, NULL, &parent_distributed_dm)); 77 PetscCallMPI(MPI_Comm_size(comm, &size)); 78 if (parent_distributed_dm) { 79 PetscCall(DMDestroy(&user->spanstats.dm)); 80 user->spanstats.dm = parent_distributed_dm; 81 } else if (size > 1) { 82 PetscCall(PetscPrintf(comm, "WARNING: Turbulent spanwise statistics: parent DM could not be distributed accross %d ranks.\n", size)); 83 } 84 } 85 { 86 PetscBool is_simplex = PETSC_FALSE; 87 PetscCall(DMPlexIsSimplex(user->spanstats.dm, &is_simplex)); 88 PetscCheck(is_simplex != PETSC_TRUE, comm, PETSC_ERR_ARG_WRONGSTATE, "Spanwise Statistics is not implemented for non-tensor product grids"); 89 } 90 91 PetscCall(PetscObjectSetName((PetscObject)user->spanstats.dm, "Spanwise_Stats")); 92 PetscCall(DMSetOptionsPrefix(user->spanstats.dm, "turbulence_spanstats_")); 93 PetscCall(DMSetFromOptions(user->spanstats.dm)); 94 PetscCall(DMViewFromOptions(user->spanstats.dm, NULL, "-dm_view")); 95 96 // Create FE space for parent DM 97 PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, user->app_ctx->degree, 1, user->app_ctx->q_extra, 1, &user->spanstats.num_comp_stats, 98 user->spanstats.dm)); 99 100 // Create Section for data 101 PetscCall(DMGetLocalSection(user->spanstats.dm, §ion)); 102 PetscCall(PetscSectionSetFieldName(section, 0, "")); 103 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY, "MeanDensity")); 104 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE, "MeanPressure")); 105 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_SQUARED, "MeanPressureSquared")); 106 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_X, "MeanPressureVelocityX")); 107 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Y, "MeanPressureVelocityY")); 108 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Z, "MeanPressureVelocityZ")); 109 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE, "MeanDensityTemperature")); 110 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_X, "MeanDensityTemperatureFluxX")); 111 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Y, "MeanDensityTemperatureFluxY")); 112 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Z, "MeanDensityTemperatureFluxZ")); 113 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_X, "MeanMomentumX")); 114 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Y, "MeanMomentumY")); 115 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Z, "MeanMomentumZ")); 116 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XX, "MeanMomentumFluxXX")); 117 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YY, "MeanMomentumFluxYY")); 118 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_ZZ, "MeanMomentumFluxZZ")); 119 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YZ, "MeanMomentumFluxYZ")); 120 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XZ, "MeanMomentumFluxXZ")); 121 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XY, "MeanMomentumFluxXY")); 122 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_X, "MeanVelocityX")); 123 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Y, "MeanVelocityY")); 124 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Z, "MeanVelocityZ")); 125 126 PetscCall(PetscLogStagePop()); 127 PetscFunctionReturn(PETSC_SUCCESS); 128 } 129 130 /** @brief Create CeedElemRestriction for collocated data in component-major order. 131 a. Sets the strides of the restriction to component-major order 132 Number of quadrature points is used from the CeedBasis, and number of elements is used from the CeedElemRestriction. 133 */ 134 static PetscErrorCode CreateElemRestrColloc_CompMajor(Ceed ceed, CeedInt num_comp, CeedBasis basis, CeedElemRestriction elem_restr_base, 135 CeedElemRestriction *elem_restr_collocated) { 136 CeedInt num_elem_qpts, loc_num_elem; 137 138 PetscFunctionBeginUser; 139 PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(basis, &num_elem_qpts)); 140 PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_base, &loc_num_elem)); 141 142 const CeedInt strides[] = {num_comp, 1, num_elem_qpts * num_comp}; 143 PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, loc_num_elem, num_elem_qpts, num_comp, num_comp * loc_num_elem * num_elem_qpts, strides, 144 elem_restr_collocated)); 145 PetscFunctionReturn(PETSC_SUCCESS); 146 } 147 148 // Get coordinates of quadrature points 149 PetscErrorCode GetQuadratureCoords(Ceed ceed, DM dm, CeedElemRestriction elem_restr_x, CeedBasis basis_x, CeedVector x_coords, Vec *Qx_coords, 150 PetscInt *total_nqpnts) { 151 CeedElemRestriction elem_restr_qx; 152 CeedQFunction qf_quad_coords; 153 CeedOperator op_quad_coords; 154 CeedInt num_comp_x; 155 CeedSize l_vec_size; 156 OperatorApplyContext op_quad_coords_ctx; 157 158 PetscFunctionBeginUser; 159 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x)); 160 PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_x, basis_x, elem_restr_x, &elem_restr_qx)); 161 PetscCallCeed(ceed, CeedElemRestrictionGetLVectorSize(elem_restr_qx, &l_vec_size)); 162 *total_nqpnts = l_vec_size / num_comp_x; 163 164 // Create QFunction 165 PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_x, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_quad_coords)); 166 167 // Create Operator 168 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_quad_coords, NULL, NULL, &op_quad_coords)); 169 PetscCallCeed(ceed, CeedOperatorSetField(op_quad_coords, "input", elem_restr_x, basis_x, x_coords)); 170 PetscCallCeed(ceed, CeedOperatorSetField(op_quad_coords, "output", elem_restr_qx, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 171 172 PetscCall(CeedOperatorCreateLocalVecs(op_quad_coords, DMReturnVecType(dm), PetscObjectComm((PetscObject)dm), NULL, Qx_coords)); 173 PetscCall(OperatorApplyContextCreate(NULL, NULL, ceed, op_quad_coords, CEED_VECTOR_NONE, NULL, NULL, NULL, &op_quad_coords_ctx)); 174 175 PetscCall(ApplyCeedOperatorLocalToLocal(NULL, *Qx_coords, op_quad_coords_ctx)); 176 177 PetscCall(OperatorApplyContextDestroy(op_quad_coords_ctx)); 178 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qx)); 179 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_quad_coords)); 180 PetscCallCeed(ceed, CeedOperatorDestroy(&op_quad_coords)); 181 PetscFunctionReturn(PETSC_SUCCESS); 182 } 183 184 PetscErrorCode SpanStatsSetupDataCreate(Ceed ceed, User user, CeedData ceed_data, ProblemData problem, SpanStatsSetupData *stats_data) { 185 DM dm = user->spanstats.dm; 186 CeedInt num_comp_x, num_comp_stats = user->spanstats.num_comp_stats; 187 Vec X_loc; 188 DMLabel domain_label = NULL; 189 PetscInt label_value = 0, height = 0, dm_field = 0; 190 191 PetscFunctionBeginUser; 192 PetscCall(PetscNew(stats_data)); 193 194 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &(*stats_data)->elem_restr_parent_stats)); 195 PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &(*stats_data)->elem_restr_parent_x)); 196 PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, problem->q_data_size_sur, 197 &(*stats_data)->elem_restr_parent_qd)); 198 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents((*stats_data)->elem_restr_parent_x, &num_comp_x)); 199 PetscCallCeed(ceed, CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_x, &(*stats_data)->x_coord, NULL)); 200 PetscCallCeed(ceed, CeedElemRestrictionCreateVector((*stats_data)->elem_restr_parent_qd, &(*stats_data)->q_data, NULL)); 201 202 { 203 DM dm_coord; 204 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 205 PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &(*stats_data)->basis_x)); 206 PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &(*stats_data)->basis_stats)); 207 } 208 209 PetscCall(CreateElemRestrColloc_CompMajor(ceed, num_comp_stats, (*stats_data)->basis_stats, (*stats_data)->elem_restr_parent_stats, 210 &(*stats_data)->elem_restr_parent_colloc)); 211 PetscCall( 212 CreateElemRestrColloc_CompMajor(ceed, num_comp_stats, ceed_data->basis_q, ceed_data->elem_restr_q, &(*stats_data)->elem_restr_child_colloc)); 213 214 { // -- Copy DM coordinates into CeedVector 215 DM cdm; 216 PetscCall(DMGetCellCoordinateDM(dm, &cdm)); 217 if (cdm) { 218 PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc)); 219 } else { 220 PetscCall(DMGetCoordinatesLocal(dm, &X_loc)); 221 } 222 } 223 PetscCall(VecScale(X_loc, problem->dm_scale)); 224 PetscCall(VecCopyPetscToCeed(X_loc, (*stats_data)->x_coord)); 225 PetscFunctionReturn(PETSC_SUCCESS); 226 } 227 228 PetscErrorCode SpanStatsSetupDataDestroy(SpanStatsSetupData data) { 229 Ceed ceed; 230 231 PetscFunctionBeginUser; 232 PetscCall(CeedElemRestrictionGetCeed(data->elem_restr_parent_x, &ceed)); 233 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_x)); 234 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_stats)); 235 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_qd)); 236 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_parent_colloc)); 237 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&data->elem_restr_child_colloc)); 238 239 PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_x)); 240 PetscCallCeed(ceed, CeedBasisDestroy(&data->basis_stats)); 241 242 PetscCallCeed(ceed, CeedVectorDestroy(&data->x_coord)); 243 PetscCallCeed(ceed, CeedVectorDestroy(&data->q_data)); 244 245 PetscCall(PetscFree(data)); 246 PetscFunctionReturn(PETSC_SUCCESS); 247 } 248 249 // Create PetscSF for child-to-parent communication 250 PetscErrorCode CreateStatsSF(Ceed ceed, CeedData ceed_data, SpanStatsSetupData stats_data, DM parentdm, DM childdm, PetscSF *statssf) { 251 PetscInt child_num_qpnts, parent_num_qpnts; 252 CeedInt num_comp_x; 253 Vec Child_qx_coords, Parent_qx_coords; 254 255 PetscFunctionBeginUser; 256 PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)childdm), statssf)); 257 258 // Assume that child and parent have the same number of components 259 PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_x, &num_comp_x)); 260 const PetscInt num_comp_sf = num_comp_x - 1; // Number of coord components used in the creation of the SF 261 262 // Get quad_coords for child and parent DM 263 PetscCall(GetQuadratureCoords(ceed, childdm, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord, &Child_qx_coords, &child_num_qpnts)); 264 PetscCall(GetQuadratureCoords(ceed, parentdm, stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord, &Parent_qx_coords, 265 &parent_num_qpnts)); 266 267 { // Remove z component of coordinates for matching 268 const PetscReal *child_quad_coords, *parent_quad_coords; 269 PetscReal *child_coords, *parent_coords; 270 271 PetscCall(VecGetArrayRead(Child_qx_coords, &child_quad_coords)); 272 PetscCall(VecGetArrayRead(Parent_qx_coords, &parent_quad_coords)); 273 274 PetscCall(PetscMalloc2(child_num_qpnts * 2, &child_coords, parent_num_qpnts * 2, &parent_coords)); 275 for (int i = 0; i < child_num_qpnts; i++) { 276 child_coords[0 + i * num_comp_sf] = child_quad_coords[0 + i * num_comp_x]; 277 child_coords[1 + i * num_comp_sf] = child_quad_coords[1 + i * num_comp_x]; 278 } 279 for (int i = 0; i < parent_num_qpnts; i++) { 280 parent_coords[0 + i * num_comp_sf] = parent_quad_coords[0 + i * num_comp_x]; 281 parent_coords[1 + i * num_comp_sf] = parent_quad_coords[1 + i * num_comp_x]; 282 } 283 PetscCall(VecRestoreArrayRead(Child_qx_coords, &child_quad_coords)); 284 PetscCall(VecRestoreArrayRead(Parent_qx_coords, &parent_quad_coords)); 285 286 PetscCall(PetscSFSetGraphFromCoordinates(*statssf, parent_num_qpnts, child_num_qpnts, num_comp_sf, 1e-12, parent_coords, child_coords)); 287 PetscCall(PetscFree2(child_coords, parent_coords)); 288 } 289 290 PetscCall(PetscSFViewFromOptions(*statssf, NULL, "-spanstats_sf_view")); 291 292 PetscCall(VecDestroy(&Child_qx_coords)); 293 PetscCall(VecDestroy(&Parent_qx_coords)); 294 PetscFunctionReturn(PETSC_SUCCESS); 295 } 296 297 // @brief Setup RHS and LHS for L^2 projection of statistics 298 PetscErrorCode SetupL2ProjectionStats(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) { 299 CeedOperator op_mass, op_setup_sur, op_proj_rhs; 300 CeedQFunction qf_mass, qf_stats_proj; 301 CeedInt q_data_size, num_comp_stats = user->spanstats.num_comp_stats; 302 MPI_Comm comm = PetscObjectComm((PetscObject)user->spanstats.dm); 303 304 PetscFunctionBeginUser; 305 // -- Create Operator for RHS of L^2 projection of statistics 306 // Simply take collocated parent data (with quadrature weight already applied) and multiply by weight function. 307 // Therefore, an Identity QF is sufficient 308 PetscCallCeed(ceed, CeedQFunctionCreateIdentity(ceed, num_comp_stats, CEED_EVAL_NONE, CEED_EVAL_INTERP, &qf_stats_proj)); 309 310 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_proj, NULL, NULL, &op_proj_rhs)); 311 PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "input", stats_data->elem_restr_parent_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 312 PetscCallCeed(ceed, CeedOperatorSetField(op_proj_rhs, "output", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 313 314 PetscCall(OperatorApplyContextCreate(NULL, user->spanstats.dm, ceed, op_proj_rhs, NULL, NULL, NULL, NULL, &user->spanstats.op_proj_rhs_ctx)); 315 PetscCall(CeedOperatorCreateLocalVecs(op_proj_rhs, DMReturnVecType(user->spanstats.dm), comm, &user->spanstats.Parent_Stats_loc, NULL)); 316 317 // -- Setup LHS of L^2 projection 318 // Get q_data for mass matrix operator 319 PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur)); 320 PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", stats_data->elem_restr_parent_x, stats_data->basis_x, CEED_VECTOR_ACTIVE)); 321 PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, stats_data->basis_x, CEED_VECTOR_NONE)); 322 PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 323 PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, stats_data->x_coord, stats_data->q_data, CEED_REQUEST_IMMEDIATE)); 324 325 // CEED Restriction 326 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size)); 327 328 // Create Mass CeedOperator 329 PetscCall(CreateMassQFunction(ceed, num_comp_stats, q_data_size, &qf_mass)); 330 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, &op_mass)); 331 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "u", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 332 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_NONE, stats_data->q_data)); 333 PetscCallCeed(ceed, CeedOperatorSetField(op_mass, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 334 335 { // Setup KSP for L^2 projection 336 Mat mat_mass; 337 KSP ksp; 338 339 PetscCall(MatCeedCreate(user->spanstats.dm, user->spanstats.dm, op_mass, NULL, &mat_mass)); 340 341 PetscCall(KSPCreate(comm, &ksp)); 342 PetscCall(KSPSetOptionsPrefix(ksp, "turbulence_spanstats_")); 343 { 344 PC pc; 345 PetscCall(KSPGetPC(ksp, &pc)); 346 PetscCall(PCSetType(pc, PCJACOBI)); 347 PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL)); 348 PetscCall(KSPSetType(ksp, KSPCG)); 349 PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL)); 350 PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 351 } 352 PetscCall(KSPSetFromOptions_WithMatCeed(ksp, mat_mass)); 353 user->spanstats.ksp = ksp; 354 PetscCall(MatDestroy(&mat_mass)); 355 } 356 357 // Cleanup 358 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass)); 359 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_proj)); 360 PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass)); 361 PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur)); 362 PetscCallCeed(ceed, CeedOperatorDestroy(&op_proj_rhs)); 363 PetscFunctionReturn(PETSC_SUCCESS); 364 } 365 366 // Create CeedOperator for statistics collection 367 PetscErrorCode CreateStatisticCollectionOperator(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data, ProblemData problem) { 368 CeedInt num_comp_stats = user->spanstats.num_comp_stats, num_comp_x = problem->dim, num_comp_q; 369 Turbulence_SpanStatsContext collect_ctx; 370 NewtonianIdealGasContext newtonian_ig_ctx; 371 CeedQFunctionContext collect_context; 372 CeedQFunction qf_stats_collect; 373 CeedOperator op_stats_collect; 374 375 PetscFunctionBeginUser; 376 PetscCallCeed(ceed, CeedBasisGetNumComponents(ceed_data->basis_q, &num_comp_q)); 377 378 // Create Operator for statistics collection 379 switch (user->phys->state_var) { 380 case STATEVAR_PRIMITIVE: 381 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &qf_stats_collect)); 382 break; 383 case STATEVAR_CONSERVATIVE: 384 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &qf_stats_collect)); 385 break; 386 default: 387 SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "No statisics collection available for chosen state variable"); 388 } 389 390 if (user->spanstats.do_mms_test) { 391 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect)); 392 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest, ChildStatsCollectionMMSTest_loc, &qf_stats_collect)); 393 } 394 395 { // Setup Collection Context 396 PetscCall(PetscNew(&collect_ctx)); 397 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfunction_context, CEED_MEM_HOST, &newtonian_ig_ctx)); 398 collect_ctx->gas = *newtonian_ig_ctx; 399 400 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &collect_context)); 401 PetscCallCeed(ceed, CeedQFunctionContextSetData(collect_context, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx)); 402 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(collect_context, CEED_MEM_HOST, FreeContextPetsc)); 403 404 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(collect_context, "solution time", 405 offsetof(struct Turbulence_SpanStatsContext_, solution_time), 1, "Current solution time")); 406 PetscCallCeed( 407 ceed, CeedQFunctionContextRegisterDouble(collect_context, "previous time", offsetof(struct Turbulence_SpanStatsContext_, previous_time), 1, 408 "Previous time statistics collection was done")); 409 410 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfunction_context, &newtonian_ig_ctx)); 411 } 412 413 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_stats_collect, collect_context)); 414 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&collect_context)); 415 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP)); 416 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q_data", problem->q_data_size_vol, CEED_EVAL_NONE)); 417 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP)); 418 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE)); 419 420 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_collect, NULL, NULL, &op_stats_collect)); 421 PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 422 PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q_data", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 423 PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord)); 424 PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "v", stats_data->elem_restr_child_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 425 426 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "solution time", &user->spanstats.solution_time_label)); 427 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "previous time", &user->spanstats.previous_time_label)); 428 429 PetscCall(OperatorApplyContextCreate(user->dm, user->spanstats.dm, user->ceed, op_stats_collect, user->q_ceed, NULL, NULL, NULL, 430 &user->spanstats.op_stats_collect_ctx)); 431 432 PetscCall(CeedOperatorCreateLocalVecs(op_stats_collect, DMReturnVecType(user->spanstats.dm), PetscObjectComm((PetscObject)user->spanstats.dm), NULL, 433 &user->spanstats.Child_Stats_loc)); 434 PetscCall(VecZeroEntries(user->spanstats.Child_Stats_loc)); 435 436 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect)); 437 PetscCallCeed(ceed, CeedOperatorDestroy(&op_stats_collect)); 438 PetscFunctionReturn(PETSC_SUCCESS); 439 } 440 441 // Creates operator for calculating error of method of manufactured solution (MMS) test 442 PetscErrorCode SetupMMSErrorChecking(Ceed ceed, User user, CeedData ceed_data, SpanStatsSetupData stats_data) { 443 CeedInt num_comp_stats = user->spanstats.num_comp_stats, num_comp_x, q_data_size; 444 CeedQFunction qf_error; 445 CeedOperator op_error; 446 CeedVector x_ceed, y_ceed; 447 448 PetscFunctionBeginUser; 449 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(stats_data->elem_restr_parent_qd, &q_data_size)); 450 PetscCallCeed(ceed, CeedBasisGetNumComponents(stats_data->basis_x, &num_comp_x)); 451 452 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest_Error, ChildStatsCollectionMMSTest_Error_loc, &qf_error)); 453 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP)); 454 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE)); 455 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP)); 456 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP)); 457 458 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error)); 459 PetscCallCeed(ceed, CeedOperatorSetField(op_error, "q", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 460 PetscCallCeed(ceed, CeedOperatorSetField(op_error, "qdata", stats_data->elem_restr_parent_qd, CEED_BASIS_NONE, stats_data->q_data)); 461 PetscCallCeed(ceed, CeedOperatorSetField(op_error, "x", stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord)); 462 PetscCallCeed(ceed, CeedOperatorSetField(op_error, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 463 464 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &x_ceed, NULL)); 465 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &y_ceed, NULL)); 466 PetscCall(OperatorApplyContextCreate(user->spanstats.dm, user->spanstats.dm, user->ceed, op_error, x_ceed, y_ceed, NULL, NULL, 467 &user->spanstats.mms_error_ctx)); 468 469 PetscCallCeed(ceed, CeedOperatorDestroy(&op_error)); 470 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_error)); 471 PetscCallCeed(ceed, CeedVectorDestroy(&x_ceed)); 472 PetscCallCeed(ceed, CeedVectorDestroy(&y_ceed)); 473 PetscFunctionReturn(PETSC_SUCCESS); 474 } 475 476 // Setup for statistics collection 477 PetscErrorCode TurbulenceStatisticsSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) { 478 SpanStatsSetupData stats_data; 479 PetscLogStage stage_stats_setup; 480 481 PetscFunctionBeginUser; 482 PetscCall(PetscLogStageGetId("Stats Setup", &stage_stats_setup)); 483 if (stage_stats_setup == -1) PetscCall(PetscLogStageRegister("Stats Setup", &stage_stats_setup)); 484 PetscCall(PetscLogStagePush(stage_stats_setup)); 485 486 // Create parent DM 487 PetscCall(CreateStatsDM(user, problem, user->app_ctx->degree)); 488 489 // Create necessary CeedObjects for setting up statistics 490 PetscCall(SpanStatsSetupDataCreate(ceed, user, ceed_data, problem, &stats_data)); 491 492 // Create SF for communicating child data back their respective parents 493 PetscCall(CreateStatsSF(ceed, ceed_data, stats_data, user->dm, user->spanstats.dm, &user->spanstats.sf)); 494 495 // Create CeedOperators for statistics collection 496 PetscCall(CreateStatisticCollectionOperator(ceed, user, ceed_data, stats_data, problem)); 497 498 // Setup KSP and Mat for L^2 projection of statistics 499 PetscCall(SetupL2ProjectionStats(ceed, user, ceed_data, stats_data)); 500 501 PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_turbulence_spanstats_mms", &user->spanstats.do_mms_test, NULL)); 502 if (user->spanstats.do_mms_test) { 503 PetscCall(SetupMMSErrorChecking(ceed, user, ceed_data, stats_data)); 504 } 505 506 { // Setup stats viewer with prefix 507 PetscViewerType viewer_type; 508 PetscCall(PetscViewerGetType(user->app_ctx->turb_spanstats_viewer, &viewer_type)); 509 PetscCall(PetscOptionsSetValue(NULL, "-ts_monitor_turbulence_spanstats_viewer_type", viewer_type)); 510 511 PetscCall(PetscViewerSetOptionsPrefix(user->app_ctx->turb_spanstats_viewer, "ts_monitor_turbulence_spanstats_")); 512 PetscCall(PetscViewerSetFromOptions(user->app_ctx->turb_spanstats_viewer)); 513 } 514 515 PetscCall(SpanStatsSetupDataDestroy(stats_data)); 516 PetscCall(PetscLogStagePop()); 517 PetscFunctionReturn(PETSC_SUCCESS); 518 } 519 520 // Collect statistics based on the solution Q 521 PetscErrorCode CollectStatistics(User user, PetscScalar solution_time, Vec Q) { 522 SpanStatsData user_stats = user->spanstats; 523 524 PetscFunctionBeginUser; 525 PetscLogStage stage_stats_collect; 526 PetscCall(PetscLogStageGetId("Stats Collect", &stage_stats_collect)); 527 if (stage_stats_collect == -1) PetscCall(PetscLogStageRegister("Stats Collect", &stage_stats_collect)); 528 PetscCall(PetscLogStagePush(stage_stats_collect)); 529 530 PetscCall(UpdateBoundaryValues(user, user->Q_loc, solution_time)); 531 PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.solution_time_label, &solution_time)); 532 PetscCall(DMGlobalToLocal(user->dm, Q, INSERT_VALUES, user->Q_loc)); 533 PetscCall(ApplyAddCeedOperatorLocalToLocal(user->Q_loc, user_stats.Child_Stats_loc, user_stats.op_stats_collect_ctx)); 534 535 PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user_stats.op_stats_collect_ctx->op, user_stats.previous_time_label, &solution_time)); 536 537 PetscCall(PetscLogStagePop()); 538 PetscFunctionReturn(PETSC_SUCCESS); 539 } 540 541 // Process the child statistics into parent statistics and project them onto stats 542 PetscErrorCode ProcessStatistics(User user, Vec stats) { 543 SpanStatsData user_stats = user->spanstats; 544 const PetscScalar *child_stats; 545 PetscScalar *parent_stats; 546 MPI_Datatype unit; 547 Vec RHS; 548 549 PetscFunctionBeginUser; 550 PetscLogStage stage_stats_process; 551 PetscCall(PetscLogStageGetId("Stats Process", &stage_stats_process)); 552 if (stage_stats_process == -1) PetscCall(PetscLogStageRegister("Stats Process", &stage_stats_process)); 553 PetscCall(PetscLogStagePush(stage_stats_process)); 554 555 PetscCall(VecZeroEntries(user_stats.Parent_Stats_loc)); 556 557 PetscCall(VecGetArrayRead(user_stats.Child_Stats_loc, &child_stats)); 558 PetscCall(VecGetArray(user_stats.Parent_Stats_loc, &parent_stats)); 559 560 if (user_stats.num_comp_stats == 1) unit = MPIU_REAL; 561 else { 562 PetscCallMPI(MPI_Type_contiguous(user_stats.num_comp_stats, MPIU_REAL, &unit)); 563 PetscCallMPI(MPI_Type_commit(&unit)); 564 } 565 566 PetscCall(PetscSFReduceBegin(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM)); 567 PetscCall(PetscSFReduceEnd(user_stats.sf, unit, child_stats, parent_stats, MPI_SUM)); 568 569 PetscCall(VecRestoreArrayRead(user_stats.Child_Stats_loc, &child_stats)); 570 PetscCall(VecRestoreArray(user_stats.Parent_Stats_loc, &parent_stats)); 571 PetscCallMPI(MPI_Type_free(&unit)); 572 573 PetscReal solution_time; 574 PetscCall(DMGetOutputSequenceNumber(user_stats.dm, NULL, &solution_time)); 575 PetscReal summing_duration = solution_time - user->app_ctx->cont_time; 576 PetscCall(VecScale(user_stats.Parent_Stats_loc, 1 / (summing_duration * user_stats.span_width))); 577 578 // L^2 projection with the parent_data 579 PetscCall(DMGetGlobalVector(user_stats.dm, &RHS)); 580 PetscCall(ApplyCeedOperatorLocalToGlobal(user_stats.Parent_Stats_loc, RHS, user_stats.op_proj_rhs_ctx)); 581 582 PetscCall(KSPSolve(user_stats.ksp, RHS, stats)); 583 584 PetscCall(DMRestoreGlobalVector(user_stats.dm, &RHS)); 585 PetscCall(PetscLogStagePop()); 586 PetscFunctionReturn(PETSC_SUCCESS); 587 } 588 589 // TSMonitor for the statistics collection and processing 590 PetscErrorCode TSMonitor_TurbulenceStatistics(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, void *ctx) { 591 User user = (User)ctx; 592 Vec stats; 593 TSConvergedReason reason; 594 PetscInt collect_interval = user->app_ctx->turb_spanstats_collect_interval, viewer_interval = user->app_ctx->turb_spanstats_viewer_interval; 595 596 PetscFunctionBeginUser; 597 PetscCall(TSGetConvergedReason(ts, &reason)); 598 // Do not collect or process on the first step of the run (ie. on the initial condition) 599 if (steps == user->app_ctx->cont_steps) PetscFunctionReturn(PETSC_SUCCESS); 600 601 PetscBool run_processing_and_viewer = (steps % viewer_interval == 0 && viewer_interval != -1) || reason != TS_CONVERGED_ITERATING; 602 603 if (steps % collect_interval == 0 || run_processing_and_viewer) { 604 PetscCall(CollectStatistics(user, solution_time, Q)); 605 606 if (run_processing_and_viewer) { 607 PetscCall(DMSetOutputSequenceNumber(user->spanstats.dm, steps, solution_time)); 608 PetscCall(DMGetGlobalVector(user->spanstats.dm, &stats)); 609 PetscCall(ProcessStatistics(user, stats)); 610 if (user->app_ctx->test_type == TESTTYPE_NONE) { 611 PetscCall(PetscViewerPushFormat(user->app_ctx->turb_spanstats_viewer, user->app_ctx->turb_spanstats_viewer_format)); 612 PetscCall(VecView(stats, user->app_ctx->turb_spanstats_viewer)); 613 PetscCall(PetscViewerPopFormat(user->app_ctx->turb_spanstats_viewer)); 614 } 615 if (user->app_ctx->test_type == TESTTYPE_TURB_SPANSTATS && reason != TS_CONVERGED_ITERATING) { 616 PetscCall(RegressionTest(user->app_ctx, stats)); 617 } 618 if (user->spanstats.do_mms_test && reason != TS_CONVERGED_ITERATING) { 619 Vec error; 620 PetscCall(VecDuplicate(stats, &error)); 621 PetscCall(ApplyCeedOperatorGlobalToGlobal(stats, error, user->spanstats.mms_error_ctx)); 622 PetscScalar error_sq = 0; 623 PetscCall(VecSum(error, &error_sq)); 624 PetscScalar l2_error = sqrt(error_sq); 625 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error)); 626 } 627 PetscCall(DMRestoreGlobalVector(user->spanstats.dm, &stats)); 628 } 629 } 630 PetscFunctionReturn(PETSC_SUCCESS); 631 } 632 633 PetscErrorCode TurbulenceStatisticsDestroy(User user, CeedData ceed_data) { 634 PetscFunctionBeginUser; 635 PetscCall(VecDestroy(&user->spanstats.Child_Stats_loc)); 636 PetscCall(VecDestroy(&user->spanstats.Parent_Stats_loc)); 637 638 PetscCall(OperatorApplyContextDestroy(user->spanstats.op_stats_collect_ctx)); 639 PetscCall(OperatorApplyContextDestroy(user->spanstats.op_proj_rhs_ctx)); 640 PetscCall(OperatorApplyContextDestroy(user->spanstats.mms_error_ctx)); 641 642 PetscCall(KSPDestroy(&user->spanstats.ksp)); 643 PetscCall(PetscSFDestroy(&user->spanstats.sf)); 644 PetscCall(DMDestroy(&user->spanstats.dm)); 645 PetscFunctionReturn(PETSC_SUCCESS); 646 } 647