xref: /honee/src/smartsim/sgs_dd_training.c (revision 8fc6ab98e527355f6c031781f7c3aa8a26105128)
1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3 
4 #include "../../qfunctions/sgs_dd_training.h"
5 
6 #include <differential_filter.h>
7 #include <navierstokes.h>
8 #include <petscdmplex.h>
9 #include <smartsim-impl.h>
10 
11 typedef struct {
12   DM                   dm_dd_training;
13   PetscInt             num_comp_dd_inputs, write_data_interval, num_filter_widths;
14   PetscScalar          filter_widths[16];
15   OperatorApplyContext op_training_data_calc_ctx;
16   DiffFilterData       diff_filter;
17   NodalProjectionData  filtered_grad_velo_proj;
18   size_t               training_data_array_dims[2];
19   PetscBool            overwrite_training_data;
20 } *SGS_DD_TrainingData;
21 
22 #define SGS_DD_TRAIN_KEY "SGS Data Driven Training"
23 
24 static PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData *sgs_dd_train) {
25   SGS_DD_TrainingData sgs_dd_train_ = *sgs_dd_train;
26 
27   PetscFunctionBeginUser;
28   if (!sgs_dd_train_) PetscFunctionReturn(PETSC_SUCCESS);
29   PetscCall(OperatorApplyContextDestroy(sgs_dd_train_->op_training_data_calc_ctx));
30   PetscCall(NodalProjectionDataDestroy(&sgs_dd_train_->filtered_grad_velo_proj));
31   PetscCall(DMDestroy(&sgs_dd_train_->dm_dd_training));
32   PetscCall(DifferentialFilterDataDestroy(&sgs_dd_train_->diff_filter));
33   PetscCall(PetscFree(sgs_dd_train_));
34   *sgs_dd_train = NULL;
35   PetscFunctionReturn(PETSC_SUCCESS);
36 }
37 
38 typedef struct {
39   CeedElemRestriction  elem_restr_grid_aniso;
40   CeedVector           grid_aniso_ceed;
41   CeedQFunctionContext sgs_dd_train_qfctx;
42 } *SGS_DD_TrainingSetupData;
43 
44 static PetscErrorCode SGS_DD_TrainingSetupDataDestroy(SGS_DD_TrainingSetupData sgs_dd_train_setup_data) {
45   Ceed ceed;
46 
47   PetscFunctionBeginUser;
48   PetscCall(CeedElemRestrictionGetCeed(sgs_dd_train_setup_data->elem_restr_grid_aniso, &ceed));
49 
50   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_train_setup_data->elem_restr_grid_aniso));
51   PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_train_setup_data->grid_aniso_ceed));
52   PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_train_setup_data->sgs_dd_train_qfctx));
53   PetscCall(PetscFree(sgs_dd_train_setup_data));
54   PetscCheck(CeedDestroy(&ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "Destroying Ceed object failed");
55   PetscFunctionReturn(PETSC_SUCCESS);
56 }
57 
58 // @brief Create DM for storing data-drive SGS model inputs
59 static PetscErrorCode SGS_DD_TrainingCreateDM(DM dm_source, DM *dm_dd_training, PetscInt degree, PetscInt q_extra, PetscInt *num_components) {
60   PetscSection section;
61 
62   PetscFunctionBeginUser;
63   *num_components = 12;
64 
65   PetscCall(DMClone(dm_source, dm_dd_training));
66   PetscCall(DMSetMatrixPreallocateSkip(*dm_dd_training, PETSC_TRUE));
67   PetscCall(PetscObjectSetName((PetscObject)*dm_dd_training, "Data-Driven SGS Training Data"));
68 
69   PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_dd_training));
70 
71   PetscCall(DMGetLocalSection(*dm_dd_training, &section));
72   PetscCall(PetscSectionSetFieldName(section, 0, "Data-Driven SGS Training Data"));
73   PetscCall(PetscSectionSetComponentName(section, 0, 0, "SGSInput1"));
74   PetscCall(PetscSectionSetComponentName(section, 0, 1, "SGSInput2"));
75   PetscCall(PetscSectionSetComponentName(section, 0, 2, "SGSInput3"));
76   PetscCall(PetscSectionSetComponentName(section, 0, 3, "SGSInput4"));
77   PetscCall(PetscSectionSetComponentName(section, 0, 4, "SGSInput5"));
78   PetscCall(PetscSectionSetComponentName(section, 0, 5, "SGSInput6"));
79   PetscCall(PetscSectionSetComponentName(section, 0, 6, "FilteredSGSXX"));
80   PetscCall(PetscSectionSetComponentName(section, 0, 7, "FilteredSGSYY"));
81   PetscCall(PetscSectionSetComponentName(section, 0, 8, "FilteredSGSZZ"));
82   PetscCall(PetscSectionSetComponentName(section, 0, 9, "FilteredSGSYZ"));
83   PetscCall(PetscSectionSetComponentName(section, 0, 10, "FilteredSGSXZ"));
84   PetscCall(PetscSectionSetComponentName(section, 0, 11, "FilteredSGSXY"));
85   PetscFunctionReturn(PETSC_SUCCESS);
86 };
87 
88 // @brief Create CeedOperator to calculate training data for data-drive SGS model at nodes
89 static PetscErrorCode SetupTrainingDataCalculation(Ceed ceed, Honee honee, ProblemData problem, SGS_DD_TrainingSetupData sgs_dd_train_setup_data) {
90   SGS_DD_TrainingData sgs_dd_train;
91   CeedQFunction       qf_sgs_dd_train;
92   CeedOperator        op_sgs_dd_train;
93   CeedInt             num_comp_grad_velo, num_comp_grid_aniso;
94   CeedVector          inv_multiplicity, filtered_fields;
95   CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs_train;
96   DMLabel             domain_label = NULL;
97   PetscInt            label_value = 0, height = 0, dm_field = 0;
98 
99   PetscFunctionBeginUser;
100   PetscCall(HoneeGetContainer(honee, SGS_DD_TRAIN_KEY, &sgs_dd_train));
101   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_train_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso));
102 
103   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_train->dm_dd_training, domain_label, label_value, height, dm_field, &elem_restr_sgs_train));
104   PetscCall(GetInverseMultiplicity(ceed, sgs_dd_train->dm_dd_training, domain_label, label_value, height, dm_field, PETSC_TRUE,
105                                    &elem_restr_inv_multiplicity, &inv_multiplicity));
106 
107   CeedElemRestriction elem_restr_filtered_state;
108   CeedInt             num_comp_filtered_state;
109   {  // -- Setup filtered velocity gradient projection
110     CeedBasis         basis_filtered_state;
111     CeedOperatorField op_field;
112     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->diff_filter->op_rhs_ctx->op, "v0", &op_field));
113     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_filtered_state, &basis_filtered_state, NULL));
114     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_state, &num_comp_filtered_state));
115     PetscCall(VelocityGradientProjectionSetup(ceed, honee, problem, STATEVAR_PRIMITIVE, elem_restr_filtered_state, basis_filtered_state,
116                                               &sgs_dd_train->filtered_grad_velo_proj));
117     PetscCallCeed(ceed, CeedBasisDestroy(&basis_filtered_state));
118     // Get velocity gradient information
119     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->filtered_grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field));
120     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo));
121     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo));
122   }
123 
124   CeedElemRestriction elem_restr_filtered_velo_prod;
125   CeedInt             num_comp_filtered_velo_prod;
126   {  // Get filtered velocity product information
127     CeedOperatorField op_field;
128     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->diff_filter->op_rhs_ctx->op, "v1", &op_field));
129     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_filtered_velo_prod));
130     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_velo_prod, &num_comp_filtered_velo_prod));
131   }
132 
133   // -- Create operator for generating training data at nodes
134   // Differential Filter only provides filtered primitive variables
135   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSGS_DDAnisotropicTrainingDataNodal_Prim,
136                                                   ComputeSGS_DDAnisotropicTrainingDataNodal_Prim_loc, &qf_sgs_dd_train));
137 
138   PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_train, sgs_dd_train_setup_data->sgs_dd_train_qfctx));
139   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "q", num_comp_filtered_state, CEED_EVAL_NONE));
140   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "velocity product", num_comp_filtered_velo_prod, CEED_EVAL_NONE));
141   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE));
142   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE));
143   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "inverse multiplicity", 1, CEED_EVAL_NONE));
144   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_train, "training data", sgs_dd_train->num_comp_dd_inputs, CEED_EVAL_NONE));
145 
146   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_filtered_state, &filtered_fields, NULL));
147   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_train, NULL, NULL, &op_sgs_dd_train));
148   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "q", elem_restr_filtered_state, CEED_BASIS_NONE, filtered_fields));
149   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "velocity product", elem_restr_filtered_velo_prod, CEED_BASIS_NONE, filtered_fields));
150   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
151   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "anisotropy tensor", sgs_dd_train_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE,
152                                            sgs_dd_train_setup_data->grid_aniso_ceed));
153   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity));
154   PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "training data", elem_restr_sgs_train, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
155 
156   PetscCall(OperatorApplyContextCreate(sgs_dd_train->filtered_grad_velo_proj->dm, sgs_dd_train->dm_dd_training, ceed, op_sgs_dd_train, NULL, NULL,
157                                        NULL, NULL, &sgs_dd_train->op_training_data_calc_ctx));
158 
159   PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity));
160   PetscCallCeed(ceed, CeedVectorDestroy(&filtered_fields));
161   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity));
162   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_filtered_state));
163   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo));
164   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_filtered_velo_prod));
165   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_train));
166   PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_train));
167   PetscFunctionReturn(PETSC_SUCCESS);
168 }
169 
170 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, Honee honee) {
171   SGS_DDTrainingContext    sgsdd_train_ctx;
172   SGS_DD_TrainingSetupData sgs_dd_train_setup_data;
173   SGS_DD_TrainingData      sgs_dd_train;
174   ProblemData              problem = honee->problem_data;
175 
176   PetscFunctionBeginUser;
177 
178   PetscCall(PetscNew(&sgsdd_train_ctx));
179   PetscCall(PetscNew(&sgs_dd_train_setup_data));
180   PetscCall(PetscNew(&sgs_dd_train));
181   PetscCall(HoneeSetContainer(honee, SGS_DD_TRAIN_KEY, sgs_dd_train, (PetscCtxDestroyFn *)SGS_DD_TrainingDataDestroy));
182   PetscCall(DifferentialFilterSetup(honee, &sgs_dd_train->diff_filter));
183 
184   sgs_dd_train->overwrite_training_data = PETSC_TRUE;
185   sgs_dd_train->write_data_interval     = 1;
186   sgs_dd_train->num_filter_widths       = sizeof(sgs_dd_train->filter_widths) / sizeof(sgs_dd_train->filter_widths[0]);
187   PetscOptionsBegin(honee->comm, NULL, "SGS Data-Driven Training Options", NULL);
188   PetscCall(PetscOptionsInt("-sgs_train_write_data_interval", "Number of timesteps between writing data into database", NULL,
189                             sgs_dd_train->write_data_interval, &sgs_dd_train->write_data_interval, NULL));
190   PetscCall(PetscOptionsBool("-sgs_train_overwrite_data", "Overwrite old training data in the database", NULL, sgs_dd_train->overwrite_training_data,
191                              &sgs_dd_train->overwrite_training_data, NULL));
192   PetscCall(PetscOptionsRealArray("-sgs_train_filter_width_scales", "Scales of each filter width put into training database", NULL,
193                                   sgs_dd_train->filter_widths, &sgs_dd_train->num_filter_widths, NULL));
194   PetscOptionsEnd();
195 
196   // -- Create DM for storing training data
197   PetscCall(SGS_DD_TrainingCreateDM(honee->dm, &sgs_dd_train->dm_dd_training, honee->app_ctx->degree, honee->app_ctx->q_extra,
198                                     &sgs_dd_train->num_comp_dd_inputs));
199 
200   {  // -- Create QFunction Context
201     NewtonianIdealGasContext gas;
202     PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &gas));
203     sgsdd_train_ctx->gas = *gas;
204     PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &gas));
205     PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &sgs_dd_train_setup_data->sgs_dd_train_qfctx));
206     PetscCallCeed(ceed, CeedQFunctionContextSetData(sgs_dd_train_setup_data->sgs_dd_train_qfctx, CEED_MEM_HOST, CEED_USE_POINTER,
207                                                     sizeof(*sgsdd_train_ctx), sgsdd_train_ctx));
208     PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_train_setup_data->sgs_dd_train_qfctx, CEED_MEM_HOST, FreeContextPetsc));
209   }
210 
211   {  // -- Send training data array info to SmartRedis database
212     PetscMPIInt  rank, num_ranks;
213     SmartSimData smartsim;
214     PetscCall(HoneeGetSmartSimData(honee, &smartsim));
215     PetscCallMPI(MPI_Comm_rank(honee->comm, &rank));
216     PetscCallMPI(MPI_Comm_size(honee->comm, &num_ranks));
217 
218     {
219       PetscSection global_section;
220       PetscInt     num_dofs, num_comps, local_min_max[2] = {0.}, global_min_max[2] = {0.};
221 
222       PetscCall(DMGetGlobalSection(sgs_dd_train->dm_dd_training, &global_section));
223       PetscCall(DMGetGlobalVectorInfo(sgs_dd_train->dm_dd_training, &num_dofs, NULL, NULL));
224       PetscCall(PetscSectionGetFieldComponents(global_section, 0, &num_comps));
225       local_min_max[0] = num_dofs;
226       PetscCall(PetscGlobalMinMaxInt(honee->comm, local_min_max, global_min_max));
227 
228       sgs_dd_train->training_data_array_dims[0] = global_min_max[0] / num_comps;
229       sgs_dd_train->training_data_array_dims[1] = num_comps;
230     }
231 
232     if (rank % smartsim->collocated_database_num_ranks == 0) {
233       {  // Communicate info on simulation size
234         const char tensor_name[]  = "sizeInfo";
235         size_t     array_info_dim = 6;
236         PetscInt64 array_info[6] = {0}, num_features = 6;
237 
238         array_info[0] = sgs_dd_train->training_data_array_dims[0];
239         array_info[1] = sgs_dd_train->training_data_array_dims[1];
240         array_info[2] = num_features;
241         array_info[3] = num_ranks;
242         array_info[4] = smartsim->collocated_database_num_ranks;
243         array_info[5] = rank;
244 
245         PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
246         PetscCallSmartRedis(
247             put_tensor(smartsim->client, tensor_name, strlen(tensor_name), array_info, &array_info_dim, 1, SRTensorTypeInt64, SRMemLayoutContiguous));
248         PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name)));
249         PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
250       }
251 
252       {  // Send array that communicates if tensors are overwritten in database
253         const char tensor_name[]       = "tensor-ow";
254         PetscInt64 tensor_overwrite[2] = {sgs_dd_train->overwrite_training_data};
255         size_t     dim_2[1]            = {2};
256 
257         PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
258         PetscCallSmartRedis(
259             put_tensor(smartsim->client, tensor_name, strlen(tensor_name), tensor_overwrite, dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous));
260         PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name)));
261         PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
262       }
263 
264       {  // Communicate number of filter widths used
265         const char tensor_name[]     = "num_filter_widths";
266         PetscInt64 num_filter_widths = sgs_dd_train->num_filter_widths;
267         size_t     dim_2             = 1;
268 
269         PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
270         PetscCallSmartRedis(
271             put_tensor(smartsim->client, tensor_name, strlen(tensor_name), &num_filter_widths, &dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous));
272         PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name)));
273         PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
274       }
275     }
276   }
277 
278   // -- Compute and store anisotropy tensor
279   PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, honee, &sgs_dd_train_setup_data->elem_restr_grid_aniso,
280                                                      &sgs_dd_train_setup_data->grid_aniso_ceed));
281 
282   // -- Create Nodal Evaluation Operator
283   PetscCall(SetupTrainingDataCalculation(ceed, honee, problem, sgs_dd_train_setup_data));
284 
285   PetscCall(SGS_DD_TrainingSetupDataDestroy(sgs_dd_train_setup_data));
286   PetscFunctionReturn(PETSC_SUCCESS);
287 }
288 
289 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx) {
290   Honee               honee = (Honee)ctx;
291   Ceed                ceed  = honee->ceed;
292   SGS_DD_TrainingData sgs_dd_train;
293   SmartSimData        smartsim;
294   Vec                 TrainingData;
295   PetscMPIInt         rank;
296 
297   PetscFunctionBeginUser;
298   PetscCall(HoneeGetSmartSimData(honee, &smartsim));
299   PetscCall(HoneeGetContainer(honee, SGS_DD_TRAIN_KEY, &sgs_dd_train));
300   PetscCallMPI(MPI_Comm_rank(honee->comm, &rank));
301 
302   if (step_num % sgs_dd_train->write_data_interval != 0) PetscFunctionReturn(PETSC_SUCCESS);
303   PetscCall(DMGetGlobalVector(sgs_dd_train->dm_dd_training, &TrainingData));
304 
305   for (PetscInt filter_index = 0; filter_index < sgs_dd_train->num_filter_widths; filter_index++) {
306     PetscCall(PetscLogEventBegin(HONEE_TrainDataCompute, 0, 0, 0, 0));
307     {  // -- Compute and assemble training data
308       Vec          FilteredVelocityGradient, FilteredFields, FilteredFields_loc;
309       PetscMemType filtered_fields_mem_type;
310       CeedVector   filtered_fields;
311 
312       {  // Set filter width for the current solve
313         double       filter_width_scaling[3];
314         CeedOperator op_mat;
315         Mat          A_mat;
316 
317         for (int j = 0; j < 3; j++) filter_width_scaling[j] = sgs_dd_train->filter_widths[filter_index];
318         PetscCall(KSPGetOperators(sgs_dd_train->diff_filter->ksp, &A_mat, NULL));
319         PetscCall(MatCeedGetCeedOperators(A_mat, &op_mat, NULL));
320         PetscCall(CeedOperatorSetContextDouble(op_mat, sgs_dd_train->diff_filter->filter_width_scaling_label, filter_width_scaling));
321       }
322 
323       PetscCall(DMGetGlobalVector(sgs_dd_train->diff_filter->dm_filter, &FilteredFields));
324       PetscCall(DMGetLocalVector(sgs_dd_train->diff_filter->dm_filter, &FilteredFields_loc));
325 
326       PetscCall(DifferentialFilterApply(honee, sgs_dd_train->diff_filter, solution_time, Q, FilteredFields));
327       PetscCall(DMGlobalToLocal(sgs_dd_train->diff_filter->dm_filter, FilteredFields, INSERT_VALUES, FilteredFields_loc));
328 
329       PetscCall(DMGetGlobalVector(sgs_dd_train->filtered_grad_velo_proj->dm, &FilteredVelocityGradient));
330       PetscCall(VelocityGradientProjectionApply(sgs_dd_train->filtered_grad_velo_proj, FilteredFields_loc, FilteredVelocityGradient));
331 
332       {
333         CeedOperatorField op_field;
334 
335         PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->op_training_data_calc_ctx->op, "q", &op_field));
336         PetscCallCeed(ceed, CeedOperatorFieldGetVector(op_field, &filtered_fields));
337       }
338 
339       PetscCall(VecPetscToCeed(FilteredFields_loc, &filtered_fields_mem_type, filtered_fields));  // filtered_fields is an implicit input
340       PetscCall(ApplyCeedOperatorGlobalToGlobal(FilteredVelocityGradient, TrainingData, sgs_dd_train->op_training_data_calc_ctx));
341       PetscCall(VecCeedToPetsc(filtered_fields, filtered_fields_mem_type, FilteredFields_loc));
342 
343       PetscCall(DMRestoreGlobalVector(sgs_dd_train->filtered_grad_velo_proj->dm, &FilteredVelocityGradient));
344       PetscCall(DMRestoreGlobalVector(sgs_dd_train->diff_filter->dm_filter, &FilteredFields));
345       PetscCall(DMRestoreLocalVector(sgs_dd_train->diff_filter->dm_filter, &FilteredFields_loc));
346       PetscCallCeed(ceed, CeedVectorDestroy(&filtered_fields));
347     }
348     PetscCall(PetscLogEventEnd(HONEE_TrainDataCompute, 0, 0, 0, 0));
349 
350     {  // -- Send training data to SmartSim
351       char   array_key[PETSC_MAX_PATH_LEN];
352       size_t array_key_len;
353 
354       if (sgs_dd_train->overwrite_training_data) {
355         PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT, smartsim->rank_id_name, filter_index));
356       } else {
357         PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT "%" PetscInt_FMT, smartsim->rank_id_name, step_num, filter_index));
358       }
359       PetscCall(PetscStrlen(array_key, &array_key_len));
360 
361       {
362         const PetscScalar *training_data;
363         PetscCall(VecGetArrayRead(TrainingData, &training_data));
364         PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Write, 0, 0, 0, 0));
365         PetscCallSmartRedis(put_tensor(smartsim->client, array_key, array_key_len, (void *)training_data, sgs_dd_train->training_data_array_dims, 2,
366                                        SRTensorTypeDouble, SRMemLayoutContiguous));
367         PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Write, 0, 0, 0, 0));
368         PetscCall(VecRestoreArrayRead(TrainingData, &training_data));
369       }
370     }
371   }
372 
373   if (rank % smartsim->collocated_database_num_ranks == 0) {
374     const char tensor_name[] = "step";
375     size_t     dim_2[1]      = {2};
376     PetscInt64 step_array[2] = {step_num, step_num};
377 
378     PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
379     PetscCallSmartRedis(
380         put_tensor(smartsim->client, tensor_name, strlen(tensor_name), step_array, dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous));
381     PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
382   }
383 
384   PetscCall(DMRestoreGlobalVector(sgs_dd_train->dm_dd_training, &TrainingData));
385   PetscFunctionReturn(PETSC_SUCCESS);
386 }
387 
388 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts) {
389   Honee        honee;
390   const char   check_run_key[]   = "check-run";
391   PetscReal    check_run[2]      = {1};
392   const size_t check_run_dims[1] = {2};
393   size_t       check_run_key_size;
394   SmartSimData smartsim;
395 
396   PetscFunctionBeginUser;
397   PetscCall(PetscStrlen(check_run_key, &check_run_key_size));
398   PetscCall(TSGetApplicationContext(ts, &honee));
399   PetscCall(HoneeGetSmartSimData(honee, &smartsim));
400 
401   PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
402   PetscCallSmartRedis(
403       unpack_tensor(smartsim->client, check_run_key, check_run_key_size, check_run, check_run_dims, 1, SRTensorTypeDouble, SRMemLayoutContiguous));
404   PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0));
405   if (check_run[0] == 0) {
406     PetscCall(PetscPrintf(honee->comm, "-- Simulation stopped by 'check-run' tensor in Redis database\n"));
407     PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER));
408   }
409 
410   PetscFunctionReturn(PETSC_SUCCESS);
411 }
412