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