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