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_model.h" 5 6 #include <petscdmplex.h> 7 8 #include <navierstokes.h> 9 #include <sgs_model_torch.h> 10 11 typedef PetscErrorCode (*SgsDDNodalStressEval)(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc); 12 typedef PetscErrorCode (*SgsDDNodalStressInference)(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx); 13 typedef struct { 14 DM dm_sgs, dm_dd_inputs, dm_dd_outputs; 15 PetscInt num_comp_sgs, num_comp_inputs, num_comp_outputs; 16 OperatorApplyContext op_nodal_evaluation_ctx, op_nodal_dd_inputs_ctx, op_nodal_dd_outputs_ctx, op_sgs_apply_ctx; 17 CeedVector sgs_nodal_ceed, grad_velo_ceed; 18 SgsDDNodalStressEval sgs_nodal_eval; 19 SgsDDNodalStressInference sgs_nodal_inference; 20 void *sgs_nodal_inference_ctx; 21 PetscErrorCode (*sgs_nodal_inference_ctx_destroy)(void *ctx); 22 } *SgsDDData; 23 24 // @brief Destroy `SgsDDData` object 25 static PetscErrorCode SgsDDDataDestroy(SgsDDData *sgs_dd_data) { 26 SgsDDData sgs_dd_data_ = *sgs_dd_data; 27 28 PetscFunctionBeginUser; 29 if (!sgs_dd_data_) PetscFunctionReturn(PETSC_SUCCESS); 30 Ceed ceed = sgs_dd_data_->op_sgs_apply_ctx->ceed; 31 32 PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data_->sgs_nodal_ceed)); 33 PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data_->grad_velo_ceed)); 34 PetscCall(OperatorApplyContextDestroy(sgs_dd_data_->op_nodal_evaluation_ctx)); 35 PetscCall(OperatorApplyContextDestroy(sgs_dd_data_->op_sgs_apply_ctx)); 36 PetscCall(OperatorApplyContextDestroy(sgs_dd_data_->op_nodal_dd_inputs_ctx)); 37 PetscCall(OperatorApplyContextDestroy(sgs_dd_data_->op_nodal_dd_outputs_ctx)); 38 PetscCall(DMDestroy(&sgs_dd_data_->dm_sgs)); 39 PetscCall(DMDestroy(&sgs_dd_data_->dm_dd_inputs)); 40 PetscCall(DMDestroy(&sgs_dd_data_->dm_dd_outputs)); 41 if (sgs_dd_data_->sgs_nodal_inference_ctx) PetscCall(sgs_dd_data_->sgs_nodal_inference_ctx_destroy(sgs_dd_data_->sgs_nodal_inference_ctx)); 42 PetscCall(PetscFree(sgs_dd_data_)); 43 *sgs_dd_data = NULL; 44 PetscFunctionReturn(PETSC_SUCCESS); 45 } 46 47 typedef struct { 48 CeedElemRestriction elem_restr_grid_aniso, elem_restr_sgs; 49 CeedVector grid_aniso_ceed; 50 CeedQFunctionContext sgsdd_qfctx, ifunction_qfctx; 51 SGSModelDDImplementation sgs_dd_model_implementation; 52 } *SgsDDSetupData; 53 54 #define GRAD_VELO_PROJ_KEY "Gradient of Velocity Projection" 55 #define SGS_DD_DATA_KEY "SGS Data Driven Data" 56 57 PetscErrorCode SgsDDSetupDataDestroy(SgsDDSetupData sgs_dd_setup_data) { 58 Ceed ceed; 59 60 PetscFunctionBeginUser; 61 PetscCall(CeedElemRestrictionGetCeed(sgs_dd_setup_data->elem_restr_sgs, &ceed)); 62 63 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_grid_aniso)); 64 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_sgs)); 65 PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_setup_data->grid_aniso_ceed)); 66 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->sgsdd_qfctx)); 67 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->ifunction_qfctx)); 68 PetscCall(PetscFree(sgs_dd_setup_data)); 69 PetscCheck(CeedDestroy(&ceed) == CEED_ERROR_SUCCESS, PETSC_COMM_SELF, PETSC_ERR_LIB, "Destroying Ceed object failed"); 70 PetscFunctionReturn(PETSC_SUCCESS); 71 } 72 73 // @brief Create DM for storing subgrid stress at nodes 74 static PetscErrorCode SgsDDCreateDM(DM dm_source, DM *dm_sgs, PetscInt degree, PetscInt q_extra, PetscInt *num_components) { 75 PetscSection section; 76 77 PetscFunctionBeginUser; 78 *num_components = 6; 79 80 PetscCall(DMClone(dm_source, dm_sgs)); 81 PetscCall(DMSetMatrixPreallocateSkip(*dm_sgs, PETSC_TRUE)); 82 PetscCall(PetscObjectSetName((PetscObject)*dm_sgs, "Subgrid Stress Projection")); 83 84 PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_sgs)); 85 86 PetscCall(DMGetLocalSection(*dm_sgs, §ion)); 87 PetscCall(PetscSectionSetFieldName(section, 0, "")); 88 PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMSubgridStressXX")); 89 PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMSubgridStressYY")); 90 PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMSubgridStressZZ")); 91 PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMSubgridStressYZ")); 92 PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMSubgridStressXZ")); 93 PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMSubgridStressXY")); 94 PetscFunctionReturn(PETSC_SUCCESS); 95 }; 96 97 // @brief Evaluate data-driven SGS using fused method 98 static PetscErrorCode SgsDDNodalStressEval_Fused(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) { 99 SgsDDData sgs_dd_data; 100 PetscMemType q_mem_type; 101 102 PetscFunctionBeginUser; 103 PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 104 PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed)); // q_ceed is an implicit input 105 106 PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, SGSNodal_loc, sgs_dd_data->op_nodal_evaluation_ctx)); 107 108 PetscCall(VecCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc)); 109 PetscFunctionReturn(PETSC_SUCCESS); 110 } 111 112 // @brief Create CeedOperator to calculate data-drive SGS at nodes using fused operator 113 static PetscErrorCode SgsDDSetupNodalEvaluation_Fused(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) { 114 SgsDDData sgs_dd_data; 115 CeedQFunction qf_sgs_dd_nodal; 116 CeedOperator op_sgs_dd_nodal; 117 CeedInt num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso; 118 PetscInt dim; 119 CeedVector inv_multiplicity; 120 CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs; 121 DMLabel domain_label = NULL; 122 PetscInt label_value = 0, height = 0, dm_field = 0; 123 NodalProjectionData grad_velo_proj; 124 125 PetscFunctionBeginUser; 126 PetscCall(DMGetDimension(honee->dm, &dim)); 127 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x)); 128 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q)); 129 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso)); 130 PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 131 PetscCall(HoneeGetContainer(honee, GRAD_VELO_PROJ_KEY, &grad_velo_proj)); 132 133 { // Get velocity gradient information 134 CeedOperatorField op_field; 135 PetscCallCeed(ceed, CeedOperatorGetFieldByName(grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field)); 136 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo)); 137 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo)); 138 } 139 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs)); 140 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL)); 141 142 PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity, 143 &inv_multiplicity)); 144 145 // -- Create operator for SGS DD model nodal evaluation 146 switch (honee->phys->state_var) { 147 case STATEVAR_PRIMITIVE: 148 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Prim, ComputeSgsDDNodal_Prim_loc, &qf_sgs_dd_nodal)); 149 break; 150 case STATEVAR_CONSERVATIVE: 151 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Conserv, ComputeSgsDDNodal_Conserv_loc, &qf_sgs_dd_nodal)); 152 break; 153 case STATEVAR_ENTROPY: 154 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Entropy, ComputeSgsDDNodal_Entropy_loc, &qf_sgs_dd_nodal)); 155 break; 156 } 157 158 // Mesh/geometry order and solution basis order may differ, therefore must interpolate 159 CeedBasis basis_x_to_q; 160 PetscCallCeed(ceed, CeedBasisCreateProjection(honee->basis_x, honee->basis_q, &basis_x_to_q)); 161 162 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx)); 163 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE)); 164 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP)); 165 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE)); 166 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 167 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE)); 168 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE)); 169 170 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal)); 171 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "q", honee->elem_restr_q, CEED_BASIS_NONE, honee->q_ceed)); 172 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "x", honee->elem_restr_x, basis_x_to_q, honee->x_coord)); 173 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 174 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 175 sgs_dd_setup_data->grid_aniso_ceed)); 176 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 177 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 178 179 PetscCall(OperatorApplyContextCreate(grad_velo_proj->dm, sgs_dd_data->dm_sgs, ceed, op_sgs_dd_nodal, NULL, sgs_dd_data->sgs_nodal_ceed, NULL, NULL, 180 &sgs_dd_data->op_nodal_evaluation_ctx)); 181 182 sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs; 183 sgs_dd_data->sgs_nodal_eval = SgsDDNodalStressEval_Fused; 184 185 PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity)); 186 PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q)); 187 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity)); 188 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo)); 189 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_nodal)); 190 PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_nodal)); 191 PetscFunctionReturn(PETSC_SUCCESS); 192 } 193 194 // @brief Setup data-driven model inference using libCEED native implementation 195 static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Ceed(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data, 196 CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs, 197 CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity, 198 void **ctx) { 199 CeedQFunction qf_sgs_dd_inference; 200 CeedOperator op_sgs_dd_inference; 201 OperatorApplyContext *op_context = (OperatorApplyContext *)ctx; 202 203 PetscFunctionBeginUser; 204 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inference, ComputeSgsDDNodal_Sequential_Inference_loc, 205 &qf_sgs_dd_inference)); 206 207 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inference, sgs_dd_setup_data->sgsdd_qfctx)); 208 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE)); 209 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "inverse multiplicity", 1, CEED_EVAL_NONE)); 210 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inference, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE)); 211 212 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inference, NULL, NULL, &op_sgs_dd_inference)); 213 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 214 PetscCallCeed(ceed, 215 CeedOperatorSetField(op_sgs_dd_inference, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 216 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 217 218 PetscCall(OperatorApplyContextCreate(sgs_dd_data->dm_dd_inputs, sgs_dd_data->dm_dd_outputs, ceed, op_sgs_dd_inference, NULL, NULL, NULL, NULL, 219 op_context)); 220 sgs_dd_data->sgs_nodal_inference_ctx_destroy = (PetscErrorCode(*)(void *))OperatorApplyContextDestroy; 221 222 PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inference)); 223 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inference)); 224 PetscFunctionReturn(PETSC_SUCCESS); 225 } 226 227 // @brief Perform data-driven model inference using libCEED native implementation 228 PetscErrorCode SgsDDNodalStressEval_Sequential_Ceed(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) { 229 OperatorApplyContext op_context = *(OperatorApplyContext *)ctx; 230 231 PetscFunctionBeginUser; 232 PetscCall(PetscLogEventBegin(HONEE_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL)); 233 PetscCall(PetscLogEventBegin(HONEE_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL)); 234 PetscCall(PetscLogGpuTimeBegin()); 235 PetscCall(ApplyCeedOperatorLocalToLocal(DD_Inputs_loc, DD_Outputs_loc, op_context)); 236 PetscCall(PetscLogGpuTimeEnd()); 237 PetscCall(PetscLogEventEnd(HONEE_SgsModelDDInference, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL)); 238 PetscCall(PetscLogEventEnd(HONEE_SgsModelDDData, DD_Inputs_loc, DD_Outputs_loc, NULL, NULL)); 239 PetscFunctionReturn(PETSC_SUCCESS); 240 } 241 242 // @brief Setup data-driven model inference using libtorch 243 static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Torch(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data, 244 CeedElemRestriction elem_restr_dd_inputs, CeedElemRestriction elem_restr_dd_outputs, 245 CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity, 246 void **ctx) { 247 const char *ceed_resource; 248 char model_path[PETSC_MAX_PATH_LEN] = ""; 249 TorchDeviceType model_device_type; 250 251 PetscFunctionBeginUser; 252 PetscCallCeed(ceed, CeedGetResource(ceed, &ceed_resource)); 253 if (strstr(ceed_resource, "/gpu/cuda")) model_device_type = TORCH_DEVICE_CUDA; 254 else if (strstr(ceed_resource, "/gpu/hip")) model_device_type = TORCH_DEVICE_HIP; 255 // On-device XPU is not working reliably currently, default to CPU inference evaluation 256 // else if (strstr(ceed_resource, "/gpu/sycl")) model_device_type = TORCH_DEVICE_XPU; 257 else model_device_type = TORCH_DEVICE_CPU; 258 PetscCall(PetscOptionsGetEnum(NULL, NULL, "-sgs_model_dd_torch_model_device", TorchDeviceTypes, (PetscEnum *)&model_device_type, NULL)); 259 PetscCall(PetscOptionsGetString(NULL, NULL, "-sgs_model_dd_torch_model_path", model_path, sizeof(model_path), NULL)); 260 261 PetscCall(LoadModel_Torch(model_path, model_device_type)); 262 PetscFunctionReturn(PETSC_SUCCESS); 263 } 264 265 // @brief Perform data-driven model inference using libtorch 266 static PetscErrorCode SgsDDNodalStressEval_Sequential_Torch(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) { 267 static PetscBool run_through = PETSC_FALSE; 268 269 PetscFunctionBeginUser; 270 if (!run_through) { 271 PetscCall(VecViewFromOptions(DD_Inputs_loc, NULL, "-dd_inputs_loc_view")); 272 } 273 PetscCall(ModelInference_Torch(DD_Inputs_loc, DD_Outputs_loc)); 274 if (!run_through) { 275 PetscCall(VecViewFromOptions(DD_Outputs_loc, NULL, "-dd_outputs_loc_view")); 276 run_through = PETSC_TRUE; 277 } 278 PetscFunctionReturn(PETSC_SUCCESS); 279 } 280 281 // @brief Evaluate data-driven SGS using sequential method 282 PetscErrorCode SgsDDNodalStressEval_Sequential(Honee honee, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) { 283 SgsDDData sgs_dd_data; 284 PetscMemType q_mem_type; 285 Vec DD_Inputs_loc, DD_Outputs_loc; 286 287 PetscFunctionBeginUser; 288 PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 289 PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc)); 290 PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc)); 291 PetscCall(VecPetscToCeed(Q_loc, &q_mem_type, honee->q_ceed)); // q_ceed is an implicit input 292 293 PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, DD_Inputs_loc, sgs_dd_data->op_nodal_dd_inputs_ctx)); 294 PetscCall(sgs_dd_data->sgs_nodal_inference(DD_Inputs_loc, DD_Outputs_loc, &sgs_dd_data->sgs_nodal_inference_ctx)); 295 PetscCall(ApplyCeedOperatorLocalToLocal(DD_Outputs_loc, SGSNodal_loc, sgs_dd_data->op_nodal_dd_outputs_ctx)); 296 297 PetscCall(VecCeedToPetsc(honee->q_ceed, q_mem_type, Q_loc)); 298 PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc)); 299 PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc)); 300 PetscFunctionReturn(PETSC_SUCCESS); 301 } 302 303 // @brief Create CeedOperator to calculate data-drive SGS at nodes using sequentially-applied operators 304 static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) { 305 SgsDDData sgs_dd_data; 306 CeedInt num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso, num_comp_eigvec = 9 + 1; 307 PetscInt dim; 308 CeedVector inv_multiplicity, eigvec; 309 NodalProjectionData grad_velo_proj; 310 CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs, elem_restr_eigvec, elem_restr_dd_inputs, 311 elem_restr_dd_outputs; 312 DMLabel domain_label = NULL; 313 PetscInt label_value = 0, height = 0, dm_field = 0; 314 315 PetscFunctionBeginUser; 316 PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 317 { // Create DMs for data-driven input and output values 318 PetscSection section; 319 PetscInt degree, q_extra; 320 { // Get degree and number of quadrature points from dm_sgs 321 PetscFE fe; 322 PetscSpace basis; 323 PetscQuadrature quadrature; 324 PetscInt num_qpnts; 325 PetscCall(DMGetField(sgs_dd_data->dm_sgs, 0, NULL, (PetscObject *)&fe)); 326 PetscCall(PetscFEGetBasisSpace(fe, &basis)); 327 PetscCall(PetscSpaceGetDegree(basis, °ree, NULL)); 328 PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 329 PetscCall(PetscQuadratureGetOrder(quadrature, &num_qpnts)); 330 q_extra = degree - num_qpnts; 331 } 332 333 PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_inputs)); 334 PetscCall(DMSetMatrixPreallocateSkip(sgs_dd_data->dm_dd_inputs, PETSC_TRUE)); 335 PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_inputs, "Data-Driven Model Inputs")); 336 PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &sgs_dd_data->num_comp_inputs, sgs_dd_data->dm_dd_inputs)); 337 PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_inputs, §ion)); 338 PetscCall(PetscSectionSetFieldName(section, 0, "")); 339 for (CeedInt i = 0; i < sgs_dd_data->num_comp_inputs; i++) { 340 char component_name[PETSC_MAX_PATH_LEN]; 341 342 PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenInput%" CeedInt_FMT, i + 1)); 343 PetscCall(PetscSectionSetComponentName(section, 0, i, component_name)); 344 } 345 346 PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_outputs)); 347 PetscCall(DMSetMatrixPreallocateSkip(sgs_dd_data->dm_dd_outputs, PETSC_TRUE)); 348 PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_outputs, "Data-Driven Model Outputs")); 349 PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, &sgs_dd_data->num_comp_outputs, sgs_dd_data->dm_dd_outputs)); 350 PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_outputs, §ion)); 351 PetscCall(PetscSectionSetFieldName(section, 0, "")); 352 for (CeedInt i = 0; i < sgs_dd_data->num_comp_outputs; i++) { 353 char component_name[PETSC_MAX_PATH_LEN]; 354 355 PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenOutput%" CeedInt_FMT, i + 1)); 356 PetscCall(PetscSectionSetComponentName(section, 0, i, component_name)); 357 } 358 } 359 360 PetscCall(DMGetDimension(honee->dm, &dim)); 361 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x)); 362 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q)); 363 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso)); 364 PetscCall(HoneeGetContainer(honee, GRAD_VELO_PROJ_KEY, &grad_velo_proj)); 365 366 { // Get velocity gradient information 367 CeedOperatorField op_field; 368 PetscCallCeed(ceed, CeedOperatorGetFieldByName(grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field)); 369 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo)); 370 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo)); 371 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_grad_velo, &sgs_dd_data->grad_velo_ceed, NULL)); 372 } 373 374 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs)); 375 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL)); 376 PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, num_comp_eigvec, 377 &elem_restr_eigvec)); 378 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_eigvec, &eigvec, NULL)); 379 380 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_inputs, domain_label, label_value, height, dm_field, &elem_restr_dd_inputs)); 381 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_outputs, domain_label, label_value, height, dm_field, &elem_restr_dd_outputs)); 382 383 PetscCall(GetInverseMultiplicity(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, PETSC_FALSE, &elem_restr_inv_multiplicity, 384 &inv_multiplicity)); 385 386 { // Create operator for data-driven input evaluation 387 CeedQFunction qf_sgs_dd_inputs; 388 CeedOperator op_sgs_dd_inputs; 389 390 switch (honee->phys->state_var) { 391 case STATEVAR_PRIMITIVE: 392 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Prim, 393 ComputeSgsDDNodal_Sequential_Inputs_Prim_loc, &qf_sgs_dd_inputs)); 394 break; 395 case STATEVAR_CONSERVATIVE: 396 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Conserv, 397 ComputeSgsDDNodal_Sequential_Inputs_Conserv_loc, &qf_sgs_dd_inputs)); 398 break; 399 case STATEVAR_ENTROPY: 400 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Entropy, 401 ComputeSgsDDNodal_Sequential_Inputs_Entropy_loc, &qf_sgs_dd_inputs)); 402 break; 403 } 404 405 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inputs, sgs_dd_setup_data->sgsdd_qfctx)); 406 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "q", num_comp_q, CEED_EVAL_NONE)); 407 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE)); 408 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 409 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "inverse multiplicity", 1, CEED_EVAL_NONE)); 410 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE)); 411 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE)); 412 413 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inputs, NULL, NULL, &op_sgs_dd_inputs)); 414 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "q", honee->elem_restr_q, CEED_BASIS_NONE, honee->q_ceed)); 415 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 416 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 417 sgs_dd_setup_data->grid_aniso_ceed)); 418 PetscCallCeed(ceed, 419 CeedOperatorSetField(op_sgs_dd_inputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 420 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec)); 421 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 422 423 PetscCall(OperatorApplyContextCreate(grad_velo_proj->dm, sgs_dd_data->dm_dd_inputs, ceed, op_sgs_dd_inputs, NULL, NULL, NULL, NULL, 424 &sgs_dd_data->op_nodal_dd_inputs_ctx)); 425 PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inputs)); 426 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inputs)); 427 } 428 429 { // Create operator for data-driven output handling 430 CeedQFunction qf_sgs_dd_outputs; 431 CeedOperator op_sgs_dd_outputs; 432 433 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Outputs, ComputeSgsDDNodal_Sequential_Outputs_loc, 434 &qf_sgs_dd_outputs)); 435 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_outputs, sgs_dd_setup_data->sgsdd_qfctx)); 436 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE)); 437 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 438 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "inverse multiplicity", 1, CEED_EVAL_NONE)); 439 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE)); 440 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_outputs, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE)); 441 442 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_outputs, NULL, NULL, &op_sgs_dd_outputs)); 443 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 444 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 445 sgs_dd_setup_data->grid_aniso_ceed)); 446 PetscCallCeed(ceed, 447 CeedOperatorSetField(op_sgs_dd_outputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 448 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec)); 449 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 450 451 PetscCall(OperatorApplyContextCreate(sgs_dd_data->dm_dd_outputs, sgs_dd_data->dm_sgs, ceed, op_sgs_dd_outputs, NULL, sgs_dd_data->sgs_nodal_ceed, 452 NULL, NULL, &sgs_dd_data->op_nodal_dd_outputs_ctx)); 453 PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_outputs)); 454 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_outputs)); 455 } 456 457 sgs_dd_data->sgs_nodal_eval = SgsDDNodalStressEval_Sequential; 458 459 if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_CEED) { 460 sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Ceed; 461 PetscCall(SgsDDSetupNodalEvaluation_Sequential_Ceed(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs, 462 elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx)); 463 } else if (sgs_dd_setup_data->sgs_dd_model_implementation == SGS_MODEL_DD_SEQENTIAL_TORCH) { 464 sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Torch; 465 PetscCall(SgsDDSetupNodalEvaluation_Sequential_Torch(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs, 466 elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx)); 467 } 468 469 sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs; 470 471 PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity)); 472 PetscCallCeed(ceed, CeedVectorDestroy(&eigvec)); 473 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity)); 474 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_eigvec)); 475 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_inputs)); 476 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_outputs)); 477 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_grad_velo)); 478 PetscFunctionReturn(PETSC_SUCCESS); 479 } 480 481 // @brief Create CeedOperator to compute SGS contribution to the residual 482 static PetscErrorCode SgsSetupNodalIFunction(Ceed ceed, Honee honee, SgsDDSetupData sgs_dd_setup_data) { 483 SgsDDData sgs_dd_data; 484 CeedInt num_comp_q, q_data_size, num_comp_x; 485 PetscInt dim; 486 CeedQFunction qf_sgs_apply; 487 CeedOperator op_sgs_apply; 488 CeedBasis basis_sgs; 489 CeedVector q_data; 490 CeedElemRestriction elem_restr_qd; 491 492 PetscFunctionBeginUser; 493 PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 494 PetscCall(DMGetDimension(honee->dm, &dim)); 495 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_q, &num_comp_q)); 496 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(honee->elem_restr_x, &num_comp_x)); 497 498 { 499 DMLabel domain_label = NULL; 500 PetscInt label_value = 0, height = 0, dm_field = 0; 501 502 PetscCall(DMPlexCeedBasisCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &basis_sgs)); 503 PetscCall(QDataGet(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, 504 &q_data, &q_data_size)); 505 } 506 507 switch (honee->phys->state_var) { 508 case STATEVAR_PRIMITIVE: 509 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Prim, IFunction_NodalSgs_Prim_loc, &qf_sgs_apply)); 510 break; 511 case STATEVAR_CONSERVATIVE: 512 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Conserv, IFunction_NodalSgs_Conserv_loc, &qf_sgs_apply)); 513 break; 514 case STATEVAR_ENTROPY: 515 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Entropy, IFunction_NodalSgs_Entropy_loc, &qf_sgs_apply)); 516 break; 517 } 518 519 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_apply, sgs_dd_setup_data->ifunction_qfctx)); 520 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "q", num_comp_q, CEED_EVAL_INTERP)); 521 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "qdata", q_data_size, CEED_EVAL_NONE)); 522 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_INTERP)); 523 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_apply, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 524 525 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_apply, NULL, NULL, &op_sgs_apply)); 526 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 527 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data)); 528 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "km_sgs", sgs_dd_setup_data->elem_restr_sgs, basis_sgs, sgs_dd_data->sgs_nodal_ceed)); 529 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE)); 530 531 PetscCall(OperatorApplyContextCreate(honee->dm, honee->dm, ceed, op_sgs_apply, honee->q_ceed, honee->g_ceed, NULL, NULL, 532 &sgs_dd_data->op_sgs_apply_ctx)); 533 534 PetscCallCeed(ceed, CeedVectorDestroy(&q_data)); 535 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd)); 536 PetscCallCeed(ceed, CeedBasisDestroy(&basis_sgs)); 537 PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_apply)); 538 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_apply)); 539 PetscFunctionReturn(PETSC_SUCCESS); 540 } 541 542 // @brief Calculate and add data-driven SGS residual to the global residual 543 PetscErrorCode SgsDDApplyIFunction(Honee honee, const Vec Q_loc, Vec G_loc) { 544 SgsDDData sgs_dd_data; 545 Vec VelocityGradient, SGSNodal_loc; 546 PetscMemType sgs_nodal_mem_type; 547 NodalProjectionData grad_velo_proj; 548 549 PetscFunctionBeginUser; 550 PetscCall(PetscLogEventBegin(HONEE_SgsModel, Q_loc, G_loc, NULL, NULL)); 551 PetscCall(HoneeGetContainer(honee, SGS_DD_DATA_KEY, &sgs_dd_data)); 552 PetscCall(HoneeGetContainer(honee, GRAD_VELO_PROJ_KEY, &grad_velo_proj)); 553 PetscCall(DMGetGlobalVector(grad_velo_proj->dm, &VelocityGradient)); 554 PetscCall(VelocityGradientProjectionApply(grad_velo_proj, Q_loc, VelocityGradient)); 555 556 // -- Compute Nodal SGS tensor 557 PetscCall(DMGetLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc)); 558 PetscCall(sgs_dd_data->sgs_nodal_eval(honee, Q_loc, VelocityGradient, SGSNodal_loc)); 559 560 // -- Compute contribution of the SGS stress 561 PetscCall(VecPetscToCeed(SGSNodal_loc, &sgs_nodal_mem_type, sgs_dd_data->sgs_nodal_ceed)); // sgs_nodal_ceed is an implicit input 562 PetscCall(ApplyAddCeedOperatorLocalToLocal(Q_loc, G_loc, sgs_dd_data->op_sgs_apply_ctx)); 563 564 // -- Return local SGS vector 565 PetscCall(VecCeedToPetsc(sgs_dd_data->sgs_nodal_ceed, sgs_nodal_mem_type, SGSNodal_loc)); 566 PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc)); 567 PetscCall(DMRestoreGlobalVector(grad_velo_proj->dm, &VelocityGradient)); 568 PetscCall(PetscLogEventEnd(HONEE_SgsModel, Q_loc, G_loc, NULL, NULL)); 569 PetscFunctionReturn(PETSC_SUCCESS); 570 } 571 572 // @brief B = A^T, A is NxM, B is MxN 573 static PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) { 574 PetscFunctionBeginUser; 575 for (PetscInt i = 0; i < N; i++) { 576 for (PetscInt j = 0; j < M; j++) { 577 B[j * N + i] = A[i * M + j]; 578 } 579 } 580 PetscFunctionReturn(PETSC_SUCCESS); 581 } 582 583 // @brief Read neural network coefficients from file and put into context struct 584 static PetscErrorCode SgsDDContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SgsDDContext *psgsdd_ctx) { 585 SgsDDContext sgsdd_ctx; 586 PetscInt num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons; 587 char file_path[PETSC_MAX_PATH_LEN]; 588 PetscScalar *temp; 589 590 PetscFunctionBeginUser; 591 { 592 SgsDDContext sgsdd_temp; 593 PetscCall(PetscNew(&sgsdd_temp)); 594 *sgsdd_temp = **psgsdd_ctx; 595 sgsdd_temp->offsets.bias1 = 0; 596 sgsdd_temp->offsets.bias2 = sgsdd_temp->offsets.bias1 + num_neurons; 597 sgsdd_temp->offsets.weight1 = sgsdd_temp->offsets.bias2 + num_neurons; 598 sgsdd_temp->offsets.weight2 = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs; 599 sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons; 600 PetscInt total_num_scalars = sgsdd_temp->offsets.out_scaling + 2 * num_outputs; 601 sgsdd_temp->total_bytes = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]); 602 PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx)); 603 *sgsdd_ctx = *sgsdd_temp; 604 PetscCall(PetscFree(sgsdd_temp)); 605 } 606 607 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat")); 608 PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1])); 609 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat")); 610 PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2])); 611 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat")); 612 PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling])); 613 614 { 615 PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp)); 616 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat")); 617 PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp)); 618 PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons)); 619 PetscCall(PetscFree(temp)); 620 } 621 { 622 PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp)); 623 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat")); 624 PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp)); 625 PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs)); 626 PetscCall(PetscFree(temp)); 627 } 628 629 PetscCall(PetscFree(*psgsdd_ctx)); 630 *psgsdd_ctx = sgsdd_ctx; 631 PetscFunctionReturn(PETSC_SUCCESS); 632 } 633 634 PetscErrorCode SgsDDSetup(Ceed ceed, Honee honee, ProblemData problem) { 635 PetscReal alpha = 0; 636 SgsDDContext sgsdd_ctx; 637 MPI_Comm comm = honee->comm; 638 char sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters"; 639 SgsDDSetupData sgs_dd_setup_data; 640 NewtonianIdealGasContext newt_ctx; 641 NodalProjectionData grad_velo_proj; 642 SgsDDData sgs_dd_data; 643 644 PetscFunctionBeginUser; 645 PetscCall(VelocityGradientProjectionSetup(ceed, honee, problem, honee->phys->state_var, honee->elem_restr_q, honee->basis_q, &grad_velo_proj)); 646 PetscCall(HoneeSetContainer(honee, GRAD_VELO_PROJ_KEY, grad_velo_proj, (PetscCtxDestroyFn *)NodalProjectionDataDestroy)); 647 648 PetscCall(PetscNew(&sgs_dd_data)); 649 sgs_dd_data->num_comp_inputs = 6; 650 sgs_dd_data->num_comp_outputs = 6; 651 PetscCall(HoneeSetContainer(honee, SGS_DD_DATA_KEY, sgs_dd_data, (PetscCtxDestroyFn *)SgsDDDataDestroy)); 652 653 PetscCall(PetscNew(&sgs_dd_setup_data)); 654 655 PetscOptionsBegin(comm, NULL, "SGS Data-Driven Model Options", NULL); 656 PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL)); 657 PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir, 658 sgs_dd_dir, sizeof(sgs_dd_dir), NULL)); 659 PetscCall(PetscOptionsDeprecated("-sgs_model_dd_use_fused", NULL, "libCEED 0.12.0", "Use -sgs_model_dd_type instead")); 660 sgs_dd_setup_data->sgs_dd_model_implementation = SGS_MODEL_DD_FUSED; 661 PetscCall(PetscOptionsEnum("-sgs_model_dd_implementation", "Data-Driven SGS model implementation", NULL, SGSModelDDImplementations, 662 (PetscEnum)sgs_dd_setup_data->sgs_dd_model_implementation, (PetscEnum *)&sgs_dd_setup_data->sgs_dd_model_implementation, 663 NULL)); 664 PetscOptionsEnd(); 665 666 PetscCall(PetscNew(&sgsdd_ctx)); 667 *sgsdd_ctx = (struct SgsDDContext_){ 668 .num_layers = 1, 669 .num_inputs = 6, 670 .num_outputs = 6, 671 .num_neurons = 20, 672 .alpha = alpha, 673 }; 674 675 PetscCall(SgsDDContextFill(comm, sgs_dd_dir, &sgsdd_ctx)); 676 677 // -- Create DM for storing SGS tensor at nodes 678 PetscCall(SgsDDCreateDM(honee->dm, &sgs_dd_data->dm_sgs, honee->app_ctx->degree, honee->app_ctx->q_extra, &sgs_dd_data->num_comp_sgs)); 679 680 PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &newt_ctx)); 681 sgsdd_ctx->newt_ctx = *newt_ctx; 682 PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &newt_ctx)); 683 PetscCallCeed(ceed, CeedQFunctionContextCreate(honee->ceed, &sgs_dd_setup_data->sgsdd_qfctx)); 684 PetscCallCeed(ceed, 685 CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx)); 686 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 687 688 PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(problem->apply_vol_ifunction.qfctx, &sgs_dd_setup_data->ifunction_qfctx)); 689 690 // -- Compute and store anisotropy tensor 691 PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, honee, &sgs_dd_setup_data->elem_restr_grid_aniso, &sgs_dd_setup_data->grid_aniso_ceed)); 692 693 // -- Create Nodal Evaluation Operator 694 switch (sgs_dd_setup_data->sgs_dd_model_implementation) { 695 case SGS_MODEL_DD_FUSED: 696 PetscCall(SgsDDSetupNodalEvaluation_Fused(ceed, honee, sgs_dd_setup_data)); 697 break; 698 case SGS_MODEL_DD_SEQENTIAL_CEED: 699 case SGS_MODEL_DD_SEQENTIAL_TORCH: 700 PetscCall(SgsDDSetupNodalEvaluation_Sequential(ceed, honee, sgs_dd_setup_data)); 701 break; 702 } 703 704 // -- Create Operator to evalutate residual of SGS stress 705 PetscCall(SgsSetupNodalIFunction(ceed, honee, sgs_dd_setup_data)); 706 707 PetscCall(SgsDDSetupDataDestroy(sgs_dd_setup_data)); 708 PetscFunctionReturn(PETSC_SUCCESS); 709 } 710