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