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