162b7942eSJames Wright // Copyright (c) 2017-2023, Lawrence Livermore National Security, LLC and other CEED contributors. 262b7942eSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 362b7942eSJames Wright // 462b7942eSJames Wright // SPDX-License-Identifier: BSD-2-Clause 562b7942eSJames Wright // 662b7942eSJames Wright // This file is part of CEED: http://github.com/ceed 762b7942eSJames Wright 862b7942eSJames Wright #include "../qfunctions/sgs_dd_model.h" 962b7942eSJames Wright 1062b7942eSJames Wright #include <petscdmplex.h> 1162b7942eSJames Wright 1262b7942eSJames Wright #include "../navierstokes.h" 1362b7942eSJames Wright 14c38c977aSJames Wright typedef struct { 15c38c977aSJames Wright CeedElemRestriction elem_restr_grid_aniso, elem_restr_sgs; 16c38c977aSJames Wright CeedVector grid_aniso_ceed; 1740816385SJames Wright CeedQFunctionContext sgsdd_qfctx, ifunction_qfctx; 18ad494f68SJames Wright } *SgsDDSetupData; 19c38c977aSJames Wright 20ad494f68SJames Wright PetscErrorCode SgsDDSetupDataDestroy(SgsDDSetupData sgs_dd_setup_data) { 21b4c37c5cSJames Wright Ceed ceed; 22ad494f68SJames Wright 23c38c977aSJames Wright PetscFunctionBeginUser; 24b4c37c5cSJames Wright PetscCall(CeedElemRestrictionGetCeed(sgs_dd_setup_data->elem_restr_sgs, &ceed)); 25ad494f68SJames Wright 26b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_grid_aniso)); 27b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&sgs_dd_setup_data->elem_restr_sgs)); 28b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_setup_data->grid_aniso_ceed)); 29b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->sgsdd_qfctx)); 3040816385SJames Wright PetscCallCeed(ceed, CeedQFunctionContextDestroy(&sgs_dd_setup_data->ifunction_qfctx)); 31c38c977aSJames Wright PetscCall(PetscFree(sgs_dd_setup_data)); 32d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 33c38c977aSJames Wright } 34c38c977aSJames Wright 35ee1455b7SJames Wright // @brief Create DM for storing subgrid stress at nodes 36ad494f68SJames Wright static PetscErrorCode SgsDDCreateDM(DM dm_source, DM *dm_sgs, PetscInt degree, PetscInt q_extra, PetscInt *num_components) { 37ee1455b7SJames Wright PetscSection section; 38ee1455b7SJames Wright 39ee1455b7SJames Wright PetscFunctionBeginUser; 40ee1455b7SJames Wright *num_components = 6; 41ee1455b7SJames Wright 42ee1455b7SJames Wright PetscCall(DMClone(dm_source, dm_sgs)); 43ee1455b7SJames Wright PetscCall(PetscObjectSetName((PetscObject)*dm_sgs, "Subgrid Stress Projection")); 44ee1455b7SJames Wright 45da4ca0cfSJames Wright PetscCall(DMSetupByOrder_FEM(PETSC_TRUE, PETSC_TRUE, degree, 1, q_extra, 1, num_components, *dm_sgs)); 46ee1455b7SJames Wright 47ee1455b7SJames Wright PetscCall(DMGetLocalSection(*dm_sgs, §ion)); 48ee1455b7SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 49ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 0, "KMSubgridStressXX")); 50ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 1, "KMSubgridStressYY")); 51ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 2, "KMSubgridStressZZ")); 52ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 3, "KMSubgridStressYZ")); 53ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 4, "KMSubgridStressXZ")); 54ee1455b7SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, 5, "KMSubgridStressXY")); 55d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 56ee1455b7SJames Wright }; 57ee1455b7SJames Wright 58ad494f68SJames Wright // @brief Evaluate data-driven SGS using fused method 59ad494f68SJames Wright static PetscErrorCode SgsDDNodalStressEval_Fused(User user, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) { 60cceb3143SJames Wright SgsDDData sgs_dd_data = user->sgs_dd_data; 61cceb3143SJames Wright PetscMemType q_mem_type; 62cceb3143SJames Wright 63cceb3143SJames Wright PetscFunctionBeginUser; 64cceb3143SJames Wright PetscCall(VecP2C(Q_loc, &q_mem_type, user->q_ceed)); // q_ceed is an implicit input 65cceb3143SJames Wright 66cceb3143SJames Wright PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, SGSNodal_loc, sgs_dd_data->op_nodal_evaluation_ctx)); 67cceb3143SJames Wright 68cceb3143SJames Wright PetscCall(VecC2P(user->q_ceed, q_mem_type, Q_loc)); 69cceb3143SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 70cceb3143SJames Wright } 71cceb3143SJames Wright 72*b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using fused operator 73ad494f68SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Fused(Ceed ceed, User user, CeedData ceed_data, SgsDDSetupData sgs_dd_setup_data) { 7442454adaSJames Wright SgsDDData sgs_dd_data = user->sgs_dd_data; 75ee1455b7SJames Wright CeedQFunction qf_multiplicity, qf_sgs_dd_nodal; 76ee1455b7SJames Wright CeedOperator op_multiplicity, op_sgs_dd_nodal; 7715c18037SJames Wright CeedInt num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso; 78defe8520SJames Wright PetscInt dim; 794b0f6111SJames Wright CeedVector multiplicity, inv_multiplicity; 80ee1455b7SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs; 8115c18037SJames Wright DMLabel domain_label = NULL; 8215c18037SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 83ee1455b7SJames Wright 84ee1455b7SJames Wright PetscFunctionBeginUser; 85ee1455b7SJames Wright PetscCall(DMGetDimension(user->dm, &dim)); 86b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x)); 87b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q)); 88b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso)); 89ee1455b7SJames Wright 90ee1455b7SJames Wright { // Get velocity gradient information 91ee1455b7SJames Wright CeedOperatorField op_field; 92b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(user->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field)); 93b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo)); 94b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo)); 95ee1455b7SJames Wright } 9615c18037SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs)); 97b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL)); 98ee1455b7SJames Wright 99ee1455b7SJames Wright // -- Create inverse multiplicity for correcting nodal assembly 100b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &multiplicity, NULL)); 101b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, multiplicity)); 10215c18037SJames Wright PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, 1, &elem_restr_inv_multiplicity)); 103b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_inv_multiplicity, &inv_multiplicity, NULL)); 104ee1455b7SJames Wright 105b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, InverseMultiplicity, InverseMultiplicity_loc, &qf_multiplicity)); 106b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_multiplicity, "multiplicity", num_comp_q, CEED_EVAL_NONE)); 107b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_multiplicity, "inverse multiplicity", 1, CEED_EVAL_NONE)); 108ee1455b7SJames Wright 109b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_multiplicity, NULL, NULL, &op_multiplicity)); 110b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetName(op_multiplicity, "SGS DD Model - Create Multiplicity Scaling")); 11158e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_multiplicity, "multiplicity", ceed_data->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 11258e1cbfdSJeremy L Thompson PetscCallCeed(ceed, 11358e1cbfdSJeremy L Thompson CeedOperatorSetField(op_multiplicity, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 114ee1455b7SJames Wright 115b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorApply(op_multiplicity, multiplicity, inv_multiplicity, CEED_REQUEST_IMMEDIATE)); 116ee1455b7SJames Wright 117ee1455b7SJames Wright // -- Create operator for SGS DD model nodal evaluation 118ee1455b7SJames Wright switch (user->phys->state_var) { 119ee1455b7SJames Wright case STATEVAR_PRIMITIVE: 120ad494f68SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Prim, ComputeSgsDDNodal_Prim_loc, &qf_sgs_dd_nodal)); 121ee1455b7SJames Wright break; 122ee1455b7SJames Wright case STATEVAR_CONSERVATIVE: 123ad494f68SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Conserv, ComputeSgsDDNodal_Conserv_loc, &qf_sgs_dd_nodal)); 124ee1455b7SJames Wright break; 125ee1455b7SJames Wright default: 126ad494f68SJames Wright SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "Data-driven SGS nodal evaluation not available for chosen state variable"); 127ee1455b7SJames Wright } 128ee1455b7SJames Wright 129ee1455b7SJames Wright // Mesh/geometry order and solution basis order may differ, therefore must interpolate 130ee1455b7SJames Wright CeedBasis basis_x_to_q; 131b4c37c5cSJames Wright PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_x_to_q)); 132ee1455b7SJames Wright 133b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_nodal, sgs_dd_setup_data->sgsdd_qfctx)); 134b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "q", num_comp_q, CEED_EVAL_NONE)); 135b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "x", num_comp_x, CEED_EVAL_INTERP)); 136b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE)); 137b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 138b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_nodal, "inverse multiplicity", 1, CEED_EVAL_NONE)); 139b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_nodal, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE)); 140ee1455b7SJames Wright 141b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_nodal, NULL, NULL, &op_sgs_dd_nodal)); 14258e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "q", ceed_data->elem_restr_q, CEED_BASIS_NONE, user->q_ceed)); 143b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "x", ceed_data->elem_restr_x, basis_x_to_q, ceed_data->x_coord)); 14458e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 14558e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 146b4c37c5cSJames Wright sgs_dd_setup_data->grid_aniso_ceed)); 14758e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 14858e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_nodal, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 149ee1455b7SJames Wright 1504b0f6111SJames Wright 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, 1514b0f6111SJames Wright NULL, &sgs_dd_data->op_nodal_evaluation_ctx)); 152ee1455b7SJames Wright 153ee1455b7SJames Wright sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs; 154ad494f68SJames Wright sgs_dd_data->sgs_nodal_eval = SgsDDNodalStressEval_Fused; 155ee1455b7SJames Wright 156b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&multiplicity)); 157b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity)); 158b4c37c5cSJames Wright PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q)); 159b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity)); 160b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_multiplicity)); 161b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_nodal)); 162b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_multiplicity)); 163b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_nodal)); 164d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 165ee1455b7SJames Wright } 166ee1455b7SJames Wright 167*b87d60b3SJames Wright // @brief Setup data-driven model inference using internal (libCEED native) implementation 168*b87d60b3SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential_Internal(Ceed ceed, SgsDDData sgs_dd_data, SgsDDSetupData sgs_dd_setup_data, 169*b87d60b3SJames Wright CeedElemRestriction elem_restr_dd_inputs, 170*b87d60b3SJames Wright CeedElemRestriction elem_restr_dd_outputs, 171*b87d60b3SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, CeedVector inv_multiplicity, 172*b87d60b3SJames Wright void **ctx) { 173*b87d60b3SJames Wright CeedQFunction qf_sgs_dd_inference; 174*b87d60b3SJames Wright CeedOperator op_sgs_dd_inference; 175*b87d60b3SJames Wright OperatorApplyContext *op_context = (OperatorApplyContext *)ctx; 176*b87d60b3SJames Wright 177*b87d60b3SJames Wright PetscFunctionBeginUser; 178*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inference, ComputeSgsDDNodal_Sequential_Inference_loc, 179*b87d60b3SJames Wright &qf_sgs_dd_inference)); 180*b87d60b3SJames Wright 181*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inference, sgs_dd_setup_data->sgsdd_qfctx)); 182*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE)); 183*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inference, "inverse multiplicity", 1, CEED_EVAL_NONE)); 184*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inference, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE)); 185*b87d60b3SJames Wright 186*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inference, NULL, NULL, &op_sgs_dd_inference)); 187*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 188*b87d60b3SJames Wright PetscCallCeed(ceed, 189*b87d60b3SJames Wright CeedOperatorSetField(op_sgs_dd_inference, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 190*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inference, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 191*b87d60b3SJames Wright 192*b87d60b3SJames Wright PetscCall(OperatorApplyContextCreate(sgs_dd_data->dm_dd_inputs, sgs_dd_data->dm_dd_outputs, ceed, op_sgs_dd_inference, NULL, NULL, NULL, NULL, 193*b87d60b3SJames Wright op_context)); 194*b87d60b3SJames Wright sgs_dd_data->sgs_nodal_inference_ctx_destroy = (PetscErrorCode(*)(void *))OperatorApplyContextDestroy; 195*b87d60b3SJames Wright 196*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inference)); 197*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inference)); 198*b87d60b3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 199*b87d60b3SJames Wright } 200*b87d60b3SJames Wright 201*b87d60b3SJames Wright // @brief Perform data-driven model inference using internal (libCEED native) implementation 202*b87d60b3SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential_Internal(Vec DD_Inputs_loc, Vec DD_Outputs_loc, void *ctx) { 203*b87d60b3SJames Wright OperatorApplyContext op_context = *(OperatorApplyContext *)ctx; 204*b87d60b3SJames Wright 205*b87d60b3SJames Wright PetscFunctionBeginUser; 206*b87d60b3SJames Wright PetscCall(ApplyCeedOperatorLocalToLocal(DD_Inputs_loc, DD_Outputs_loc, op_context)); 207*b87d60b3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 208*b87d60b3SJames Wright } 209*b87d60b3SJames Wright 210*b87d60b3SJames Wright // @brief Evaluate data-driven SGS using sequential method 211*b87d60b3SJames Wright PetscErrorCode SgsDDNodalStressEval_Sequential(User user, Vec Q_loc, Vec VelocityGradient, Vec SGSNodal_loc) { 212*b87d60b3SJames Wright SgsDDData sgs_dd_data = user->sgs_dd_data; 213*b87d60b3SJames Wright PetscMemType q_mem_type; 214*b87d60b3SJames Wright Vec DD_Inputs_loc, DD_Outputs_loc; 215*b87d60b3SJames Wright 216*b87d60b3SJames Wright PetscFunctionBeginUser; 217*b87d60b3SJames Wright PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc)); 218*b87d60b3SJames Wright PetscCall(DMGetLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc)); 219*b87d60b3SJames Wright PetscCall(VecP2C(Q_loc, &q_mem_type, user->q_ceed)); // q_ceed is an implicit input 220*b87d60b3SJames Wright 221*b87d60b3SJames Wright PetscCall(ApplyCeedOperatorGlobalToLocal(VelocityGradient, DD_Inputs_loc, sgs_dd_data->op_nodal_dd_inputs_ctx)); 222*b87d60b3SJames Wright PetscCall(sgs_dd_data->sgs_nodal_inference(DD_Inputs_loc, DD_Outputs_loc, &sgs_dd_data->sgs_nodal_inference_ctx)); 223*b87d60b3SJames Wright PetscCall(ApplyCeedOperatorLocalToLocal(DD_Outputs_loc, SGSNodal_loc, sgs_dd_data->op_nodal_dd_outputs_ctx)); 224*b87d60b3SJames Wright 225*b87d60b3SJames Wright PetscCall(VecC2P(user->q_ceed, q_mem_type, Q_loc)); 226*b87d60b3SJames Wright PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_inputs, &DD_Inputs_loc)); 227*b87d60b3SJames Wright PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_dd_outputs, &DD_Outputs_loc)); 228*b87d60b3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 229*b87d60b3SJames Wright } 230*b87d60b3SJames Wright 231*b87d60b3SJames Wright // @brief Create CeedOperator to calculate data-drive SGS at nodes using sequentially-applied operators 232*b87d60b3SJames Wright static PetscErrorCode SgsDDSetupNodalEvaluation_Sequential(Ceed ceed, User user, CeedData ceed_data, SgsDDSetupData sgs_dd_setup_data) { 233*b87d60b3SJames Wright SgsDDData sgs_dd_data = user->sgs_dd_data; 234*b87d60b3SJames Wright CeedInt num_comp_q, num_comp_grad_velo, num_comp_x, num_comp_grid_aniso, num_comp_eigvec = 9 + 1; 235*b87d60b3SJames Wright PetscInt dim; 236*b87d60b3SJames Wright CeedVector inv_multiplicity, eigvec; 237*b87d60b3SJames Wright CeedElemRestriction elem_restr_inv_multiplicity, elem_restr_grad_velo, elem_restr_sgs, elem_restr_eigvec, elem_restr_dd_inputs, 238*b87d60b3SJames Wright elem_restr_dd_outputs; 239*b87d60b3SJames Wright DMLabel domain_label = NULL; 240*b87d60b3SJames Wright PetscInt label_value = 0, height = 0, dm_field = 0; 241*b87d60b3SJames Wright 242*b87d60b3SJames Wright PetscFunctionBeginUser; 243*b87d60b3SJames Wright { // Create DMs for data-driven input and output values 244*b87d60b3SJames Wright PetscSection section; 245*b87d60b3SJames Wright PetscInt degree, q_extra; 246*b87d60b3SJames Wright { // Get degree and number of quadrature points from dm_sgs 247*b87d60b3SJames Wright PetscFE fe; 248*b87d60b3SJames Wright PetscSpace basis; 249*b87d60b3SJames Wright PetscQuadrature quadrature; 250*b87d60b3SJames Wright PetscInt num_qpnts; 251*b87d60b3SJames Wright PetscCall(DMGetField(sgs_dd_data->dm_sgs, 0, NULL, (PetscObject *)&fe)); 252*b87d60b3SJames Wright PetscCall(PetscFEGetBasisSpace(fe, &basis)); 253*b87d60b3SJames Wright PetscCall(PetscSpaceGetDegree(basis, °ree, NULL)); 254*b87d60b3SJames Wright PetscCall(PetscFEGetQuadrature(fe, &quadrature)); 255*b87d60b3SJames Wright PetscCall(PetscQuadratureGetOrder(quadrature, &num_qpnts)); 256*b87d60b3SJames Wright q_extra = degree - num_qpnts; 257*b87d60b3SJames Wright } 258*b87d60b3SJames Wright 259*b87d60b3SJames Wright PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_inputs)); 260*b87d60b3SJames Wright PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_inputs, "Data-Driven Model Inputs")); 261*b87d60b3SJames Wright 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*b87d60b3SJames Wright PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_inputs, §ion)); 263*b87d60b3SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 264*b87d60b3SJames Wright for (CeedInt i = 0; i < sgs_dd_data->num_comp_inputs; i++) { 265*b87d60b3SJames Wright char component_name[PETSC_MAX_PATH_LEN]; 266*b87d60b3SJames Wright 267*b87d60b3SJames Wright PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenInput%" CeedInt_FMT, i + 1)); 268*b87d60b3SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, i, component_name)); 269*b87d60b3SJames Wright } 270*b87d60b3SJames Wright 271*b87d60b3SJames Wright PetscCall(DMClone(sgs_dd_data->dm_sgs, &sgs_dd_data->dm_dd_outputs)); 272*b87d60b3SJames Wright PetscCall(PetscObjectSetName((PetscObject)sgs_dd_data->dm_dd_outputs, "Data-Driven Model Outputs")); 273*b87d60b3SJames Wright 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*b87d60b3SJames Wright PetscCall(DMGetLocalSection(sgs_dd_data->dm_dd_outputs, §ion)); 275*b87d60b3SJames Wright PetscCall(PetscSectionSetFieldName(section, 0, "")); 276*b87d60b3SJames Wright for (CeedInt i = 0; i < sgs_dd_data->num_comp_outputs; i++) { 277*b87d60b3SJames Wright char component_name[PETSC_MAX_PATH_LEN]; 278*b87d60b3SJames Wright 279*b87d60b3SJames Wright PetscCall(PetscSNPrintf(component_name, sizeof component_name, "DataDrivenOutput%" CeedInt_FMT, i + 1)); 280*b87d60b3SJames Wright PetscCall(PetscSectionSetComponentName(section, 0, i, component_name)); 281*b87d60b3SJames Wright } 282*b87d60b3SJames Wright } 283*b87d60b3SJames Wright 284*b87d60b3SJames Wright PetscCall(DMGetDimension(user->dm, &dim)); 285*b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x)); 286*b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q)); 287*b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(sgs_dd_setup_data->elem_restr_grid_aniso, &num_comp_grid_aniso)); 288*b87d60b3SJames Wright 289*b87d60b3SJames Wright { // Get velocity gradient information 290*b87d60b3SJames Wright CeedOperatorField op_field; 291*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorGetFieldByName(user->grad_velo_proj->l2_rhs_ctx->op, "velocity gradient", &op_field)); 292*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(op_field, &elem_restr_grad_velo)); 293*b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_grad_velo, &num_comp_grad_velo)); 294*b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_grad_velo, &sgs_dd_data->grad_velo_ceed, NULL)); 295*b87d60b3SJames Wright } 296*b87d60b3SJames Wright 297*b87d60b3SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, dm_field, &elem_restr_sgs)); 298*b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_sgs, &sgs_dd_data->sgs_nodal_ceed, NULL)); 299*b87d60b3SJames Wright PetscCall( 300*b87d60b3SJames Wright DMPlexCeedElemRestrictionCollocatedCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, num_comp_eigvec, &elem_restr_eigvec)); 301*b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_eigvec, &eigvec, NULL)); 302*b87d60b3SJames Wright 303*b87d60b3SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_inputs, domain_label, label_value, height, dm_field, &elem_restr_dd_inputs)); 304*b87d60b3SJames Wright PetscCall(DMPlexCeedElemRestrictionCreate(ceed, sgs_dd_data->dm_dd_outputs, domain_label, label_value, height, dm_field, &elem_restr_dd_outputs)); 305*b87d60b3SJames Wright 306*b87d60b3SJames Wright { // Create inverse multiplicity for correcting nodal assembly 307*b87d60b3SJames Wright CeedQFunction qf_multiplicity; 308*b87d60b3SJames Wright CeedOperator op_multiplicity; 309*b87d60b3SJames Wright CeedVector multiplicity; 310*b87d60b3SJames Wright 311*b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &multiplicity, NULL)); 312*b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(ceed_data->elem_restr_q, multiplicity)); 313*b87d60b3SJames Wright PetscCall( 314*b87d60b3SJames Wright DMPlexCeedElemRestrictionCollocatedCreate(ceed, sgs_dd_data->dm_sgs, domain_label, label_value, height, 1, &elem_restr_inv_multiplicity)); 315*b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_inv_multiplicity, &inv_multiplicity, NULL)); 316*b87d60b3SJames Wright 317*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, InverseMultiplicity, InverseMultiplicity_loc, &qf_multiplicity)); 318*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_multiplicity, "multiplicity", num_comp_q, CEED_EVAL_NONE)); 319*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_multiplicity, "inverse multiplicity", 1, CEED_EVAL_NONE)); 320*b87d60b3SJames Wright 321*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_multiplicity, NULL, NULL, &op_multiplicity)); 322*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetName(op_multiplicity, "SGS DD Model - Create Multiplicity Scaling")); 323*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_multiplicity, "multiplicity", ceed_data->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 324*b87d60b3SJames Wright PetscCallCeed(ceed, 325*b87d60b3SJames Wright CeedOperatorSetField(op_multiplicity, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 326*b87d60b3SJames Wright 327*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorApply(op_multiplicity, multiplicity, inv_multiplicity, CEED_REQUEST_IMMEDIATE)); 328*b87d60b3SJames Wright 329*b87d60b3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&multiplicity)); 330*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_multiplicity)); 331*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_multiplicity)); 332*b87d60b3SJames Wright } 333*b87d60b3SJames Wright 334*b87d60b3SJames Wright { // Create operator for data-driven input evaluation 335*b87d60b3SJames Wright CeedQFunction qf_sgs_dd_inputs; 336*b87d60b3SJames Wright CeedOperator op_sgs_dd_inputs; 337*b87d60b3SJames Wright 338*b87d60b3SJames Wright switch (user->phys->state_var) { 339*b87d60b3SJames Wright case STATEVAR_PRIMITIVE: 340*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Prim, 341*b87d60b3SJames Wright ComputeSgsDDNodal_Sequential_Inputs_Prim_loc, &qf_sgs_dd_inputs)); 342*b87d60b3SJames Wright break; 343*b87d60b3SJames Wright case STATEVAR_CONSERVATIVE: 344*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Inputs_Conserv, 345*b87d60b3SJames Wright ComputeSgsDDNodal_Sequential_Inputs_Conserv_loc, &qf_sgs_dd_inputs)); 346*b87d60b3SJames Wright break; 347*b87d60b3SJames Wright default: 348*b87d60b3SJames Wright SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, 349*b87d60b3SJames Wright "Data-driven SGS nodal input evaluation not available for chosen state variable"); 350*b87d60b3SJames Wright } 351*b87d60b3SJames Wright 352*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_inputs, sgs_dd_setup_data->sgsdd_qfctx)); 353*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "q", num_comp_q, CEED_EVAL_NONE)); 354*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "gradient velocity", num_comp_grad_velo, CEED_EVAL_NONE)); 355*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 356*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_inputs, "inverse multiplicity", 1, CEED_EVAL_NONE)); 357*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE)); 358*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_inputs, "model inputs", sgs_dd_data->num_comp_inputs, CEED_EVAL_NONE)); 359*b87d60b3SJames Wright 360*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_inputs, NULL, NULL, &op_sgs_dd_inputs)); 361*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "q", ceed_data->elem_restr_q, CEED_BASIS_NONE, user->q_ceed)); 362*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "gradient velocity", elem_restr_grad_velo, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 363*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 364*b87d60b3SJames Wright sgs_dd_setup_data->grid_aniso_ceed)); 365*b87d60b3SJames Wright PetscCallCeed(ceed, 366*b87d60b3SJames Wright CeedOperatorSetField(op_sgs_dd_inputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 367*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec)); 368*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_inputs, "model inputs", elem_restr_dd_inputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 369*b87d60b3SJames Wright 370*b87d60b3SJames Wright PetscCall(OperatorApplyContextCreate(user->grad_velo_proj->dm, sgs_dd_data->dm_dd_inputs, ceed, op_sgs_dd_inputs, NULL, NULL, NULL, NULL, 371*b87d60b3SJames Wright &sgs_dd_data->op_nodal_dd_inputs_ctx)); 372*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_inputs)); 373*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_inputs)); 374*b87d60b3SJames Wright } 375*b87d60b3SJames Wright 376*b87d60b3SJames Wright { // Create operator for data-driven output handling 377*b87d60b3SJames Wright CeedQFunction qf_sgs_dd_outputs; 378*b87d60b3SJames Wright CeedOperator op_sgs_dd_outputs; 379*b87d60b3SJames Wright 380*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, ComputeSgsDDNodal_Sequential_Outputs, ComputeSgsDDNodal_Sequential_Outputs_loc, 381*b87d60b3SJames Wright &qf_sgs_dd_outputs)); 382*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_dd_outputs, sgs_dd_setup_data->sgsdd_qfctx)); 383*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "model outputs", sgs_dd_data->num_comp_outputs, CEED_EVAL_NONE)); 384*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "anisotropy tensor", num_comp_grid_aniso, CEED_EVAL_NONE)); 385*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "inverse multiplicity", 1, CEED_EVAL_NONE)); 386*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_dd_outputs, "eigenvectors", num_comp_eigvec, CEED_EVAL_NONE)); 387*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_dd_outputs, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_NONE)); 388*b87d60b3SJames Wright 389*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_dd_outputs, NULL, NULL, &op_sgs_dd_outputs)); 390*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "model outputs", elem_restr_dd_outputs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 391*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "anisotropy tensor", sgs_dd_setup_data->elem_restr_grid_aniso, CEED_BASIS_NONE, 392*b87d60b3SJames Wright sgs_dd_setup_data->grid_aniso_ceed)); 393*b87d60b3SJames Wright PetscCallCeed(ceed, 394*b87d60b3SJames Wright CeedOperatorSetField(op_sgs_dd_outputs, "inverse multiplicity", elem_restr_inv_multiplicity, CEED_BASIS_NONE, inv_multiplicity)); 395*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "eigenvectors", elem_restr_eigvec, CEED_BASIS_NONE, eigvec)); 396*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_dd_outputs, "km_sgs", elem_restr_sgs, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 397*b87d60b3SJames Wright 398*b87d60b3SJames Wright 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*b87d60b3SJames Wright NULL, NULL, &sgs_dd_data->op_nodal_dd_outputs_ctx)); 400*b87d60b3SJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_dd_outputs)); 401*b87d60b3SJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_dd_outputs)); 402*b87d60b3SJames Wright } 403*b87d60b3SJames Wright 404*b87d60b3SJames Wright sgs_dd_data->sgs_nodal_inference = SgsDDNodalStressEval_Sequential_Internal; 405*b87d60b3SJames Wright sgs_dd_data->sgs_nodal_eval = SgsDDNodalStressEval_Sequential; 406*b87d60b3SJames Wright PetscCall(SgsDDSetupNodalEvaluation_Sequential_Internal(ceed, sgs_dd_data, sgs_dd_setup_data, elem_restr_dd_inputs, elem_restr_dd_outputs, 407*b87d60b3SJames Wright elem_restr_inv_multiplicity, inv_multiplicity, &sgs_dd_data->sgs_nodal_inference_ctx)); 408*b87d60b3SJames Wright 409*b87d60b3SJames Wright sgs_dd_setup_data->elem_restr_sgs = elem_restr_sgs; 410*b87d60b3SJames Wright 411*b87d60b3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&inv_multiplicity)); 412*b87d60b3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&eigvec)); 413*b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_inv_multiplicity)); 414*b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_eigvec)); 415*b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_inputs)); 416*b87d60b3SJames Wright PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dd_outputs)); 417*b87d60b3SJames Wright PetscFunctionReturn(PETSC_SUCCESS); 418*b87d60b3SJames Wright } 419*b87d60b3SJames Wright 4209c678832SJames Wright // @brief Create CeedOperator to compute SGS contribution to the residual 421ad494f68SJames Wright static PetscErrorCode SgsSetupNodalIFunction(Ceed ceed, User user, CeedData ceed_data, SgsDDSetupData sgs_dd_setup_data) { 42242454adaSJames Wright SgsDDData sgs_dd_data = user->sgs_dd_data; 42367263decSKenneth E. Jansen CeedInt num_comp_q, num_comp_qd, num_comp_x; 424defe8520SJames Wright PetscInt dim; 4259c678832SJames Wright CeedQFunction qf_sgs_apply; 4269c678832SJames Wright CeedOperator op_sgs_apply; 4279c678832SJames Wright CeedBasis basis_sgs; 4289c678832SJames Wright 4299c678832SJames Wright PetscFunctionBeginUser; 4309c678832SJames Wright PetscCall(DMGetDimension(user->dm, &dim)); 431b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_q, &num_comp_q)); 432b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_qd_i, &num_comp_qd)); 433b4c37c5cSJames Wright PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(ceed_data->elem_restr_x, &num_comp_x)); 4349c678832SJames Wright 43567263decSKenneth E. Jansen PetscCall(CreateBasisFromPlex(ceed, sgs_dd_data->dm_sgs, 0, 0, 0, 0, &basis_sgs)); 4369c678832SJames Wright 4379c678832SJames Wright switch (user->phys->state_var) { 4389c678832SJames Wright case STATEVAR_PRIMITIVE: 43942454adaSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Prim, IFunction_NodalSgs_Prim_loc, &qf_sgs_apply)); 4409c678832SJames Wright break; 4419c678832SJames Wright case STATEVAR_CONSERVATIVE: 44242454adaSJames Wright PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, IFunction_NodalSgs_Conserv, IFunction_NodalSgs_Conserv_loc, &qf_sgs_apply)); 4439c678832SJames Wright break; 4449c678832SJames Wright default: 4459c678832SJames Wright SETERRQ(PetscObjectComm((PetscObject)user->dm), PETSC_ERR_SUP, "Nodal SGS evaluation not available for chosen state variable"); 4469c678832SJames Wright } 4479c678832SJames Wright 44840816385SJames Wright PetscCallCeed(ceed, CeedQFunctionSetContext(qf_sgs_apply, sgs_dd_setup_data->ifunction_qfctx)); 449b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "q", num_comp_q, CEED_EVAL_INTERP)); 450b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "qdata", num_comp_qd, CEED_EVAL_NONE)); 451b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddInput(qf_sgs_apply, "km_sgs", sgs_dd_data->num_comp_sgs, CEED_EVAL_INTERP)); 452b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_sgs_apply, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD)); 4539c678832SJames Wright 454b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_sgs_apply, NULL, NULL, &op_sgs_apply)); 455b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 45658e1cbfdSJeremy L Thompson PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data)); 457b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "km_sgs", sgs_dd_setup_data->elem_restr_sgs, basis_sgs, sgs_dd_data->sgs_nodal_ceed)); 458b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorSetField(op_sgs_apply, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE)); 4599c678832SJames Wright 4609c678832SJames Wright PetscCall( 4619c678832SJames Wright OperatorApplyContextCreate(user->dm, user->dm, ceed, op_sgs_apply, user->q_ceed, user->g_ceed, NULL, NULL, &sgs_dd_data->op_sgs_apply_ctx)); 4629c678832SJames Wright 463b4c37c5cSJames Wright PetscCallCeed(ceed, CeedOperatorDestroy(&op_sgs_apply)); 464b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_sgs_apply)); 465d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4669c678832SJames Wright } 4679c678832SJames Wright 4689c678832SJames Wright // @brief Calculate and add data-driven SGS residual to the global residual 469ad494f68SJames Wright PetscErrorCode SgsDDApplyIFunction(User user, const Vec Q_loc, Vec G_loc) { 47042454adaSJames Wright SgsDDData sgs_dd_data = user->sgs_dd_data; 4719c678832SJames Wright Vec VelocityGradient, SGSNodal_loc; 472cceb3143SJames Wright PetscMemType sgs_nodal_mem_type; 4739c678832SJames Wright 4749c678832SJames Wright PetscFunctionBeginUser; 4759c678832SJames Wright PetscCall(DMGetGlobalVector(user->grad_velo_proj->dm, &VelocityGradient)); 476f6ac214eSJames Wright PetscCall(VelocityGradientProjectionApply(user->grad_velo_proj, Q_loc, VelocityGradient)); 4779c678832SJames Wright 4789c678832SJames Wright // -- Compute Nodal SGS tensor 4799c678832SJames Wright PetscCall(DMGetLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc)); 480cceb3143SJames Wright PetscCall(sgs_dd_data->sgs_nodal_eval(user, Q_loc, VelocityGradient, SGSNodal_loc)); 4819c678832SJames Wright 4829c678832SJames Wright // -- Compute contribution of the SGS stress 483cceb3143SJames Wright PetscCall(VecP2C(SGSNodal_loc, &sgs_nodal_mem_type, sgs_dd_data->sgs_nodal_ceed)); // sgs_nodal_ceed is an implicit input 4849c678832SJames Wright PetscCall(ApplyAddCeedOperatorLocalToLocal(Q_loc, G_loc, sgs_dd_data->op_sgs_apply_ctx)); 4859c678832SJames Wright 4869c678832SJames Wright // -- Return local SGS vector 4879c678832SJames Wright PetscCall(VecC2P(sgs_dd_data->sgs_nodal_ceed, sgs_nodal_mem_type, SGSNodal_loc)); 4889c678832SJames Wright PetscCall(DMRestoreLocalVector(sgs_dd_data->dm_sgs, &SGSNodal_loc)); 4899c678832SJames Wright PetscCall(DMRestoreGlobalVector(user->grad_velo_proj->dm, &VelocityGradient)); 490d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 4919c678832SJames Wright } 4929c678832SJames Wright 49362b7942eSJames Wright // @brief B = A^T, A is NxM, B is MxN 494cceb3143SJames Wright static PetscErrorCode TransposeMatrix(const PetscScalar *A, PetscScalar *B, const PetscInt N, const PetscInt M) { 49562b7942eSJames Wright PetscFunctionBeginUser; 49662b7942eSJames Wright for (PetscInt i = 0; i < N; i++) { 49762b7942eSJames Wright for (PetscInt j = 0; j < M; j++) { 49862b7942eSJames Wright B[j * N + i] = A[i * M + j]; 49962b7942eSJames Wright } 50062b7942eSJames Wright } 501d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 50262b7942eSJames Wright } 50362b7942eSJames Wright 50462b7942eSJames Wright // @brief Read neural network coefficients from file and put into context struct 505ad494f68SJames Wright static PetscErrorCode SgsDDContextFill(MPI_Comm comm, char data_dir[PETSC_MAX_PATH_LEN], SgsDDContext *psgsdd_ctx) { 506ad494f68SJames Wright SgsDDContext sgsdd_ctx; 50762b7942eSJames Wright PetscInt num_inputs = (*psgsdd_ctx)->num_inputs, num_outputs = (*psgsdd_ctx)->num_outputs, num_neurons = (*psgsdd_ctx)->num_neurons; 50862b7942eSJames Wright char file_path[PETSC_MAX_PATH_LEN]; 50962b7942eSJames Wright PetscScalar *temp; 51062b7942eSJames Wright 51162b7942eSJames Wright PetscFunctionBeginUser; 51262b7942eSJames Wright { 513ad494f68SJames Wright SgsDDContext sgsdd_temp; 51462b7942eSJames Wright PetscCall(PetscNew(&sgsdd_temp)); 51562b7942eSJames Wright *sgsdd_temp = **psgsdd_ctx; 51662b7942eSJames Wright sgsdd_temp->offsets.bias1 = 0; 51762b7942eSJames Wright sgsdd_temp->offsets.bias2 = sgsdd_temp->offsets.bias1 + num_neurons; 51862b7942eSJames Wright sgsdd_temp->offsets.weight1 = sgsdd_temp->offsets.bias2 + num_neurons; 51962b7942eSJames Wright sgsdd_temp->offsets.weight2 = sgsdd_temp->offsets.weight1 + num_neurons * num_inputs; 52062b7942eSJames Wright sgsdd_temp->offsets.out_scaling = sgsdd_temp->offsets.weight2 + num_inputs * num_neurons; 52162b7942eSJames Wright PetscInt total_num_scalars = sgsdd_temp->offsets.out_scaling + 2 * num_outputs; 52262b7942eSJames Wright sgsdd_temp->total_bytes = sizeof(*sgsdd_ctx) + total_num_scalars * sizeof(sgsdd_ctx->data[0]); 52362b7942eSJames Wright PetscCall(PetscMalloc(sgsdd_temp->total_bytes, &sgsdd_ctx)); 52462b7942eSJames Wright *sgsdd_ctx = *sgsdd_temp; 52562b7942eSJames Wright PetscCall(PetscFree(sgsdd_temp)); 52662b7942eSJames Wright } 52762b7942eSJames Wright 52862b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b1.dat")); 52942454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias1])); 53062b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "b2.dat")); 53142454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.bias2])); 53262b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "OutScaling.dat")); 53342454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, &sgsdd_ctx->data[sgsdd_ctx->offsets.out_scaling])); 53462b7942eSJames Wright 53562b7942eSJames Wright { 53662b7942eSJames Wright PetscCall(PetscMalloc1(num_inputs * num_neurons, &temp)); 53762b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w1.dat")); 53842454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp)); 53962b7942eSJames Wright PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight1], num_inputs, num_neurons)); 54062b7942eSJames Wright PetscCall(PetscFree(temp)); 54162b7942eSJames Wright } 54262b7942eSJames Wright { 54362b7942eSJames Wright PetscCall(PetscMalloc1(num_outputs * num_neurons, &temp)); 54462b7942eSJames Wright PetscCall(PetscSNPrintf(file_path, sizeof file_path, "%s/%s", data_dir, "w2.dat")); 54542454adaSJames Wright PetscCall(PhastaDatFileReadToArrayReal(comm, file_path, temp)); 54662b7942eSJames Wright PetscCall(TransposeMatrix(temp, &sgsdd_ctx->data[sgsdd_ctx->offsets.weight2], num_neurons, num_outputs)); 54762b7942eSJames Wright PetscCall(PetscFree(temp)); 54862b7942eSJames Wright } 54962b7942eSJames Wright 55062b7942eSJames Wright PetscCall(PetscFree(*psgsdd_ctx)); 55162b7942eSJames Wright *psgsdd_ctx = sgsdd_ctx; 552d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 55362b7942eSJames Wright } 55462b7942eSJames Wright 555ad494f68SJames Wright PetscErrorCode SgsDDSetup(Ceed ceed, User user, CeedData ceed_data, ProblemData *problem) { 556ee1455b7SJames Wright PetscReal alpha = 0; 557ad494f68SJames Wright SgsDDContext sgsdd_ctx; 55862b7942eSJames Wright MPI_Comm comm = user->comm; 559ee1455b7SJames Wright char sgs_dd_dir[PETSC_MAX_PATH_LEN] = "./dd_sgs_parameters"; 560ad494f68SJames Wright SgsDDSetupData sgs_dd_setup_data; 561*b87d60b3SJames Wright PetscBool use_fused; 562ee1455b7SJames Wright NewtonianIdealGasContext gas; 56362b7942eSJames Wright 56406f41313SJames Wright PetscFunctionBeginUser; 565f6ac214eSJames Wright PetscCall(VelocityGradientProjectionSetup(ceed, user, ceed_data, problem, user->phys->state_var, ceed_data->elem_restr_q, ceed_data->basis_q, 566f6ac214eSJames Wright &user->grad_velo_proj)); 5679ab09d51SJames Wright 568*b87d60b3SJames Wright PetscCall(PetscNew(&user->sgs_dd_data)); 569*b87d60b3SJames Wright user->sgs_dd_data->num_comp_inputs = 6; 570*b87d60b3SJames Wright user->sgs_dd_data->num_comp_outputs = 6; 57162b7942eSJames Wright 572*b87d60b3SJames Wright use_fused = PETSC_TRUE; 5734b0f6111SJames Wright PetscOptionsBegin(comm, NULL, "SGS Data-Driven Model Options", NULL); 57462b7942eSJames Wright PetscCall(PetscOptionsReal("-sgs_model_dd_leakyrelu_alpha", "Slope parameter for Leaky ReLU activation function", NULL, alpha, &alpha, NULL)); 57562b7942eSJames Wright PetscCall(PetscOptionsString("-sgs_model_dd_parameter_dir", "Path to directory with model parameters (weights, biases, etc.)", NULL, sgs_dd_dir, 57662b7942eSJames Wright sgs_dd_dir, sizeof(sgs_dd_dir), NULL)); 577*b87d60b3SJames Wright PetscCall( 578*b87d60b3SJames Wright PetscOptionsBool("-sgs_model_dd_use_fused", "Use the fused SGS DD model evaluation instead of sequential", NULL, use_fused, &use_fused, NULL)); 57962b7942eSJames Wright PetscOptionsEnd(); 58062b7942eSJames Wright 581*b87d60b3SJames Wright PetscCall(PetscNew(&sgsdd_ctx)); 5824b0f6111SJames Wright sgsdd_ctx->num_layers = 1; 58362b7942eSJames Wright sgsdd_ctx->num_inputs = 6; 58462b7942eSJames Wright sgsdd_ctx->num_outputs = 6; 58562b7942eSJames Wright sgsdd_ctx->num_neurons = 20; 58662b7942eSJames Wright sgsdd_ctx->alpha = alpha; 58762b7942eSJames Wright 588ad494f68SJames Wright PetscCall(SgsDDContextFill(comm, sgs_dd_dir, &sgsdd_ctx)); 58962b7942eSJames Wright 590ee1455b7SJames Wright // -- Create DM for storing SGS tensor at nodes 591ad494f68SJames Wright 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)); 592ee1455b7SJames Wright 593c38c977aSJames Wright PetscCall(PetscNew(&sgs_dd_setup_data)); 594c38c977aSJames Wright 595b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas)); 596ee1455b7SJames Wright sgsdd_ctx->gas = *gas; 597b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas)); 598b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextCreate(user->ceed, &sgs_dd_setup_data->sgsdd_qfctx)); 599b4c37c5cSJames Wright PetscCallCeed(ceed, 600b4c37c5cSJames Wright CeedQFunctionContextSetData(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, CEED_USE_POINTER, sgsdd_ctx->total_bytes, sgsdd_ctx)); 601b4c37c5cSJames Wright PetscCallCeed(ceed, CeedQFunctionContextSetDataDestroy(sgs_dd_setup_data->sgsdd_qfctx, CEED_MEM_HOST, FreeContextPetsc)); 602ee1455b7SJames Wright 60340816385SJames Wright PetscCallCeed(ceed, CeedQFunctionContextReferenceCopy(problem->apply_vol_ifunction.qfunction_context, &sgs_dd_setup_data->ifunction_qfctx)); 60440816385SJames Wright 605c38c977aSJames Wright // -- Compute and store anisotropy tensor 606c38c977aSJames Wright PetscCall(GridAnisotropyTensorProjectionSetupApply(ceed, user, ceed_data, &sgs_dd_setup_data->elem_restr_grid_aniso, 607c38c977aSJames Wright &sgs_dd_setup_data->grid_aniso_ceed)); 608c38c977aSJames Wright 609ee1455b7SJames Wright // -- Create Nodal Evaluation Operator 610*b87d60b3SJames Wright if (use_fused) PetscCall(SgsDDSetupNodalEvaluation_Fused(ceed, user, ceed_data, sgs_dd_setup_data)); 611*b87d60b3SJames Wright else PetscCall(SgsDDSetupNodalEvaluation_Sequential(ceed, user, ceed_data, sgs_dd_setup_data)); 612ee1455b7SJames Wright 6139c678832SJames Wright // -- Create Operator to evalutate residual of SGS stress 614ad494f68SJames Wright PetscCall(SgsSetupNodalIFunction(ceed, user, ceed_data, sgs_dd_setup_data)); 6159c678832SJames Wright 616ad494f68SJames Wright PetscCall(SgsDDSetupDataDestroy(sgs_dd_setup_data)); 617d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 61862b7942eSJames Wright } 619ee1455b7SJames Wright 62042454adaSJames Wright PetscErrorCode SgsDDDataDestroy(SgsDDData sgs_dd_data) { 621ee1455b7SJames Wright PetscFunctionBeginUser; 622d949ddfcSJames Wright if (!sgs_dd_data) PetscFunctionReturn(PETSC_SUCCESS); 623b4c37c5cSJames Wright Ceed ceed = sgs_dd_data->op_sgs_apply_ctx->ceed; 624ee1455b7SJames Wright 625b4c37c5cSJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->sgs_nodal_ceed)); 626*b87d60b3SJames Wright PetscCallCeed(ceed, CeedVectorDestroy(&sgs_dd_data->grad_velo_ceed)); 627ee1455b7SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_evaluation_ctx)); 62841edf198SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_sgs_apply_ctx)); 629*b87d60b3SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_dd_inputs_ctx)); 630*b87d60b3SJames Wright PetscCall(OperatorApplyContextDestroy(sgs_dd_data->op_nodal_dd_outputs_ctx)); 631ee1455b7SJames Wright PetscCall(DMDestroy(&sgs_dd_data->dm_sgs)); 632*b87d60b3SJames Wright PetscCall(DMDestroy(&sgs_dd_data->dm_dd_inputs)); 633*b87d60b3SJames Wright PetscCall(DMDestroy(&sgs_dd_data->dm_dd_outputs)); 634*b87d60b3SJames Wright if (sgs_dd_data->sgs_nodal_inference_ctx) PetscCall(sgs_dd_data->sgs_nodal_inference_ctx_destroy(sgs_dd_data->sgs_nodal_inference_ctx)); 635ee1455b7SJames Wright PetscCall(PetscFree(sgs_dd_data)); 636d949ddfcSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 637ee1455b7SJames Wright } 638