1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3aa0b7f76SJames Wright 4aa0b7f76SJames Wright #include "../../qfunctions/sgs_dd_training.h" 5aa0b7f76SJames Wright 68d78d7c8SJames Wright #include <differential_filter.h> 7149fb536SJames Wright #include <navierstokes.h> 88d78d7c8SJames Wright #include <petscdmplex.h> 9149fb536SJames Wright #include <smartsim.h> 10aa0b7f76SJames Wright 11aa0b7f76SJames Wright typedef struct { 1239169b57SJames Wright DM dm_dd_training; 1339169b57SJames Wright PetscInt num_comp_dd_inputs, write_data_interval, num_filter_widths; 1439169b57SJames Wright PetscScalar filter_widths[16]; 1539169b57SJames Wright OperatorApplyContext op_training_data_calc_ctx; 16cb8a476cSJames Wright DiffFilterData diff_filter; 1739169b57SJames Wright NodalProjectionData filtered_grad_velo_proj; 1839169b57SJames Wright size_t training_data_array_dims[2]; 1939169b57SJames Wright PetscBool overwrite_training_data; 2039169b57SJames Wright } *SGS_DD_TrainingData; 2139169b57SJames Wright 2239169b57SJames Wright #define SGS_DD_TRAIN_KEY "SGS Data Driven Training" 2339169b57SJames Wright 2439169b57SJames Wright PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData *sgs_dd_train) { 2539169b57SJames Wright SGS_DD_TrainingData sgs_dd_train_ = *sgs_dd_train; 2639169b57SJames Wright 2739169b57SJames Wright PetscFunctionBeginUser; 2839169b57SJames Wright if (!sgs_dd_train_) PetscFunctionReturn(PETSC_SUCCESS); 2939169b57SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_train_->op_training_data_calc_ctx)); 3039169b57SJames Wright PetscCall(NodalProjectionDataDestroy(&sgs_dd_train_->filtered_grad_velo_proj)); 3139169b57SJames Wright PetscCall(DMDestroy(&sgs_dd_train_->dm_dd_training)); 32cb8a476cSJames Wright PetscCall(DifferentialFilterDataDestroy(&sgs_dd_train_->diff_filter)); 335206a5a0SJames Wright PetscCall(PetscFree(sgs_dd_train_)); 345206a5a0SJames Wright *sgs_dd_train = NULL; 3539169b57SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 3639169b57SJames Wright } 3739169b57SJames Wright 3839169b57SJames Wright typedef struct { 39aa0b7f76SJames Wright CeedElemRestriction elem_restr_grid_aniso; 40aa0b7f76SJames Wright CeedVector grid_aniso_ceed; 41aa0b7f76SJames Wright CeedQFunctionContext sgs_dd_train_qfctx; 42aa0b7f76SJames Wright } *SGS_DD_TrainingSetupData; 43aa0b7f76SJames Wright 44aa0b7f76SJames Wright static PetscErrorCode SGS_DD_TrainingSetupDataDestroy(SGS_DD_TrainingSetupData sgs_dd_train_setup_data) { 45aa0b7f76SJames Wright Ceed ceed; 46aa0b7f76SJames Wright 47aa0b7f76SJames Wright PetscFunctionBeginUser; 48aa0b7f76SJames Wright PetscCall(CeedElemRestrictionGetCeed(sgs_dd_train_setup_data->elem_restr_grid_aniso, &ceed)); 49aa0b7f76SJames Wright 50aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_train_setup_data->elem_restr_grid_aniso)); 51aa0b7f76SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_train_setup_data->grid_aniso_ceed)); 52aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_train_setup_data->sgs_dd_train_qfctx)); 53aa0b7f76SJames Wright PetscCall(PetscFree(sgs_dd_train_setup_data)); 54519781aeSJames Wright PetscCheck(CeedDestroy(&ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "Destroying Ceed object failed"); 55aa0b7f76SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 56aa0b7f76SJames Wright } 57aa0b7f76SJames Wright 58aa0b7f76SJames Wright // @brief Create DM for storing data-drive SGS model inputs 59aa0b7f76SJames Wright static PetscErrorCode SGS_DD_TrainingCreateDM(DM dm_source, DM *dm_dd_training, PetscInt degree, PetscInt q_extra, PetscInt *num_components) { 60aa0b7f76SJames Wright PetscSection section; 61aa0b7f76SJames Wright 62aa0b7f76SJames Wright PetscFunctionBeginUser; 63aa0b7f76SJames Wright *num_components = 12; 64aa0b7f76SJames Wright 65aa0b7f76SJames Wright PetscCall(DMClone(dm_source, dm_dd_training)); 660dee9b8eSJames Wright PetscCall(DMSetMatrixPreallocateSkip(*dm_dd_training, PETSC_TRUE)); 67aa0b7f76SJames Wright PetscCall(PetscObjectSetName((PetscObject)*dm_dd_training, "Data-Driven SGS Training Data")); 68aa0b7f76SJames Wright 69aa0b7f76SJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_dd_training)); 70aa0b7f76SJames Wright 71aa0b7f76SJames Wright PetscCall(DMGetLocalSection(*dm_dd_training, §ion)); 72aa0b7f76SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "Data-Driven SGS Training Data")); 73aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "SGSInput1")); 74aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "SGSInput2")); 75aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "SGSInput3")); 76aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "SGSInput4")); 77aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "SGSInput5")); 78aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 5, "SGSInput6")); 79aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 6, "FilteredSGSXX")); 80aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 7, "FilteredSGSYY")); 81aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 8, "FilteredSGSZZ")); 82aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 9, "FilteredSGSYZ")); 83aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 10, "FilteredSGSXZ")); 84aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 11, "FilteredSGSXY")); 85aa0b7f76SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 86aa0b7f76SJames Wright }; 87aa0b7f76SJames Wright 88aa0b7f76SJames Wright // @brief Create CeedOperator to calculate training data for data-drive SGS model at nodes 89e3663b90SJames Wright static PetscErrorCode SetupTrainingDataCalculation(Ceed ceed, Honee honee, ProblemData problem, SGS_DD_TrainingSetupData sgs_dd_train_setup_data) { 9039169b57SJames Wright SGS_DD_TrainingData sgs_dd_train; 91aa0b7f76SJames Wright CeedQFunction qf_sgs_dd_train; 92aa0b7f76SJames Wright CeedOperator op_sgs_dd_train; 93aa0b7f76SJames Wright CeedInt num_comp_grad_velo, num_comp_grid_aniso; 94aa0b7f76SJames Wright CeedVector inv_multiplicity, filtered_fields; 95aa0b7f76SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs_train; 96aa0b7f76SJames Wright DMLabel domain_label = NULL; 97aa0b7f76SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 98aa0b7f76SJames Wright 99aa0b7f76SJames Wright PetscFunctionBeginUser; 1000c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, SGS_DD_TRAIN_KEY, &sgs_dd_train)); 101aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_train_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso)); 102aa0b7f76SJames Wright 103aa0b7f76SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_train->dm_dd_training, domain_label, label_value, height, dm_field, &elem_restr_sgs_train)); 1045930f037SJames Wright PetscCall(GetInverseMultiplicity(ceed, sgs_dd_train->dm_dd_training, domain_label, label_value, height, dm_field, PETSC_TRUE, 1055930f037SJames Wright &elem_restr_inv_multiplicity, &inv_multiplicity)); 106aa0b7f76SJames Wright 107aa0b7f76SJames Wright CeedElemRestriction elem_restr_filtered_state; 108aa0b7f76SJames Wright CeedInt num_comp_filtered_state; 109aa0b7f76SJames Wright { // -- Setup filtered velocity gradient projection 110aa0b7f76SJames Wright CeedBasis basis_filtered_state; 111aa0b7f76SJames Wright CeedOperatorField op_field; 112cb8a476cSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->diff_filter->op_rhs_ctx->op, "v0", &op_field)); 11301e19bfaSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_filtered_state, &basis_filtered_state, NULL)); 114aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_state, &num_comp_filtered_state)); 115e3663b90SJames Wright PetscCall(VelocityGradientProjectionSetup(ceed, honee, problem, STATEVAR_PRIMITIVE, elem_restr_filtered_state, basis_filtered_state, 116aa0b7f76SJames Wright &sgs_dd_train->filtered_grad_velo_proj)); 117fff85bd3SJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_filtered_state)); 118aa0b7f76SJames Wright // Get velocity gradient information 119aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->filtered_grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field)); 120aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo)); 121aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo)); 122aa0b7f76SJames Wright } 123aa0b7f76SJames Wright 124aa0b7f76SJames Wright CeedElemRestriction elem_restr_filtered_velo_prod; 125aa0b7f76SJames Wright CeedInt num_comp_filtered_velo_prod; 126aa0b7f76SJames Wright { // Get filtered velocity product information 127aa0b7f76SJames Wright CeedOperatorField op_field; 128cb8a476cSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->diff_filter->op_rhs_ctx->op, "v1", &op_field)); 129aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_filtered_velo_prod)); 130aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_velo_prod, &num_comp_filtered_velo_prod)); 131aa0b7f76SJames Wright } 132aa0b7f76SJames Wright 133aa0b7f76SJames Wright // -- Create operator for generating training data at nodes 134aa0b7f76SJames Wright // Differential Filter only provides filtered primitive variables 135aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSGS_DDAnisotropicTrainingDataNodal_Prim, 136aa0b7f76SJames Wright ComputeSGS_DDAnisotropicTrainingDataNodal_Prim_loc, &qf_sgs_dd_train)); 137aa0b7f76SJames Wright 138aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_train, sgs_dd_train_setup_data->sgs_dd_train_qfctx)); 139aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "q", num_comp_filtered_state, CEED_EVAL_NONE)); 140aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "velocity product", num_comp_filtered_velo_prod, CEED_EVAL_NONE)); 141aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE)); 142aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 143aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "inverse multiplicity", 1, CEED_EVAL_NONE)); 144aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_train, "training data", sgs_dd_train->num_comp_dd_inputs, CEED_EVAL_NONE)); 145aa0b7f76SJames Wright 146aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_filtered_state, &filtered_fields, NULL)); 147aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_train, NULL, NULL, &op_sgs_dd_train)); 148c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "q", elem_restr_filtered_state, CEED_BASIS_NONE, filtered_fields)); 149c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "velocity product", elem_restr_filtered_velo_prod, CEED_BASIS_NONE, filtered_fields)); 150c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 151c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "anisotropy tensor", sgs_dd_train_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 152c6271fa9SJeremy L Thompson sgs_dd_train_setup_data->grid_aniso_ceed)); 153c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 154c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "training data", elem_restr_sgs_train, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 155aa0b7f76SJames Wright 156aa0b7f76SJames Wright PetscCall(OperatorApplyContextCreate(sgs_dd_train->filtered_grad_velo_proj->dm, sgs_dd_train->dm_dd_training, ceed, op_sgs_dd_train, NULL, NULL, 157aa0b7f76SJames Wright NULL, NULL, &sgs_dd_train->op_training_data_calc_ctx)); 158aa0b7f76SJames Wright 159aa0b7f76SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity)); 160aa0b7f76SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&filtered_fields)); 161aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity)); 162fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_filtered_state)); 163fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo)); 164fff85bd3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_filtered_velo_prod)); 165aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_train)); 166aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_train)); 167aa0b7f76SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 168aa0b7f76SJames Wright } 169aa0b7f76SJames Wright 170e3663b90SJames Wright PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, Honee honee, ProblemData problem) { 1717f3a2123SJames Wright SGS_DDTrainingContext sgsdd_train_ctx; 172aa0b7f76SJames Wright SGS_DD_TrainingSetupData sgs_dd_train_setup_data; 17339169b57SJames Wright SGS_DD_TrainingData sgs_dd_train; 174aa0b7f76SJames Wright 175aa0b7f76SJames Wright PetscFunctionBeginUser; 176aa0b7f76SJames Wright 1777f3a2123SJames Wright PetscCall(PetscNew(&sgsdd_train_ctx)); 178aa0b7f76SJames Wright PetscCall(PetscNew(&sgs_dd_train_setup_data)); 17939169b57SJames Wright PetscCall(PetscNew(&sgs_dd_train)); 1800c70a8bcSJames Wright PetscCall(HoneeSetContainer(honee, SGS_DD_TRAIN_KEY, sgs_dd_train, (PetscCtxDestroyFn *)SGS_DD_TrainingDataDestroy)); 181cb8a476cSJames Wright PetscCall(DifferentialFilterSetup(honee, &sgs_dd_train->diff_filter)); 182aa0b7f76SJames Wright 183aa0b7f76SJames Wright sgs_dd_train->overwrite_training_data = PETSC_TRUE; 184aa0b7f76SJames Wright sgs_dd_train->write_data_interval = 1; 1856ea7c1aeSJames Wright sgs_dd_train->num_filter_widths = sizeof(sgs_dd_train->filter_widths) / sizeof(sgs_dd_train->filter_widths[0]); 1860c373b74SJames Wright PetscOptionsBegin(honee->comm, NULL, "SGS Data-Driven Training Options", NULL); 187aa0b7f76SJames Wright PetscCall(PetscOptionsInt("-sgs_train_write_data_interval", "Number of timesteps between writing data into database", NULL, 188aa0b7f76SJames Wright sgs_dd_train->write_data_interval, &sgs_dd_train->write_data_interval, NULL)); 189aa0b7f76SJames Wright PetscCall(PetscOptionsBool("-sgs_train_overwrite_data", "Overwrite old training data in the database", NULL, sgs_dd_train->overwrite_training_data, 190aa0b7f76SJames Wright &sgs_dd_train->overwrite_training_data, NULL)); 1916ea7c1aeSJames Wright PetscCall(PetscOptionsRealArray("-sgs_train_filter_width_scales", "Scales of each filter width put into training database", NULL, 1926ea7c1aeSJames Wright sgs_dd_train->filter_widths, &sgs_dd_train->num_filter_widths, NULL)); 193aa0b7f76SJames Wright PetscOptionsEnd(); 194aa0b7f76SJames Wright 195aa0b7f76SJames Wright // -- Create DM for storing training data 1960c373b74SJames Wright PetscCall(SGS_DD_TrainingCreateDM(honee->dm, &sgs_dd_train->dm_dd_training, honee->app_ctx->degree, honee->app_ctx->q_extra, 197aa0b7f76SJames Wright &sgs_dd_train->num_comp_dd_inputs)); 198aa0b7f76SJames Wright 199aa0b7f76SJames Wright { // -- Create QFunction Context 200aa0b7f76SJames Wright NewtonianIdealGasContext gas; 201e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &gas)); 2027f3a2123SJames Wright sgsdd_train_ctx->gas = *gas; 203e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &gas)); 2040c373b74SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &sgs_dd_train_setup_data->sgs_dd_train_qfctx)); 205aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(sgs_dd_train_setup_data->sgs_dd_train_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, 2067f3a2123SJames Wright sizeof(*sgsdd_train_ctx), sgsdd_train_ctx)); 207aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_train_setup_data->sgs_dd_train_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 208aa0b7f76SJames Wright } 209aa0b7f76SJames Wright 210aa0b7f76SJames Wright { // -- Send training data array info to SmartRedis database 211aa0b7f76SJames Wright PetscMPIInt rank, num_ranks; 2127ebeccb9SJames Wright SmartSimData smartsim; 2137ebeccb9SJames Wright PetscCall(HoneeGetSmartSimData(honee, &smartsim)); 2140c373b74SJames Wright PetscCallMPI(MPI_Comm_rank(honee->comm, &rank)); 2150c373b74SJames Wright PetscCallMPI(MPI_Comm_size(honee->comm, &num_ranks)); 216aa0b7f76SJames Wright 217aa0b7f76SJames Wright { 218aa0b7f76SJames Wright PetscSection global_section; 2197ff16c02SJames Wright PetscInt num_dofs, num_comps, local_min_max[2] = {0.}, global_min_max[2] = {0.}; 2207ff16c02SJames Wright 221aa0b7f76SJames Wright PetscCall(DMGetGlobalSection(sgs_dd_train->dm_dd_training, &global_section)); 222aa0b7f76SJames Wright PetscCall(DMGetGlobalVectorInfo(sgs_dd_train->dm_dd_training, &num_dofs, NULL, NULL)); 223aa0b7f76SJames Wright PetscCall(PetscSectionGetFieldComponents(global_section, 0, &num_comps)); 2247ff16c02SJames Wright local_min_max[0] = num_dofs; 2250c373b74SJames Wright PetscCall(PetscGlobalMinMaxInt(honee->comm, local_min_max, global_min_max)); 2267ff16c02SJames Wright 2277ff16c02SJames Wright sgs_dd_train->training_data_array_dims[0] = global_min_max[0] / num_comps; 228aa0b7f76SJames Wright sgs_dd_train->training_data_array_dims[1] = num_comps; 229aa0b7f76SJames Wright } 230aa0b7f76SJames Wright 231aa0b7f76SJames Wright if (rank % smartsim->collocated_database_num_ranks == 0) { 2324fa1625aSJames Wright { // Communicate info on simulation size 2334fa1625aSJames Wright const char tensor_name[] = "sizeInfo"; 234aa0b7f76SJames Wright size_t array_info_dim = 6; 235f4632befSRiccardo Balin PetscInt64 array_info[6] = {0}, num_features = 6; 236aa0b7f76SJames Wright 237aa0b7f76SJames Wright array_info[0] = sgs_dd_train->training_data_array_dims[0]; 238aa0b7f76SJames Wright array_info[1] = sgs_dd_train->training_data_array_dims[1]; 239aa0b7f76SJames Wright array_info[2] = num_features; 240aa0b7f76SJames Wright array_info[3] = num_ranks; 241aa0b7f76SJames Wright array_info[4] = smartsim->collocated_database_num_ranks; 242aa0b7f76SJames Wright array_info[5] = rank; 243aa0b7f76SJames Wright 244ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0)); 24543e9749fSJames Wright PetscCallSmartRedis( 2464fa1625aSJames Wright put_tensor(smartsim->client, tensor_name, strlen(tensor_name), array_info, &array_info_dim, 1, SRTensorTypeInt64, SRMemLayoutContiguous)); 2474fa1625aSJames Wright PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name))); 248ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0)); 2494fa1625aSJames Wright } 250aa0b7f76SJames Wright 2514fa1625aSJames Wright { // Send array that communicates if tensors are overwritten in database 2524fa1625aSJames Wright const char tensor_name[] = "tensor-ow"; 253f4632befSRiccardo Balin PetscInt64 tensor_overwrite[2] = {sgs_dd_train->overwrite_training_data}; 254aa0b7f76SJames Wright size_t dim_2[1] = {2}; 2554fa1625aSJames Wright 256ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0)); 25743e9749fSJames Wright PetscCallSmartRedis( 2584fa1625aSJames Wright put_tensor(smartsim->client, tensor_name, strlen(tensor_name), tensor_overwrite, dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous)); 2594fa1625aSJames Wright PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name))); 260ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0)); 261aa0b7f76SJames Wright } 2626ea7c1aeSJames Wright 2636ea7c1aeSJames Wright { // Communicate number of filter widths used 2646ea7c1aeSJames Wright const char tensor_name[] = "num_filter_widths"; 2656ea7c1aeSJames Wright PetscInt64 num_filter_widths = sgs_dd_train->num_filter_widths; 2666ea7c1aeSJames Wright size_t dim_2 = 1; 2676ea7c1aeSJames Wright 268ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0)); 26943e9749fSJames Wright PetscCallSmartRedis( 2706ea7c1aeSJames Wright put_tensor(smartsim->client, tensor_name, strlen(tensor_name), &num_filter_widths, &dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous)); 2716ea7c1aeSJames Wright PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name))); 272ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0)); 2736ea7c1aeSJames Wright } 274aa0b7f76SJames Wright } 2754fa1625aSJames Wright } 276aa0b7f76SJames Wright 277aa0b7f76SJames Wright // -- Compute and store anisotropy tensor 278e3663b90SJames Wright PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, honee, &sgs_dd_train_setup_data->elem_restr_grid_aniso, 279aa0b7f76SJames Wright &sgs_dd_train_setup_data->grid_aniso_ceed)); 280aa0b7f76SJames Wright 281aa0b7f76SJames Wright // -- Create Nodal Evaluation Operator 282e3663b90SJames Wright PetscCall(SetupTrainingDataCalculation(ceed, honee, problem, sgs_dd_train_setup_data)); 283aa0b7f76SJames Wright 284aa0b7f76SJames Wright PetscCall(SGS_DD_TrainingSetupDataDestroy(sgs_dd_train_setup_data)); 285aa0b7f76SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 286aa0b7f76SJames Wright } 287aa0b7f76SJames Wright 288aa0b7f76SJames Wright PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx) { 2890c373b74SJames Wright Honee honee = (Honee)ctx; 2900c373b74SJames Wright Ceed ceed = honee->ceed; 29139169b57SJames Wright SGS_DD_TrainingData sgs_dd_train; 2927ebeccb9SJames Wright SmartSimData smartsim; 293aa0b7f76SJames Wright Vec TrainingData; 294ee81c423SRiccardo Balin PetscMPIInt rank; 295aa0b7f76SJames Wright 296aa0b7f76SJames Wright PetscFunctionBeginUser; 2977ebeccb9SJames Wright PetscCall(HoneeGetSmartSimData(honee, &smartsim)); 2980c70a8bcSJames Wright PetscCall(HoneeGetContainer(honee, SGS_DD_TRAIN_KEY, &sgs_dd_train)); 2990c373b74SJames Wright PetscCallMPI(MPI_Comm_rank(honee->comm, &rank)); 300ee81c423SRiccardo Balin 301aa0b7f76SJames Wright if (step_num % sgs_dd_train->write_data_interval != 0) PetscFunctionReturn(PETSC_SUCCESS); 302aa0b7f76SJames Wright PetscCall(DMGetGlobalVector(sgs_dd_train->dm_dd_training, &TrainingData)); 303aa0b7f76SJames Wright 3046ea7c1aeSJames Wright for (PetscInt filter_index = 0; filter_index < sgs_dd_train->num_filter_widths; filter_index++) { 305ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_TrainDataCompute, 0, 0, 0, 0)); 306aa0b7f76SJames Wright { // -- Compute and assemble training data 307aa0b7f76SJames Wright Vec FilteredVelocityGradient, FilteredFields, FilteredFields_loc; 308aa0b7f76SJames Wright PetscMemType filtered_fields_mem_type; 309aa0b7f76SJames Wright CeedVector filtered_fields; 310aa0b7f76SJames Wright 3116ea7c1aeSJames Wright { // Set filter width for the current solve 3126ea7c1aeSJames Wright double filter_width_scaling[3]; 3136ea7c1aeSJames Wright CeedOperator op_mat; 3146ea7c1aeSJames Wright Mat A_mat; 3156ea7c1aeSJames Wright 3166ea7c1aeSJames Wright for (int j = 0; j < 3; j++) filter_width_scaling[j] = sgs_dd_train->filter_widths[filter_index]; 317cb8a476cSJames Wright PetscCall(KSPGetOperators(sgs_dd_train->diff_filter->ksp, &A_mat, NULL)); 3186ea7c1aeSJames Wright PetscCall(MatCeedGetCeedOperators(A_mat, &op_mat, NULL)); 319cb8a476cSJames Wright PetscCall(CeedOperatorSetContextDouble(op_mat, sgs_dd_train->diff_filter->filter_width_scaling_label, filter_width_scaling)); 3206ea7c1aeSJames Wright } 3216ea7c1aeSJames Wright 322cb8a476cSJames Wright PetscCall(DMGetGlobalVector(sgs_dd_train->diff_filter->dm_filter, &FilteredFields)); 323cb8a476cSJames Wright PetscCall(DMGetLocalVector(sgs_dd_train->diff_filter->dm_filter, &FilteredFields_loc)); 324aa0b7f76SJames Wright 325cb8a476cSJames Wright PetscCall(DifferentialFilterApply(honee, sgs_dd_train->diff_filter, solution_time, Q, FilteredFields)); 326cb8a476cSJames Wright PetscCall(DMGlobalToLocal(sgs_dd_train->diff_filter->dm_filter, FilteredFields, INSERT_VALUES, FilteredFields_loc)); 327aa0b7f76SJames Wright 328aa0b7f76SJames Wright PetscCall(DMGetGlobalVector(sgs_dd_train->filtered_grad_velo_proj->dm, &FilteredVelocityGradient)); 329aa0b7f76SJames Wright PetscCall(VelocityGradientProjectionApply(sgs_dd_train->filtered_grad_velo_proj, FilteredFields_loc, FilteredVelocityGradient)); 330aa0b7f76SJames Wright 331aa0b7f76SJames Wright { 332aa0b7f76SJames Wright CeedOperatorField op_field; 3336ea7c1aeSJames Wright 334aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->op_training_data_calc_ctx->op, "q", &op_field)); 335aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetVector(op_field, &filtered_fields)); 336aa0b7f76SJames Wright } 3376ea7c1aeSJames Wright 338a7dac1d5SJames Wright PetscCall(VecPetscToCeed(FilteredFields_loc, &filtered_fields_mem_type, filtered_fields)); // filtered_fields is an implicit input 339aa0b7f76SJames Wright PetscCall(ApplyCeedOperatorGlobalToGlobal(FilteredVelocityGradient, TrainingData, sgs_dd_train->op_training_data_calc_ctx)); 340a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(filtered_fields, filtered_fields_mem_type, FilteredFields_loc)); 341aa0b7f76SJames Wright 342aa0b7f76SJames Wright PetscCall(DMRestoreGlobalVector(sgs_dd_train->filtered_grad_velo_proj->dm, &FilteredVelocityGradient)); 343cb8a476cSJames Wright PetscCall(DMRestoreGlobalVector(sgs_dd_train->diff_filter->dm_filter, &FilteredFields)); 344cb8a476cSJames Wright PetscCall(DMRestoreLocalVector(sgs_dd_train->diff_filter->dm_filter, &FilteredFields_loc)); 345fff85bd3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&filtered_fields)); 346aa0b7f76SJames Wright } 347ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_TrainDataCompute, 0, 0, 0, 0)); 348aa0b7f76SJames Wright 349aa0b7f76SJames Wright { // -- Send training data to SmartSim 350aa0b7f76SJames Wright char array_key[PETSC_MAX_PATH_LEN]; 351aa0b7f76SJames Wright size_t array_key_len; 352aa0b7f76SJames Wright 353aa0b7f76SJames Wright if (sgs_dd_train->overwrite_training_data) { 3546ea7c1aeSJames Wright PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT, smartsim->rank_id_name, filter_index)); 355aa0b7f76SJames Wright } else { 3566ea7c1aeSJames Wright PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT "%" PetscInt_FMT, smartsim->rank_id_name, step_num, filter_index)); 357aa0b7f76SJames Wright } 358aa0b7f76SJames Wright PetscCall(PetscStrlen(array_key, &array_key_len)); 359aa0b7f76SJames Wright 360aa0b7f76SJames Wright { 361aa0b7f76SJames Wright const PetscScalar *training_data; 362aa0b7f76SJames Wright PetscCall(VecGetArrayRead(TrainingData, &training_data)); 363*e171caa6SJames Wright PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Write, 0, 0, 0, 0)); 36443e9749fSJames Wright PetscCallSmartRedis(put_tensor(smartsim->client, array_key, array_key_len, (void *)training_data, sgs_dd_train->training_data_array_dims, 2, 365aa0b7f76SJames Wright SRTensorTypeDouble, SRMemLayoutContiguous)); 366*e171caa6SJames Wright PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Write, 0, 0, 0, 0)); 367aa0b7f76SJames Wright PetscCall(VecRestoreArrayRead(TrainingData, &training_data)); 368aa0b7f76SJames Wright } 369ee81c423SRiccardo Balin } 370ee81c423SRiccardo Balin } 371aa0b7f76SJames Wright 372aa0b7f76SJames Wright if (rank % smartsim->collocated_database_num_ranks == 0) { 3736ea7c1aeSJames Wright const char tensor_name[] = "step"; 374aa0b7f76SJames Wright size_t dim_2[1] = {2}; 375f4632befSRiccardo Balin PetscInt64 step_array[2] = {step_num, step_num}; 3766ea7c1aeSJames Wright 377ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0)); 37843e9749fSJames Wright PetscCallSmartRedis( 3796ea7c1aeSJames Wright put_tensor(smartsim->client, tensor_name, strlen(tensor_name), step_array, dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous)); 380ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0)); 381aa0b7f76SJames Wright } 382aa0b7f76SJames Wright 38339169b57SJames Wright PetscCall(DMRestoreGlobalVector(sgs_dd_train->dm_dd_training, &TrainingData)); 384aa0b7f76SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 385aa0b7f76SJames Wright } 386aa0b7f76SJames Wright 387632a41e1SJames Wright PetscErrorCode TSPostStep_SGS_DD_Training(TS ts) { 3880c373b74SJames Wright Honee honee; 389632a41e1SJames Wright const char check_run_key[] = "check-run"; 390632a41e1SJames Wright PetscReal check_run[2] = {1}; 391632a41e1SJames Wright const size_t check_run_dims[1] = {2}; 392632a41e1SJames Wright size_t check_run_key_size; 3937ebeccb9SJames Wright SmartSimData smartsim; 394632a41e1SJames Wright 395632a41e1SJames Wright PetscFunctionBeginUser; 396632a41e1SJames Wright PetscCall(PetscStrlen(check_run_key, &check_run_key_size)); 3970c373b74SJames Wright PetscCall(TSGetApplicationContext(ts, &honee)); 3987ebeccb9SJames Wright PetscCall(HoneeGetSmartSimData(honee, &smartsim)); 399632a41e1SJames Wright 400ea615d4cSJames Wright PetscCall(PetscLogEventBegin(HONEE_SmartRedis_Meta, 0, 0, 0, 0)); 40143e9749fSJames Wright PetscCallSmartRedis( 402632a41e1SJames Wright unpack_tensor(smartsim->client, check_run_key, check_run_key_size, check_run, check_run_dims, 1, SRTensorTypeDouble, SRMemLayoutContiguous)); 403ea615d4cSJames Wright PetscCall(PetscLogEventEnd(HONEE_SmartRedis_Meta, 0, 0, 0, 0)); 404632a41e1SJames Wright if (check_run[0] == 0) { 4050c373b74SJames Wright PetscCall(PetscPrintf(honee->comm, "-- Simulation stopped by 'check-run' tensor in Redis database\n")); 406632a41e1SJames Wright PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER)); 407632a41e1SJames Wright } 408632a41e1SJames Wright 409632a41e1SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 410632a41e1SJames Wright } 411