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