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(DMSetMatrixPreallocateSkip(*dm_dd_training, PETSC_TRUE)); 39 PetscCall(PetscObjectSetName((PetscObject)*dm_dd_training, "Data-Driven SGS Training Data")); 40 41 PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_dd_training)); 42 43 PetscCall(DMGetLocalSection(*dm_dd_training, §ion)); 44 PetscCall(PetscSectionSetFieldName(section, 0, "Data-Driven SGS Training Data")); 45 PetscCall(PetscSectionSetComponentName(section, 0, 0, "SGSInput1")); 46 PetscCall(PetscSectionSetComponentName(section, 0, 1, "SGSInput2")); 47 PetscCall(PetscSectionSetComponentName(section, 0, 2, "SGSInput3")); 48 PetscCall(PetscSectionSetComponentName(section, 0, 3, "SGSInput4")); 49 PetscCall(PetscSectionSetComponentName(section, 0, 4, "SGSInput5")); 50 PetscCall(PetscSectionSetComponentName(section, 0, 5, "SGSInput6")); 51 PetscCall(PetscSectionSetComponentName(section, 0, 6, "FilteredSGSXX")); 52 PetscCall(PetscSectionSetComponentName(section, 0, 7, "FilteredSGSYY")); 53 PetscCall(PetscSectionSetComponentName(section, 0, 8, "FilteredSGSZZ")); 54 PetscCall(PetscSectionSetComponentName(section, 0, 9, "FilteredSGSYZ")); 55 PetscCall(PetscSectionSetComponentName(section, 0, 10, "FilteredSGSXZ")); 56 PetscCall(PetscSectionSetComponentName(section, 0, 11, "FilteredSGSXY")); 57 PetscFunctionReturn(PETSC_SUCCESS); 58 }; 59 60 // @brief Create CeedOperator to calculate training data for data-drive SGS model at nodes 61 static PetscErrorCode SetupTrainingDataCalculation(Ceed ceed, Honee honee, CeedData ceed_data, ProblemData problem, 62 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, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_filtered_state)); 86 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_state, &num_comp_filtered_state)); 87 PetscCallCeed(ceed, CeedOperatorFieldGetBasis(op_field, &basis_filtered_state)); 88 PetscCall(VelocityGradientProjectionSetup(ceed, honee, ceed_data, problem, STATEVAR_PRIMITIVE, elem_restr_filtered_state, basis_filtered_state, 89 &sgs_dd_train->filtered_grad_velo_proj)); 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, CeedQFunctionDestroy(&qf_sgs_dd_train)); 135 PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_train)); 136 PetscFunctionReturn(PETSC_SUCCESS); 137 } 138 139 PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, Honee honee, CeedData ceed_data, ProblemData problem) { 140 SGS_DDTrainingContext sgsdd_train_qfctx; 141 SGS_DD_TrainingSetupData sgs_dd_train_setup_data; 142 143 PetscFunctionBeginUser; 144 if (!honee->diff_filter) PetscCall(DifferentialFilterSetup(ceed, honee, ceed_data, problem)); 145 if (!honee->smartsim) PetscCall(SmartSimSetup(honee)); 146 147 PetscCall(PetscNew(&sgsdd_train_qfctx)); 148 PetscCall(PetscNew(&sgs_dd_train_setup_data)); 149 PetscCall(PetscNew(&honee->sgs_dd_train)); 150 SGS_DD_TrainingData sgs_dd_train = honee->sgs_dd_train; 151 152 sgs_dd_train->overwrite_training_data = PETSC_TRUE; 153 sgs_dd_train->write_data_interval = 1; 154 sgs_dd_train->num_filter_widths = sizeof(sgs_dd_train->filter_widths) / sizeof(sgs_dd_train->filter_widths[0]); 155 PetscOptionsBegin(honee->comm, NULL, "SGS Data-Driven Training Options", NULL); 156 PetscCall(PetscOptionsInt("-sgs_train_write_data_interval", "Number of timesteps between writing data into database", NULL, 157 sgs_dd_train->write_data_interval, &sgs_dd_train->write_data_interval, NULL)); 158 PetscCall(PetscOptionsBool("-sgs_train_overwrite_data", "Overwrite old training data in the database", NULL, sgs_dd_train->overwrite_training_data, 159 &sgs_dd_train->overwrite_training_data, NULL)); 160 PetscCall(PetscOptionsRealArray("-sgs_train_filter_width_scales", "Scales of each filter width put into training database", NULL, 161 sgs_dd_train->filter_widths, &sgs_dd_train->num_filter_widths, NULL)); 162 PetscOptionsEnd(); 163 164 // -- Create DM for storing training data 165 PetscCall(SGS_DD_TrainingCreateDM(honee->dm, &sgs_dd_train->dm_dd_training, honee->app_ctx->degree, honee->app_ctx->q_extra, 166 &sgs_dd_train->num_comp_dd_inputs)); 167 168 { // -- Create QFunction Context 169 NewtonianIdealGasContext gas; 170 PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &gas)); 171 sgsdd_train_qfctx->gas = *gas; 172 PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &gas)); 173 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &sgs_dd_train_setup_data->sgs_dd_train_qfctx)); 174 PetscCallCeed(ceed, CeedQFunctionContextSetData(sgs_dd_train_setup_data->sgs_dd_train_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, 175 sizeof(*sgsdd_train_qfctx), sgsdd_train_qfctx)); 176 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_train_setup_data->sgs_dd_train_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 177 } 178 179 { // -- Send training data array info to SmartRedis database 180 PetscMPIInt rank, num_ranks; 181 SmartSimData smartsim = honee->smartsim; 182 PetscCallMPI(MPI_Comm_rank(honee->comm, &rank)); 183 PetscCallMPI(MPI_Comm_size(honee->comm, &num_ranks)); 184 185 { 186 PetscSection global_section; 187 PetscInt num_dofs, num_comps, local_min_max[2] = {0.}, global_min_max[2] = {0.}; 188 189 PetscCall(DMGetGlobalSection(sgs_dd_train->dm_dd_training, &global_section)); 190 PetscCall(DMGetGlobalVectorInfo(sgs_dd_train->dm_dd_training, &num_dofs, NULL, NULL)); 191 PetscCall(PetscSectionGetFieldComponents(global_section, 0, &num_comps)); 192 local_min_max[0] = num_dofs; 193 PetscCall(PetscGlobalMinMaxInt(honee->comm, local_min_max, global_min_max)); 194 195 sgs_dd_train->training_data_array_dims[0] = global_min_max[0] / num_comps; 196 sgs_dd_train->training_data_array_dims[1] = num_comps; 197 } 198 199 if (rank % smartsim->collocated_database_num_ranks == 0) { 200 { // Communicate info on simulation size 201 const char tensor_name[] = "sizeInfo"; 202 size_t array_info_dim = 6; 203 PetscInt64 array_info[6] = {0}, num_features = 6; 204 205 array_info[0] = sgs_dd_train->training_data_array_dims[0]; 206 array_info[1] = sgs_dd_train->training_data_array_dims[1]; 207 array_info[2] = num_features; 208 array_info[3] = num_ranks; 209 array_info[4] = smartsim->collocated_database_num_ranks; 210 array_info[5] = rank; 211 212 PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 213 PetscCallSmartRedis( 214 put_tensor(smartsim->client, tensor_name, strlen(tensor_name), array_info, &array_info_dim, 1, SRTensorTypeInt64, SRMemLayoutContiguous)); 215 PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name))); 216 PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 217 } 218 219 { // Send array that communicates if tensors are overwritten in database 220 const char tensor_name[] = "tensor-ow"; 221 PetscInt64 tensor_overwrite[2] = {sgs_dd_train->overwrite_training_data}; 222 size_t dim_2[1] = {2}; 223 224 PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 225 PetscCallSmartRedis( 226 put_tensor(smartsim->client, tensor_name, strlen(tensor_name), tensor_overwrite, dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous)); 227 PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name))); 228 PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 229 } 230 231 { // Communicate number of filter widths used 232 const char tensor_name[] = "num_filter_widths"; 233 PetscInt64 num_filter_widths = sgs_dd_train->num_filter_widths; 234 size_t dim_2 = 1; 235 236 PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 237 PetscCallSmartRedis( 238 put_tensor(smartsim->client, tensor_name, strlen(tensor_name), &num_filter_widths, &dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous)); 239 PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name))); 240 PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 241 } 242 } 243 } 244 245 // -- Compute and store anisotropy tensor 246 PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, honee, ceed_data, &sgs_dd_train_setup_data->elem_restr_grid_aniso, 247 &sgs_dd_train_setup_data->grid_aniso_ceed)); 248 249 // -- Create Nodal Evaluation Operator 250 PetscCall(SetupTrainingDataCalculation(ceed, honee, ceed_data, problem, sgs_dd_train_setup_data)); 251 252 PetscCall(SGS_DD_TrainingSetupDataDestroy(sgs_dd_train_setup_data)); 253 PetscFunctionReturn(PETSC_SUCCESS); 254 } 255 256 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx) { 257 Honee honee = (Honee)ctx; 258 Ceed ceed = honee->ceed; 259 SGS_DD_TrainingData sgs_dd_train = honee->sgs_dd_train; 260 SmartSimData smartsim = honee->smartsim; 261 Vec TrainingData; 262 PetscMPIInt rank; 263 264 PetscFunctionBeginUser; 265 266 PetscCallMPI(MPI_Comm_rank(honee->comm, &rank)); 267 268 if (step_num % sgs_dd_train->write_data_interval != 0) PetscFunctionReturn(PETSC_SUCCESS); 269 PetscCall(DMGetGlobalVector(sgs_dd_train->dm_dd_training, &TrainingData)); 270 271 for (PetscInt filter_index = 0; filter_index < sgs_dd_train->num_filter_widths; filter_index++) { 272 PetscCall(PetscLogEventBegin(FLUIDS_TrainDataCompute, 0, 0, 0, 0)); 273 { // -- Compute and assemble training data 274 Vec FilteredVelocityGradient, FilteredFields, FilteredFields_loc; 275 PetscMemType filtered_fields_mem_type; 276 CeedVector filtered_fields; 277 278 { // Set filter width for the current solve 279 double filter_width_scaling[3]; 280 CeedOperator op_mat; 281 Mat A_mat; 282 283 for (int j = 0; j < 3; j++) filter_width_scaling[j] = sgs_dd_train->filter_widths[filter_index]; 284 PetscCall(KSPGetOperators(honee->diff_filter->ksp, &A_mat, NULL)); 285 PetscCall(MatCeedGetCeedOperators(A_mat, &op_mat, NULL)); 286 PetscCall(CeedOperatorSetContextDouble(op_mat, honee->diff_filter->filter_width_scaling_label, filter_width_scaling)); 287 } 288 289 PetscCall(DMGetGlobalVector(honee->diff_filter->dm_filter, &FilteredFields)); 290 PetscCall(DMGetLocalVector(honee->diff_filter->dm_filter, &FilteredFields_loc)); 291 292 PetscCall(DifferentialFilterApply(honee, solution_time, Q, FilteredFields)); 293 PetscCall(DMGlobalToLocal(honee->diff_filter->dm_filter, FilteredFields, INSERT_VALUES, FilteredFields_loc)); 294 295 PetscCall(DMGetGlobalVector(sgs_dd_train->filtered_grad_velo_proj->dm, &FilteredVelocityGradient)); 296 PetscCall(VelocityGradientProjectionApply(sgs_dd_train->filtered_grad_velo_proj, FilteredFields_loc, FilteredVelocityGradient)); 297 298 { 299 CeedOperatorField op_field; 300 301 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->op_training_data_calc_ctx->op, "q", &op_field)); 302 PetscCallCeed(ceed, CeedOperatorFieldGetVector(op_field, &filtered_fields)); 303 } 304 305 PetscCall(VecPetscToCeed(FilteredFields_loc, &filtered_fields_mem_type, filtered_fields)); // filtered_fields is an implicit input 306 PetscCall(ApplyCeedOperatorGlobalToGlobal(FilteredVelocityGradient, TrainingData, sgs_dd_train->op_training_data_calc_ctx)); 307 PetscCall(VecCeedToPetsc(filtered_fields, filtered_fields_mem_type, FilteredFields_loc)); 308 309 PetscCall(DMRestoreGlobalVector(sgs_dd_train->filtered_grad_velo_proj->dm, &FilteredVelocityGradient)); 310 PetscCall(DMRestoreGlobalVector(honee->diff_filter->dm_filter, &FilteredFields)); 311 PetscCall(DMRestoreLocalVector(honee->diff_filter->dm_filter, &FilteredFields_loc)); 312 } 313 PetscCall(PetscLogEventEnd(FLUIDS_TrainDataCompute, 0, 0, 0, 0)); 314 315 { // -- Send training data to SmartSim 316 char array_key[PETSC_MAX_PATH_LEN]; 317 size_t array_key_len; 318 319 if (sgs_dd_train->overwrite_training_data) { 320 PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT, smartsim->rank_id_name, filter_index)); 321 } else { 322 PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT "%" PetscInt_FMT, smartsim->rank_id_name, step_num, filter_index)); 323 } 324 PetscCall(PetscStrlen(array_key, &array_key_len)); 325 326 { 327 const PetscScalar *training_data; 328 PetscCall(VecGetArrayRead(TrainingData, &training_data)); 329 PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Train, 0, 0, 0, 0)); 330 PetscCallSmartRedis(put_tensor(smartsim->client, array_key, array_key_len, (void *)training_data, sgs_dd_train->training_data_array_dims, 2, 331 SRTensorTypeDouble, SRMemLayoutContiguous)); 332 PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Train, 0, 0, 0, 0)); 333 PetscCall(VecRestoreArrayRead(TrainingData, &training_data)); 334 } 335 } 336 } 337 338 if (rank % smartsim->collocated_database_num_ranks == 0) { 339 const char tensor_name[] = "step"; 340 size_t dim_2[1] = {2}; 341 PetscInt64 step_array[2] = {step_num, step_num}; 342 343 PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 344 PetscCallSmartRedis( 345 put_tensor(smartsim->client, tensor_name, strlen(tensor_name), step_array, dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous)); 346 PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 347 } 348 349 PetscCall(DMRestoreGlobalVector(honee->sgs_dd_train->dm_dd_training, &TrainingData)); 350 PetscFunctionReturn(PETSC_SUCCESS); 351 } 352 353 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts) { 354 Honee honee; 355 const char check_run_key[] = "check-run"; 356 PetscReal check_run[2] = {1}; 357 const size_t check_run_dims[1] = {2}; 358 size_t check_run_key_size; 359 360 PetscFunctionBeginUser; 361 PetscCall(PetscStrlen(check_run_key, &check_run_key_size)); 362 PetscCall(TSGetApplicationContext(ts, &honee)); 363 SmartSimData smartsim = honee->smartsim; 364 365 PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 366 PetscCallSmartRedis( 367 unpack_tensor(smartsim->client, check_run_key, check_run_key_size, check_run, check_run_dims, 1, SRTensorTypeDouble, SRMemLayoutContiguous)); 368 PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 369 if (check_run[0] == 0) { 370 PetscCall(PetscPrintf(honee->comm, "-- Simulation stopped by 'check-run' tensor in Redis database\n")); 371 PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER)); 372 } 373 374 PetscFunctionReturn(PETSC_SUCCESS); 375 } 376 377 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train) { 378 PetscFunctionBeginUser; 379 if (!sgs_dd_train) PetscFunctionReturn(PETSC_SUCCESS); 380 381 PetscCall(OperatorApplyContextDestroy(sgs_dd_train->op_training_data_calc_ctx)); 382 PetscCall(NodalProjectionDataDestroy(sgs_dd_train->filtered_grad_velo_proj)); 383 PetscCall(DMDestroy(&sgs_dd_train->dm_dd_training)); 384 PetscCall(PetscFree(sgs_dd_train)); 385 386 PetscFunctionReturn(PETSC_SUCCESS); 387 } 388