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