1b30619f6SJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2b30619f6SJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3b30619f6SJames Wright /// @file 4b30619f6SJames Wright /// Functions for setting up and performing spanwise-statistics collection of CFL and Peclet number 5b30619f6SJames Wright 6b30619f6SJames Wright #include "../qfunctions/spanstats/cflpe.h" 7b30619f6SJames Wright #include "../qfunctions/advection_types.h" 8b30619f6SJames Wright 9b30619f6SJames Wright #include <ceed.h> 10b30619f6SJames Wright #include <petscdmplex.h> 11b30619f6SJames Wright #include <petscsf.h> 12b30619f6SJames Wright #include <spanstats.h> 13b30619f6SJames Wright 14b30619f6SJames Wright #include <navierstokes.h> 15b30619f6SJames Wright #include <petsc_ops.h> 16b30619f6SJames Wright 17b30619f6SJames Wright // @brief Create CeedOperator for statistics collection 18b30619f6SJames Wright static PetscErrorCode CreateStatisticCollectionOperator(Honee honee, SpanStatsCtx spanstats, SpanStatsSetupData stats_data, ProblemData problem) { 19b30619f6SJames Wright Ceed ceed = honee->ceed; 20b30619f6SJames Wright CeedInt num_comp_stats = spanstats->num_comp_stats, num_comp_x, num_comp_q, q_data_size; 21*e3db12f8SJames Wright PetscInt dim; 22b30619f6SJames Wright CflPe_SpanStatsContext collect_ctx; 23b30619f6SJames Wright CeedQFunctionContext collect_qfctx; 24b30619f6SJames Wright CeedQFunction qf_stats_collect = NULL; 25b30619f6SJames Wright CeedOperator op_stats_collect; 26b30619f6SJames Wright CeedVector q_data; 27b30619f6SJames Wright CeedElemRestriction elem_restr_qd; 28b30619f6SJames Wright 29b30619f6SJames Wright PetscFunctionBeginUser; 30b30619f6SJames Wright PetscCall(DMGetDimension(honee->dm, &dim)); 31b30619f6SJames Wright num_comp_x = dim; 32b30619f6SJames Wright PetscCallCeed(ceed, CeedBasisGetNumComponents(honee->basis_q, &num_comp_q)); 33b30619f6SJames Wright 34b30619f6SJames Wright // Create Operator for statistics collection 35b30619f6SJames Wright switch (dim) { 36b30619f6SJames Wright case 2: 37b30619f6SJames Wright switch (honee->phys->state_var) { 38b30619f6SJames Wright case STATEVAR_PRIMITIVE: 39b30619f6SJames Wright PetscCallCeed(ceed, 40b30619f6SJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_2D_Prim, ChildStatsCollection_2D_Prim_loc, &qf_stats_collect)); 41b30619f6SJames Wright break; 42b30619f6SJames Wright case STATEVAR_CONSERVATIVE: 439eadbee4SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_2D_Conserv, ChildStatsCollection_2D_Conserv_loc, 449eadbee4SJames Wright &qf_stats_collect)); 45b30619f6SJames Wright break; 46b30619f6SJames Wright case STATEVAR_ENTROPY: 479eadbee4SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_2D_Entropy, ChildStatsCollection_2D_Entropy_loc, 489eadbee4SJames Wright &qf_stats_collect)); 49b30619f6SJames Wright break; 50b30619f6SJames Wright } 51b30619f6SJames Wright break; 52b30619f6SJames Wright case 3: 53b30619f6SJames Wright switch (honee->phys->state_var) { 54b30619f6SJames Wright case STATEVAR_PRIMITIVE: 55b30619f6SJames Wright PetscCallCeed(ceed, 56b30619f6SJames Wright CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_3D_Prim, ChildStatsCollection_3D_Prim_loc, &qf_stats_collect)); 57b30619f6SJames Wright break; 58b30619f6SJames Wright case STATEVAR_CONSERVATIVE: 599eadbee4SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_3D_Conserv, ChildStatsCollection_3D_Conserv_loc, 609eadbee4SJames Wright &qf_stats_collect)); 61b30619f6SJames Wright break; 62b30619f6SJames Wright case STATEVAR_ENTROPY: 639eadbee4SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_3D_Entropy, ChildStatsCollection_3D_Entropy_loc, 649eadbee4SJames Wright &qf_stats_collect)); 65b30619f6SJames Wright break; 66b30619f6SJames Wright } 67b30619f6SJames Wright break; 68b30619f6SJames Wright } 69b30619f6SJames Wright PetscCheck(qf_stats_collect, PetscObjectComm((PetscObject)honee->dm), PETSC_ERR_SUP, 70b30619f6SJames Wright "Could not create CFL/Pe spanstats collection QFunction for dim %" PetscInt_FMT " and state variable %s", dim, 71b30619f6SJames Wright StateVariables[honee->phys->state_var]); 72b30619f6SJames Wright 73*e3db12f8SJames Wright PetscCall(QDataGet(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, 74*e3db12f8SJames Wright &q_data, &q_data_size)); 75b30619f6SJames Wright { // Setup Collection Context 76b30619f6SJames Wright PetscBool is_advection; 77b30619f6SJames Wright 78b30619f6SJames Wright PetscCall(PetscNew(&collect_ctx)); 79b30619f6SJames Wright PetscCall(PetscStrcmp("advection", honee->app_ctx->problem_name, &is_advection)); 80b30619f6SJames Wright if (is_advection) { 81b30619f6SJames Wright AdvectionContext advection_ctx; 82b30619f6SJames Wright 83b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &advection_ctx)); 84cde3d787SJames Wright collect_ctx->newt_ctx = (struct NewtonianIdealGasContext_){0}; 85b30619f6SJames Wright collect_ctx->diffusion_coeff = advection_ctx->diffusion_coeff; 86b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &advection_ctx)); 87b30619f6SJames Wright } else { 88b30619f6SJames Wright NewtonianIdealGasContext newtonian_ig_ctx; 89b30619f6SJames Wright 90b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx)); 91cde3d787SJames Wright collect_ctx->newt_ctx = *newtonian_ig_ctx; 92cde3d787SJames Wright collect_ctx->diffusion_coeff = newtonian_ig_ctx->gas.mu; 93b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx)); 94b30619f6SJames Wright } 95b30619f6SJames Wright 96b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &collect_qfctx)); 97b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(collect_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx)); 98b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(collect_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 99b30619f6SJames Wright 100b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(collect_qfctx, "solution time", offsetof(struct CflPe_SpanStatsContext_, solution_time), 1, 101b30619f6SJames Wright "Current solution time")); 102b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(collect_qfctx, "previous time", offsetof(struct CflPe_SpanStatsContext_, previous_time), 1, 103b30619f6SJames Wright "Previous time statistics collection was done")); 104b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(collect_qfctx, "timestep", offsetof(struct CflPe_SpanStatsContext_, timestep), 1, 105b30619f6SJames Wright "Current timestep taken by TS")); 106b30619f6SJames Wright } 107b30619f6SJames Wright 108b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_stats_collect, collect_qfctx)); 109b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&collect_qfctx)); 110b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP)); 111b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q_data", q_data_size, CEED_EVAL_NONE)); 112b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP)); 113b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE)); 114b30619f6SJames Wright 115b30619f6SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_collect, NULL, NULL, &op_stats_collect)); 116b30619f6SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 117b30619f6SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q_data", elem_restr_qd, CEED_BASIS_NONE, q_data)); 118b30619f6SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "x", honee->elem_restr_x, honee->basis_x, honee->x_coord)); 119b30619f6SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "v", stats_data->elem_restr_child_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 120b30619f6SJames Wright 121b30619f6SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "solution time", &spanstats->solution_time_label)); 122b30619f6SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "previous time", &spanstats->previous_time_label)); 123b30619f6SJames Wright PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "timestep", &spanstats->timestep_label)); 124b30619f6SJames Wright 125b30619f6SJames Wright PetscCall(OperatorApplyContextCreate(honee->dm, spanstats->dm, honee->ceed, op_stats_collect, honee->q_ceed, NULL, NULL, NULL, 126b30619f6SJames Wright &spanstats->op_stats_collect_ctx)); 127b30619f6SJames Wright 128b30619f6SJames Wright PetscCall(CeedOperatorCreateLocalVecs(op_stats_collect, DMReturnVecType(spanstats->dm), PETSC_COMM_SELF, NULL, &spanstats->Child_Stats_loc)); 129b30619f6SJames Wright PetscCall(VecZeroEntries(spanstats->Child_Stats_loc)); 130b30619f6SJames Wright 131b30619f6SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 132b30619f6SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 133b30619f6SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect)); 134b30619f6SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_stats_collect)); 135b30619f6SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 136b30619f6SJames Wright } 137b30619f6SJames Wright 138b30619f6SJames Wright // @brief Setup for statistics collection of CFL and Peclet number 139b30619f6SJames Wright PetscErrorCode SpanwiseStatisticsSetup_CflPe(TS ts, PetscViewerAndFormat *ctx) { 140b30619f6SJames Wright Honee honee; 141b30619f6SJames Wright SpanStatsSetupData stats_setup_data; 142b30619f6SJames Wright SpanStatsCtx spanstats; 143b30619f6SJames Wright const char *prefix = "ts_monitor_spanstats_cflpe_"; 144b30619f6SJames Wright PetscBool is_ascii; 145b30619f6SJames Wright 146b30619f6SJames Wright PetscFunctionBeginUser; 147b30619f6SJames Wright PetscCall(TSGetApplicationContext(ts, &honee)); 148b30619f6SJames Wright PetscCall(PetscObjectTypeCompare((PetscObject)ctx->viewer, PETSCVIEWERASCII, &is_ascii)); 149b30619f6SJames Wright PetscCheck(!is_ascii || honee->app_ctx->test_type != TESTTYPE_NONE, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, 150b30619f6SJames Wright "Turbulence spanwise statistics does not support ASCII viewers"); 151b30619f6SJames Wright 152b30619f6SJames Wright PetscCall(SpanwiseStatisticsSetupInitialize(honee, honee->app_ctx->degree, 6, prefix, &stats_setup_data, &spanstats)); 153b30619f6SJames Wright 154b30619f6SJames Wright { // Create Section for data 155b30619f6SJames Wright PetscSection section; 156b30619f6SJames Wright 157b30619f6SJames Wright PetscCall(DMGetLocalSection(spanstats->dm, §ion)); 158b30619f6SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 159b30619f6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "MeanCFL")); 160b30619f6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "MeanCFLSquared")); 161b30619f6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "MeanCFLCubed")); 162b30619f6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "MeanPe")); 163b30619f6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "MeanPeSquared")); 164b30619f6SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 5, "MeanPeCubed")); 165b30619f6SJames Wright } 166b30619f6SJames Wright 167b30619f6SJames Wright // Create CeedOperators for statistics collection 168b30619f6SJames Wright PetscCall(CreateStatisticCollectionOperator(honee, spanstats, stats_setup_data, honee->problem_data)); 169b30619f6SJames Wright 170b30619f6SJames Wright ctx->data = spanstats; 1715206a5a0SJames Wright ctx->data_destroy = (PetscCtxDestroyFn *)SpanStatsCtxDestroy; 172b30619f6SJames Wright 173b30619f6SJames Wright PetscCall(SpanwiseStatisticsSetupFinalize(ts, honee, spanstats, ctx, &stats_setup_data)); 174b30619f6SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 175b30619f6SJames Wright } 176b30619f6SJames Wright 177b30619f6SJames Wright // @brief TSMonitor for collection and processing of CFL and Peclet number statistics 178b30619f6SJames Wright PetscErrorCode TSMonitor_SpanwiseStatisticsCflPe(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx) { 179b30619f6SJames Wright Honee honee; 180b30619f6SJames Wright Vec stats; 181b30619f6SJames Wright TSConvergedReason reason; 182b30619f6SJames Wright SpanStatsCtx spanstats = ctx->data; 183b30619f6SJames Wright PetscInt collect_interval = spanstats->collect_interval, viewer_interval = ctx->view_interval; 184b30619f6SJames Wright PetscReal timestep; 185b30619f6SJames Wright 186b30619f6SJames Wright PetscFunctionBeginUser; 187b30619f6SJames Wright PetscCall(TSGetApplicationContext(ts, &honee)); 188b30619f6SJames Wright PetscCall(TSGetConvergedReason(ts, &reason)); 189b30619f6SJames Wright // Do not collect or process on the first step of the run (ie. on the initial condition) 190b30619f6SJames Wright if (steps == honee->app_ctx->cont_steps) PetscFunctionReturn(PETSC_SUCCESS); 191b30619f6SJames Wright 192b30619f6SJames Wright PetscBool run_processing_and_viewer = (steps % viewer_interval == 0 && viewer_interval != -1) || reason != TS_CONVERGED_ITERATING; 193b30619f6SJames Wright 194b30619f6SJames Wright if (steps % collect_interval == 0 || run_processing_and_viewer) { 195b30619f6SJames Wright PetscCall(TSGetTimeStep(ts, ×tep)); 196b30619f6SJames Wright PetscCallCeed(honee->ceed, CeedOperatorSetContextDouble(spanstats->op_stats_collect_ctx->op, spanstats->timestep_label, ×tep)); 197b30619f6SJames Wright PetscCall(SpanwiseStatisticsCollect(honee, spanstats, solution_time, Q)); 198b30619f6SJames Wright 199b30619f6SJames Wright if (run_processing_and_viewer) { 200b30619f6SJames Wright PetscCall(DMSetOutputSequenceNumber(spanstats->dm, steps, solution_time)); 201b30619f6SJames Wright PetscCall(DMGetGlobalVector(spanstats->dm, &stats)); 202b30619f6SJames Wright PetscCall(SpanwiseStatisticsProcess(honee, spanstats, stats)); 203b30619f6SJames Wright 2040aab7249SJames Wright if (honee->app_ctx->test_type == TESTTYPE_NONE) { 205b30619f6SJames Wright PetscCall(PetscViewerPushFormat(ctx->viewer, ctx->format)); 206b30619f6SJames Wright PetscCall(VecView(stats, ctx->viewer)); 207b30619f6SJames Wright PetscCall(PetscViewerPopFormat(ctx->viewer)); 208b30619f6SJames Wright { 209b30619f6SJames Wright PetscInt tab_level; 210b30619f6SJames Wright PetscViewer viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ts)); 211b30619f6SJames Wright CeedScalar second = honee->units->second; 212b30619f6SJames Wright const char *filename; 213b30619f6SJames Wright 2140aab7249SJames Wright PetscCall(PetscViewerFileGetName(ctx->viewer, &filename)); 215b30619f6SJames Wright PetscCall(PetscObjectGetTabLevel((PetscObject)ts, &tab_level)); 216b30619f6SJames Wright PetscCall(PetscViewerASCIIAddTab(viewer, tab_level + 1)); 217b30619f6SJames Wright if (filename) { 218b30619f6SJames Wright PetscCall(PetscViewerASCIIPrintf( 2190aab7249SJames Wright viewer, "Spanwise cflpe statistics written to %s for time span (%0.12e,%0.12e] and step span [%" PetscInt_FMT ",%" PetscInt_FMT "]\n", 220b30619f6SJames Wright filename, spanstats->initial_solution_time / second, solution_time / second, spanstats->initial_solution_step, steps)); 221b30619f6SJames Wright } else { 222b30619f6SJames Wright PetscCall(PetscViewerASCIIPrintf( 223b30619f6SJames Wright viewer, "Spanwise statistics (%s) file written for time span (%0.12e,%0.12e] and step span [%" PetscInt_FMT ",%" PetscInt_FMT "]\n", 224b30619f6SJames Wright spanstats->prefix, spanstats->initial_solution_time / second, solution_time / second, spanstats->initial_solution_step, steps)); 225b30619f6SJames Wright } 226b30619f6SJames Wright PetscCall(PetscViewerASCIISubtractTab(viewer, tab_level + 1)); 227b30619f6SJames Wright } 2280aab7249SJames Wright } else if (honee->app_ctx->test_type == TESTTYPE_SPANSTATS && reason != TS_CONVERGED_ITERATING) { 2290aab7249SJames Wright PetscCall(RegressionTest(honee->app_ctx, stats)); 2300aab7249SJames Wright } 231b30619f6SJames Wright 232b30619f6SJames Wright PetscCall(DMRestoreGlobalVector(spanstats->dm, &stats)); 233b30619f6SJames Wright } 234b30619f6SJames Wright } 235b30619f6SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 236b30619f6SJames Wright } 237