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