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