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