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