xref: /honee/src/strong_boundary_conditions.c (revision 21ba7ba458ca7fbd5416be1f0b0586f2885778ca)
1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
22eb7bf1fSJames Wright // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
32eb7bf1fSJames Wright //
42eb7bf1fSJames Wright // SPDX-License-Identifier: BSD-2-Clause
52eb7bf1fSJames Wright //
62eb7bf1fSJames Wright // This file is part of CEED:  http://github.com/ceed
72eb7bf1fSJames Wright 
82eb7bf1fSJames Wright #include "../qfunctions/strong_boundary_conditions.h"
92eb7bf1fSJames Wright 
102eb7bf1fSJames Wright #include <ceed.h>
112eb7bf1fSJames Wright #include <petscdmplex.h>
122eb7bf1fSJames Wright 
13149fb536SJames Wright #include <navierstokes.h>
142eb7bf1fSJames Wright #include "../problems/stg_shur14.h"
152eb7bf1fSJames Wright 
16991aef52SJames Wright PetscErrorCode SetupStrongSTG_Ceed(Ceed ceed, CeedData ceed_data, DM dm, ProblemData problem, SimpleBC bc, Physics phys, CeedOperator op_strong_bc) {
17*21ba7ba4SJames Wright   CeedInt             num_comp_x = problem->dim, num_comp_q = 5, stg_data_size = 1, dXdx_size;
186f188493SJames Wright   CeedVector          multiplicity, x_stored, scale_stored, stg_data, dXdx;
19*21ba7ba4SJames Wright   CeedBasis           basis_x_to_q_face;
20*21ba7ba4SJames Wright   CeedElemRestriction elem_restr_x_face, elem_restr_q_face, elem_restr_x_stored, elem_restr_scale, elem_restr_stgdata, elem_restr_dXdx;
212eb7bf1fSJames Wright   CeedQFunction       qf_setup, qf_strongbc, qf_stgdata;
226f188493SJames Wright   CeedOperator        op_setup, op_strong_bc_sub, op_stgdata;
232eb7bf1fSJames Wright   DMLabel             domain_label;
24*21ba7ba4SJames Wright   const PetscInt      dm_field = 0, height_face = 1, height_cell = 0;
25*21ba7ba4SJames Wright   PetscInt dim;
26*21ba7ba4SJames Wright   MPI_Comm            comm = PetscObjectComm((PetscObject)dm);
272eb7bf1fSJames Wright 
282eb7bf1fSJames Wright   PetscFunctionBeginUser;
292eb7bf1fSJames Wright   PetscCall(DMGetLabel(dm, "Face Sets", &domain_label));
30*21ba7ba4SJames Wright   PetscCall(DMGetDimension(dm, &dim));
31*21ba7ba4SJames Wright   dXdx_size = num_comp_x * (dim - height_cell);
322eb7bf1fSJames Wright 
33866f9b4aSJames Wright   {  // Basis
34*21ba7ba4SJames Wright     CeedBasis basis_x_face, basis_q_face;
35866f9b4aSJames Wright     DM        dm_coord;
36866f9b4aSJames Wright 
37866f9b4aSJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
38866f9b4aSJames Wright     DMLabel  label       = NULL;
39866f9b4aSJames Wright     PetscInt label_value = 0;
40*21ba7ba4SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm, label, label_value, height_face, dm_field, &basis_q_face));
41*21ba7ba4SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, label_value, height_face, dm_field, &basis_x_face));
42866f9b4aSJames Wright 
43*21ba7ba4SJames Wright     PetscCallCeed(ceed, CeedBasisCreateProjection(basis_x_face, basis_q_face, &basis_x_to_q_face));
44866f9b4aSJames Wright 
45*21ba7ba4SJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_q_face));
46*21ba7ba4SJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_face));
47866f9b4aSJames Wright   }
482eb7bf1fSJames Wright 
492eb7bf1fSJames Wright   // Setup QFunction
50b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, SetupStrongBC, SetupStrongBC_loc, &qf_setup));
51b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "x", num_comp_x, CEED_EVAL_INTERP));
52*21ba7ba4SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "dxdX", dXdx_size, CEED_EVAL_GRAD));
53b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup, "multiplicity", num_comp_q, CEED_EVAL_NONE));
54b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "x stored", num_comp_x, CEED_EVAL_NONE));
55b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "scale", 1, CEED_EVAL_NONE));
566f188493SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup, "dXdx", dXdx_size, CEED_EVAL_NONE));
572eb7bf1fSJames Wright 
5851a85e6cSJames Wright   // Setup STG QFunctions
5942454adaSJames Wright   PetscCall(SetupStrongStg_PreProcessing(ceed, problem, num_comp_x, stg_data_size, dXdx_size, &qf_stgdata));
6051a85e6cSJames Wright   PetscCall(SetupStrongStg_QF(ceed, problem, num_comp_x, num_comp_q, stg_data_size, dXdx_size, &qf_strongbc));
612eb7bf1fSJames Wright 
622eb7bf1fSJames Wright   // Compute contribution on each boundary face
632eb7bf1fSJames Wright   for (CeedInt i = 0; i < bc->num_inflow; i++) {
64*21ba7ba4SJames Wright     PetscInt num_face_orientations;  // Ignored for now
65*21ba7ba4SJames Wright     char    *face_orientation_label_name;
66*21ba7ba4SJames Wright     PetscCall(DMPlexCreateFaceLabel(dm, bc->inflows[i], &num_face_orientations, &face_orientation_label_name));
6767263decSKenneth E. Jansen 
68*21ba7ba4SJames Wright     IS              is_orientation_values;
69*21ba7ba4SJames Wright     DMLabel         face_orientation_label;
70*21ba7ba4SJames Wright     PetscInt        num_orientations_local;
71*21ba7ba4SJames Wright     PetscInt        minmax_orientation_values[2];
72*21ba7ba4SJames Wright     const PetscInt *orientations_local;
73*21ba7ba4SJames Wright     PetscCall(DMGetLabel(dm, face_orientation_label_name, &face_orientation_label));
74*21ba7ba4SJames Wright     PetscCall(DMLabelGetValueIS(face_orientation_label, &is_orientation_values));
75*21ba7ba4SJames Wright     PetscCall(ISSort(is_orientation_values));
76*21ba7ba4SJames Wright     PetscCall(ISGetLocalSize(is_orientation_values, &num_orientations_local));
77*21ba7ba4SJames Wright     PetscCall(ISGetIndices(is_orientation_values, &orientations_local));
78*21ba7ba4SJames Wright     {
79*21ba7ba4SJames Wright       const PetscInt minmax_orientation_values_loc[2] = {orientations_local[0], orientations_local[num_orientations_local - 1]};
80*21ba7ba4SJames Wright 
81*21ba7ba4SJames Wright       PetscCall(PetscGlobalMinMaxInt(comm, minmax_orientation_values_loc, minmax_orientation_values));
82*21ba7ba4SJames Wright     }
83*21ba7ba4SJames Wright 
84*21ba7ba4SJames Wright     for (PetscInt orientation = minmax_orientation_values[0]; orientation <= minmax_orientation_values[1]; orientation++) {
85*21ba7ba4SJames Wright       CeedBasis basis_x_to_q_cell;
86*21ba7ba4SJames Wright       CeedElemRestriction elem_restr_x_cell;
87*21ba7ba4SJames Wright       {  // Skip orientation if no process has this face with this orientation
88*21ba7ba4SJames Wright         PetscInt location_local, location_global = -1;
89*21ba7ba4SJames Wright 
90*21ba7ba4SJames Wright         PetscCall(PetscFindInt(orientation, num_orientations_local, orientations_local, &location_local));
91*21ba7ba4SJames Wright         PetscCall(MPIU_Allreduce(&location_local, &location_global, 1, MPIU_INT, MPI_MAX, comm));
92*21ba7ba4SJames Wright         if (location_global < 0) continue;
93*21ba7ba4SJames Wright       }
94*21ba7ba4SJames Wright 
95*21ba7ba4SJames Wright       {
96*21ba7ba4SJames Wright         CeedBasis basis_x_cell_to_face, basis_q_face;
97*21ba7ba4SJames Wright 
98*21ba7ba4SJames Wright         PetscCall(DMPlexCeedBasisCellToFaceCoordinateCreate(ceed, dm, face_orientation_label, orientation, orientation, &basis_x_cell_to_face));
99*21ba7ba4SJames Wright         PetscCall(CreateBasisFromPlex(ceed, dm, face_orientation_label, orientation, height_face, dm_field, &basis_q_face));
100*21ba7ba4SJames Wright         PetscCallCeed(ceed, CeedBasisCreateProjection(basis_x_cell_to_face, basis_q_face, &basis_x_to_q_cell));
101*21ba7ba4SJames Wright         PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_cell_to_face));
102*21ba7ba4SJames Wright         PetscCallCeed(ceed, CeedBasisDestroy(&basis_q_face));
103*21ba7ba4SJames Wright       }
104*21ba7ba4SJames Wright 
105*21ba7ba4SJames Wright       PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, bc->inflows[i], height_face, dm_field, &elem_restr_q_face));
106*21ba7ba4SJames Wright       PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, bc->inflows[i], height_face, &elem_restr_x_face));
107*21ba7ba4SJames Wright       PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, face_orientation_label, orientation, height_cell, &elem_restr_x_cell));
108*21ba7ba4SJames Wright       PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_q_face, &multiplicity, NULL));
109*21ba7ba4SJames Wright       PetscCallCeed(ceed, CeedElemRestrictionGetMultiplicity(elem_restr_q_face, multiplicity));
110*21ba7ba4SJames Wright 
111*21ba7ba4SJames Wright       PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, dm, domain_label, bc->inflows[i], height_face, num_comp_x, &elem_restr_x_stored));
112*21ba7ba4SJames Wright       PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, dm, domain_label, bc->inflows[i], height_face, 1, &elem_restr_scale));
113*21ba7ba4SJames Wright       PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, dm, domain_label, bc->inflows[i], height_face, stg_data_size, &elem_restr_stgdata));
114*21ba7ba4SJames Wright       PetscCall(DMPlexCeedElemRestrictionCollocatedCreate(ceed, dm, domain_label, bc->inflows[i], height_face, dXdx_size, &elem_restr_dXdx));
115b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_x_stored, &x_stored, NULL));
116b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_scale, &scale_stored, NULL));
117b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_stgdata, &stg_data, NULL));
1186f188493SJames Wright       PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_dXdx, &dXdx, NULL));
1192eb7bf1fSJames Wright 
1202eb7bf1fSJames Wright       // -- Setup Operator
121b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup, NULL, NULL, &op_setup));
1226f188493SJames Wright       PetscCallCeed(ceed, CeedOperatorSetName(op_setup, "Precomputed data for strong boundary conditions"));
123*21ba7ba4SJames Wright       PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "x", elem_restr_x_face, basis_x_to_q_face, CEED_VECTOR_ACTIVE));
124*21ba7ba4SJames Wright       PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "dxdX", elem_restr_x_cell, basis_x_to_q_cell, CEED_VECTOR_ACTIVE));
125*21ba7ba4SJames Wright       PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "multiplicity", elem_restr_q_face, CEED_BASIS_NONE, multiplicity));
12658e1cbfdSJeremy L Thompson       PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "x stored", elem_restr_x_stored, CEED_BASIS_NONE, x_stored));
12758e1cbfdSJeremy L Thompson       PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "scale", elem_restr_scale, CEED_BASIS_NONE, scale_stored));
12858e1cbfdSJeremy L Thompson       PetscCallCeed(ceed, CeedOperatorSetField(op_setup, "dXdx", elem_restr_dXdx, CEED_BASIS_NONE, dXdx));
1292eb7bf1fSJames Wright 
1302eb7bf1fSJames Wright       // -- Compute geometric factors
131b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedOperatorApply(op_setup, ceed_data->x_coord, CEED_VECTOR_NONE, CEED_REQUEST_IMMEDIATE));
1322eb7bf1fSJames Wright 
1332eb7bf1fSJames Wright       // -- Compute STGData
134b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_stgdata, NULL, NULL, &op_stgdata));
13558e1cbfdSJeremy L Thompson       PetscCallCeed(ceed, CeedOperatorSetField(op_stgdata, "dXdx", elem_restr_dXdx, CEED_BASIS_NONE, dXdx));
13658e1cbfdSJeremy L Thompson       PetscCallCeed(ceed, CeedOperatorSetField(op_stgdata, "x", elem_restr_x_stored, CEED_BASIS_NONE, x_stored));
13758e1cbfdSJeremy L Thompson       PetscCallCeed(ceed, CeedOperatorSetField(op_stgdata, "stg data", elem_restr_stgdata, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
1382eb7bf1fSJames Wright 
139b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedOperatorApply(op_stgdata, CEED_VECTOR_NONE, stg_data, CEED_REQUEST_IMMEDIATE));
1402eb7bf1fSJames Wright 
141*21ba7ba4SJames Wright       // -- Setup BC QFunctions
142b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_strongbc, NULL, NULL, &op_strong_bc_sub));
143b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedOperatorSetName(op_strong_bc_sub, "Strong STG"));
1442eb7bf1fSJames Wright 
14558e1cbfdSJeremy L Thompson       PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "dXdx", elem_restr_dXdx, CEED_BASIS_NONE, dXdx));
14658e1cbfdSJeremy L Thompson       PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "x", elem_restr_x_stored, CEED_BASIS_NONE, x_stored));
14758e1cbfdSJeremy L Thompson       PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "scale", elem_restr_scale, CEED_BASIS_NONE, scale_stored));
14858e1cbfdSJeremy L Thompson       PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "stg data", elem_restr_stgdata, CEED_BASIS_NONE, stg_data));
149*21ba7ba4SJames Wright       PetscCallCeed(ceed, CeedOperatorSetField(op_strong_bc_sub, "q", elem_restr_q_face, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
1502eb7bf1fSJames Wright 
1512eb7bf1fSJames Wright       // -- Add to composite operator
152b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_strong_bc, op_strong_bc_sub));
1532eb7bf1fSJames Wright 
154b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedVectorDestroy(&multiplicity));
155b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedVectorDestroy(&x_stored));
156b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedVectorDestroy(&scale_stored));
157b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedVectorDestroy(&stg_data));
1586f188493SJames Wright       PetscCallCeed(ceed, CeedVectorDestroy(&dXdx));
159*21ba7ba4SJames Wright       PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q_cell));
160*21ba7ba4SJames Wright       PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_face));
161*21ba7ba4SJames Wright       PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_cell));
162*21ba7ba4SJames Wright       PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_face));
163b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_stored));
164b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_scale));
165b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_stgdata));
1666f188493SJames Wright       PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_dXdx));
167b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedOperatorDestroy(&op_strong_bc_sub));
168b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup));
169b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedOperatorDestroy(&op_stgdata));
1702eb7bf1fSJames Wright     }
171*21ba7ba4SJames Wright   }
1722eb7bf1fSJames Wright 
173b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_strong_bc, "solution time", &phys->stg_solution_time_label));
1742eb7bf1fSJames Wright 
175*21ba7ba4SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_to_q_face));
17651a85e6cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_strongbc));
17751a85e6cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_stgdata));
178b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup));
179d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1802eb7bf1fSJames Wright }
1812eb7bf1fSJames Wright 
1822eb7bf1fSJames Wright PetscErrorCode DMPlexInsertBoundaryValues_StrongBCCeed(DM dm, PetscBool insert_essential, Vec Q_loc, PetscReal time, Vec face_geom_FVM,
1832eb7bf1fSJames Wright                                                        Vec cell_geom_FVM, Vec grad_FVM) {
1842eb7bf1fSJames Wright   Vec  boundary_mask;
1852eb7bf1fSJames Wright   User user;
1862eb7bf1fSJames Wright 
1872eb7bf1fSJames Wright   PetscFunctionBeginUser;
1882eb7bf1fSJames Wright   PetscCall(DMGetApplicationContext(dm, &user));
1892eb7bf1fSJames Wright 
1902eb7bf1fSJames Wright   if (user->phys->stg_solution_time_label) {
191b4c37c5cSJames Wright     PetscCallCeed(user->ceed, CeedOperatorSetContextDouble(user->op_strong_bc_ctx->op, user->phys->stg_solution_time_label, &time));
1922eb7bf1fSJames Wright   }
1932eb7bf1fSJames Wright 
1942eb7bf1fSJames Wright   // Mask Strong BC entries
1952eb7bf1fSJames Wright   PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
1962eb7bf1fSJames Wright   PetscCall(VecPointwiseMult(Q_loc, Q_loc, boundary_mask));
1972eb7bf1fSJames Wright   PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
1982eb7bf1fSJames Wright 
1992eb7bf1fSJames Wright   PetscCall(ApplyAddCeedOperatorLocalToLocal(NULL, Q_loc, user->op_strong_bc_ctx));
200d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2012eb7bf1fSJames Wright }
2022eb7bf1fSJames Wright 
203991aef52SJames Wright PetscErrorCode SetupStrongBC_Ceed(Ceed ceed, CeedData ceed_data, DM dm, User user, ProblemData problem, SimpleBC bc) {
2042eb7bf1fSJames Wright   CeedOperator op_strong_bc;
2052eb7bf1fSJames Wright 
2062eb7bf1fSJames Wright   PetscFunctionBeginUser;
2072eb7bf1fSJames Wright   {
2082eb7bf1fSJames Wright     Vec boundary_mask, global_vec;
2092eb7bf1fSJames Wright 
2102eb7bf1fSJames Wright     PetscCall(DMGetNamedLocalVector(dm, "boundary mask", &boundary_mask));
2112eb7bf1fSJames Wright     PetscCall(DMGetGlobalVector(dm, &global_vec));
2122eb7bf1fSJames Wright     PetscCall(VecZeroEntries(boundary_mask));
2132eb7bf1fSJames Wright     PetscCall(VecSet(global_vec, 1.0));
2142eb7bf1fSJames Wright     PetscCall(DMGlobalToLocal(dm, global_vec, INSERT_VALUES, boundary_mask));
2152eb7bf1fSJames Wright     PetscCall(DMRestoreNamedLocalVector(dm, "boundary mask", &boundary_mask));
2162eb7bf1fSJames Wright     PetscCall(DMRestoreGlobalVector(dm, &global_vec));
2172eb7bf1fSJames Wright   }
2182eb7bf1fSJames Wright 
219b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_strong_bc));
2202eb7bf1fSJames Wright   {
2212eb7bf1fSJames Wright     PetscBool use_strongstg = PETSC_FALSE;
2222eb7bf1fSJames Wright     PetscCall(PetscOptionsGetBool(NULL, NULL, "-stg_strong", &use_strongstg, NULL));
2232eb7bf1fSJames Wright 
2242eb7bf1fSJames Wright     if (use_strongstg) {
2256f188493SJames Wright       PetscCall(SetupStrongSTG_Ceed(ceed, ceed_data, dm, problem, bc, user->phys, op_strong_bc));
2262eb7bf1fSJames Wright     }
2272eb7bf1fSJames Wright   }
2282eb7bf1fSJames Wright 
2292eb7bf1fSJames Wright   PetscCall(OperatorApplyContextCreate(NULL, NULL, ceed, op_strong_bc, CEED_VECTOR_NONE, NULL, NULL, NULL, &user->op_strong_bc_ctx));
2302eb7bf1fSJames Wright 
2312eb7bf1fSJames Wright   PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexInsertBoundaryValues_C", DMPlexInsertBoundaryValues_StrongBCCeed));
232d7b7c37aSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_strong_bc));
233d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2342eb7bf1fSJames Wright }
235