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, Honee honee, DM dm, ProblemData problem, Physics phys, CeedOperator op_strong_bc) { 13 CeedInt num_comp_x, 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 MPI_Comm comm = honee->comm; 23 24 PetscFunctionBeginUser; 25 PetscCall(DMGetLabel(dm, "Face Sets", &domain_label)); 26 PetscCall(DMGetDimension(dm, &dim)); 27 num_comp_x = dim; 28 dXdx_size = num_comp_x * (dim - height_cell); 29 30 { // Basis 31 CeedBasis basis_x_face, basis_q_face; 32 DM dm_coord; 33 34 PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 35 DMLabel label = NULL; 36 PetscInt label_value = 0; 37 PetscCall(CreateBasisFromPlex(ceed, dm, label, label_value, height_face, dm_field, &basis_q_face)); 38 PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, label_value, height_face, dm_field, &basis_x_face)); 39 40 PetscCallCeed(ceed, CeedBasisCreateProjection(basis_x_face, basis_q_face, &basis_x_to_q_face)); 41 42 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q_face)); 43 PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_face)); 44 } 45 46 // Setup QFunction 47 PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, SetupStrongBC, SetupStrongBC_loc, &qf_setup)); 48 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "x", num_comp_x, CEED_EVAL_INTERP)); 49 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "dxdX", dXdx_size, CEED_EVAL_GRAD)); 50 PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "multiplicity", num_comp_q, CEED_EVAL_NONE)); 51 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "x stored", num_comp_x, CEED_EVAL_NONE)); 52 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "scale", 1, CEED_EVAL_NONE)); 53 PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "dXdx", dXdx_size, CEED_EVAL_NONE)); 54 55 // Setup STG QFunctions 56 PetscCall(SetupStrongStg_PreProcessing(ceed, problem, num_comp_x, stg_data_size, dXdx_size, &qf_stgdata)); 57 PetscCall(SetupStrongStg_QF(ceed, problem, num_comp_x, num_comp_q, stg_data_size, dXdx_size, &qf_strongbc)); 58 59 // Get face label values for the inflow BC 60 PetscInt num_inflow_faces; 61 const PetscInt *inflow_faces; 62 { 63 BCDefinition bc_def = NULL; 64 for (PetscCount b = 0; b < problem->num_bc_defs; b++) { 65 BCDefinition bc_def_ = problem->bc_defs[b]; 66 const char *name; 67 68 PetscCall(BCDefinitionGetInfo(bc_def_, &name, NULL, NULL)); 69 if (!strcmp(name, "inflow")) { 70 bc_def = bc_def_; 71 break; 72 } 73 } 74 PetscCheck(bc_def, comm, PETSC_ERR_SUP, "Could not find BCDefintion named 'inflow' for STG setup."); 75 PetscCall(BCDefinitionGetInfo(bc_def, NULL, &num_inflow_faces, &inflow_faces)); 76 } 77 78 // Compute contribution on each boundary face 79 for (CeedInt i = 0; i < num_inflow_faces; i++) { 80 DMLabel face_orientation_label; 81 PetscInt num_orientations_values, *orientation_values; 82 83 { 84 char *face_orientation_label_name; 85 86 PetscCall(DMPlexCreateFaceLabel(dm, inflow_faces[i], &face_orientation_label_name)); 87 PetscCall(DMGetLabel(dm, face_orientation_label_name, &face_orientation_label)); 88 PetscCall(PetscFree(face_orientation_label_name)); 89 } 90 PetscCall(DMLabelCreateGlobalValueArray(dm, face_orientation_label, &num_orientations_values, &orientation_values)); 91 for (PetscInt o = 0; o < num_orientations_values; o++) { 92 CeedBasis basis_x_to_q_cell; 93 CeedElemRestriction elem_restr_x_cell; 94 PetscInt orientation = orientation_values[o]; 95 96 { 97 CeedBasis basis_x_cell_to_face, basis_q_face; 98 99 PetscCall(DMPlexCeedBasisCellToFaceCoordinateCreate(ceed, dm, face_orientation_label, orientation, orientation, &basis_x_cell_to_face)); 100 PetscCall(CreateBasisFromPlex(ceed, dm, face_orientation_label, orientation, height_face, dm_field, &basis_q_face)); 101 PetscCallCeed(ceed, CeedBasisCreateProjection(basis_x_cell_to_face, basis_q_face, &basis_x_to_q_cell)); 102 PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_cell_to_face)); 103 PetscCallCeed(ceed, CeedBasisDestroy(&basis_q_face)); 104 } 105 106 PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, face_orientation_label, orientation, height_face, dm_field, &elem_restr_q_face)); 107 PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, face_orientation_label, orientation, height_face, &elem_restr_x_face)); 108 PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, face_orientation_label, orientation, height_cell, &elem_restr_x_cell)); 109 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_q_face, &multiplicity, NULL)); 110 PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(elem_restr_q_face, multiplicity)); 111 112 PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, dm, face_orientation_label, orientation, height_face, num_comp_x, 113 &elem_restr_x_stored)); 114 PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, dm, face_orientation_label, orientation, height_face, 1, &elem_restr_scale)); 115 PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, dm, face_orientation_label, orientation, height_face, stg_data_size, 116 &elem_restr_stgdata)); 117 PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, dm, face_orientation_label, orientation, height_face, dXdx_size, &elem_restr_dXdx)); 118 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_x_stored, &x_stored, NULL)); 119 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_scale, &scale_stored, NULL)); 120 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_stgdata, &stg_data, NULL)); 121 PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_dXdx, &dXdx, NULL)); 122 123 // -- Setup Operator 124 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup)); 125 PetscCallCeed(ceed, CeedOperatorSetName(op_setup, "Precomputed data for strong boundary conditions")); 126 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "x", elem_restr_x_face, basis_x_to_q_face, CEED_VECTOR_ACTIVE)); 127 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "dxdX", elem_restr_x_cell, basis_x_to_q_cell, CEED_VECTOR_ACTIVE)); 128 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "multiplicity", elem_restr_q_face, CEED_BASIS_NONE, multiplicity)); 129 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "x stored", elem_restr_x_stored, CEED_BASIS_NONE, x_stored)); 130 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "scale", elem_restr_scale, CEED_BASIS_NONE, scale_stored)); 131 PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "dXdx", elem_restr_dXdx, CEED_BASIS_NONE, dXdx)); 132 133 // -- Compute geometric factors 134 PetscCallCeed(ceed, CeedOperatorApply(op_setup, honee->x_coord, CEED_VECTOR_NONE, CEED_REQUEST_IMMEDIATE)); 135 136 // -- Compute STGData 137 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stgdata, NULL, NULL, &op_stgdata)); 138 PetscCallCeed(ceed, CeedOperatorSetField(op_stgdata, "dXdx", elem_restr_dXdx, CEED_BASIS_NONE, dXdx)); 139 PetscCallCeed(ceed, CeedOperatorSetField(op_stgdata, "x", elem_restr_x_stored, CEED_BASIS_NONE, x_stored)); 140 PetscCallCeed(ceed, CeedOperatorSetField(op_stgdata, "stg data", elem_restr_stgdata, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 141 142 PetscCallCeed(ceed, CeedOperatorApply(op_stgdata, CEED_VECTOR_NONE, stg_data, CEED_REQUEST_IMMEDIATE)); 143 144 // -- Setup BC QFunctions 145 PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_strongbc, NULL, NULL, &op_strong_bc_sub)); 146 PetscCallCeed(ceed, CeedOperatorSetName(op_strong_bc_sub, "Strong STG")); 147 148 PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "dXdx", elem_restr_dXdx, CEED_BASIS_NONE, dXdx)); 149 PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "x", elem_restr_x_stored, CEED_BASIS_NONE, x_stored)); 150 PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "scale", elem_restr_scale, CEED_BASIS_NONE, scale_stored)); 151 PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "stg data", elem_restr_stgdata, CEED_BASIS_NONE, stg_data)); 152 PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "q", elem_restr_q_face, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE)); 153 154 // -- Add to composite operator 155 PetscCallCeed(ceed, CeedOperatorCompositeAddSub(op_strong_bc, op_strong_bc_sub)); 156 157 PetscCallCeed(ceed, CeedVectorDestroy(&multiplicity)); 158 PetscCallCeed(ceed, CeedVectorDestroy(&x_stored)); 159 PetscCallCeed(ceed, CeedVectorDestroy(&scale_stored)); 160 PetscCallCeed(ceed, CeedVectorDestroy(&stg_data)); 161 PetscCallCeed(ceed, CeedVectorDestroy(&dXdx)); 162 PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q_cell)); 163 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_face)); 164 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_cell)); 165 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_face)); 166 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_stored)); 167 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_scale)); 168 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_stgdata)); 169 PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dXdx)); 170 PetscCallCeed(ceed, CeedOperatorDestroy(&op_strong_bc_sub)); 171 PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup)); 172 PetscCallCeed(ceed, CeedOperatorDestroy(&op_stgdata)); 173 } 174 PetscCall(PetscFree(orientation_values)); 175 } 176 177 PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_strong_bc, "solution time", &phys->stg_solution_time_label)); 178 179 PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q_face)); 180 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_strongbc)); 181 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stgdata)); 182 PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup)); 183 PetscFunctionReturn(PETSC_SUCCESS); 184 } 185 186 PetscErrorCode DMPlexInsertBoundaryValues_StrongBCCeed(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, 187 Vec cell_geom_FVM, Vec grad_FVM) { 188 Vec boundary_mask; 189 Honee honee; 190 191 PetscFunctionBeginUser; 192 PetscCall(PetscLogEventBegin(HONEE_StrongBCInsert, dm, Q_loc, 0, 0)); 193 PetscCall(DMGetApplicationContext(dm, &honee)); 194 195 if (honee->phys->stg_solution_time_label) { 196 PetscCallCeed(honee->ceed, CeedOperatorSetContextDouble(honee->op_strong_bc_ctx->op, honee->phys->stg_solution_time_label, &time)); 197 } 198 199 // Mask Strong BC entries 200 PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 201 PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask)); 202 PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 203 204 PetscCall(PetscLogEventBegin(HONEE_StrongBCCeed, dm, Q_loc, 0, 0)); 205 PetscCall(ApplyAddCeedOperatorLocalToLocal(NULL, Q_loc, honee->op_strong_bc_ctx)); 206 PetscCall(PetscLogEventEnd(HONEE_StrongBCCeed, dm, Q_loc, 0, 0)); 207 PetscCall(PetscLogEventEnd(HONEE_StrongBCInsert, dm, Q_loc, 0, 0)); 208 PetscFunctionReturn(PETSC_SUCCESS); 209 } 210 211 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, DM dm, Honee honee, ProblemData problem) { 212 CeedOperator op_strong_bc; 213 214 PetscFunctionBeginUser; 215 { 216 Vec boundary_mask, global_vec; 217 218 PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 219 PetscCall(DMGetGlobalVector(dm, &global_vec)); 220 PetscCall(VecZeroEntries(boundary_mask)); 221 PetscCall(VecSet(global_vec, 1.0)); 222 PetscCall(DMGlobalToLocal(dm, global_vec, INSERT_VALUES, boundary_mask)); 223 PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 224 PetscCall(DMRestoreGlobalVector(dm, &global_vec)); 225 } 226 227 PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, &op_strong_bc)); 228 { 229 PetscBool use_strongstg = PETSC_FALSE; 230 PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL)); 231 232 if (use_strongstg) { 233 PetscCall(SetupStrongSTG_Ceed(ceed, honee, dm, problem, honee->phys, op_strong_bc)); 234 } 235 } 236 237 PetscCall(OperatorApplyContextCreate(NULL, NULL, ceed, op_strong_bc, CEED_VECTOR_NONE, NULL, NULL, NULL, &honee->op_strong_bc_ctx)); 238 239 PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_StrongBCCeed)); 240 PetscCallCeed(ceed, CeedOperatorDestroy(&op_strong_bc)); 241 PetscFunctionReturn(PETSC_SUCCESS); 242 } 243