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