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