1*2eb7bf1fSJames Wright // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors. 2*2eb7bf1fSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details. 3*2eb7bf1fSJames Wright // 4*2eb7bf1fSJames Wright // SPDX-License-Identifier: BSD-2-Clause 5*2eb7bf1fSJames Wright // 6*2eb7bf1fSJames Wright // This file is part of CEED: http://github.com/ceed 7*2eb7bf1fSJames Wright 8*2eb7bf1fSJames Wright #include "../qfunctions/strong_boundary_conditions.h" 9*2eb7bf1fSJames Wright 10*2eb7bf1fSJames Wright #include <ceed.h> 11*2eb7bf1fSJames Wright #include <petscdmplex.h> 12*2eb7bf1fSJames Wright 13*2eb7bf1fSJames Wright #include "../navierstokes.h" 14*2eb7bf1fSJames Wright #include "../problems/stg_shur14.h" 15*2eb7bf1fSJames Wright 16*2eb7bf1fSJames Wright PetscErrorCode SetupStrongSTG_Ceed(Ceed ceed, CeedData ceed_data, DM dm, AppCtx app_ctx, ProblemData *problem, SimpleBC bc, Physics phys, 17*2eb7bf1fSJames Wright CeedInt Q_sur, CeedInt q_data_size_sur, CeedOperator op_strong_bc) { 18*2eb7bf1fSJames Wright CeedInt num_comp_x = problem->dim, num_comp_q = 5, num_elem, elem_size, stg_data_size = 1; 19*2eb7bf1fSJames Wright CeedVector multiplicity, x_stored, scale_stored, q_data_sur, stg_data; 20*2eb7bf1fSJames Wright CeedBasis basis_x_to_q_sur; 21*2eb7bf1fSJames Wright CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_x_stored, elem_restr_scale, elem_restr_qd_sur, elem_restr_stgdata; 22*2eb7bf1fSJames Wright CeedQFunction qf_setup, qf_strongbc, qf_stgdata; 23*2eb7bf1fSJames Wright CeedOperator op_setup, op_strong_bc_sub, op_setup_sur, op_stgdata; 24*2eb7bf1fSJames Wright DMLabel domain_label; 25*2eb7bf1fSJames Wright 26*2eb7bf1fSJames Wright PetscFunctionBeginUser; 27*2eb7bf1fSJames Wright PetscCall(DMGetLabel(dm, "Face Sets", &domain_label)); 28*2eb7bf1fSJames Wright 29*2eb7bf1fSJames Wright // Basis 30*2eb7bf1fSJames Wright CeedInt height = 1; 31*2eb7bf1fSJames Wright PetscCall(CeedBasisCreateProjection(ceed_data->basis_x_sur, ceed_data->basis_q_sur, &basis_x_to_q_sur)); 32*2eb7bf1fSJames Wright 33*2eb7bf1fSJames Wright // Setup QFunction 34*2eb7bf1fSJames Wright CeedQFunctionCreateInterior(ceed, 1, SetupStrongBC, SetupStrongBC_loc, &qf_setup); 35*2eb7bf1fSJames Wright CeedQFunctionAddInput(qf_setup, "x", num_comp_x, CEED_EVAL_INTERP); 36*2eb7bf1fSJames Wright CeedQFunctionAddInput(qf_setup, "multiplicity", num_comp_q, CEED_EVAL_NONE); 37*2eb7bf1fSJames Wright CeedQFunctionAddOutput(qf_setup, "x stored", num_comp_x, CEED_EVAL_NONE); 38*2eb7bf1fSJames Wright CeedQFunctionAddOutput(qf_setup, "scale", 1, CEED_EVAL_NONE); 39*2eb7bf1fSJames Wright 40*2eb7bf1fSJames Wright // Setup STG Setup QFunction 41*2eb7bf1fSJames Wright PetscCall(SetupStrongSTG_PreProcessing(ceed, problem, num_comp_x, stg_data_size, q_data_size_sur, &qf_stgdata)); 42*2eb7bf1fSJames Wright 43*2eb7bf1fSJames Wright // Compute contribution on each boundary face 44*2eb7bf1fSJames Wright for (CeedInt i = 0; i < bc->num_inflow; i++) { 45*2eb7bf1fSJames Wright // -- Restrictions 46*2eb7bf1fSJames Wright PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, bc->inflows[i], Q_sur, q_data_size_sur, &elem_restr_q_sur, &elem_restr_x_sur, 47*2eb7bf1fSJames Wright &elem_restr_qd_sur)); 48*2eb7bf1fSJames Wright CeedElemRestrictionCreateVector(elem_restr_q_sur, &multiplicity, NULL); 49*2eb7bf1fSJames Wright CeedElemRestrictionGetMultiplicity(elem_restr_q_sur, multiplicity); 50*2eb7bf1fSJames Wright CeedElemRestrictionGetNumElements(elem_restr_q_sur, &num_elem); 51*2eb7bf1fSJames Wright CeedElemRestrictionGetElementSize(elem_restr_q_sur, &elem_size); 52*2eb7bf1fSJames Wright CeedElemRestrictionCreateStrided(ceed, num_elem, elem_size, num_comp_x, num_elem * elem_size * num_comp_x, CEED_STRIDES_BACKEND, 53*2eb7bf1fSJames Wright &elem_restr_x_stored); 54*2eb7bf1fSJames Wright CeedElemRestrictionCreateVector(elem_restr_x_stored, &x_stored, NULL); 55*2eb7bf1fSJames Wright 56*2eb7bf1fSJames Wright CeedElemRestrictionCreateStrided(ceed, num_elem, elem_size, 1, num_elem * elem_size, CEED_STRIDES_BACKEND, &elem_restr_scale); 57*2eb7bf1fSJames Wright CeedElemRestrictionCreateVector(elem_restr_scale, &scale_stored, NULL); 58*2eb7bf1fSJames Wright 59*2eb7bf1fSJames Wright CeedElemRestrictionCreateStrided(ceed, num_elem, elem_size, stg_data_size, num_elem * elem_size, CEED_STRIDES_BACKEND, &elem_restr_stgdata); 60*2eb7bf1fSJames Wright CeedElemRestrictionCreateVector(elem_restr_stgdata, &stg_data, NULL); 61*2eb7bf1fSJames Wright 62*2eb7bf1fSJames Wright CeedVectorCreate(ceed, q_data_size_sur * num_elem * elem_size, &q_data_sur); 63*2eb7bf1fSJames Wright 64*2eb7bf1fSJames Wright // -- Setup Operator 65*2eb7bf1fSJames Wright CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup); 66*2eb7bf1fSJames Wright CeedOperatorSetName(op_setup, "surface geometric data"); 67*2eb7bf1fSJames Wright CeedOperatorSetField(op_setup, "x", elem_restr_x_sur, basis_x_to_q_sur, CEED_VECTOR_ACTIVE); 68*2eb7bf1fSJames Wright CeedOperatorSetField(op_setup, "multiplicity", elem_restr_q_sur, CEED_BASIS_COLLOCATED, multiplicity); 69*2eb7bf1fSJames Wright CeedOperatorSetField(op_setup, "x stored", elem_restr_x_stored, CEED_BASIS_COLLOCATED, x_stored); 70*2eb7bf1fSJames Wright CeedOperatorSetField(op_setup, "scale", elem_restr_scale, CEED_BASIS_COLLOCATED, scale_stored); 71*2eb7bf1fSJames Wright 72*2eb7bf1fSJames Wright // -- Compute geometric factors 73*2eb7bf1fSJames Wright CeedOperatorApply(op_setup, ceed_data->x_coord, CEED_VECTOR_NONE, CEED_REQUEST_IMMEDIATE); 74*2eb7bf1fSJames Wright 75*2eb7bf1fSJames Wright // -- Compute QData for the surface 76*2eb7bf1fSJames Wright CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur); 77*2eb7bf1fSJames Wright CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, ceed_data->basis_xc_sur, CEED_VECTOR_ACTIVE); 78*2eb7bf1fSJames Wright CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_xc_sur, CEED_VECTOR_NONE); 79*2eb7bf1fSJames Wright CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_sur, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 80*2eb7bf1fSJames Wright 81*2eb7bf1fSJames Wright CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE); 82*2eb7bf1fSJames Wright 83*2eb7bf1fSJames Wright // -- Compute STGData 84*2eb7bf1fSJames Wright CeedOperatorCreate(ceed, qf_stgdata, NULL, NULL, &op_stgdata); 85*2eb7bf1fSJames Wright CeedOperatorSetField(op_stgdata, "surface qdata", elem_restr_qd_sur, CEED_BASIS_COLLOCATED, q_data_sur); 86*2eb7bf1fSJames Wright CeedOperatorSetField(op_stgdata, "x", elem_restr_x_stored, CEED_BASIS_COLLOCATED, x_stored); 87*2eb7bf1fSJames Wright CeedOperatorSetField(op_stgdata, "stg data", elem_restr_stgdata, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 88*2eb7bf1fSJames Wright CeedOperatorSetNumQuadraturePoints(op_stgdata, elem_size); 89*2eb7bf1fSJames Wright 90*2eb7bf1fSJames Wright CeedOperatorApply(op_stgdata, CEED_VECTOR_NONE, stg_data, CEED_REQUEST_IMMEDIATE); 91*2eb7bf1fSJames Wright 92*2eb7bf1fSJames Wright // -- Setup BC QFunctions 93*2eb7bf1fSJames Wright SetupStrongSTG_QF(ceed, problem, num_comp_x, num_comp_q, stg_data_size, q_data_size_sur, &qf_strongbc); 94*2eb7bf1fSJames Wright CeedOperatorCreate(ceed, qf_strongbc, NULL, NULL, &op_strong_bc_sub); 95*2eb7bf1fSJames Wright CeedOperatorSetName(op_strong_bc_sub, "Strong STG"); 96*2eb7bf1fSJames Wright 97*2eb7bf1fSJames Wright CeedOperatorSetField(op_strong_bc_sub, "surface qdata", elem_restr_qd_sur, CEED_BASIS_COLLOCATED, q_data_sur); 98*2eb7bf1fSJames Wright CeedOperatorSetField(op_strong_bc_sub, "x", elem_restr_x_stored, CEED_BASIS_COLLOCATED, x_stored); 99*2eb7bf1fSJames Wright CeedOperatorSetField(op_strong_bc_sub, "scale", elem_restr_scale, CEED_BASIS_COLLOCATED, scale_stored); 100*2eb7bf1fSJames Wright CeedOperatorSetField(op_strong_bc_sub, "stg data", elem_restr_stgdata, CEED_BASIS_COLLOCATED, stg_data); 101*2eb7bf1fSJames Wright CeedOperatorSetField(op_strong_bc_sub, "q", elem_restr_q_sur, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 102*2eb7bf1fSJames Wright CeedOperatorSetNumQuadraturePoints(op_strong_bc_sub, elem_size); 103*2eb7bf1fSJames Wright 104*2eb7bf1fSJames Wright // -- Add to composite operator 105*2eb7bf1fSJames Wright CeedCompositeOperatorAddSub(op_strong_bc, op_strong_bc_sub); 106*2eb7bf1fSJames Wright 107*2eb7bf1fSJames Wright CeedVectorDestroy(&q_data_sur); 108*2eb7bf1fSJames Wright CeedVectorDestroy(&multiplicity); 109*2eb7bf1fSJames Wright CeedVectorDestroy(&x_stored); 110*2eb7bf1fSJames Wright CeedVectorDestroy(&scale_stored); 111*2eb7bf1fSJames Wright CeedVectorDestroy(&stg_data); 112*2eb7bf1fSJames Wright CeedElemRestrictionDestroy(&elem_restr_x_sur); 113*2eb7bf1fSJames Wright CeedElemRestrictionDestroy(&elem_restr_q_sur); 114*2eb7bf1fSJames Wright CeedElemRestrictionDestroy(&elem_restr_qd_sur); 115*2eb7bf1fSJames Wright CeedElemRestrictionDestroy(&elem_restr_x_stored); 116*2eb7bf1fSJames Wright CeedElemRestrictionDestroy(&elem_restr_scale); 117*2eb7bf1fSJames Wright CeedElemRestrictionDestroy(&elem_restr_stgdata); 118*2eb7bf1fSJames Wright CeedQFunctionDestroy(&qf_strongbc); 119*2eb7bf1fSJames Wright CeedQFunctionDestroy(&qf_stgdata); 120*2eb7bf1fSJames Wright CeedOperatorDestroy(&op_setup_sur); 121*2eb7bf1fSJames Wright CeedOperatorDestroy(&op_strong_bc_sub); 122*2eb7bf1fSJames Wright CeedOperatorDestroy(&op_setup); 123*2eb7bf1fSJames Wright CeedOperatorDestroy(&op_stgdata); 124*2eb7bf1fSJames Wright } 125*2eb7bf1fSJames Wright 126*2eb7bf1fSJames Wright CeedOperatorGetContextFieldLabel(op_strong_bc, "solution time", &phys->stg_solution_time_label); 127*2eb7bf1fSJames Wright 128*2eb7bf1fSJames Wright CeedBasisDestroy(&basis_x_to_q_sur); 129*2eb7bf1fSJames Wright CeedQFunctionDestroy(&qf_setup); 130*2eb7bf1fSJames Wright 131*2eb7bf1fSJames Wright PetscFunctionReturn(0); 132*2eb7bf1fSJames Wright } 133*2eb7bf1fSJames Wright 134*2eb7bf1fSJames Wright PetscErrorCode DMPlexInsertBoundaryValues_StrongBCCeed(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM, 135*2eb7bf1fSJames Wright Vec cell_geom_FVM, Vec grad_FVM) { 136*2eb7bf1fSJames Wright Vec boundary_mask; 137*2eb7bf1fSJames Wright User user; 138*2eb7bf1fSJames Wright 139*2eb7bf1fSJames Wright PetscFunctionBeginUser; 140*2eb7bf1fSJames Wright PetscCall(DMGetApplicationContext(dm, &user)); 141*2eb7bf1fSJames Wright 142*2eb7bf1fSJames Wright if (user->phys->stg_solution_time_label) { 143*2eb7bf1fSJames Wright CeedOperatorSetContextDouble(user->op_strong_bc_ctx->op, user->phys->stg_solution_time_label, &time); 144*2eb7bf1fSJames Wright } 145*2eb7bf1fSJames Wright 146*2eb7bf1fSJames Wright // Mask Strong BC entries 147*2eb7bf1fSJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 148*2eb7bf1fSJames Wright PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask)); 149*2eb7bf1fSJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 150*2eb7bf1fSJames Wright 151*2eb7bf1fSJames Wright PetscCall(ApplyAddCeedOperatorLocalToLocal(NULL, Q_loc, user->op_strong_bc_ctx)); 152*2eb7bf1fSJames Wright 153*2eb7bf1fSJames Wright PetscFunctionReturn(0); 154*2eb7bf1fSJames Wright } 155*2eb7bf1fSJames Wright 156*2eb7bf1fSJames Wright PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc, CeedInt Q_sur, 157*2eb7bf1fSJames Wright CeedInt q_data_size_sur) { 158*2eb7bf1fSJames Wright CeedOperator op_strong_bc; 159*2eb7bf1fSJames Wright 160*2eb7bf1fSJames Wright PetscFunctionBeginUser; 161*2eb7bf1fSJames Wright { 162*2eb7bf1fSJames Wright Vec boundary_mask, global_vec; 163*2eb7bf1fSJames Wright 164*2eb7bf1fSJames Wright PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask)); 165*2eb7bf1fSJames Wright PetscCall(DMGetGlobalVector(dm, &global_vec)); 166*2eb7bf1fSJames Wright PetscCall(VecZeroEntries(boundary_mask)); 167*2eb7bf1fSJames Wright PetscCall(VecSet(global_vec, 1.0)); 168*2eb7bf1fSJames Wright PetscCall(DMGlobalToLocal(dm, global_vec, INSERT_VALUES, boundary_mask)); 169*2eb7bf1fSJames Wright PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask)); 170*2eb7bf1fSJames Wright PetscCall(DMRestoreGlobalVector(dm, &global_vec)); 171*2eb7bf1fSJames Wright } 172*2eb7bf1fSJames Wright 173*2eb7bf1fSJames Wright CeedCompositeOperatorCreate(ceed, &op_strong_bc); 174*2eb7bf1fSJames Wright { 175*2eb7bf1fSJames Wright PetscBool use_strongstg = PETSC_FALSE; 176*2eb7bf1fSJames Wright PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL)); 177*2eb7bf1fSJames Wright 178*2eb7bf1fSJames Wright if (use_strongstg) { 179*2eb7bf1fSJames Wright PetscCall(SetupStrongSTG_Ceed(ceed, ceed_data, dm, app_ctx, problem, bc, user->phys, Q_sur, q_data_size_sur, op_strong_bc)); 180*2eb7bf1fSJames Wright } 181*2eb7bf1fSJames Wright } 182*2eb7bf1fSJames Wright 183*2eb7bf1fSJames Wright PetscCall(OperatorApplyContextCreate(NULL, NULL, ceed, op_strong_bc, CEED_VECTOR_NONE, NULL, NULL, NULL, &user->op_strong_bc_ctx)); 184*2eb7bf1fSJames Wright 185*2eb7bf1fSJames Wright PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_StrongBCCeed)); 186*2eb7bf1fSJames Wright PetscFunctionReturn(0); 187*2eb7bf1fSJames Wright } 188