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