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 turbulence quantities 5 6 #include "../../qfunctions/spanstats/turbulence.h" 7 8 #include <ceed.h> 9 #include <petscdmplex.h> 10 #include <petscsf.h> 11 #include <spanstats.h> 12 13 #include <navierstokes.h> 14 #include <petsc_ops.h> 15 16 // @brief Create CeedOperator for statistics collection 17 static PetscErrorCode CreateStatisticCollectionOperator(Honee honee, SpanStatsCtx spanstats, SpanStatsSetupData stats_data, ProblemData problem) { 18 Ceed ceed = honee->ceed; 19 CeedInt num_comp_stats = spanstats->num_comp_stats, num_comp_x, num_comp_q, q_data_size; 20 PetscInt dim; 21 Turbulence_SpanStatsContext collect_ctx; 22 NewtonianIdealGasContext newtonian_ig_ctx; 23 CeedQFunctionContext collect_qfctx; 24 CeedQFunction qf_stats_collect; 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 (honee->phys->state_var) { 36 case STATEVAR_PRIMITIVE: 37 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Prim, ChildStatsCollection_Prim_loc, &qf_stats_collect)); 38 break; 39 case STATEVAR_CONSERVATIVE: 40 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Conserv, ChildStatsCollection_Conserv_loc, &qf_stats_collect)); 41 break; 42 case STATEVAR_ENTROPY: 43 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollection_Entropy, ChildStatsCollection_Entropy_loc, &qf_stats_collect)); 44 break; 45 } 46 47 if (spanstats->do_mms_test) { 48 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect)); 49 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest, ChildStatsCollectionMMSTest_loc, &qf_stats_collect)); 50 } 51 52 PetscCall(QDataGet(ceed, honee->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, 53 &q_data, &q_data_size)); 54 55 { // Setup Collection Context 56 PetscCall(PetscNew(&collect_ctx)); 57 PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx)); 58 collect_ctx->newt_ctx = *newtonian_ig_ctx; 59 60 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &collect_qfctx)); 61 PetscCallCeed(ceed, CeedQFunctionContextSetData(collect_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx)); 62 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(collect_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 63 64 PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(collect_qfctx, "solution time", 65 offsetof(struct Turbulence_SpanStatsContext_, solution_time), 1, "Current solution time")); 66 PetscCallCeed(ceed, 67 CeedQFunctionContextRegisterDouble(collect_qfctx, "previous time", offsetof(struct Turbulence_SpanStatsContext_, previous_time), 1, 68 "Previous time statistics collection was done")); 69 70 PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx)); 71 } 72 73 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_stats_collect, collect_qfctx)); 74 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&collect_qfctx)); 75 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP)); 76 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q_data", q_data_size, CEED_EVAL_NONE)); 77 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP)); 78 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE)); 79 80 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_collect, NULL, NULL, &op_stats_collect)); 81 PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 82 PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q_data", elem_restr_qd, CEED_BASIS_NONE, q_data)); 83 PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "x", honee->elem_restr_x, honee->basis_x, honee->x_coord)); 84 PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "v", stats_data->elem_restr_child_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 85 86 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "solution time", &spanstats->solution_time_label)); 87 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "previous time", &spanstats->previous_time_label)); 88 89 PetscCall(OperatorApplyContextCreate(honee->dm, spanstats->dm, honee->ceed, op_stats_collect, honee->q_ceed, NULL, NULL, NULL, 90 &spanstats->op_stats_collect_ctx)); 91 92 PetscCall(CeedOperatorCreateLocalVecs(op_stats_collect, DMReturnVecType(spanstats->dm), PETSC_COMM_SELF, NULL, &spanstats->Child_Stats_loc)); 93 PetscCall(VecZeroEntries(spanstats->Child_Stats_loc)); 94 95 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 96 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 97 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect)); 98 PetscCallCeed(ceed, CeedOperatorDestroy(&op_stats_collect)); 99 PetscFunctionReturn(PETSC_SUCCESS); 100 } 101 102 // @brief Creates operator for calculating error of method of manufactured solution (MMS) test 103 static PetscErrorCode SetupMMSErrorChecking(Honee honee, SpanStatsCtx spanstats, SpanStatsSetupData stats_data) { 104 Ceed ceed = honee->ceed; 105 CeedInt num_comp_stats = spanstats->num_comp_stats, num_comp_x, q_data_size; 106 CeedQFunction qf_error; 107 CeedOperator op_error; 108 CeedVector x_ceed, y_ceed; 109 CeedVector q_data; 110 CeedElemRestriction elem_restr_parent_qd; 111 112 PetscFunctionBeginUser; 113 PetscCall(QDataGet(ceed, spanstats->dm, DMLABEL_DEFAULT, DMLABEL_DEFAULT_VALUE, stats_data->elem_restr_parent_x, stats_data->basis_x, 114 stats_data->x_coord, &elem_restr_parent_qd, &q_data, &q_data_size)); 115 PetscCallCeed(ceed, CeedBasisGetNumComponents(stats_data->basis_x, &num_comp_x)); 116 117 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ChildStatsCollectionMMSTest_Error, ChildStatsCollectionMMSTest_Error_loc, &qf_error)); 118 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "q", num_comp_stats, CEED_EVAL_INTERP)); 119 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE)); 120 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_error, "x", num_comp_x, CEED_EVAL_INTERP)); 121 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_error, "v", num_comp_stats, CEED_EVAL_INTERP)); 122 123 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error)); 124 PetscCallCeed(ceed, CeedOperatorSetField(op_error, "q", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 125 PetscCallCeed(ceed, CeedOperatorSetField(op_error, "qdata", elem_restr_parent_qd, CEED_BASIS_NONE, q_data)); 126 PetscCallCeed(ceed, CeedOperatorSetField(op_error, "x", stats_data->elem_restr_parent_x, stats_data->basis_x, stats_data->x_coord)); 127 PetscCallCeed(ceed, CeedOperatorSetField(op_error, "v", stats_data->elem_restr_parent_stats, stats_data->basis_stats, CEED_VECTOR_ACTIVE)); 128 129 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &x_ceed, NULL)); 130 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(stats_data->elem_restr_parent_stats, &y_ceed, NULL)); 131 PetscCall(OperatorApplyContextCreate(spanstats->dm, spanstats->dm, honee->ceed, op_error, x_ceed, y_ceed, NULL, NULL, &spanstats->mms_error_ctx)); 132 133 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 134 PetscCallCeed(ceed, CeedVectorDestroy(&x_ceed)); 135 PetscCallCeed(ceed, CeedVectorDestroy(&y_ceed)); 136 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_parent_qd)); 137 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_error)); 138 PetscCallCeed(ceed, CeedOperatorDestroy(&op_error)); 139 PetscFunctionReturn(PETSC_SUCCESS); 140 } 141 142 // @brief Setup for statistics collection for turbulence quantities 143 PetscErrorCode SpanwiseStatisticsSetup_Turbulence(TS ts, PetscViewerAndFormat *ctx) { 144 Honee honee; 145 SpanStatsSetupData stats_setup_data; 146 SpanStatsCtx spanstats; 147 const char *prefix = "ts_monitor_spanstats_turbulence_"; 148 PetscBool is_ascii; 149 150 PetscFunctionBeginUser; 151 PetscCall(TSGetApplicationContext(ts, &honee)); 152 PetscCall(PetscObjectTypeCompare((PetscObject)ctx->viewer, PETSCVIEWERASCII, &is_ascii)); 153 PetscCheck(!is_ascii || honee->app_ctx->test_type != TESTTYPE_NONE, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, 154 "Turbulence spanwise statistics does not support ASCII viewers"); 155 156 PetscCall(SpanwiseStatisticsSetupInitialize(honee, honee->app_ctx->degree, TURB_NUM_COMPONENTS, prefix, &stats_setup_data, &spanstats)); 157 158 { // Create Section for data 159 PetscSection section; 160 161 PetscCall(DMGetLocalSection(spanstats->dm, §ion)); 162 PetscCall(PetscSectionSetFieldName(section, 0, "")); 163 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY, "MeanDensity")); 164 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE, "MeanPressure")); 165 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_SQUARED, "MeanPressureSquared")); 166 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_X, "MeanPressureVelocityX")); 167 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Y, "MeanPressureVelocityY")); 168 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_PRESSURE_VELOCITY_Z, "MeanPressureVelocityZ")); 169 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE, "MeanDensityTemperature")); 170 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_X, "MeanDensityTemperatureFluxX")); 171 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Y, "MeanDensityTemperatureFluxY")); 172 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_DENSITY_TEMPERATURE_FLUX_Z, "MeanDensityTemperatureFluxZ")); 173 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_X, "MeanMomentumX")); 174 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Y, "MeanMomentumY")); 175 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUM_Z, "MeanMomentumZ")); 176 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XX, "MeanMomentumFluxXX")); 177 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YY, "MeanMomentumFluxYY")); 178 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_ZZ, "MeanMomentumFluxZZ")); 179 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_YZ, "MeanMomentumFluxYZ")); 180 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XZ, "MeanMomentumFluxXZ")); 181 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_MOMENTUMFLUX_XY, "MeanMomentumFluxXY")); 182 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_X, "MeanVelocityX")); 183 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Y, "MeanVelocityY")); 184 PetscCall(PetscSectionSetComponentName(section, 0, TURB_MEAN_VELOCITY_Z, "MeanVelocityZ")); 185 } 186 187 // Create CeedOperators for statistics collection 188 PetscCall(CreateStatisticCollectionOperator(honee, spanstats, stats_setup_data, honee->problem_data)); 189 190 PetscCall(PetscOptionsGetBool(NULL, NULL, "-ts_monitor_spanstats_turbulence_mms", &spanstats->do_mms_test, NULL)); 191 if (spanstats->do_mms_test) { 192 PetscCall(SetupMMSErrorChecking(honee, spanstats, stats_setup_data)); 193 } 194 195 ctx->data = spanstats; 196 ctx->data_destroy = (PetscCtxDestroyFn *)SpanStatsCtxDestroy; 197 198 PetscCall(SpanwiseStatisticsSetupFinalize(ts, honee, spanstats, ctx, &stats_setup_data)); 199 PetscFunctionReturn(PETSC_SUCCESS); 200 } 201 202 // @brief TSMonitor for collection and processing of turbulence statistics 203 PetscErrorCode TSMonitor_SpanwiseStatisticsTurbulence(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx) { 204 Honee honee; 205 Vec stats; 206 TSConvergedReason reason; 207 SpanStatsCtx spanstats = ctx->data; 208 PetscInt collect_interval = spanstats->collect_interval, viewer_interval = ctx->view_interval; 209 210 PetscFunctionBeginUser; 211 PetscCall(TSGetApplicationContext(ts, &honee)); 212 PetscCall(TSGetConvergedReason(ts, &reason)); 213 // Do not collect or process on the first step of the run (ie. on the initial condition) 214 if (steps == honee->app_ctx->cont_steps) PetscFunctionReturn(PETSC_SUCCESS); 215 216 PetscBool run_processing_and_viewer = (steps % viewer_interval == 0 && viewer_interval != -1) || reason != TS_CONVERGED_ITERATING; 217 218 if (steps % collect_interval == 0 || run_processing_and_viewer) { 219 PetscCall(SpanwiseStatisticsCollect(honee, spanstats, solution_time, Q)); 220 221 if (run_processing_and_viewer) { 222 PetscCall(DMSetOutputSequenceNumber(spanstats->dm, steps, solution_time)); 223 PetscCall(DMGetGlobalVector(spanstats->dm, &stats)); 224 PetscCall(SpanwiseStatisticsProcess(honee, spanstats, stats)); 225 if (honee->app_ctx->test_type == TESTTYPE_NONE) { 226 PetscCall(PetscViewerPushFormat(ctx->viewer, ctx->format)); 227 PetscCall(VecView(stats, ctx->viewer)); 228 PetscCall(PetscViewerPopFormat(ctx->viewer)); 229 { 230 PetscInt tab_level; 231 PetscViewer viewer = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)ts)); 232 CeedScalar second = honee->units->second; 233 const char *filename; 234 235 PetscCall(PetscViewerFileGetName(ctx->viewer, &filename)); 236 PetscCall(PetscObjectGetTabLevel((PetscObject)ts, &tab_level)); 237 PetscCall(PetscViewerASCIIAddTab(viewer, tab_level + 1)); 238 if (filename) { 239 PetscCall(PetscViewerASCIIPrintf(viewer, 240 "Spanwise turbulent statistics written to %s for time span (%0.12e,%0.12e] and step span [%" PetscInt_FMT 241 ",%" PetscInt_FMT "]\n", 242 filename, spanstats->initial_solution_time / second, solution_time / second, 243 spanstats->initial_solution_step, steps)); 244 } else { 245 PetscCall(PetscViewerASCIIPrintf( 246 viewer, "Spanwise statistics (%s) file written for time span (%0.12e,%0.12e] and step span [%" PetscInt_FMT ",%" PetscInt_FMT "]\n", 247 spanstats->prefix, spanstats->initial_solution_time / second, solution_time / second, spanstats->initial_solution_step, steps)); 248 } 249 PetscCall(PetscViewerASCIISubtractTab(viewer, tab_level + 1)); 250 } 251 } 252 if (honee->app_ctx->test_type == TESTTYPE_TURB_SPANSTATS && reason != TS_CONVERGED_ITERATING) { 253 PetscCall(RegressionTest(honee->app_ctx, stats)); 254 } 255 if (spanstats->do_mms_test && reason != TS_CONVERGED_ITERATING) { 256 Vec error; 257 PetscCall(VecDuplicate(stats, &error)); 258 PetscCall(ApplyCeedOperatorGlobalToGlobal(stats, error, spanstats->mms_error_ctx)); 259 PetscScalar error_sq = 0; 260 PetscCall(VecSum(error, &error_sq)); 261 PetscScalar l2_error = sqrt(error_sq); 262 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "l2 error: %.5e\n", l2_error)); 263 } 264 PetscCall(DMRestoreGlobalVector(spanstats->dm, &stats)); 265 } 266 } 267 PetscFunctionReturn(PETSC_SUCCESS); 268 } 269