1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 #include "../../qfunctions/sgs_dd_training.h" 9 10 #include <petscdmplex.h> 11 12 #include "../../include/smartsim.h" 13 #include "../../navierstokes.h" 14 15 typedef struct { 16 CeedElemRestriction elem_restr_grid_aniso; 17 CeedVector grid_aniso_ceed; 18 CeedQFunctionContext sgs_dd_train_qfctx; 19 } *SGS_DD_TrainingSetupData; 20 21 static PetscErrorCode SGS_DD_TrainingSetupDataDestroy(SGS_DD_TrainingSetupData sgs_dd_train_setup_data) { 22 Ceed ceed; 23 24 PetscFunctionBeginUser; 25 PetscCall(CeedElemRestrictionGetCeed(sgs_dd_train_setup_data->elem_restr_grid_aniso, &ceed)); 26 27 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_train_setup_data->elem_restr_grid_aniso)); 28 PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_train_setup_data->grid_aniso_ceed)); 29 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_train_setup_data->sgs_dd_train_qfctx)); 30 PetscCall(PetscFree(sgs_dd_train_setup_data)); 31 PetscFunctionReturn(PETSC_SUCCESS); 32 } 33 34 // @brief Create DM for storing data-drive SGS model inputs 35 static PetscErrorCode SGS_DD_TrainingCreateDM(DM dm_source, DM *dm_dd_training, PetscInt degree, PetscInt q_extra, PetscInt *num_components) { 36 PetscSection section; 37 38 PetscFunctionBeginUser; 39 *num_components = 12; 40 41 PetscCall(DMClone(dm_source, dm_dd_training)); 42 PetscCall(PetscObjectSetName((PetscObject)*dm_dd_training, "Data-Driven SGS Training Data")); 43 44 PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_dd_training)); 45 46 PetscCall(DMGetLocalSection(*dm_dd_training, §ion)); 47 PetscCall(PetscSectionSetFieldName(section, 0, "Data-Driven SGS Training Data")); 48 PetscCall(PetscSectionSetComponentName(section, 0, 0, "SGSInput1")); 49 PetscCall(PetscSectionSetComponentName(section, 0, 1, "SGSInput2")); 50 PetscCall(PetscSectionSetComponentName(section, 0, 2, "SGSInput3")); 51 PetscCall(PetscSectionSetComponentName(section, 0, 3, "SGSInput4")); 52 PetscCall(PetscSectionSetComponentName(section, 0, 4, "SGSInput5")); 53 PetscCall(PetscSectionSetComponentName(section, 0, 5, "SGSInput6")); 54 PetscCall(PetscSectionSetComponentName(section, 0, 6, "FilteredSGSXX")); 55 PetscCall(PetscSectionSetComponentName(section, 0, 7, "FilteredSGSYY")); 56 PetscCall(PetscSectionSetComponentName(section, 0, 8, "FilteredSGSZZ")); 57 PetscCall(PetscSectionSetComponentName(section, 0, 9, "FilteredSGSYZ")); 58 PetscCall(PetscSectionSetComponentName(section, 0, 10, "FilteredSGSXZ")); 59 PetscCall(PetscSectionSetComponentName(section, 0, 11, "FilteredSGSXY")); 60 PetscFunctionReturn(PETSC_SUCCESS); 61 }; 62 63 // @brief Create CeedOperator to calculate training data for data-drive SGS model at nodes 64 static PetscErrorCode SetupTrainingDataCalculation(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem, 65 SGS_DD_TrainingSetupData sgs_dd_train_setup_data) { 66 SGS_DD_TrainingData sgs_dd_train = user->sgs_dd_train; 67 CeedQFunction qf_sgs_dd_train; 68 CeedOperator op_sgs_dd_train; 69 CeedInt num_comp_grad_velo, num_comp_grid_aniso; 70 CeedVector inv_multiplicity, filtered_fields; 71 CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs_train; 72 DMLabel domain_label = NULL; 73 PetscInt label_value = 0, height = 0, dm_field = 0; 74 75 PetscFunctionBeginUser; 76 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_train_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso)); 77 78 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_train->dm_dd_training, domain_label, label_value, height, dm_field, &elem_restr_sgs_train)); 79 PetscCall(GetInverseMultiplicity(ceed, sgs_dd_train->dm_dd_training, domain_label, label_value, height, dm_field, PETSC_TRUE, 80 &elem_restr_inv_multiplicity, &inv_multiplicity)); 81 82 CeedElemRestriction elem_restr_filtered_state; 83 CeedInt num_comp_filtered_state; 84 { // -- Setup filtered velocity gradient projection 85 CeedBasis basis_filtered_state; 86 CeedOperatorField op_field; 87 PetscCallCeed(ceed, CeedOperatorGetFieldByName(user->diff_filter->op_rhs_ctx->op, "v0", &op_field)); 88 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_filtered_state)); 89 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_state, &num_comp_filtered_state)); 90 PetscCallCeed(ceed, CeedOperatorFieldGetBasis(op_field, &basis_filtered_state)); 91 PetscCall(VelocityGradientProjectionSetup(ceed, user, ceed_data, problem, STATEVAR_PRIMITIVE, elem_restr_filtered_state, basis_filtered_state, 92 &sgs_dd_train->filtered_grad_velo_proj)); 93 // Get velocity gradient information 94 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->filtered_grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field)); 95 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo)); 96 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo)); 97 } 98 99 CeedElemRestriction elem_restr_filtered_velo_prod; 100 CeedInt num_comp_filtered_velo_prod; 101 { // Get filtered velocity product information 102 CeedOperatorField op_field; 103 PetscCallCeed(ceed, CeedOperatorGetFieldByName(user->diff_filter->op_rhs_ctx->op, "v1", &op_field)); 104 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_filtered_velo_prod)); 105 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_velo_prod, &num_comp_filtered_velo_prod)); 106 } 107 108 // -- Create operator for generating training data at nodes 109 // Differential Filter only provides filtered primitive variables 110 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSGS_DDAnisotropicTrainingDataNodal_Prim, 111 ComputeSGS_DDAnisotropicTrainingDataNodal_Prim_loc, &qf_sgs_dd_train)); 112 113 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_train, sgs_dd_train_setup_data->sgs_dd_train_qfctx)); 114 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "q", num_comp_filtered_state, CEED_EVAL_NONE)); 115 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "velocity product", num_comp_filtered_velo_prod, CEED_EVAL_NONE)); 116 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE)); 117 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 118 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "inverse multiplicity", 1, CEED_EVAL_NONE)); 119 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_train, "training data", sgs_dd_train->num_comp_dd_inputs, CEED_EVAL_NONE)); 120 121 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_filtered_state, &filtered_fields, NULL)); 122 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_train, NULL, NULL, &op_sgs_dd_train)); 123 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "q", elem_restr_filtered_state, CEED_BASIS_NONE, filtered_fields)); 124 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "velocity product", elem_restr_filtered_velo_prod, CEED_BASIS_NONE, filtered_fields)); 125 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 126 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "anisotropy tensor", sgs_dd_train_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 127 sgs_dd_train_setup_data->grid_aniso_ceed)); 128 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 129 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "training data", elem_restr_sgs_train, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 130 131 PetscCall(OperatorApplyContextCreate(sgs_dd_train->filtered_grad_velo_proj->dm, sgs_dd_train->dm_dd_training, ceed, op_sgs_dd_train, NULL, NULL, 132 NULL, NULL, &sgs_dd_train->op_training_data_calc_ctx)); 133 134 PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity)); 135 PetscCallCeed(ceed, CeedVectorDestroy(&filtered_fields)); 136 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity)); 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, User user, CeedData ceed_data, ProblemData *problem) { 143 SGS_DDTrainingContext sgsdd_train_qfctx; 144 SGS_DD_TrainingSetupData sgs_dd_train_setup_data; 145 146 PetscFunctionBeginUser; 147 if (!user->diff_filter) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem)); 148 if (!user->smartsim) PetscCall(SmartSimSetup(user)); 149 150 PetscCall(PetscNew(&sgsdd_train_qfctx)); 151 PetscCall(PetscNew(&sgs_dd_train_setup_data)); 152 PetscCall(PetscNew(&user->sgs_dd_train)); 153 SGS_DD_TrainingData sgs_dd_train = user->sgs_dd_train; 154 155 sgs_dd_train->overwrite_training_data = PETSC_TRUE; 156 sgs_dd_train->write_data_interval = 1; 157 PetscOptionsBegin(user->comm, NULL, "SGS Data-Driven Training Options", NULL); 158 PetscCall(PetscOptionsInt("-sgs_train_write_data_interval", "Number of timesteps between writing data into database", NULL, 159 sgs_dd_train->write_data_interval, &sgs_dd_train->write_data_interval, NULL)); 160 PetscCall(PetscOptionsBool("-sgs_train_overwrite_data", "Overwrite old training data in the database", NULL, sgs_dd_train->overwrite_training_data, 161 &sgs_dd_train->overwrite_training_data, NULL)); 162 PetscOptionsEnd(); 163 164 // -- Create DM for storing training data 165 PetscCall(SGS_DD_TrainingCreateDM(user->dm, &sgs_dd_train->dm_dd_training, user->app_ctx->degree, user->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.qfunction_context, CEED_MEM_HOST, &gas)); 171 sgsdd_train_qfctx->gas = *gas; 172 PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas)); 173 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->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 = user->smartsim; 182 PetscCallMPI(MPI_Comm_rank(user->comm, &rank)); 183 PetscCallMPI(MPI_Comm_size(user->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(user->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 PetscSmartRedisCall( 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 PetscSmartRedisCall( 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 } 232 233 // -- Compute and store anisotropy tensor 234 PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, user, ceed_data, &sgs_dd_train_setup_data->elem_restr_grid_aniso, 235 &sgs_dd_train_setup_data->grid_aniso_ceed)); 236 237 // -- Create Nodal Evaluation Operator 238 PetscCall(SetupTrainingDataCalculation(ceed, user, ceed_data, problem, sgs_dd_train_setup_data)); 239 240 PetscCall(SGS_DD_TrainingSetupDataDestroy(sgs_dd_train_setup_data)); 241 PetscFunctionReturn(PETSC_SUCCESS); 242 } 243 244 PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx) { 245 User user = (User)ctx; 246 Ceed ceed = user->ceed; 247 SGS_DD_TrainingData sgs_dd_train = user->sgs_dd_train; 248 SmartSimData smartsim = user->smartsim; 249 Vec TrainingData; 250 251 PetscFunctionBeginUser; 252 if (step_num % sgs_dd_train->write_data_interval != 0) PetscFunctionReturn(PETSC_SUCCESS); 253 PetscCall(DMGetGlobalVector(sgs_dd_train->dm_dd_training, &TrainingData)); 254 255 PetscCall(PetscLogEventBegin(FLUIDS_TrainDataCompute, 0, 0, 0, 0)); 256 { // -- Compute and assemble training data 257 Vec FilteredVelocityGradient, FilteredFields, FilteredFields_loc; 258 PetscMemType filtered_fields_mem_type; 259 CeedVector filtered_fields; 260 261 PetscCall(DMGetGlobalVector(user->diff_filter->dm_filter, &FilteredFields)); 262 PetscCall(DMGetLocalVector(user->diff_filter->dm_filter, &FilteredFields_loc)); 263 264 PetscCall(DifferentialFilterApply(user, solution_time, Q, FilteredFields)); 265 PetscCall(DMGlobalToLocal(user->diff_filter->dm_filter, FilteredFields, INSERT_VALUES, FilteredFields_loc)); 266 267 PetscCall(DMGetGlobalVector(sgs_dd_train->filtered_grad_velo_proj->dm, &FilteredVelocityGradient)); 268 PetscCall(VelocityGradientProjectionApply(sgs_dd_train->filtered_grad_velo_proj, FilteredFields_loc, FilteredVelocityGradient)); 269 270 { 271 CeedOperatorField op_field; 272 PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->op_training_data_calc_ctx->op, "q", &op_field)); 273 PetscCallCeed(ceed, CeedOperatorFieldGetVector(op_field, &filtered_fields)); 274 } 275 PetscCall(VecPetscToCeed(FilteredFields_loc, &filtered_fields_mem_type, filtered_fields)); // filtered_fields is an implicit input 276 277 PetscCall(ApplyCeedOperatorGlobalToGlobal(FilteredVelocityGradient, TrainingData, sgs_dd_train->op_training_data_calc_ctx)); 278 279 PetscCall(VecCeedToPetsc(filtered_fields, filtered_fields_mem_type, FilteredFields_loc)); 280 281 PetscCall(DMRestoreGlobalVector(sgs_dd_train->filtered_grad_velo_proj->dm, &FilteredVelocityGradient)); 282 PetscCall(DMRestoreGlobalVector(user->diff_filter->dm_filter, &FilteredFields)); 283 PetscCall(DMRestoreLocalVector(user->diff_filter->dm_filter, &FilteredFields_loc)); 284 } 285 PetscCall(PetscLogEventEnd(FLUIDS_TrainDataCompute, 0, 0, 0, 0)); 286 287 { // -- Send training data to SmartSim 288 char array_key[PETSC_MAX_PATH_LEN]; 289 size_t array_key_len; 290 PetscMPIInt rank; 291 292 PetscCallMPI(MPI_Comm_rank(user->comm, &rank)); 293 294 if (sgs_dd_train->overwrite_training_data) { 295 PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s", smartsim->rank_id_name)); 296 } else { 297 PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT, smartsim->rank_id_name, step_num)); 298 } 299 PetscCall(PetscStrlen(array_key, &array_key_len)); 300 301 { 302 const PetscScalar *training_data; 303 PetscCall(VecGetArrayRead(TrainingData, &training_data)); 304 PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Train, 0, 0, 0, 0)); 305 PetscSmartRedisCall(put_tensor(smartsim->client, array_key, array_key_len, (void *)training_data, sgs_dd_train->training_data_array_dims, 2, 306 SRTensorTypeDouble, SRMemLayoutContiguous)); 307 PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Train, 0, 0, 0, 0)); 308 PetscCall(VecRestoreArrayRead(TrainingData, &training_data)); 309 } 310 311 if (rank % smartsim->collocated_database_num_ranks == 0) { 312 size_t dim_2[1] = {2}; 313 PetscInt64 step_array[2] = {step_num, step_num}; 314 PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 315 PetscSmartRedisCall(put_tensor(smartsim->client, "step", 4, step_array, dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous)); 316 PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 317 } 318 } 319 320 PetscCall(DMRestoreGlobalVector(user->sgs_dd_train->dm_dd_training, &TrainingData)); 321 PetscFunctionReturn(PETSC_SUCCESS); 322 } 323 324 PetscErrorCode TSPostStep_SGS_DD_Training(TS ts) { 325 User user; 326 const char check_run_key[] = "check-run"; 327 PetscReal check_run[2] = {1}; 328 const size_t check_run_dims[1] = {2}; 329 size_t check_run_key_size; 330 331 PetscFunctionBeginUser; 332 PetscCall(PetscStrlen(check_run_key, &check_run_key_size)); 333 PetscCall(TSGetApplicationContext(ts, &user)); 334 SmartSimData smartsim = user->smartsim; 335 336 PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 337 PetscSmartRedisCall( 338 unpack_tensor(smartsim->client, check_run_key, check_run_key_size, check_run, check_run_dims, 1, SRTensorTypeDouble, SRMemLayoutContiguous)); 339 PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 340 if (check_run[0] == 0) { 341 PetscCall(PetscPrintf(user->comm, "-- Simulation stopped by 'check-run' tensor in Redis database\n")); 342 PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER)); 343 } 344 345 PetscFunctionReturn(PETSC_SUCCESS); 346 } 347 348 PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train) { 349 PetscFunctionBeginUser; 350 if (!sgs_dd_train) PetscFunctionReturn(PETSC_SUCCESS); 351 352 PetscCall(OperatorApplyContextDestroy(sgs_dd_train->op_training_data_calc_ctx)); 353 PetscCall(NodalProjectionDataDestroy(sgs_dd_train->filtered_grad_velo_proj)); 354 PetscCall(DMDestroy(&sgs_dd_train->dm_dd_training)); 355 PetscCall(PetscFree(sgs_dd_train)); 356 357 PetscFunctionReturn(PETSC_SUCCESS); 358 } 359