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