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