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 PetscCheck(CeedDestroy(&ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "Destroying Ceed object failed"); 28 PetscFunctionReturn(PETSC_SUCCESS); 29 } 30 31 // @brief Create DM for storing data-drive SGS model inputs 32 static PetscErrorCode SGS_DD_TrainingCreateDM(DM dm_source, DM *dm_dd_training, PetscInt degree, PetscInt q_extra, PetscInt *num_components) { 33 PetscSection section; 34 35 PetscFunctionBeginUser; 36 *num_components = 12; 37 38 PetscCall(DMClone(dm_source, dm_dd_training)); 39 PetscCall(DMSetMatrixPreallocateSkip(*dm_dd_training, PETSC_TRUE)); 40 PetscCall(PetscObjectSetName((PetscObject)*dm_dd_training, "Data-Driven SGS Training Data")); 41 42 PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_dd_training)); 43 44 PetscCall(DMGetLocalSection(*dm_dd_training, §ion)); 45 PetscCall(PetscSectionSetFieldName(section, 0, "Data-Driven SGS Training Data")); 46 PetscCall(PetscSectionSetComponentName(section, 0, 0, "SGSInput1")); 47 PetscCall(PetscSectionSetComponentName(section, 0, 1, "SGSInput2")); 48 PetscCall(PetscSectionSetComponentName(section, 0, 2, "SGSInput3")); 49 PetscCall(PetscSectionSetComponentName(section, 0, 3, "SGSInput4")); 50 PetscCall(PetscSectionSetComponentName(section, 0, 4, "SGSInput5")); 51 PetscCall(PetscSectionSetComponentName(section, 0, 5, "SGSInput6")); 52 PetscCall(PetscSectionSetComponentName(section, 0, 6, "FilteredSGSXX")); 53 PetscCall(PetscSectionSetComponentName(section, 0, 7, "FilteredSGSYY")); 54 PetscCall(PetscSectionSetComponentName(section, 0, 8, "FilteredSGSZZ")); 55 PetscCall(PetscSectionSetComponentName(section, 0, 9, "FilteredSGSYZ")); 56 PetscCall(PetscSectionSetComponentName(section, 0, 10, "FilteredSGSXZ")); 57 PetscCall(PetscSectionSetComponentName(section, 0, 11, "FilteredSGSXY")); 58 PetscFunctionReturn(PETSC_SUCCESS); 59 }; 60 61 // @brief Create CeedOperator to calculate training data for data-drive SGS model at nodes 62 static PetscErrorCode SetupTrainingDataCalculation(Ceed ceed, Honee honee, ProblemData problem, SGS_DD_TrainingSetupData sgs_dd_train_setup_data) { 63 SGS_DD_TrainingData sgs_dd_train = honee->sgs_dd_train; 64 CeedQFunction qf_sgs_dd_train; 65 CeedOperator op_sgs_dd_train; 66 CeedInt num_comp_grad_velo, num_comp_grid_aniso; 67 CeedVector inv_multiplicity, filtered_fields; 68 CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs_train; 69 DMLabel domain_label = NULL; 70 PetscInt label_value = 0, height = 0, dm_field = 0; 71 72 PetscFunctionBeginUser; 73 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_train_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso)); 74 75 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_train->dm_dd_training, domain_label, label_value, height, dm_field, &elem_restr_sgs_train)); 76 PetscCall(GetInverseMultiplicity(ceed, sgs_dd_train->dm_dd_training, domain_label, label_value, height, dm_field, PETSC_TRUE, 77 &elem_restr_inv_multiplicity, &inv_multiplicity)); 78 79 CeedElemRestriction elem_restr_filtered_state; 80 CeedInt num_comp_filtered_state; 81 { // -- Setup filtered velocity gradient projection 82 CeedBasis basis_filtered_state; 83 CeedOperatorField op_field; 84 PetscCallCeed(ceed, CeedOperatorGetFieldByName(honee->diff_filter->op_rhs_ctx->op, "v0", &op_field)); 85 PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_filtered_state, &basis_filtered_state, NULL)); 86 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_state, &num_comp_filtered_state)); 87 PetscCall(VelocityGradientProjectionSetup(ceed, honee, problem, STATEVAR_PRIMITIVE, elem_restr_filtered_state, basis_filtered_state, 88 &sgs_dd_train->filtered_grad_velo_proj)); 89 PetscCallCeed(ceed, CeedBasisDestroy(&basis_filtered_state)); 90 // Get velocity gradient information 91 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->filtered_grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field)); 92 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo)); 93 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo)); 94 } 95 96 CeedElemRestriction elem_restr_filtered_velo_prod; 97 CeedInt num_comp_filtered_velo_prod; 98 { // Get filtered velocity product information 99 CeedOperatorField op_field; 100 PetscCallCeed(ceed, CeedOperatorGetFieldByName(honee->diff_filter->op_rhs_ctx->op, "v1", &op_field)); 101 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_filtered_velo_prod)); 102 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_velo_prod, &num_comp_filtered_velo_prod)); 103 } 104 105 // -- Create operator for generating training data at nodes 106 // Differential Filter only provides filtered primitive variables 107 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSGS_DDAnisotropicTrainingDataNodal_Prim, 108 ComputeSGS_DDAnisotropicTrainingDataNodal_Prim_loc, &qf_sgs_dd_train)); 109 110 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_train, sgs_dd_train_setup_data->sgs_dd_train_qfctx)); 111 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "q", num_comp_filtered_state, CEED_EVAL_NONE)); 112 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "velocity product", num_comp_filtered_velo_prod, CEED_EVAL_NONE)); 113 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE)); 114 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 115 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "inverse multiplicity", 1, CEED_EVAL_NONE)); 116 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_train, "training data", sgs_dd_train->num_comp_dd_inputs, CEED_EVAL_NONE)); 117 118 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_filtered_state, &filtered_fields, NULL)); 119 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_train, NULL, NULL, &op_sgs_dd_train)); 120 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "q", elem_restr_filtered_state, CEED_BASIS_NONE, filtered_fields)); 121 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "velocity product", elem_restr_filtered_velo_prod, CEED_BASIS_NONE, filtered_fields)); 122 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 123 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "anisotropy tensor", sgs_dd_train_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 124 sgs_dd_train_setup_data->grid_aniso_ceed)); 125 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 126 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "training data", elem_restr_sgs_train, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 127 128 PetscCall(OperatorApplyContextCreate(sgs_dd_train->filtered_grad_velo_proj->dm, sgs_dd_train->dm_dd_training, ceed, op_sgs_dd_train, NULL, NULL, 129 NULL, NULL, &sgs_dd_train->op_training_data_calc_ctx)); 130 131 PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity)); 132 PetscCallCeed(ceed, CeedVectorDestroy(&filtered_fields)); 133 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity)); 134 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_filtered_state)); 135 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo)); 136 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_filtered_velo_prod)); 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, Honee honee, ProblemData problem) { 143 SGS_DDTrainingContext sgsdd_train_qfctx; 144 SGS_DD_TrainingSetupData sgs_dd_train_setup_data; 145 146 PetscFunctionBeginUser; 147 if (!honee->diff_filter) PetscCall(DifferentialFilterSetup(ceed, honee, problem)); 148 if (!honee->smartsim) PetscCall(SmartSimSetup(honee)); 149 150 PetscCall(PetscNew(&sgsdd_train_qfctx)); 151 PetscCall(PetscNew(&sgs_dd_train_setup_data)); 152 PetscCall(PetscNew(&honee->sgs_dd_train)); 153 SGS_DD_TrainingData sgs_dd_train = honee->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(honee->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(honee->dm, &sgs_dd_train->dm_dd_training, honee->app_ctx->degree, honee->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.qfctx, CEED_MEM_HOST, &gas)); 174 sgsdd_train_qfctx->gas = *gas; 175 PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &gas)); 176 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->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 = honee->smartsim; 185 PetscCallMPI(MPI_Comm_rank(honee->comm, &rank)); 186 PetscCallMPI(MPI_Comm_size(honee->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(honee->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(HONEE_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(HONEE_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(HONEE_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(HONEE_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(HONEE_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(HONEE_SmartRedis_Meta, 0, 0, 0, 0)); 244 } 245 } 246 } 247 248 // -- Compute and store anisotropy tensor 249 PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, honee, &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, honee, 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 Honee honee = (Honee)ctx; 261 Ceed ceed = honee->ceed; 262 SGS_DD_TrainingData sgs_dd_train = honee->sgs_dd_train; 263 SmartSimData smartsim = honee->smartsim; 264 Vec TrainingData; 265 PetscMPIInt rank; 266 267 PetscFunctionBeginUser; 268 269 PetscCallMPI(MPI_Comm_rank(honee->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(HONEE_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(honee->diff_filter->ksp, &A_mat, NULL)); 288 PetscCall(MatCeedGetCeedOperators(A_mat, &op_mat, NULL)); 289 PetscCall(CeedOperatorSetContextDouble(op_mat, honee->diff_filter->filter_width_scaling_label, filter_width_scaling)); 290 } 291 292 PetscCall(DMGetGlobalVector(honee->diff_filter->dm_filter, &FilteredFields)); 293 PetscCall(DMGetLocalVector(honee->diff_filter->dm_filter, &FilteredFields_loc)); 294 295 PetscCall(DifferentialFilterApply(honee, solution_time, Q, FilteredFields)); 296 PetscCall(DMGlobalToLocal(honee->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(honee->diff_filter->dm_filter, &FilteredFields)); 314 PetscCall(DMRestoreLocalVector(honee->diff_filter->dm_filter, &FilteredFields_loc)); 315 PetscCallCeed(ceed, CeedVectorDestroy(&filtered_fields)); 316 } 317 PetscCall(PetscLogEventEnd(HONEE_TrainDataCompute, 0, 0, 0, 0)); 318 319 { // -- Send training data to SmartSim 320 char array_key[PETSC_MAX_PATH_LEN]; 321 size_t array_key_len; 322 323 if (sgs_dd_train->overwrite_training_data) { 324 PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT, smartsim->rank_id_name, filter_index)); 325 } else { 326 PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT "%" PetscInt_FMT, smartsim->rank_id_name, step_num, filter_index)); 327 } 328 PetscCall(PetscStrlen(array_key, &array_key_len)); 329 330 { 331 const PetscScalar *training_data; 332 PetscCall(VecGetArrayRead(TrainingData, &training_data)); 333 PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Train, 0, 0, 0, 0)); 334 PetscCallSmartRedis(put_tensor(smartsim->client, array_key, array_key_len, (void *)training_data, sgs_dd_train->training_data_array_dims, 2, 335 SRTensorTypeDouble, SRMemLayoutContiguous)); 336 PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Train, 0, 0, 0, 0)); 337 PetscCall(VecRestoreArrayRead(TrainingData, &training_data)); 338 } 339 } 340 } 341 342 if (rank % smartsim->collocated_database_num_ranks == 0) { 343 const char tensor_name[] = "step"; 344 size_t dim_2[1] = {2}; 345 PetscInt64 step_array[2] = {step_num, step_num}; 346 347 PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0)); 348 PetscCallSmartRedis( 349 put_tensor(smartsim->client, tensor_name, strlen(tensor_name), step_array, dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous)); 350 PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0)); 351 } 352 353 PetscCall(DMRestoreGlobalVector(honee->sgs_dd_train->dm_dd_training, &TrainingData)); 354 PetscFunctionReturn(PETSC_SUCCESS); 355 } 356 357 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts) { 358 Honee honee; 359 const char check_run_key[] = "check-run"; 360 PetscReal check_run[2] = {1}; 361 const size_t check_run_dims[1] = {2}; 362 size_t check_run_key_size; 363 364 PetscFunctionBeginUser; 365 PetscCall(PetscStrlen(check_run_key, &check_run_key_size)); 366 PetscCall(TSGetApplicationContext(ts, &honee)); 367 SmartSimData smartsim = honee->smartsim; 368 369 PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0)); 370 PetscCallSmartRedis( 371 unpack_tensor(smartsim->client, check_run_key, check_run_key_size, check_run, check_run_dims, 1, SRTensorTypeDouble, SRMemLayoutContiguous)); 372 PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0)); 373 if (check_run[0] == 0) { 374 PetscCall(PetscPrintf(honee->comm, "-- Simulation stopped by 'check-run' tensor in Redis database\n")); 375 PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER)); 376 } 377 378 PetscFunctionReturn(PETSC_SUCCESS); 379 } 380 381 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train) { 382 PetscFunctionBeginUser; 383 if (!sgs_dd_train) PetscFunctionReturn(PETSC_SUCCESS); 384 385 PetscCall(OperatorApplyContextDestroy(sgs_dd_train->op_training_data_calc_ctx)); 386 PetscCall(NodalProjectionDataDestroy(sgs_dd_train->filtered_grad_velo_proj)); 387 PetscCall(DMDestroy(&sgs_dd_train->dm_dd_training)); 388 PetscCall(PetscFree(sgs_dd_train)); 389 390 PetscFunctionReturn(PETSC_SUCCESS); 391 } 392