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 6aa0b7f76SJames Wright #include <petscdmplex.h> 7aa0b7f76SJames Wright 8149fb536SJames Wright #include <navierstokes.h> 9149fb536SJames Wright #include <smartsim.h> 10aa0b7f76SJames Wright 11aa0b7f76SJames Wright typedef struct { 12aa0b7f76SJames Wright CeedElemRestriction elem_restr_grid_aniso; 13aa0b7f76SJames Wright CeedVector grid_aniso_ceed; 14aa0b7f76SJames Wright CeedQFunctionContext sgs_dd_train_qfctx; 15aa0b7f76SJames Wright } *SGS_DD_TrainingSetupData; 16aa0b7f76SJames Wright 17aa0b7f76SJames Wright static PetscErrorCode SGS_DD_TrainingSetupDataDestroy(SGS_DD_TrainingSetupData sgs_dd_train_setup_data) { 18aa0b7f76SJames Wright Ceed ceed; 19aa0b7f76SJames Wright 20aa0b7f76SJames Wright PetscFunctionBeginUser; 21aa0b7f76SJames Wright PetscCall(CeedElemRestrictionGetCeed(sgs_dd_train_setup_data->elem_restr_grid_aniso, &ceed)); 22aa0b7f76SJames Wright 23aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_train_setup_data->elem_restr_grid_aniso)); 24aa0b7f76SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_train_setup_data->grid_aniso_ceed)); 25aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_train_setup_data->sgs_dd_train_qfctx)); 26aa0b7f76SJames Wright PetscCall(PetscFree(sgs_dd_train_setup_data)); 27aa0b7f76SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 28aa0b7f76SJames Wright } 29aa0b7f76SJames Wright 30aa0b7f76SJames Wright // @brief Create DM for storing data-drive SGS model inputs 31aa0b7f76SJames Wright static PetscErrorCode SGS_DD_TrainingCreateDM(DM dm_source, DM *dm_dd_training, PetscInt degree, PetscInt q_extra, PetscInt *num_components) { 32aa0b7f76SJames Wright PetscSection section; 33aa0b7f76SJames Wright 34aa0b7f76SJames Wright PetscFunctionBeginUser; 35aa0b7f76SJames Wright *num_components = 12; 36aa0b7f76SJames Wright 37aa0b7f76SJames Wright PetscCall(DMClone(dm_source, dm_dd_training)); 38*0dee9b8eSJames Wright PetscCall(DMSetMatrixPreallocateSkip(*dm_dd_training, PETSC_TRUE)); 39aa0b7f76SJames Wright PetscCall(PetscObjectSetName((PetscObject)*dm_dd_training, "Data-Driven SGS Training Data")); 40aa0b7f76SJames Wright 41aa0b7f76SJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_dd_training)); 42aa0b7f76SJames Wright 43aa0b7f76SJames Wright PetscCall(DMGetLocalSection(*dm_dd_training, §ion)); 44aa0b7f76SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "Data-Driven SGS Training Data")); 45aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "SGSInput1")); 46aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "SGSInput2")); 47aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "SGSInput3")); 48aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "SGSInput4")); 49aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "SGSInput5")); 50aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 5, "SGSInput6")); 51aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 6, "FilteredSGSXX")); 52aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 7, "FilteredSGSYY")); 53aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 8, "FilteredSGSZZ")); 54aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 9, "FilteredSGSYZ")); 55aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 10, "FilteredSGSXZ")); 56aa0b7f76SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 11, "FilteredSGSXY")); 57aa0b7f76SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 58aa0b7f76SJames Wright }; 59aa0b7f76SJames Wright 60aa0b7f76SJames Wright // @brief Create CeedOperator to calculate training data for data-drive SGS model at nodes 61991aef52SJames Wright static PetscErrorCode SetupTrainingDataCalculation(Ceed ceed, User user, CeedData ceed_data, ProblemData problem, 62aa0b7f76SJames Wright SGS_DD_TrainingSetupData sgs_dd_train_setup_data) { 63aa0b7f76SJames Wright SGS_DD_TrainingData sgs_dd_train = user->sgs_dd_train; 64aa0b7f76SJames Wright CeedQFunction qf_sgs_dd_train; 65aa0b7f76SJames Wright CeedOperator op_sgs_dd_train; 66aa0b7f76SJames Wright CeedInt num_comp_grad_velo, num_comp_grid_aniso; 67aa0b7f76SJames Wright CeedVector inv_multiplicity, filtered_fields; 68aa0b7f76SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs_train; 69aa0b7f76SJames Wright DMLabel domain_label = NULL; 70aa0b7f76SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 71aa0b7f76SJames Wright 72aa0b7f76SJames Wright PetscFunctionBeginUser; 73aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_train_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso)); 74aa0b7f76SJames Wright 75aa0b7f76SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_train->dm_dd_training, domain_label, label_value, height, dm_field, &elem_restr_sgs_train)); 765930f037SJames Wright PetscCall(GetInverseMultiplicity(ceed, sgs_dd_train->dm_dd_training, domain_label, label_value, height, dm_field, PETSC_TRUE, 775930f037SJames Wright &elem_restr_inv_multiplicity, &inv_multiplicity)); 78aa0b7f76SJames Wright 79aa0b7f76SJames Wright CeedElemRestriction elem_restr_filtered_state; 80aa0b7f76SJames Wright CeedInt num_comp_filtered_state; 81aa0b7f76SJames Wright { // -- Setup filtered velocity gradient projection 82aa0b7f76SJames Wright CeedBasis basis_filtered_state; 83aa0b7f76SJames Wright CeedOperatorField op_field; 84aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(user->diff_filter->op_rhs_ctx->op, "v0", &op_field)); 85aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_filtered_state)); 86aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_state, &num_comp_filtered_state)); 87aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetBasis(op_field, &basis_filtered_state)); 88aa0b7f76SJames Wright PetscCall(VelocityGradientProjectionSetup(ceed, user, ceed_data, problem, STATEVAR_PRIMITIVE, elem_restr_filtered_state, basis_filtered_state, 89aa0b7f76SJames Wright &sgs_dd_train->filtered_grad_velo_proj)); 90aa0b7f76SJames Wright // Get velocity gradient information 91aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->filtered_grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field)); 92aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo)); 93aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo)); 94aa0b7f76SJames Wright } 95aa0b7f76SJames Wright 96aa0b7f76SJames Wright CeedElemRestriction elem_restr_filtered_velo_prod; 97aa0b7f76SJames Wright CeedInt num_comp_filtered_velo_prod; 98aa0b7f76SJames Wright { // Get filtered velocity product information 99aa0b7f76SJames Wright CeedOperatorField op_field; 100aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(user->diff_filter->op_rhs_ctx->op, "v1", &op_field)); 101aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_filtered_velo_prod)); 102aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_filtered_velo_prod, &num_comp_filtered_velo_prod)); 103aa0b7f76SJames Wright } 104aa0b7f76SJames Wright 105aa0b7f76SJames Wright // -- Create operator for generating training data at nodes 106aa0b7f76SJames Wright // Differential Filter only provides filtered primitive variables 107aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSGS_DDAnisotropicTrainingDataNodal_Prim, 108aa0b7f76SJames Wright ComputeSGS_DDAnisotropicTrainingDataNodal_Prim_loc, &qf_sgs_dd_train)); 109aa0b7f76SJames Wright 110aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_train, sgs_dd_train_setup_data->sgs_dd_train_qfctx)); 111aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "q", num_comp_filtered_state, CEED_EVAL_NONE)); 112aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "velocity product", num_comp_filtered_velo_prod, CEED_EVAL_NONE)); 113aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE)); 114aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 115aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_train, "inverse multiplicity", 1, CEED_EVAL_NONE)); 116aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_train, "training data", sgs_dd_train->num_comp_dd_inputs, CEED_EVAL_NONE)); 117aa0b7f76SJames Wright 118aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_filtered_state, &filtered_fields, NULL)); 119aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_train, NULL, NULL, &op_sgs_dd_train)); 120c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "q", elem_restr_filtered_state, CEED_BASIS_NONE, filtered_fields)); 121c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "velocity product", elem_restr_filtered_velo_prod, CEED_BASIS_NONE, filtered_fields)); 122c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 123c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "anisotropy tensor", sgs_dd_train_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 124c6271fa9SJeremy L Thompson sgs_dd_train_setup_data->grid_aniso_ceed)); 125c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 126c6271fa9SJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_train, "training data", elem_restr_sgs_train, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 127aa0b7f76SJames Wright 128aa0b7f76SJames Wright PetscCall(OperatorApplyContextCreate(sgs_dd_train->filtered_grad_velo_proj->dm, sgs_dd_train->dm_dd_training, ceed, op_sgs_dd_train, NULL, NULL, 129aa0b7f76SJames Wright NULL, NULL, &sgs_dd_train->op_training_data_calc_ctx)); 130aa0b7f76SJames Wright 131aa0b7f76SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity)); 132aa0b7f76SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&filtered_fields)); 133aa0b7f76SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity)); 134aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_train)); 135aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_train)); 136aa0b7f76SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 137aa0b7f76SJames Wright } 138aa0b7f76SJames Wright 139991aef52SJames Wright PetscErrorCode SGS_DD_TrainingSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData problem) { 140aa0b7f76SJames Wright SGS_DDTrainingContext sgsdd_train_qfctx; 141aa0b7f76SJames Wright SGS_DD_TrainingSetupData sgs_dd_train_setup_data; 142aa0b7f76SJames Wright 143aa0b7f76SJames Wright PetscFunctionBeginUser; 144aa0b7f76SJames Wright if (!user->diff_filter) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem)); 145aa0b7f76SJames Wright if (!user->smartsim) PetscCall(SmartSimSetup(user)); 146aa0b7f76SJames Wright 147aa0b7f76SJames Wright PetscCall(PetscNew(&sgsdd_train_qfctx)); 148aa0b7f76SJames Wright PetscCall(PetscNew(&sgs_dd_train_setup_data)); 149aa0b7f76SJames Wright PetscCall(PetscNew(&user->sgs_dd_train)); 150aa0b7f76SJames Wright SGS_DD_TrainingData sgs_dd_train = user->sgs_dd_train; 151aa0b7f76SJames Wright 152aa0b7f76SJames Wright sgs_dd_train->overwrite_training_data = PETSC_TRUE; 153aa0b7f76SJames Wright sgs_dd_train->write_data_interval = 1; 1546ea7c1aeSJames Wright sgs_dd_train->num_filter_widths = sizeof(sgs_dd_train->filter_widths) / sizeof(sgs_dd_train->filter_widths[0]); 155aa0b7f76SJames Wright PetscOptionsBegin(user->comm, NULL, "SGS Data-Driven Training Options", NULL); 156aa0b7f76SJames Wright PetscCall(PetscOptionsInt("-sgs_train_write_data_interval", "Number of timesteps between writing data into database", NULL, 157aa0b7f76SJames Wright sgs_dd_train->write_data_interval, &sgs_dd_train->write_data_interval, NULL)); 158aa0b7f76SJames Wright PetscCall(PetscOptionsBool("-sgs_train_overwrite_data", "Overwrite old training data in the database", NULL, sgs_dd_train->overwrite_training_data, 159aa0b7f76SJames Wright &sgs_dd_train->overwrite_training_data, NULL)); 1606ea7c1aeSJames Wright PetscCall(PetscOptionsRealArray("-sgs_train_filter_width_scales", "Scales of each filter width put into training database", NULL, 1616ea7c1aeSJames Wright sgs_dd_train->filter_widths, &sgs_dd_train->num_filter_widths, NULL)); 162aa0b7f76SJames Wright PetscOptionsEnd(); 163aa0b7f76SJames Wright 164aa0b7f76SJames Wright // -- Create DM for storing training data 165aa0b7f76SJames Wright PetscCall(SGS_DD_TrainingCreateDM(user->dm, &sgs_dd_train->dm_dd_training, user->app_ctx->degree, user->app_ctx->q_extra, 166aa0b7f76SJames Wright &sgs_dd_train->num_comp_dd_inputs)); 167aa0b7f76SJames Wright 168aa0b7f76SJames Wright { // -- Create QFunction Context 169aa0b7f76SJames Wright NewtonianIdealGasContext gas; 170e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &gas)); 171aa0b7f76SJames Wright sgsdd_train_qfctx->gas = *gas; 172e07531f7SJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &gas)); 173aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &sgs_dd_train_setup_data->sgs_dd_train_qfctx)); 174aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetData(sgs_dd_train_setup_data->sgs_dd_train_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, 175aa0b7f76SJames Wright sizeof(*sgsdd_train_qfctx), sgsdd_train_qfctx)); 176aa0b7f76SJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_train_setup_data->sgs_dd_train_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 177aa0b7f76SJames Wright } 178aa0b7f76SJames Wright 179aa0b7f76SJames Wright { // -- Send training data array info to SmartRedis database 180aa0b7f76SJames Wright PetscMPIInt rank, num_ranks; 181aa0b7f76SJames Wright SmartSimData smartsim = user->smartsim; 182aa0b7f76SJames Wright PetscCallMPI(MPI_Comm_rank(user->comm, &rank)); 183aa0b7f76SJames Wright PetscCallMPI(MPI_Comm_size(user->comm, &num_ranks)); 184aa0b7f76SJames Wright 185aa0b7f76SJames Wright { 186aa0b7f76SJames Wright PetscSection global_section; 1877ff16c02SJames Wright PetscInt num_dofs, num_comps, local_min_max[2] = {0.}, global_min_max[2] = {0.}; 1887ff16c02SJames Wright 189aa0b7f76SJames Wright PetscCall(DMGetGlobalSection(sgs_dd_train->dm_dd_training, &global_section)); 190aa0b7f76SJames Wright PetscCall(DMGetGlobalVectorInfo(sgs_dd_train->dm_dd_training, &num_dofs, NULL, NULL)); 191aa0b7f76SJames Wright PetscCall(PetscSectionGetFieldComponents(global_section, 0, &num_comps)); 1927ff16c02SJames Wright local_min_max[0] = num_dofs; 1937ff16c02SJames Wright PetscCall(PetscGlobalMinMaxInt(user->comm, local_min_max, global_min_max)); 1947ff16c02SJames Wright 1957ff16c02SJames Wright sgs_dd_train->training_data_array_dims[0] = global_min_max[0] / num_comps; 196aa0b7f76SJames Wright sgs_dd_train->training_data_array_dims[1] = num_comps; 197aa0b7f76SJames Wright } 198aa0b7f76SJames Wright 199aa0b7f76SJames Wright if (rank % smartsim->collocated_database_num_ranks == 0) { 2004fa1625aSJames Wright { // Communicate info on simulation size 2014fa1625aSJames Wright const char tensor_name[] = "sizeInfo"; 202aa0b7f76SJames Wright size_t array_info_dim = 6; 203f4632befSRiccardo Balin PetscInt64 array_info[6] = {0}, num_features = 6; 204aa0b7f76SJames Wright 205aa0b7f76SJames Wright array_info[0] = sgs_dd_train->training_data_array_dims[0]; 206aa0b7f76SJames Wright array_info[1] = sgs_dd_train->training_data_array_dims[1]; 207aa0b7f76SJames Wright array_info[2] = num_features; 208aa0b7f76SJames Wright array_info[3] = num_ranks; 209aa0b7f76SJames Wright array_info[4] = smartsim->collocated_database_num_ranks; 210aa0b7f76SJames Wright array_info[5] = rank; 211aa0b7f76SJames Wright 212ad2e713eSRiccardo Balin PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 21343e9749fSJames Wright PetscCallSmartRedis( 2144fa1625aSJames Wright put_tensor(smartsim->client, tensor_name, strlen(tensor_name), array_info, &array_info_dim, 1, SRTensorTypeInt64, SRMemLayoutContiguous)); 2154fa1625aSJames Wright PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name))); 216ad2e713eSRiccardo Balin PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 2174fa1625aSJames Wright } 218aa0b7f76SJames Wright 2194fa1625aSJames Wright { // Send array that communicates if tensors are overwritten in database 2204fa1625aSJames Wright const char tensor_name[] = "tensor-ow"; 221f4632befSRiccardo Balin PetscInt64 tensor_overwrite[2] = {sgs_dd_train->overwrite_training_data}; 222aa0b7f76SJames Wright size_t dim_2[1] = {2}; 2234fa1625aSJames Wright 224ad2e713eSRiccardo Balin PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 22543e9749fSJames Wright PetscCallSmartRedis( 2264fa1625aSJames Wright put_tensor(smartsim->client, tensor_name, strlen(tensor_name), tensor_overwrite, dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous)); 2274fa1625aSJames Wright PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name))); 228ad2e713eSRiccardo Balin PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 229aa0b7f76SJames Wright } 2306ea7c1aeSJames Wright 2316ea7c1aeSJames Wright { // Communicate number of filter widths used 2326ea7c1aeSJames Wright const char tensor_name[] = "num_filter_widths"; 2336ea7c1aeSJames Wright PetscInt64 num_filter_widths = sgs_dd_train->num_filter_widths; 2346ea7c1aeSJames Wright size_t dim_2 = 1; 2356ea7c1aeSJames Wright 2366ea7c1aeSJames Wright PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 23743e9749fSJames Wright PetscCallSmartRedis( 2386ea7c1aeSJames Wright put_tensor(smartsim->client, tensor_name, strlen(tensor_name), &num_filter_widths, &dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous)); 2396ea7c1aeSJames Wright PetscCall(SmartRedisVerifyPutTensor(smartsim->client, tensor_name, strlen(tensor_name))); 2406ea7c1aeSJames Wright PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 2416ea7c1aeSJames Wright } 242aa0b7f76SJames Wright } 2434fa1625aSJames Wright } 244aa0b7f76SJames Wright 245aa0b7f76SJames Wright // -- Compute and store anisotropy tensor 246aa0b7f76SJames Wright PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, user, ceed_data, &sgs_dd_train_setup_data->elem_restr_grid_aniso, 247aa0b7f76SJames Wright &sgs_dd_train_setup_data->grid_aniso_ceed)); 248aa0b7f76SJames Wright 249aa0b7f76SJames Wright // -- Create Nodal Evaluation Operator 250aa0b7f76SJames Wright PetscCall(SetupTrainingDataCalculation(ceed, user, ceed_data, problem, sgs_dd_train_setup_data)); 251aa0b7f76SJames Wright 252aa0b7f76SJames Wright PetscCall(SGS_DD_TrainingSetupDataDestroy(sgs_dd_train_setup_data)); 253aa0b7f76SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 254aa0b7f76SJames Wright } 255aa0b7f76SJames Wright 256aa0b7f76SJames Wright PetscErrorCode TSMonitor_SGS_DD_Training(TS ts, PetscInt step_num, PetscReal solution_time, Vec Q, void *ctx) { 257aa0b7f76SJames Wright User user = (User)ctx; 258aa0b7f76SJames Wright Ceed ceed = user->ceed; 259aa0b7f76SJames Wright SGS_DD_TrainingData sgs_dd_train = user->sgs_dd_train; 260aa0b7f76SJames Wright SmartSimData smartsim = user->smartsim; 261aa0b7f76SJames Wright Vec TrainingData; 262ee81c423SRiccardo Balin PetscMPIInt rank; 263aa0b7f76SJames Wright 264aa0b7f76SJames Wright PetscFunctionBeginUser; 265ee81c423SRiccardo Balin 266ee81c423SRiccardo Balin PetscCallMPI(MPI_Comm_rank(user->comm, &rank)); 267ee81c423SRiccardo Balin 268aa0b7f76SJames Wright if (step_num % sgs_dd_train->write_data_interval != 0) PetscFunctionReturn(PETSC_SUCCESS); 269aa0b7f76SJames Wright PetscCall(DMGetGlobalVector(sgs_dd_train->dm_dd_training, &TrainingData)); 270aa0b7f76SJames Wright 2716ea7c1aeSJames Wright for (PetscInt filter_index = 0; filter_index < sgs_dd_train->num_filter_widths; filter_index++) { 272ad2e713eSRiccardo Balin PetscCall(PetscLogEventBegin(FLUIDS_TrainDataCompute, 0, 0, 0, 0)); 273aa0b7f76SJames Wright { // -- Compute and assemble training data 274aa0b7f76SJames Wright Vec FilteredVelocityGradient, FilteredFields, FilteredFields_loc; 275aa0b7f76SJames Wright PetscMemType filtered_fields_mem_type; 276aa0b7f76SJames Wright CeedVector filtered_fields; 277aa0b7f76SJames Wright 2786ea7c1aeSJames Wright { // Set filter width for the current solve 2796ea7c1aeSJames Wright double filter_width_scaling[3]; 2806ea7c1aeSJames Wright CeedOperator op_mat; 2816ea7c1aeSJames Wright Mat A_mat; 2826ea7c1aeSJames Wright 2836ea7c1aeSJames Wright for (int j = 0; j < 3; j++) filter_width_scaling[j] = sgs_dd_train->filter_widths[filter_index]; 2846ea7c1aeSJames Wright PetscCall(KSPGetOperators(user->diff_filter->ksp, &A_mat, NULL)); 2856ea7c1aeSJames Wright PetscCall(MatCeedGetCeedOperators(A_mat, &op_mat, NULL)); 2866ea7c1aeSJames Wright PetscCall(CeedOperatorSetContextDouble(op_mat, user->diff_filter->filter_width_scaling_label, filter_width_scaling)); 2876ea7c1aeSJames Wright } 2886ea7c1aeSJames Wright 289aa0b7f76SJames Wright PetscCall(DMGetGlobalVector(user->diff_filter->dm_filter, &FilteredFields)); 290aa0b7f76SJames Wright PetscCall(DMGetLocalVector(user->diff_filter->dm_filter, &FilteredFields_loc)); 291aa0b7f76SJames Wright 292aa0b7f76SJames Wright PetscCall(DifferentialFilterApply(user, solution_time, Q, FilteredFields)); 293aa0b7f76SJames Wright PetscCall(DMGlobalToLocal(user->diff_filter->dm_filter, FilteredFields, INSERT_VALUES, FilteredFields_loc)); 294aa0b7f76SJames Wright 295aa0b7f76SJames Wright PetscCall(DMGetGlobalVector(sgs_dd_train->filtered_grad_velo_proj->dm, &FilteredVelocityGradient)); 296aa0b7f76SJames Wright PetscCall(VelocityGradientProjectionApply(sgs_dd_train->filtered_grad_velo_proj, FilteredFields_loc, FilteredVelocityGradient)); 297aa0b7f76SJames Wright 298aa0b7f76SJames Wright { 299aa0b7f76SJames Wright CeedOperatorField op_field; 3006ea7c1aeSJames Wright 301aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(sgs_dd_train->op_training_data_calc_ctx->op, "q", &op_field)); 302aa0b7f76SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetVector(op_field, &filtered_fields)); 303aa0b7f76SJames Wright } 3046ea7c1aeSJames Wright 305a7dac1d5SJames Wright PetscCall(VecPetscToCeed(FilteredFields_loc, &filtered_fields_mem_type, filtered_fields)); // filtered_fields is an implicit input 306aa0b7f76SJames Wright PetscCall(ApplyCeedOperatorGlobalToGlobal(FilteredVelocityGradient, TrainingData, sgs_dd_train->op_training_data_calc_ctx)); 307a7dac1d5SJames Wright PetscCall(VecCeedToPetsc(filtered_fields, filtered_fields_mem_type, FilteredFields_loc)); 308aa0b7f76SJames Wright 309aa0b7f76SJames Wright PetscCall(DMRestoreGlobalVector(sgs_dd_train->filtered_grad_velo_proj->dm, &FilteredVelocityGradient)); 310aa0b7f76SJames Wright PetscCall(DMRestoreGlobalVector(user->diff_filter->dm_filter, &FilteredFields)); 311aa0b7f76SJames Wright PetscCall(DMRestoreLocalVector(user->diff_filter->dm_filter, &FilteredFields_loc)); 312aa0b7f76SJames Wright } 313ad2e713eSRiccardo Balin PetscCall(PetscLogEventEnd(FLUIDS_TrainDataCompute, 0, 0, 0, 0)); 314aa0b7f76SJames Wright 315aa0b7f76SJames Wright { // -- Send training data to SmartSim 316aa0b7f76SJames Wright char array_key[PETSC_MAX_PATH_LEN]; 317aa0b7f76SJames Wright size_t array_key_len; 318aa0b7f76SJames Wright 319aa0b7f76SJames Wright if (sgs_dd_train->overwrite_training_data) { 3206ea7c1aeSJames Wright PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT, smartsim->rank_id_name, filter_index)); 321aa0b7f76SJames Wright } else { 3226ea7c1aeSJames Wright PetscCall(PetscSNPrintf(array_key, sizeof array_key, "%s.%" PetscInt_FMT "%" PetscInt_FMT, smartsim->rank_id_name, step_num, filter_index)); 323aa0b7f76SJames Wright } 324aa0b7f76SJames Wright PetscCall(PetscStrlen(array_key, &array_key_len)); 325aa0b7f76SJames Wright 326aa0b7f76SJames Wright { 327aa0b7f76SJames Wright const PetscScalar *training_data; 328aa0b7f76SJames Wright PetscCall(VecGetArrayRead(TrainingData, &training_data)); 329ad2e713eSRiccardo Balin PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Train, 0, 0, 0, 0)); 33043e9749fSJames Wright PetscCallSmartRedis(put_tensor(smartsim->client, array_key, array_key_len, (void *)training_data, sgs_dd_train->training_data_array_dims, 2, 331aa0b7f76SJames Wright SRTensorTypeDouble, SRMemLayoutContiguous)); 332ad2e713eSRiccardo Balin PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Train, 0, 0, 0, 0)); 333aa0b7f76SJames Wright PetscCall(VecRestoreArrayRead(TrainingData, &training_data)); 334aa0b7f76SJames Wright } 335ee81c423SRiccardo Balin } 336ee81c423SRiccardo Balin } 337aa0b7f76SJames Wright 338aa0b7f76SJames Wright if (rank % smartsim->collocated_database_num_ranks == 0) { 3396ea7c1aeSJames Wright const char tensor_name[] = "step"; 340aa0b7f76SJames Wright size_t dim_2[1] = {2}; 341f4632befSRiccardo Balin PetscInt64 step_array[2] = {step_num, step_num}; 3426ea7c1aeSJames Wright 343ad2e713eSRiccardo Balin PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 34443e9749fSJames Wright PetscCallSmartRedis( 3456ea7c1aeSJames Wright put_tensor(smartsim->client, tensor_name, strlen(tensor_name), step_array, dim_2, 1, SRTensorTypeInt64, SRMemLayoutContiguous)); 346ad2e713eSRiccardo Balin PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 347aa0b7f76SJames Wright } 348aa0b7f76SJames Wright 349aa0b7f76SJames Wright PetscCall(DMRestoreGlobalVector(user->sgs_dd_train->dm_dd_training, &TrainingData)); 350aa0b7f76SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 351aa0b7f76SJames Wright } 352aa0b7f76SJames Wright 353632a41e1SJames Wright PetscErrorCode TSPostStep_SGS_DD_Training(TS ts) { 354632a41e1SJames Wright User user; 355632a41e1SJames Wright const char check_run_key[] = "check-run"; 356632a41e1SJames Wright PetscReal check_run[2] = {1}; 357632a41e1SJames Wright const size_t check_run_dims[1] = {2}; 358632a41e1SJames Wright size_t check_run_key_size; 359632a41e1SJames Wright 360632a41e1SJames Wright PetscFunctionBeginUser; 361632a41e1SJames Wright PetscCall(PetscStrlen(check_run_key, &check_run_key_size)); 362632a41e1SJames Wright PetscCall(TSGetApplicationContext(ts, &user)); 363632a41e1SJames Wright SmartSimData smartsim = user->smartsim; 364632a41e1SJames Wright 365ad2e713eSRiccardo Balin PetscCall(PetscLogEventBegin(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 36643e9749fSJames Wright PetscCallSmartRedis( 367632a41e1SJames Wright unpack_tensor(smartsim->client, check_run_key, check_run_key_size, check_run, check_run_dims, 1, SRTensorTypeDouble, SRMemLayoutContiguous)); 368ad2e713eSRiccardo Balin PetscCall(PetscLogEventEnd(FLUIDS_SmartRedis_Meta, 0, 0, 0, 0)); 369632a41e1SJames Wright if (check_run[0] == 0) { 370632a41e1SJames Wright PetscCall(PetscPrintf(user->comm, "-- Simulation stopped by 'check-run' tensor in Redis database\n")); 371632a41e1SJames Wright PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER)); 372632a41e1SJames Wright } 373632a41e1SJames Wright 374632a41e1SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 375632a41e1SJames Wright } 376632a41e1SJames Wright 377aa0b7f76SJames Wright PetscErrorCode SGS_DD_TrainingDataDestroy(SGS_DD_TrainingData sgs_dd_train) { 378aa0b7f76SJames Wright PetscFunctionBeginUser; 379aa0b7f76SJames Wright if (!sgs_dd_train) PetscFunctionReturn(PETSC_SUCCESS); 380aa0b7f76SJames Wright 381aa0b7f76SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_train->op_training_data_calc_ctx)); 382aa0b7f76SJames Wright PetscCall(NodalProjectionDataDestroy(sgs_dd_train->filtered_grad_velo_proj)); 383aa0b7f76SJames Wright PetscCall(DMDestroy(&sgs_dd_train->dm_dd_training)); 384aa0b7f76SJames Wright PetscCall(PetscFree(sgs_dd_train)); 385aa0b7f76SJames Wright 386aa0b7f76SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 387aa0b7f76SJames Wright } 388