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