1 // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors. 2 // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause 3 4 #include "../qfunctions/strong_boundary_conditions.h" 5 6 #include <ceed.h> 7 #include <petscdmplex.h> 8 9 #include <navierstokes.h> 10 #include "../problems/stg_shur14.h" 11 12 PetscErrorCode SetupStrongSTG_Ceed(Ceed ceed, CeedData ceed_data, DM dm, ProblemData problem, SimpleBC bc, Physics phys, CeedOperator op_strong_bc) { 13 CeedInt num_comp_x = problem->dim, num_comp_q = 5, stg_data_size = 1, dXdx_size; 14 CeedVector multiplicity, x_stored, scale_stored, stg_data, dXdx; 15 CeedBasis basis_x_to_q_face; 16 CeedElemRestriction elem_restr_x_face, elem_restr_q_face, elem_restr_x_stored, elem_restr_scale, elem_restr_stgdata, elem_restr_dXdx; 17 CeedQFunction qf_setup, qf_strongbc, qf_stgdata; 18 CeedOperator op_setup, op_strong_bc_sub, op_stgdata; 19 DMLabel domain_label; 20 const PetscInt dm_field = 0, height_face = 1, height_cell = 0; 21 PetscInt dim; 22 23 PetscFunctionBeginUser; 24 PetscCall(DMGetLabel(dm, "Face Sets", &domain_label)); 25 PetscCall(DMGetDimension(dm, &dim)); 26 dXdx_size = num_comp_x * (dim - height_cell); 27 28 { // Basis 29 CeedBasis basis_x_face, basis_q_face; 30 DM dm_coord; 31 32 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 33 DMLabel label = NULL; 34 PetscInt label_value = 0; 35 PetscCall(CreateBasisFromPlex(ceed, dm, label, label_value, height_face, dm_field, &basis_q_face)); 36 PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, label_value, height_face, dm_field, &basis_x_face)); 37 38 PetscCallCeed(ceed, CeedBasisCreateProjection(basis_x_face, basis_q_face, &basis_x_to_q_face)); 39 40 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q_face)); 41 PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_face)); 42 } 43 44 // Setup QFunction 45 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, SetupStrongBC, SetupStrongBC_loc, &qf_setup)); 46 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "x", num_comp_x, CEED_EVAL_INTERP)); 47 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "dxdX", dXdx_size, CEED_EVAL_GRAD)); 48 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "multiplicity", num_comp_q, CEED_EVAL_NONE)); 49 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "x stored", num_comp_x, CEED_EVAL_NONE)); 50 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "scale", 1, CEED_EVAL_NONE)); 51 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "dXdx", dXdx_size, CEED_EVAL_NONE)); 52 53 // Setup STG QFunctions 54 PetscCall(SetupStrongStg_PreProcessing(ceed, problem, num_comp_x, stg_data_size, dXdx_size, &qf_stgdata)); 55 PetscCall(SetupStrongStg_QF(ceed, problem, num_comp_x, num_comp_q, stg_data_size, dXdx_size, &qf_strongbc)); 56 57 // Compute contribution on each boundary face 58 for (CeedInt i = 0; i < bc->num_inflow; i++) { 59 DMLabel face_orientation_label; 60 PetscInt num_orientations_values, *orientation_values; 61 62 { 63 char *face_orientation_label_name; 64 65 PetscCall(DMPlexCreateFaceLabel(dm, bc->inflows[i], &face_orientation_label_name)); 66 PetscCall(DMGetLabel(dm, face_orientation_label_name, &face_orientation_label)); 67 PetscCall(PetscFree(face_orientation_label_name)); 68 } 69 PetscCall(DMLabelCreateGlobalValueArray(dm, face_orientation_label, &num_orientations_values, &orientation_values)); 70 for (PetscInt o = 0; o < num_orientations_values; o++) { 71 CeedBasis basis_x_to_q_cell; 72 CeedElemRestriction elem_restr_x_cell; 73 PetscInt orientation = orientation_values[o]; 74 75 { 76 CeedBasis basis_x_cell_to_face, basis_q_face; 77 78 PetscCall(DMPlexCeedBasisCellToFaceCoordinateCreate(ceed, dm, face_orientation_label, orientation, orientation, &basis_x_cell_to_face)); 79 PetscCall(CreateBasisFromPlex(ceed, dm, face_orientation_label, orientation, height_face, dm_field, &basis_q_face)); 80 PetscCallCeed(ceed, CeedBasisCreateProjection(basis_x_cell_to_face, basis_q_face, &basis_x_to_q_cell)); 81 PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_cell_to_face)); 82 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q_face)); 83 } 84 85 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, bc->inflows[i], height_face, dm_field, &elem_restr_q_face)); 86 PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, bc->inflows[i], height_face, &elem_restr_x_face)); 87 PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, face_orientation_label, orientation, height_cell, &elem_restr_x_cell)); 88 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_q_face, &multiplicity, NULL)); 89 PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(elem_restr_q_face, multiplicity)); 90 91 PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, dm, domain_label, bc->inflows[i], height_face, num_comp_x, &elem_restr_x_stored)); 92 PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, dm, domain_label, bc->inflows[i], height_face, 1, &elem_restr_scale)); 93 PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, dm, domain_label, bc->inflows[i], height_face, stg_data_size, &elem_restr_stgdata)); 94 PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, dm, domain_label, bc->inflows[i], height_face, dXdx_size, &elem_restr_dXdx)); 95 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_x_stored, &x_stored, NULL)); 96 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_scale, &scale_stored, NULL)); 97 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_stgdata, &stg_data, NULL)); 98 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_dXdx, &dXdx, NULL)); 99 100 // -- Setup Operator 101 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup)); 102 PetscCallCeed(ceed, CeedOperatorSetName(op_setup, "Precomputed data for strong boundary conditions")); 103 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "x", elem_restr_x_face, basis_x_to_q_face, CEED_VECTOR_ACTIVE)); 104 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "dxdX", elem_restr_x_cell, basis_x_to_q_cell, CEED_VECTOR_ACTIVE)); 105 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "multiplicity", elem_restr_q_face, CEED_BASIS_NONE, multiplicity)); 106 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "x stored", elem_restr_x_stored, CEED_BASIS_NONE, x_stored)); 107 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "scale", elem_restr_scale, CEED_BASIS_NONE, scale_stored)); 108 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "dXdx", elem_restr_dXdx, CEED_BASIS_NONE, dXdx)); 109 110 // -- Compute geometric factors 111 PetscCallCeed(ceed, CeedOperatorApply(op_setup, ceed_data->x_coord, CEED_VECTOR_NONE, CEED_REQUEST_IMMEDIATE)); 112 113 // -- Compute STGData 114 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stgdata, NULL, NULL, &op_stgdata)); 115 PetscCallCeed(ceed, CeedOperatorSetField(op_stgdata, "dXdx", elem_restr_dXdx, CEED_BASIS_NONE, dXdx)); 116 PetscCallCeed(ceed, CeedOperatorSetField(op_stgdata, "x", elem_restr_x_stored, CEED_BASIS_NONE, x_stored)); 117 PetscCallCeed(ceed, CeedOperatorSetField(op_stgdata, "stg data", elem_restr_stgdata, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 118 119 PetscCallCeed(ceed, CeedOperatorApply(op_stgdata, CEED_VECTOR_NONE, stg_data, CEED_REQUEST_IMMEDIATE)); 120 121 // -- Setup BC QFunctions 122 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_strongbc, NULL, NULL, &op_strong_bc_sub)); 123 PetscCallCeed(ceed, CeedOperatorSetName(op_strong_bc_sub, "Strong STG")); 124 125 PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "dXdx", elem_restr_dXdx, CEED_BASIS_NONE, dXdx)); 126 PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "x", elem_restr_x_stored, CEED_BASIS_NONE, x_stored)); 127 PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "scale", elem_restr_scale, CEED_BASIS_NONE, scale_stored)); 128 PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "stg data", elem_restr_stgdata, CEED_BASIS_NONE, stg_data)); 129 PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "q", elem_restr_q_face, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 130 131 // -- Add to composite operator 132 PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_strong_bc, op_strong_bc_sub)); 133 134 PetscCallCeed(ceed, CeedVectorDestroy(&multiplicity)); 135 PetscCallCeed(ceed, CeedVectorDestroy(&x_stored)); 136 PetscCallCeed(ceed, CeedVectorDestroy(&scale_stored)); 137 PetscCallCeed(ceed, CeedVectorDestroy(&stg_data)); 138 PetscCallCeed(ceed, CeedVectorDestroy(&dXdx)); 139 PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q_cell)); 140 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_face)); 141 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_cell)); 142 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_face)); 143 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_stored)); 144 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_scale)); 145 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_stgdata)); 146 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dXdx)); 147 PetscCallCeed(ceed, CeedOperatorDestroy(&op_strong_bc_sub)); 148 PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup)); 149 PetscCallCeed(ceed, CeedOperatorDestroy(&op_stgdata)); 150 } 151 PetscCall(PetscFree(orientation_values)); 152 } 153 154 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_strong_bc, "solution time", &phys->stg_solution_time_label)); 155 156 PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q_face)); 157 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_strongbc)); 158 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stgdata)); 159 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup)); 160 PetscFunctionReturn(PETSC_SUCCESS); 161 } 162 163 PetscErrorCode DMPlexInsertBoundaryValues_StrongBCCeed(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, 164 Vec cell_geom_FVM, Vec grad_FVM) { 165 Vec boundary_mask; 166 User user; 167 168 PetscFunctionBeginUser; 169 PetscCall(DMGetApplicationContext(dm, &user)); 170 171 if (user->phys->stg_solution_time_label) { 172 PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user->op_strong_bc_ctx->op, user->phys->stg_solution_time_label, &time)); 173 } 174 175 // Mask Strong BC entries 176 PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 177 PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask)); 178 PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 179 180 PetscCall(ApplyAddCeedOperatorLocalToLocal(NULL, Q_loc, user->op_strong_bc_ctx)); 181 PetscFunctionReturn(PETSC_SUCCESS); 182 } 183 184 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, ProblemData problem, SimpleBC bc) { 185 CeedOperator op_strong_bc; 186 187 PetscFunctionBeginUser; 188 { 189 Vec boundary_mask, global_vec; 190 191 PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 192 PetscCall(DMGetGlobalVector(dm, &global_vec)); 193 PetscCall(VecZeroEntries(boundary_mask)); 194 PetscCall(VecSet(global_vec, 1.0)); 195 PetscCall(DMGlobalToLocal(dm, global_vec, INSERT_VALUES, boundary_mask)); 196 PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 197 PetscCall(DMRestoreGlobalVector(dm, &global_vec)); 198 } 199 200 PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_strong_bc)); 201 { 202 PetscBool use_strongstg = PETSC_FALSE; 203 PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL)); 204 205 if (use_strongstg) { 206 PetscCall(SetupStrongSTG_Ceed(ceed, ceed_data, dm, problem, bc, user->phys, op_strong_bc)); 207 } 208 } 209 210 PetscCall(OperatorApplyContextCreate(NULL, NULL, ceed, op_strong_bc, CEED_VECTOR_NONE, NULL, NULL, NULL, &user->op_strong_bc_ctx)); 211 212 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_StrongBCCeed)); 213 PetscCallCeed(ceed, CeedOperatorDestroy(&op_strong_bc)); 214 PetscFunctionReturn(PETSC_SUCCESS); 215 } 216