1 // Copyright (c) 2017-2024, 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/strong_boundary_conditions.h" 9 10 #include <ceed.h> 11 #include <petscdmplex.h> 12 13 #include "../navierstokes.h" 14 #include "../problems/stg_shur14.h" 15 16 PetscErrorCode SetupStrongSTG_Ceed(Ceed ceed, CeedData ceed_data, DM dm, ProblemData *problem, SimpleBC bc, Physics phys, CeedOperator op_strong_bc) { 17 CeedInt num_comp_x = problem->dim, num_comp_q = 5, stg_data_size = 1, dim_boundary = 2, dXdx_size = num_comp_x * dim_boundary; 18 CeedVector multiplicity, x_stored, scale_stored, stg_data, dXdx; 19 CeedBasis basis_x_to_q_sur; 20 CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_x_stored, elem_restr_scale, elem_restr_stgdata, elem_restr_dXdx; 21 CeedQFunction qf_setup, qf_strongbc, qf_stgdata; 22 CeedOperator op_setup, op_strong_bc_sub, op_stgdata; 23 DMLabel domain_label; 24 PetscInt dm_field = 0, height = 1; 25 26 PetscFunctionBeginUser; 27 PetscCall(DMGetLabel(dm, "Face Sets", &domain_label)); 28 29 // Basis 30 PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x_sur, ceed_data->basis_q_sur, &basis_x_to_q_sur)); 31 32 // Setup QFunction 33 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, SetupStrongBC, SetupStrongBC_loc, &qf_setup)); 34 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "x", num_comp_x, CEED_EVAL_INTERP)); 35 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "dxdX", num_comp_x * dim_boundary, CEED_EVAL_GRAD)); 36 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "multiplicity", num_comp_q, CEED_EVAL_NONE)); 37 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "x stored", num_comp_x, CEED_EVAL_NONE)); 38 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "scale", 1, CEED_EVAL_NONE)); 39 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "dXdx", dXdx_size, CEED_EVAL_NONE)); 40 41 // Setup STG Setup QFunction 42 PetscCall(SetupStrongStg_PreProcessing(ceed, problem, num_comp_x, stg_data_size, dXdx_size, &qf_stgdata)); 43 44 // Compute contribution on each boundary face 45 for (CeedInt i = 0; i < bc->num_inflow; i++) { 46 // -- Restrictions 47 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, bc->inflows[i], height, dm_field, &elem_restr_q_sur)); 48 PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, bc->inflows[i], height, &elem_restr_x_sur)); 49 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_q_sur, &multiplicity, NULL)); 50 PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(elem_restr_q_sur, multiplicity)); 51 52 PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, dm, domain_label, bc->inflows[i], height, num_comp_x, &elem_restr_x_stored)); 53 PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, dm, domain_label, bc->inflows[i], height, 1, &elem_restr_scale)); 54 PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, dm, domain_label, bc->inflows[i], height, stg_data_size, &elem_restr_stgdata)); 55 PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, dm, domain_label, bc->inflows[i], height, dXdx_size, &elem_restr_dXdx)); 56 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_x_stored, &x_stored, NULL)); 57 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_scale, &scale_stored, NULL)); 58 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_stgdata, &stg_data, NULL)); 59 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_dXdx, &dXdx, NULL)); 60 61 // -- Setup Operator 62 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup)); 63 PetscCallCeed(ceed, CeedOperatorSetName(op_setup, "Precomputed data for strong boundary conditions")); 64 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "x", elem_restr_x_sur, basis_x_to_q_sur, CEED_VECTOR_ACTIVE)); 65 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "dxdX", elem_restr_x_sur, basis_x_to_q_sur, CEED_VECTOR_ACTIVE)); 66 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "multiplicity", elem_restr_q_sur, CEED_BASIS_NONE, multiplicity)); 67 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "x stored", elem_restr_x_stored, CEED_BASIS_NONE, x_stored)); 68 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "scale", elem_restr_scale, CEED_BASIS_NONE, scale_stored)); 69 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "dXdx", elem_restr_dXdx, CEED_BASIS_NONE, dXdx)); 70 71 // -- Compute geometric factors 72 PetscCallCeed(ceed, CeedOperatorApply(op_setup, ceed_data->x_coord, CEED_VECTOR_NONE, CEED_REQUEST_IMMEDIATE)); 73 74 // -- Compute STGData 75 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stgdata, NULL, NULL, &op_stgdata)); 76 PetscCallCeed(ceed, CeedOperatorSetField(op_stgdata, "dXdx", elem_restr_dXdx, CEED_BASIS_NONE, dXdx)); 77 PetscCallCeed(ceed, CeedOperatorSetField(op_stgdata, "x", elem_restr_x_stored, CEED_BASIS_NONE, x_stored)); 78 PetscCallCeed(ceed, CeedOperatorSetField(op_stgdata, "stg data", elem_restr_stgdata, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 79 80 PetscCallCeed(ceed, CeedOperatorApply(op_stgdata, CEED_VECTOR_NONE, stg_data, CEED_REQUEST_IMMEDIATE)); 81 82 // -- Setup BC QFunctions 83 PetscCall(SetupStrongStg_QF(ceed, problem, num_comp_x, num_comp_q, stg_data_size, dXdx_size, &qf_strongbc)); 84 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_strongbc, NULL, NULL, &op_strong_bc_sub)); 85 PetscCallCeed(ceed, CeedOperatorSetName(op_strong_bc_sub, "Strong STG")); 86 87 PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "dXdx", elem_restr_dXdx, CEED_BASIS_NONE, dXdx)); 88 PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "x", elem_restr_x_stored, CEED_BASIS_NONE, x_stored)); 89 PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "scale", elem_restr_scale, CEED_BASIS_NONE, scale_stored)); 90 PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "stg data", elem_restr_stgdata, CEED_BASIS_NONE, stg_data)); 91 PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "q", elem_restr_q_sur, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 92 93 // -- Add to composite operator 94 PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_strong_bc, op_strong_bc_sub)); 95 96 PetscCallCeed(ceed, CeedVectorDestroy(&multiplicity)); 97 PetscCallCeed(ceed, CeedVectorDestroy(&x_stored)); 98 PetscCallCeed(ceed, CeedVectorDestroy(&scale_stored)); 99 PetscCallCeed(ceed, CeedVectorDestroy(&stg_data)); 100 PetscCallCeed(ceed, CeedVectorDestroy(&dXdx)); 101 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur)); 102 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur)); 103 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_stored)); 104 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_scale)); 105 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_stgdata)); 106 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dXdx)); 107 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_strongbc)); 108 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stgdata)); 109 PetscCallCeed(ceed, CeedOperatorDestroy(&op_strong_bc_sub)); 110 PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup)); 111 PetscCallCeed(ceed, CeedOperatorDestroy(&op_stgdata)); 112 } 113 114 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_strong_bc, "solution time", &phys->stg_solution_time_label)); 115 116 PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q_sur)); 117 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup)); 118 PetscFunctionReturn(PETSC_SUCCESS); 119 } 120 121 PetscErrorCode DMPlexInsertBoundaryValues_StrongBCCeed(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, 122 Vec cell_geom_FVM, Vec grad_FVM) { 123 Vec boundary_mask; 124 User user; 125 126 PetscFunctionBeginUser; 127 PetscCall(DMGetApplicationContext(dm, &user)); 128 129 if (user->phys->stg_solution_time_label) { 130 PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user->op_strong_bc_ctx->op, user->phys->stg_solution_time_label, &time)); 131 } 132 133 // Mask Strong BC entries 134 PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 135 PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask)); 136 PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 137 138 PetscCall(ApplyAddCeedOperatorLocalToLocal(NULL, Q_loc, user->op_strong_bc_ctx)); 139 PetscFunctionReturn(PETSC_SUCCESS); 140 } 141 142 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, ProblemData *problem, SimpleBC bc) { 143 CeedOperator op_strong_bc; 144 145 PetscFunctionBeginUser; 146 { 147 Vec boundary_mask, global_vec; 148 149 PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 150 PetscCall(DMGetGlobalVector(dm, &global_vec)); 151 PetscCall(VecZeroEntries(boundary_mask)); 152 PetscCall(VecSet(global_vec, 1.0)); 153 PetscCall(DMGlobalToLocal(dm, global_vec, INSERT_VALUES, boundary_mask)); 154 PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 155 PetscCall(DMRestoreGlobalVector(dm, &global_vec)); 156 } 157 158 PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_strong_bc)); 159 { 160 PetscBool use_strongstg = PETSC_FALSE; 161 PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL)); 162 163 if (use_strongstg) { 164 PetscCall(SetupStrongSTG_Ceed(ceed, ceed_data, dm, problem, bc, user->phys, op_strong_bc)); 165 } 166 } 167 168 PetscCall(OperatorApplyContextCreate(NULL, NULL, ceed, op_strong_bc, CEED_VECTOR_NONE, NULL, NULL, NULL, &user->op_strong_bc_ctx)); 169 170 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_StrongBCCeed)); 171 PetscFunctionReturn(PETSC_SUCCESS); 172 } 173