1 // Copyright (c) 2017-2023, Lawrence Livermore National Security, LLC and other CEED contributors. 2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3 // 4 // SPDX-License-Identifier: BSD-2-Clause 5 // 6 // This file is part of CEED: http://github.com/ceed 7 8 #include "../qfunctions/sgs_dd_model.h" 9 10 #include <petscdmplex.h> 11 12 #include "../navierstokes.h" 13 14 typedef struct { 15 CeedElemRestriction elem_restr_grid_aniso, elem_restr_sgs; 16 CeedVector grid_aniso_ceed; 17 CeedQFunctionContext sgsdd_qfctx; 18 } *SgsDDModelSetupData; 19 20 PetscErrorCode SgsDDModelSetupDataDestroy(SgsDDModelSetupData sgs_dd_setup_data) { 21 Ceed ceed; 22 PetscFunctionBeginUser; 23 PetscCall(CeedElemRestrictionGetCeed(sgs_dd_setup_data->elem_restr_sgs, &ceed)); 24 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_grid_aniso)); 25 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_sgs)); 26 PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_setup_data->grid_aniso_ceed)); 27 PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->sgsdd_qfctx)); 28 29 PetscCall(PetscFree(sgs_dd_setup_data)); 30 PetscFunctionReturn(PETSC_SUCCESS); 31 } 32 33 // @brief Create DM for storing subgrid stress at nodes 34 PetscErrorCode SgsDDModelCreateDM(DM dm_source, DM *dm_sgs, PetscInt degree, PetscInt q_extra, PetscInt *num_components) { 35 PetscSection section; 36 37 PetscFunctionBeginUser; 38 *num_components = 6; 39 40 PetscCall(DMClone(dm_source, dm_sgs)); 41 PetscCall(PetscObjectSetName((PetscObject)*dm_sgs, "Subgrid Stress Projection")); 42 43 PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_sgs)); 44 45 PetscCall(DMGetLocalSection(*dm_sgs, §ion)); 46 PetscCall(PetscSectionSetFieldName(section, 0, "")); 47 PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMSubgridStressXX")); 48 PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMSubgridStressYY")); 49 PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMSubgridStressZZ")); 50 PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMSubgridStressYZ")); 51 PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMSubgridStressXZ")); 52 PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMSubgridStressXY")); 53 PetscFunctionReturn(PETSC_SUCCESS); 54 }; 55 56 // @brief Create CeedOperator to calculate data-drive SGS at nodes 57 PetscErrorCode SgsDDModelSetupNodalEvaluation(Ceed ceed, User user, CeedData ceed_data, SgsDDModelSetupData sgs_dd_setup_data) { 58 SgsDDData sgs_dd_data = user->sgs_dd_data; 59 CeedQFunction qf_multiplicity, qf_sgs_dd_nodal; 60 CeedOperator op_multiplicity, op_sgs_dd_nodal; 61 CeedInt num_elem, elem_size, num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso; 62 PetscInt dim; 63 CeedVector multiplicity, inv_multiplicity; 64 CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs; 65 66 PetscFunctionBeginUser; 67 PetscCall(DMGetDimension(user->dm, &dim)); 68 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x)); 69 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q)); 70 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso)); 71 PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(ceed_data->elem_restr_q, &num_elem)); 72 PetscCallCeed(ceed, CeedElemRestrictionGetElementSize(ceed_data->elem_restr_q, &elem_size)); 73 74 { // Get velocity gradient information 75 CeedOperatorField op_field; 76 PetscCallCeed(ceed, CeedOperatorGetFieldByName(user->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field)); 77 PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo)); 78 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo)); 79 } 80 PetscCall(GetRestrictionForDomain(ceed, sgs_dd_data->dm_sgs, 0, 0, 0, 0, -1, 0, &elem_restr_sgs, NULL, NULL)); 81 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL)); 82 83 // -- Create inverse multiplicity for correcting nodal assembly 84 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &multiplicity, NULL)); 85 PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, multiplicity)); 86 PetscCallCeed( 87 ceed, CeedElemRestrictionCreateStrided(ceed, num_elem, elem_size, 1, num_elem * elem_size, CEED_STRIDES_BACKEND, &elem_restr_inv_multiplicity)); 88 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_inv_multiplicity, &inv_multiplicity, NULL)); 89 90 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, InverseMultiplicity, InverseMultiplicity_loc, &qf_multiplicity)); 91 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_multiplicity, "multiplicity", num_comp_q, CEED_EVAL_NONE)); 92 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_multiplicity, "inverse multiplicity", 1, CEED_EVAL_NONE)); 93 94 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_multiplicity, NULL, NULL, &op_multiplicity)); 95 PetscCallCeed(ceed, CeedOperatorSetName(op_multiplicity, "SGS DD Model - Create Multiplicity Scaling")); 96 PetscCallCeed(ceed, CeedOperatorSetField(op_multiplicity, "multiplicity", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE)); 97 PetscCallCeed( 98 ceed, CeedOperatorSetField(op_multiplicity, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE)); 99 100 PetscCallCeed(ceed, CeedOperatorApply(op_multiplicity, multiplicity, inv_multiplicity, CEED_REQUEST_IMMEDIATE)); 101 102 // -- Create operator for SGS DD model nodal evaluation 103 switch (user->phys->state_var) { 104 case STATEVAR_PRIMITIVE: 105 PetscCallCeed(ceed, 106 CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDAnisotropicNodal_Prim, ComputeSgsDDAnisotropicNodal_Prim_loc, &qf_sgs_dd_nodal)); 107 break; 108 case STATEVAR_CONSERVATIVE: 109 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDAnisotropicNodal_Conserv, ComputeSgsDDAnisotropicNodal_Conserv_loc, 110 &qf_sgs_dd_nodal)); 111 break; 112 default: 113 SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, 114 "Anisotropic data-driven SGS nodal evaluation not available for chosen state variable"); 115 } 116 117 // Mesh/geometry order and solution basis order may differ, therefore must interpolate 118 CeedBasis basis_x_to_q; 119 PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_x_to_q)); 120 121 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx)); 122 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE)); 123 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP)); 124 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE)); 125 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 126 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE)); 127 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE)); 128 129 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal)); 130 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "q", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, user->q_ceed)); 131 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "x", ceed_data->elem_restr_x, basis_x_to_q, ceed_data->x_coord)); 132 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE)); 133 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_COLLOCATED, 134 sgs_dd_setup_data->grid_aniso_ceed)); 135 PetscCallCeed(ceed, 136 CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_COLLOCATED, inv_multiplicity)); 137 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE)); 138 139 PetscCall(OperatorApplyContextCreate(user->grad_velo_proj->dm, sgs_dd_data->dm_sgs, ceed, op_sgs_dd_nodal, NULL, sgs_dd_data->sgs_nodal_ceed, NULL, 140 NULL, &sgs_dd_data->op_nodal_evaluation_ctx)); 141 142 sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs; 143 144 PetscCallCeed(ceed, CeedVectorDestroy(&multiplicity)); 145 PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity)); 146 PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q)); 147 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity)); 148 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_multiplicity)); 149 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_nodal)); 150 PetscCallCeed(ceed, CeedOperatorDestroy(&op_multiplicity)); 151 PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_nodal)); 152 PetscFunctionReturn(PETSC_SUCCESS); 153 } 154 155 // @brief Create CeedOperator to compute SGS contribution to the residual 156 PetscErrorCode SgsModelSetupNodalIFunction(Ceed ceed, User user, CeedData ceed_data, SgsDDModelSetupData sgs_dd_setup_data) { 157 SgsDDData sgs_dd_data = user->sgs_dd_data; 158 CeedInt num_comp_q, num_comp_qd, num_comp_x; 159 PetscInt dim; 160 CeedQFunction qf_sgs_apply; 161 CeedOperator op_sgs_apply; 162 CeedBasis basis_sgs; 163 164 PetscFunctionBeginUser; 165 PetscCall(DMGetDimension(user->dm, &dim)); 166 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q)); 167 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &num_comp_qd)); 168 PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x)); 169 170 PetscCall(CreateBasisFromPlex(ceed, sgs_dd_data->dm_sgs, 0, 0, 0, 0, &basis_sgs)); 171 172 switch (user->phys->state_var) { 173 case STATEVAR_PRIMITIVE: 174 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Prim, IFunction_NodalSgs_Prim_loc, &qf_sgs_apply)); 175 break; 176 case STATEVAR_CONSERVATIVE: 177 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Conserv, IFunction_NodalSgs_Conserv_loc, &qf_sgs_apply)); 178 break; 179 default: 180 SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "Nodal SGS evaluation not available for chosen state variable"); 181 } 182 183 PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_apply, sgs_dd_setup_data->sgsdd_qfctx)); 184 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "q", num_comp_q, CEED_EVAL_INTERP)); 185 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "qdata", num_comp_qd, CEED_EVAL_NONE)); 186 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "x", num_comp_x, CEED_EVAL_INTERP)); 187 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_INTERP)); 188 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_apply, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 189 190 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_apply, NULL, NULL, &op_sgs_apply)); 191 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 192 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data)); 193 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord)); 194 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "km_sgs", sgs_dd_setup_data->elem_restr_sgs, basis_sgs, sgs_dd_data->sgs_nodal_ceed)); 195 PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 196 197 PetscCall( 198 OperatorApplyContextCreate(user->dm, user->dm, ceed, op_sgs_apply, user->q_ceed, user->g_ceed, NULL, NULL, &sgs_dd_data->op_sgs_apply_ctx)); 199 200 PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_apply)); 201 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_apply)); 202 PetscFunctionReturn(PETSC_SUCCESS); 203 } 204 205 // @brief Calculate and add data-driven SGS residual to the global residual 206 PetscErrorCode SgsDDModelApplyIFunction(User user, const Vec Q_loc, Vec G_loc) { 207 SgsDDData sgs_dd_data = user->sgs_dd_data; 208 Vec VelocityGradient, SGSNodal_loc; 209 PetscMemType sgs_nodal_mem_type, q_mem_type; 210 211 PetscFunctionBeginUser; 212 PetscCall(DMGetGlobalVector(user->grad_velo_proj->dm, &VelocityGradient)); 213 PetscCall(VelocityGradientProjectionApply(user, Q_loc, VelocityGradient)); 214 215 // -- Compute Nodal SGS tensor 216 PetscCall(DMGetLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc)); 217 PetscCall(VecP2C(Q_loc, &q_mem_type, user->q_ceed)); // q_ceed is an implicit input 218 219 PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, SGSNodal_loc, sgs_dd_data->op_nodal_evaluation_ctx)); 220 221 PetscCall(VecC2P(user->q_ceed, q_mem_type, Q_loc)); 222 PetscCall(VecP2C(SGSNodal_loc, &sgs_nodal_mem_type, sgs_dd_data->sgs_nodal_ceed)); // sgs_nodal_ceed is an implicit input 223 224 // -- Compute contribution of the SGS stress 225 PetscCall(ApplyAddCeedOperatorLocalToLocal(Q_loc, G_loc, sgs_dd_data->op_sgs_apply_ctx)); 226 227 // -- Return local SGS vector 228 PetscCall(VecC2P(sgs_dd_data->sgs_nodal_ceed, sgs_nodal_mem_type, SGSNodal_loc)); 229 PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc)); 230 PetscCall(DMRestoreGlobalVector(user->grad_velo_proj->dm, &VelocityGradient)); 231 PetscFunctionReturn(PETSC_SUCCESS); 232 } 233 234 // @brief B = A^T, A is NxM, B is MxN 235 PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) { 236 PetscFunctionBeginUser; 237 for (PetscInt i = 0; i < N; i++) { 238 for (PetscInt j = 0; j < M; j++) { 239 B[j * N + i] = A[i * M + j]; 240 } 241 } 242 PetscFunctionReturn(PETSC_SUCCESS); 243 } 244 245 // @brief Read neural network coefficients from file and put into context struct 246 PetscErrorCode SgsDDModelContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SgsDDModelContext *psgsdd_ctx) { 247 SgsDDModelContext sgsdd_ctx; 248 PetscInt num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons; 249 char file_path[PETSC_MAX_PATH_LEN]; 250 PetscScalar *temp; 251 252 PetscFunctionBeginUser; 253 { 254 SgsDDModelContext sgsdd_temp; 255 PetscCall(PetscNew(&sgsdd_temp)); 256 *sgsdd_temp = **psgsdd_ctx; 257 sgsdd_temp->offsets.bias1 = 0; 258 sgsdd_temp->offsets.bias2 = sgsdd_temp->offsets.bias1 + num_neurons; 259 sgsdd_temp->offsets.weight1 = sgsdd_temp->offsets.bias2 + num_neurons; 260 sgsdd_temp->offsets.weight2 = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs; 261 sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons; 262 PetscInt total_num_scalars = sgsdd_temp->offsets.out_scaling + 2 * num_outputs; 263 sgsdd_temp->total_bytes = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]); 264 PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx)); 265 *sgsdd_ctx = *sgsdd_temp; 266 PetscCall(PetscFree(sgsdd_temp)); 267 } 268 269 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat")); 270 PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1])); 271 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat")); 272 PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2])); 273 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat")); 274 PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling])); 275 276 { 277 PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp)); 278 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat")); 279 PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp)); 280 PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons)); 281 PetscCall(PetscFree(temp)); 282 } 283 { 284 PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp)); 285 PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat")); 286 PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp)); 287 PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs)); 288 PetscCall(PetscFree(temp)); 289 } 290 291 PetscCall(PetscFree(*psgsdd_ctx)); 292 *psgsdd_ctx = sgsdd_ctx; 293 PetscFunctionReturn(PETSC_SUCCESS); 294 } 295 296 PetscErrorCode SgsDDModelSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) { 297 PetscReal alpha = 0; 298 SgsDDModelContext sgsdd_ctx; 299 MPI_Comm comm = user->comm; 300 char sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters"; 301 SgsDDModelSetupData sgs_dd_setup_data; 302 NewtonianIdealGasContext gas; 303 304 PetscFunctionBeginUser; 305 PetscCall(VelocityGradientProjectionSetup(ceed, user, ceed_data, problem)); 306 307 PetscCall(PetscNew(&sgsdd_ctx)); 308 309 PetscOptionsBegin(comm, NULL, "SGS Data-Driven Model Options", NULL); 310 PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL)); 311 PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir, 312 sgs_dd_dir, sizeof(sgs_dd_dir), NULL)); 313 PetscOptionsEnd(); 314 315 sgsdd_ctx->num_layers = 1; 316 sgsdd_ctx->num_inputs = 6; 317 sgsdd_ctx->num_outputs = 6; 318 sgsdd_ctx->num_neurons = 20; 319 sgsdd_ctx->alpha = alpha; 320 321 PetscCall(SgsDDModelContextFill(comm, sgs_dd_dir, &sgsdd_ctx)); 322 323 // -- Create DM for storing SGS tensor at nodes 324 PetscCall(PetscNew(&user->sgs_dd_data)); 325 PetscCall( 326 SgsDDModelCreateDM(user->dm, &user->sgs_dd_data->dm_sgs, user->app_ctx->degree, user->app_ctx->q_extra, &user->sgs_dd_data->num_comp_sgs)); 327 328 PetscCall(PetscNew(&sgs_dd_setup_data)); 329 330 PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas)); 331 sgsdd_ctx->gas = *gas; 332 PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas)); 333 PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &sgs_dd_setup_data->sgsdd_qfctx)); 334 PetscCallCeed(ceed, 335 CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx)); 336 PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 337 338 // -- Compute and store anisotropy tensor 339 PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, user, ceed_data, &sgs_dd_setup_data->elem_restr_grid_aniso, 340 &sgs_dd_setup_data->grid_aniso_ceed)); 341 342 // -- Create Nodal Evaluation Operator 343 PetscCall(SgsDDModelSetupNodalEvaluation(ceed, user, ceed_data, sgs_dd_setup_data)); 344 345 // -- Create Operator to evalutate residual of SGS stress 346 PetscCall(SgsModelSetupNodalIFunction(ceed, user, ceed_data, sgs_dd_setup_data)); 347 348 PetscCall(SgsDDModelSetupDataDestroy(sgs_dd_setup_data)); 349 PetscFunctionReturn(PETSC_SUCCESS); 350 } 351 352 PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data) { 353 PetscFunctionBeginUser; 354 if (!sgs_dd_data) PetscFunctionReturn(PETSC_SUCCESS); 355 Ceed ceed = sgs_dd_data->op_sgs_apply_ctx->ceed; 356 357 PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->sgs_nodal_ceed)); 358 PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_evaluation_ctx)); 359 PetscCall(DMDestroy(&sgs_dd_data->dm_sgs)); 360 PetscCall(PetscFree(sgs_dd_data)); 361 PetscFunctionReturn(PETSC_SUCCESS); 362 } 363