xref: /honee/src/spanstats/cflpe.c (revision c4bd9656417ac488adba016fe5c2e786414c4baf)
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 
84       PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &advection_ctx));
85       collect_ctx->newt_ctx        = (struct NewtonianIdealGasContext_){0};
86       collect_ctx->diffusion_coeff = advection_ctx->diffusion_coeff;
87       PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &advection_ctx));
88     } else {
89       NewtonianIdealGasContext newtonian_ig_ctx;
90 
91       PetscCallCeed(ceed, CeedQFunctionContextGetData(problem->apply_vol_rhs.qfctx, CEED_MEM_HOST, &newtonian_ig_ctx));
92       collect_ctx->newt_ctx        = *newtonian_ig_ctx;
93       collect_ctx->diffusion_coeff = newtonian_ig_ctx->gas.mu;
94       PetscCallCeed(ceed, CeedQFunctionContextRestoreData(problem->apply_vol_rhs.qfctx, &newtonian_ig_ctx));
95     }
96 
97     PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &collect_qfctx));
98     PetscCallCeed(ceed, CeedQFunctionContextSetData(collect_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sizeof(*collect_ctx), collect_ctx));
99     PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(collect_qfctx, CEED_MEM_HOST, FreeContextPetsc));
100 
101     PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(collect_qfctx, "solution time", offsetof(struct CflPe_SpanStatsContext_, solution_time), 1,
102                                                            "Current solution time"));
103     PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(collect_qfctx, "previous time", offsetof(struct CflPe_SpanStatsContext_, previous_time), 1,
104                                                            "Previous time statistics collection was done"));
105     PetscCallCeed(ceed, CeedQFunctionContextRegisterDouble(collect_qfctx, "timestep", offsetof(struct CflPe_SpanStatsContext_, timestep), 1,
106                                                            "Current timestep taken by TS"));
107   }
108 
109   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_stats_collect, collect_qfctx));
110   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&collect_qfctx));
111   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q", num_comp_q, CEED_EVAL_INTERP));
112   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "q_data", q_data_size, CEED_EVAL_NONE));
113   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_stats_collect, "x", num_comp_x, CEED_EVAL_INTERP));
114   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_stats_collect, "v", num_comp_stats, CEED_EVAL_NONE));
115 
116   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stats_collect, NULL, NULL, &op_stats_collect));
117   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
118   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "q_data", elem_restr_qd, CEED_BASIS_NONE, q_data));
119   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "x", honee->elem_restr_x, honee->basis_x, honee->x_coord));
120   PetscCallCeed(ceed, CeedOperatorSetField(op_stats_collect, "v", stats_data->elem_restr_child_colloc, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
121 
122   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "solution time", &spanstats->solution_time_label));
123   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "previous time", &spanstats->previous_time_label));
124   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_stats_collect, "timestep", &spanstats->timestep_label));
125 
126   PetscCall(OperatorApplyContextCreate(honee->dm, spanstats->dm, honee->ceed, op_stats_collect, honee->q_ceed, NULL, NULL, NULL,
127                                        &spanstats->op_stats_collect_ctx));
128 
129   PetscCall(CeedOperatorCreateLocalVecs(op_stats_collect, DMReturnVecType(spanstats->dm), PETSC_COMM_SELF, NULL, &spanstats->Child_Stats_loc));
130   PetscCall(VecZeroEntries(spanstats->Child_Stats_loc));
131 
132   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
133   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
134   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stats_collect));
135   PetscCallCeed(ceed, CeedOperatorDestroy(&op_stats_collect));
136   PetscFunctionReturn(PETSC_SUCCESS);
137 }
138 
139 // @brief Setup for statistics collection of CFL and Peclet number
140 PetscErrorCode SpanwiseStatisticsSetup_CflPe(TS ts, PetscViewerAndFormat *ctx) {
141   Honee              honee;
142   SpanStatsSetupData stats_setup_data;
143   SpanStatsCtx       spanstats;
144   const char        *prefix = "ts_monitor_spanstats_cflpe_";
145   PetscBool          is_ascii;
146 
147   PetscFunctionBeginUser;
148   PetscCall(TSGetApplicationContext(ts, &honee));
149   PetscCall(PetscObjectTypeCompare((PetscObject)ctx->viewer, PETSCVIEWERASCII, &is_ascii));
150   PetscCheck(!is_ascii || honee->app_ctx->test_type != TESTTYPE_NONE, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP,
151              "Turbulence spanwise statistics does not support ASCII viewers");
152 
153   PetscCall(SpanwiseStatisticsSetupInitialize(honee, honee->app_ctx->degree, 6, prefix, &stats_setup_data, &spanstats));
154 
155   {  // Create Section for data
156     PetscSection section;
157 
158     PetscCall(DMGetLocalSection(spanstats->dm, &section));
159     PetscCall(PetscSectionSetFieldName(section, 0, ""));
160     PetscCall(PetscSectionSetComponentName(section, 0, 0, "MeanCFL"));
161     PetscCall(PetscSectionSetComponentName(section, 0, 1, "MeanCFLSquared"));
162     PetscCall(PetscSectionSetComponentName(section, 0, 2, "MeanCFLCubed"));
163     PetscCall(PetscSectionSetComponentName(section, 0, 3, "MeanPe"));
164     PetscCall(PetscSectionSetComponentName(section, 0, 4, "MeanPeSquared"));
165     PetscCall(PetscSectionSetComponentName(section, 0, 5, "MeanPeCubed"));
166   }
167 
168   // Create CeedOperators for statistics collection
169   PetscCall(CreateStatisticCollectionOperator(honee, spanstats, stats_setup_data, honee->problem_data));
170 
171   ctx->data         = spanstats;
172   ctx->data_destroy = (PetscCtxDestroyFn *)SpanStatsCtxDestroy;
173 
174   PetscCall(SpanwiseStatisticsSetupFinalize(ts, honee, spanstats, ctx, &stats_setup_data));
175   PetscFunctionReturn(PETSC_SUCCESS);
176 }
177 
178 // @brief TSMonitor for collection and processing of CFL and Peclet number statistics
179 PetscErrorCode TSMonitor_SpanwiseStatisticsCflPe(TS ts, PetscInt steps, PetscReal solution_time, Vec Q, PetscViewerAndFormat *ctx) {
180   Honee             honee;
181   Vec               stats;
182   TSConvergedReason reason;
183   SpanStatsCtx      spanstats        = ctx->data;
184   PetscInt          collect_interval = spanstats->collect_interval, viewer_interval = ctx->view_interval;
185   PetscReal         timestep;
186 
187   PetscFunctionBeginUser;
188   PetscCall(TSGetApplicationContext(ts, &honee));
189   PetscCall(TSGetConvergedReason(ts, &reason));
190   // Do not collect or process on the first step of the run (ie. on the initial condition)
191   if (steps == honee->app_ctx->cont_steps) PetscFunctionReturn(PETSC_SUCCESS);
192 
193   PetscBool run_processing_and_viewer = (steps % viewer_interval == 0 && viewer_interval != -1) || reason != TS_CONVERGED_ITERATING;
194 
195   if (steps % collect_interval == 0 || run_processing_and_viewer) {
196     PetscCall(TSGetTimeStep(ts, &timestep));
197     PetscCallCeed(honee->ceed, CeedOperatorSetContextDouble(spanstats->op_stats_collect_ctx->op, spanstats->timestep_label, &timestep));
198     PetscCall(SpanwiseStatisticsCollect(honee, spanstats, solution_time, Q));
199 
200     if (run_processing_and_viewer) {
201       PetscCall(DMSetOutputSequenceNumber(spanstats->dm, steps, solution_time));
202       PetscCall(DMGetGlobalVector(spanstats->dm, &stats));
203       PetscCall(SpanwiseStatisticsProcess(honee, spanstats, stats));
204 
205       if (honee->app_ctx->test_type == TESTTYPE_NONE) {
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, "Spanwise cflpe statistics written to %s for time span (%0.12e,%0.12e] and step span [%" PetscInt_FMT ",%" PetscInt_FMT "]\n",
221                 filename, spanstats->initial_solution_time / second, solution_time / second, spanstats->initial_solution_step, steps));
222           } else {
223             PetscCall(PetscViewerASCIIPrintf(
224                 viewer, "Spanwise statistics (%s) file written for time span (%0.12e,%0.12e] and step span [%" PetscInt_FMT ",%" PetscInt_FMT "]\n",
225                 spanstats->prefix, spanstats->initial_solution_time / second, solution_time / second, spanstats->initial_solution_step, steps));
226           }
227           PetscCall(PetscViewerASCIISubtractTab(viewer, tab_level + 1));
228         }
229       } else if (honee->app_ctx->test_type == TESTTYPE_SPANSTATS && reason != TS_CONVERGED_ITERATING) {
230         PetscCall(RegressionTest(honee->app_ctx, stats));
231       }
232 
233       PetscCall(DMRestoreGlobalVector(spanstats->dm, &stats));
234     }
235   }
236   PetscFunctionReturn(PETSC_SUCCESS);
237 }
238