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