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