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