xref: /libCEED/examples/fluids/src/strong_boundary_conditions.c (revision d109957bbe678d426fc2302b20ca8bf13761f74a)
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, CeedInt q_data_size_sur,
17                                    CeedOperator op_strong_bc) {
18   CeedInt             num_comp_x = problem->dim, num_comp_q = 5, num_elem, elem_size, stg_data_size = 1;
19   CeedVector          multiplicity, x_stored, scale_stored, q_data_sur, stg_data;
20   CeedBasis           basis_x_to_q_sur;
21   CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_x_stored, elem_restr_scale, elem_restr_qd_sur, elem_restr_stgdata;
22   CeedQFunction       qf_setup, qf_strongbc, qf_stgdata;
23   CeedOperator        op_setup, op_strong_bc_sub, op_setup_sur, op_stgdata;
24   DMLabel             domain_label;
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, "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 
43   // Setup STG Setup QFunction
44   PetscCall(SetupStrongSTG_PreProcessing(ceed, problem, num_comp_x, stg_data_size, q_data_size_sur, &qf_stgdata));
45 
46   // Compute contribution on each boundary face
47   for (CeedInt i = 0; i < bc->num_inflow; i++) {
48     // -- Restrictions
49     PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, bc->inflows[i], 0, -1, -1, &elem_restr_q_sur, &elem_restr_x_sur, NULL));
50     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_q_sur, &multiplicity, NULL));
51     PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(elem_restr_q_sur, multiplicity));
52     PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_q_sur, &num_elem));
53     PetscCallCeed(ceed, CeedElemRestrictionGetElementSize(elem_restr_q_sur, &elem_size));
54     PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, bc->inflows[i], 0, elem_size, q_data_size_sur, NULL, NULL, &elem_restr_qd_sur));
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, CEED_STRIDES_BACKEND,
65                                                          &elem_restr_stgdata));
66     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_stgdata, &stg_data, NULL));
67 
68     PetscCallCeed(ceed, CeedVectorCreate(ceed, q_data_size_sur * num_elem * elem_size, &q_data_sur));
69 
70     // -- Setup Operator
71     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup));
72     PetscCallCeed(ceed, CeedOperatorSetName(op_setup, "surface geometric data"));
73     PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "x", elem_restr_x_sur, basis_x_to_q_sur, CEED_VECTOR_ACTIVE));
74     PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "multiplicity", elem_restr_q_sur, CEED_BASIS_COLLOCATED, multiplicity));
75     PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "x stored", elem_restr_x_stored, CEED_BASIS_COLLOCATED, x_stored));
76     PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "scale", elem_restr_scale, CEED_BASIS_COLLOCATED, scale_stored));
77 
78     // -- Compute geometric factors
79     PetscCallCeed(ceed, CeedOperatorApply(op_setup, ceed_data->x_coord, CEED_VECTOR_NONE, CEED_REQUEST_IMMEDIATE));
80 
81     // -- Compute QData for the surface
82     PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur));
83     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, ceed_data->basis_xc_sur, CEED_VECTOR_ACTIVE));
84     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_xc_sur, CEED_VECTOR_NONE));
85     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_sur, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
86 
87     PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE));
88 
89     // -- Compute STGData
90     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stgdata, NULL, NULL, &op_stgdata));
91     PetscCallCeed(ceed, CeedOperatorSetField(op_stgdata, "surface qdata", elem_restr_qd_sur, CEED_BASIS_COLLOCATED, q_data_sur));
92     PetscCallCeed(ceed, CeedOperatorSetField(op_stgdata, "x", elem_restr_x_stored, CEED_BASIS_COLLOCATED, x_stored));
93     PetscCallCeed(ceed, CeedOperatorSetField(op_stgdata, "stg data", elem_restr_stgdata, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
94 
95     PetscCallCeed(ceed, CeedOperatorApply(op_stgdata, CEED_VECTOR_NONE, stg_data, CEED_REQUEST_IMMEDIATE));
96 
97     // -- Setup BC QFunctions
98     SetupStrongSTG_QF(ceed, problem, num_comp_x, num_comp_q, stg_data_size, q_data_size_sur, &qf_strongbc);
99     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_strongbc, NULL, NULL, &op_strong_bc_sub));
100     PetscCallCeed(ceed, CeedOperatorSetName(op_strong_bc_sub, "Strong STG"));
101 
102     PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "surface qdata", elem_restr_qd_sur, CEED_BASIS_COLLOCATED, q_data_sur));
103     PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "x", elem_restr_x_stored, CEED_BASIS_COLLOCATED, x_stored));
104     PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "scale", elem_restr_scale, CEED_BASIS_COLLOCATED, scale_stored));
105     PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "stg data", elem_restr_stgdata, CEED_BASIS_COLLOCATED, stg_data));
106     PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "q", elem_restr_q_sur, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
107 
108     // -- Add to composite operator
109     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_strong_bc, op_strong_bc_sub));
110 
111     PetscCallCeed(ceed, CeedVectorDestroy(&q_data_sur));
112     PetscCallCeed(ceed, CeedVectorDestroy(&multiplicity));
113     PetscCallCeed(ceed, CeedVectorDestroy(&x_stored));
114     PetscCallCeed(ceed, CeedVectorDestroy(&scale_stored));
115     PetscCallCeed(ceed, CeedVectorDestroy(&stg_data));
116     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur));
117     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur));
118     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_sur));
119     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_stored));
120     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_scale));
121     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_stgdata));
122     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_strongbc));
123     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stgdata));
124     PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur));
125     PetscCallCeed(ceed, CeedOperatorDestroy(&op_strong_bc_sub));
126     PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup));
127     PetscCallCeed(ceed, CeedOperatorDestroy(&op_stgdata));
128   }
129 
130   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_strong_bc, "solution time", &phys->stg_solution_time_label));
131 
132   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q_sur));
133   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup));
134 
135   PetscFunctionReturn(PETSC_SUCCESS);
136 }
137 
138 PetscErrorCode DMPlexInsertBoundaryValues_StrongBCCeed(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM,
139                                                        Vec cell_geom_FVM, Vec grad_FVM) {
140   Vec  boundary_mask;
141   User user;
142 
143   PetscFunctionBeginUser;
144   PetscCall(DMGetApplicationContext(dm, &user));
145 
146   if (user->phys->stg_solution_time_label) {
147     PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user->op_strong_bc_ctx->op, user->phys->stg_solution_time_label, &time));
148   }
149 
150   // Mask Strong BC entries
151   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
152   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
153   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
154 
155   PetscCall(ApplyAddCeedOperatorLocalToLocal(NULL, Q_loc, user->op_strong_bc_ctx));
156 
157   PetscFunctionReturn(PETSC_SUCCESS);
158 }
159 
160 PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, ProblemData *problem, SimpleBC bc, CeedInt q_data_size_sur) {
161   CeedOperator op_strong_bc;
162 
163   PetscFunctionBeginUser;
164   {
165     Vec boundary_mask, global_vec;
166 
167     PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
168     PetscCall(DMGetGlobalVector(dm, &global_vec));
169     PetscCall(VecZeroEntries(boundary_mask));
170     PetscCall(VecSet(global_vec, 1.0));
171     PetscCall(DMGlobalToLocal(dm, global_vec, INSERT_VALUES, boundary_mask));
172     PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
173     PetscCall(DMRestoreGlobalVector(dm, &global_vec));
174   }
175 
176   PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_strong_bc));
177   {
178     PetscBool use_strongstg = PETSC_FALSE;
179     PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
180 
181     if (use_strongstg) {
182       PetscCall(SetupStrongSTG_Ceed(ceed, ceed_data, dm, problem, bc, user->phys, q_data_size_sur, op_strong_bc));
183     }
184   }
185 
186   PetscCall(OperatorApplyContextCreate(NULL, NULL, ceed, op_strong_bc, CEED_VECTOR_NONE, NULL, NULL, NULL, &user->op_strong_bc_ctx));
187 
188   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_StrongBCCeed));
189   PetscFunctionReturn(PETSC_SUCCESS);
190 }
191