xref: /libCEED/examples/fluids/src/setuplibceed.c (revision fd8f24fa4ef1026f43db0b90b0c68f0c3af48272)
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 /// @file
9 /// Setup libCEED for Navier-Stokes example using PETSc
10 
11 #include <ceed.h>
12 #include <petscdmplex.h>
13 
14 #include <petscds.h>
15 #include "../navierstokes.h"
16 
17 // Utility function to create local CEED restriction
18 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt label_value, PetscInt dm_field,
19                                          CeedElemRestriction *elem_restr) {
20   PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets_petsc;
21   CeedInt *elem_restr_offsets_ceed;
22 
23   PetscFunctionBeginUser;
24   PetscCall(
25       DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets_petsc));
26 
27   PetscCall(IntArrayP2C(num_elem * elem_size, &elem_restr_offsets_petsc, &elem_restr_offsets_ceed));
28   PetscCallCeed(ceed, CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES,
29                                                 elem_restr_offsets_ceed, elem_restr));
30   PetscCall(PetscFree(elem_restr_offsets_ceed));
31 
32   PetscFunctionReturn(PETSC_SUCCESS);
33 }
34 
35 // Utility function to get Ceed Restriction for each domain
36 PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, PetscInt label_value, PetscInt dm_field, CeedInt Q,
37                                        CeedInt q_data_size, CeedElemRestriction *elem_restr_q, CeedElemRestriction *elem_restr_x,
38                                        CeedElemRestriction *elem_restr_qd_i) {
39   DM                  dm_coord;
40   CeedInt             loc_num_elem;
41   PetscInt            dim;
42   CeedElemRestriction elem_restr_tmp;
43   PetscFunctionBeginUser;
44 
45   PetscCall(DMGetDimension(dm, &dim));
46   dim -= height;
47   PetscCall(CreateRestrictionFromPlex(ceed, dm, height, domain_label, label_value, dm_field, &elem_restr_tmp));
48   if (elem_restr_q) *elem_restr_q = elem_restr_tmp;
49   if (elem_restr_x) {
50     PetscCall(DMGetCellCoordinateDM(dm, &dm_coord));
51     if (!dm_coord) {
52       PetscCall(DMGetCoordinateDM(dm, &dm_coord));
53     }
54     PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL));
55     PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, height, domain_label, label_value, dm_field, elem_restr_x));
56   }
57   if (elem_restr_qd_i) {
58     PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_tmp, &loc_num_elem));
59     PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, loc_num_elem, Q, q_data_size, q_data_size * loc_num_elem * Q, CEED_STRIDES_BACKEND,
60                                                          elem_restr_qd_i));
61   }
62   if (!elem_restr_q) PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_tmp));
63   PetscFunctionReturn(PETSC_SUCCESS);
64 }
65 
66 PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height, CeedInt Q_sur,
67                                 CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian,
68                                 CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) {
69   CeedVector          q_data_sur, jac_data_sur;
70   CeedOperator        op_setup_sur, op_apply_bc, op_apply_bc_jacobian = NULL;
71   CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur;
72   CeedInt             num_qpts_sur;
73   PetscFunctionBeginUser;
74 
75   // --- Get number of quadrature points for the boundaries
76   PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(ceed_data->basis_q_sur, &num_qpts_sur));
77 
78   // ---- CEED Restriction
79   PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, label_value, 0, num_qpts_sur, q_data_size_sur, &elem_restr_q_sur,
80                                     &elem_restr_x_sur, &elem_restr_qd_i_sur));
81   if (jac_data_size_sur > 0) {
82     // State-dependent data will be passed from residual to Jacobian. This will be collocated.
83     PetscCall(
84         GetRestrictionForDomain(ceed, dm, height, domain_label, label_value, 0, num_qpts_sur, jac_data_size_sur, NULL, NULL, &elem_restr_jd_i_sur));
85     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL));
86   } else {
87     elem_restr_jd_i_sur = NULL;
88     jac_data_sur        = NULL;
89   }
90 
91   // ---- CEED Vector
92   CeedInt loc_num_elem_sur;
93   PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_q_sur, &loc_num_elem_sur));
94   PetscCallCeed(ceed, CeedVectorCreate(ceed, q_data_size_sur * loc_num_elem_sur * num_qpts_sur, &q_data_sur));
95 
96   // ---- CEED Operator
97   // ----- CEED Operator for Setup (geometric factors)
98   PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur));
99   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, ceed_data->basis_x_sur, CEED_VECTOR_ACTIVE));
100   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x_sur, CEED_VECTOR_NONE));
101   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
102 
103   // ----- CEED Operator for Physics
104   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc));
105   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
106   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
107   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur));
108   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord));
109   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
110   if (elem_restr_jd_i_sur)
111     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur));
112 
113   if (qf_apply_bc_jacobian) {
114     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian));
115     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
116     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
117     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur));
118     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord));
119     PetscCallCeed(ceed,
120                   CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur));
121     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
122   }
123 
124   // ----- Apply CEED operator for Setup
125   PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE));
126 
127   // ----- Apply Sub-Operator for Physics
128   PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply, op_apply_bc));
129   if (op_apply_bc_jacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_bc_jacobian));
130 
131   // ----- Cleanup
132   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_sur));
133   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data_sur));
134   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur));
135   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur));
136   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i_sur));
137   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i_sur));
138   PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur));
139   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc));
140   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc_jacobian));
141 
142   PetscFunctionReturn(PETSC_SUCCESS);
143 }
144 
145 // Utility function to create CEED Composite Operator for the entire domain
146 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol,
147                                        CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur,
148                                        CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) {
149   DMLabel domain_label;
150   PetscFunctionBeginUser;
151 
152   // Create Composite Operaters
153   PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_apply));
154   if (op_apply_ijacobian) PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_apply_ijacobian));
155 
156   // --Apply Sub-Operator for the volume
157   PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply, op_apply_vol));
158   if (op_apply_ijacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_ijacobian_vol));
159 
160   // -- Create Sub-Operator for in/outflow BCs
161   if (phys->has_neumann || 1) {
162     // --- Setup
163     PetscCall(DMGetLabel(dm, "Face Sets", &domain_label));
164 
165     // --- Create Sub-Operator for inflow boundaries
166     for (CeedInt i = 0; i < bc->num_inflow; i++) {
167       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->inflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
168                                  ceed_data->qf_apply_inflow, ceed_data->qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian));
169     }
170     // --- Create Sub-Operator for outflow boundaries
171     for (CeedInt i = 0; i < bc->num_outflow; i++) {
172       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
173                                  ceed_data->qf_apply_outflow, ceed_data->qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian));
174     }
175     // --- Create Sub-Operator for freestream boundaries
176     for (CeedInt i = 0; i < bc->num_freestream; i++) {
177       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
178                                  ceed_data->qf_apply_freestream, ceed_data->qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian));
179     }
180   }
181 
182   // ----- Get Context Labels for Operator
183   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(*op_apply, "solution time", &phys->solution_time_label));
184   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(*op_apply, "timestep size", &phys->timestep_size_label));
185 
186   PetscFunctionReturn(PETSC_SUCCESS);
187 }
188 
189 PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur,
190                                  PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian,
191                                  CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) {
192   PetscFunctionBeginUser;
193 
194   if (apply_bc.qfunction) {
195     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc));
196     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context));
197     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP));
198     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD));
199     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
200     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP));
201     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP));
202     if (jac_data_size_sur) PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
203   }
204   if (apply_bc_jacobian.qfunction) {
205     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian));
206     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context));
207     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP));
208     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD));
209     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
210     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP));
211     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
212     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP));
213   }
214   PetscFunctionReturn(PETSC_SUCCESS);
215 }
216 
217 // -----------------------------------------------------------------------------
218 // Convert DM field to DS field
219 // -----------------------------------------------------------------------------
220 PetscErrorCode DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, PetscInt *ds_field) {
221   PetscDS         ds;
222   IS              field_is;
223   const PetscInt *fields;
224   PetscInt        num_fields;
225 
226   PetscFunctionBeginUser;
227 
228   // Translate dm_field to ds_field
229   PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds, NULL));
230   PetscCall(ISGetIndices(field_is, &fields));
231   PetscCall(ISGetSize(field_is, &num_fields));
232   for (PetscInt i = 0; i < num_fields; i++) {
233     if (dm_field == fields[i]) {
234       *ds_field = i;
235       break;
236     }
237   }
238   PetscCall(ISRestoreIndices(field_is, &fields));
239 
240   if (*ds_field == -1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field);
241 
242   PetscFunctionReturn(PETSC_SUCCESS);
243 }
244 
245 // -----------------------------------------------------------------------------
246 // Utility function - convert from DMPolytopeType to CeedElemTopology
247 // -----------------------------------------------------------------------------
248 CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) {
249   switch (cell_type) {
250     case DM_POLYTOPE_TRIANGLE:
251       return CEED_TOPOLOGY_TRIANGLE;
252     case DM_POLYTOPE_QUADRILATERAL:
253       return CEED_TOPOLOGY_QUAD;
254     case DM_POLYTOPE_TETRAHEDRON:
255       return CEED_TOPOLOGY_TET;
256     case DM_POLYTOPE_HEXAHEDRON:
257       return CEED_TOPOLOGY_HEX;
258     default:
259       return 0;
260   }
261 }
262 
263 // -----------------------------------------------------------------------------
264 // Create libCEED Basis from PetscTabulation
265 // -----------------------------------------------------------------------------
266 PetscErrorCode BasisCreateFromTabulation(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt face, PetscFE fe,
267                                          PetscTabulation basis_tabulation, PetscQuadrature quadrature, CeedBasis *basis) {
268   PetscInt           first_point;
269   PetscInt           ids[1] = {label_value};
270   DMLabel            depth_label;
271   DMPolytopeType     cell_type;
272   CeedElemTopology   elem_topo;
273   PetscScalar       *q_points, *interp, *grad;
274   const PetscScalar *q_weights;
275   PetscDualSpace     dual_space;
276   PetscInt           num_dual_basis_vectors;
277   PetscInt           dim, num_comp, P, Q;
278 
279   PetscFunctionBeginUser;
280 
281   // General basis information
282   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
283   PetscCall(PetscFEGetNumComponents(fe, &num_comp));
284   PetscCall(PetscFEGetDualSpace(fe, &dual_space));
285   PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
286   P = num_dual_basis_vectors / num_comp;
287 
288   // Use depth label if no domain label present
289   if (!domain_label) {
290     PetscInt depth;
291 
292     PetscCall(DMPlexGetDepth(dm, &depth));
293     PetscCall(DMPlexGetDepthLabel(dm, &depth_label));
294     ids[0] = depth - height;
295   }
296 
297   // Get cell interp, grad, and quadrature data
298   PetscCall(DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label, 1, ids, height, &first_point, NULL));
299   PetscCall(DMPlexGetCellType(dm, first_point, &cell_type));
300   elem_topo = ElemTopologyP2C(cell_type);
301   if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMPlex topology not supported");
302   {
303     size_t             q_points_size;
304     const PetscScalar *q_points_petsc;
305     PetscInt           q_dim;
306 
307     PetscCall(PetscQuadratureGetData(quadrature, &q_dim, NULL, &Q, &q_points_petsc, &q_weights));
308     q_points_size = Q * dim * sizeof(CeedScalar);
309     PetscCall(PetscCalloc(q_points_size, &q_points));
310     for (PetscInt q = 0; q < Q; q++) {
311       for (PetscInt d = 0; d < q_dim; d++) q_points[q * dim + d] = q_points_petsc[q * q_dim + d];
312     }
313   }
314 
315   {  // Convert to libCEED orientation
316     PetscBool       is_simplex  = PETSC_FALSE;
317     IS              permutation = NULL;
318     const PetscInt *permutation_indices;
319 
320     PetscCall(DMPlexIsSimplex(dm, &is_simplex));
321     if (!is_simplex) {
322       PetscSection section;
323 
324       // -- Get permutation
325       PetscCall(DMGetLocalSection(dm, &section));
326       PetscCall(PetscSectionGetClosurePermutation(section, (PetscObject)dm, dim, num_comp * P, &permutation));
327       PetscCall(ISGetIndices(permutation, &permutation_indices));
328     }
329 
330     // -- Copy interp, grad matrices
331     PetscCall(PetscCalloc(P * Q * sizeof(CeedScalar), &interp));
332     PetscCall(PetscCalloc(P * Q * dim * sizeof(CeedScalar), &grad));
333     const CeedInt c = 0;
334     for (CeedInt q = 0; q < Q; q++) {
335       for (CeedInt p_ceed = 0; p_ceed < P; p_ceed++) {
336         CeedInt p_petsc = is_simplex ? (p_ceed * num_comp) : permutation_indices[p_ceed * num_comp];
337 
338         interp[q * P + p_ceed] = basis_tabulation->T[0][((face * Q + q) * P * num_comp + p_petsc) * num_comp + c];
339         for (CeedInt d = 0; d < dim; d++) {
340           grad[(d * Q + q) * P + p_ceed] = basis_tabulation->T[1][(((face * Q + q) * P * num_comp + p_petsc) * num_comp + c) * dim + d];
341         }
342       }
343     }
344 
345     // -- Cleanup
346     if (permutation) PetscCall(ISRestoreIndices(permutation, &permutation_indices));
347     PetscCall(ISDestroy(&permutation));
348   }
349 
350   // Finally, create libCEED basis
351   PetscCallCeed(ceed, CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points, q_weights, basis));
352   PetscCall(PetscFree(q_points));
353   PetscCall(PetscFree(interp));
354   PetscCall(PetscFree(grad));
355 
356   PetscFunctionReturn(PETSC_SUCCESS);
357 }
358 
359 // -----------------------------------------------------------------------------
360 // Get CEED Basis from DMPlex
361 // -----------------------------------------------------------------------------
362 PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis) {
363   PetscDS         ds;
364   PetscFE         fe;
365   PetscQuadrature quadrature;
366   PetscBool       is_simplex = PETSC_TRUE;
367   PetscInt        ds_field   = -1;
368 
369   PetscFunctionBeginUser;
370 
371   // Get element information
372   PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL));
373   PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field));
374   PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe));
375   PetscCall(PetscFEGetHeightSubspace(fe, height, &fe));
376   PetscCall(PetscFEGetQuadrature(fe, &quadrature));
377 
378   // Check if simplex or tensor-product mesh
379   PetscCall(DMPlexIsSimplex(dm, &is_simplex));
380 
381   // Build libCEED basis
382   if (is_simplex) {
383     PetscTabulation basis_tabulation;
384     PetscInt        num_derivatives = 1, face = 0;
385 
386     PetscCall(PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation));
387     PetscCall(BasisCreateFromTabulation(ceed, dm, domain_label, label_value, height, face, fe, basis_tabulation, quadrature, basis));
388   } else {
389     PetscDualSpace dual_space;
390     PetscInt       num_dual_basis_vectors;
391     PetscInt       dim, num_comp, P, Q;
392 
393     PetscCall(PetscFEGetSpatialDimension(fe, &dim));
394     PetscCall(PetscFEGetNumComponents(fe, &num_comp));
395     PetscCall(PetscFEGetDualSpace(fe, &dual_space));
396     PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
397     P = num_dual_basis_vectors / num_comp;
398     PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL));
399 
400     CeedInt P_1d = (CeedInt)round(pow(P, 1.0 / dim));
401     CeedInt Q_1d = (CeedInt)round(pow(Q, 1.0 / dim));
402 
403     PetscCallCeed(ceed, CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d, CEED_GAUSS, basis));
404   }
405 
406   PetscFunctionReturn(PETSC_SUCCESS);
407 }
408 
409 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc) {
410   PetscFunctionBeginUser;
411 
412   // *****************************************************************************
413   // Set up CEED objects for the interior domain (volume)
414   // *****************************************************************************
415   const PetscInt num_comp_q = 5;
416   const CeedInt  dim = problem->dim, num_comp_x = problem->dim, q_data_size_vol = problem->q_data_size_vol, jac_data_size_vol = num_comp_q + 6 + 3;
417   CeedElemRestriction elem_restr_jd_i;
418   CeedVector          jac_data;
419   CeedInt             num_qpts;
420 
421   // -----------------------------------------------------------------------------
422   // CEED Bases
423   // -----------------------------------------------------------------------------
424   DM dm_coord;
425   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
426 
427   PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, &ceed_data->basis_q));
428   PetscCall(CreateBasisFromPlex(ceed, dm_coord, 0, 0, 0, 0, &ceed_data->basis_x));
429   PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &ceed_data->basis_xc));
430   PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(ceed_data->basis_q, &num_qpts));
431 
432   // -----------------------------------------------------------------------------
433   // CEED Restrictions
434   // -----------------------------------------------------------------------------
435   // -- Create restriction
436   PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, 0, num_qpts, q_data_size_vol, &ceed_data->elem_restr_q, &ceed_data->elem_restr_x,
437                                     &ceed_data->elem_restr_qd_i));
438 
439   PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, 0, num_qpts, jac_data_size_vol, NULL, NULL, &elem_restr_jd_i));
440   // -- Create E vectors
441   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL));
442   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL));
443   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL));
444 
445   // -----------------------------------------------------------------------------
446   // CEED QFunctions
447   // -----------------------------------------------------------------------------
448   // -- Create QFunction for quadrature data
449   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, problem->setup_vol.qfunction_loc, &ceed_data->qf_setup_vol));
450   if (problem->setup_vol.qfunction_context) {
451     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_setup_vol, problem->setup_vol.qfunction_context));
452   }
453   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_vol, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
454   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT));
455   PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_setup_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
456 
457   // -- Create QFunction for ICs
458   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &ceed_data->qf_ics));
459   PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_ics, problem->ics.qfunction_context));
460   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ics, "x", num_comp_x, CEED_EVAL_INTERP));
461   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
462   PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_ics, "q0", num_comp_q, CEED_EVAL_NONE));
463 
464   // -- Create QFunction for RHS
465   if (problem->apply_vol_rhs.qfunction) {
466     PetscCallCeed(
467         ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &ceed_data->qf_rhs_vol));
468     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_rhs_vol, problem->apply_vol_rhs.qfunction_context));
469     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP));
470     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
471     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
472     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP));
473     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP));
474     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
475   }
476 
477   // -- Create QFunction for IFunction
478   if (problem->apply_vol_ifunction.qfunction) {
479     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc,
480                                                     &ceed_data->qf_ifunction_vol));
481     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context));
482     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP));
483     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
484     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP));
485     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
486     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP));
487     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP));
488     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
489     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
490   }
491 
492   CeedQFunction qf_ijacobian_vol = NULL;
493   if (problem->apply_vol_ijacobian.qfunction) {
494     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc,
495                                                     &qf_ijacobian_vol));
496     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context));
497     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP));
498     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD));
499     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
500     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "x", num_comp_x, CEED_EVAL_INTERP));
501     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
502     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP));
503     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
504   }
505 
506   // ---------------------------------------------------------------------------
507   // Element coordinates
508   // ---------------------------------------------------------------------------
509   // -- Create CEED vector
510   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL));
511 
512   // -- Copy PETSc vector in CEED vector
513   Vec X_loc;
514   {
515     DM cdm;
516     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
517     if (cdm) {
518       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
519     } else {
520       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
521     }
522   }
523   PetscCall(VecScale(X_loc, problem->dm_scale));
524   PetscCall(VecCopyP2C(X_loc, ceed_data->x_coord));
525 
526   // -----------------------------------------------------------------------------
527   // CEED vectors
528   // -----------------------------------------------------------------------------
529   // -- Create CEED vector for geometric data
530   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_qd_i, &ceed_data->q_data, NULL));
531   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL));
532 
533   // -----------------------------------------------------------------------------
534   // CEED Operators
535   // -----------------------------------------------------------------------------
536   // -- Create CEED operator for quadrature data
537   PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_vol, NULL, NULL, &ceed_data->op_setup_vol));
538   PetscCallCeed(ceed, CeedOperatorSetField(ceed_data->op_setup_vol, "dx", ceed_data->elem_restr_x, ceed_data->basis_x, CEED_VECTOR_ACTIVE));
539   PetscCallCeed(ceed, CeedOperatorSetField(ceed_data->op_setup_vol, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE));
540   PetscCallCeed(ceed, CeedOperatorSetField(ceed_data->op_setup_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
541 
542   // -- Create CEED operator for ICs
543   CeedOperator op_ics;
544   PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_ics, NULL, NULL, &op_ics));
545   PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE));
546   PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE));
547   PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
548   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label));
549   PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx));
550   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics));
551 
552   // Create CEED operator for RHS
553   if (ceed_data->qf_rhs_vol) {
554     CeedOperator op;
555     PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_rhs_vol, NULL, NULL, &op));
556     PetscCallCeed(ceed, CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
557     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
558     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data));
559     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
560     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
561     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
562     user->op_rhs_vol = op;
563   }
564 
565   // -- CEED operator for IFunction
566   if (ceed_data->qf_ifunction_vol) {
567     CeedOperator op;
568     PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_ifunction_vol, NULL, NULL, &op));
569     PetscCallCeed(ceed, CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
570     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
571     PetscCallCeed(ceed, CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed));
572     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data));
573     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
574     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
575     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
576     PetscCallCeed(ceed, CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data));
577 
578     user->op_ifunction_vol = op;
579   }
580 
581   CeedOperator op_ijacobian_vol = NULL;
582   if (qf_ijacobian_vol) {
583     CeedOperator op;
584     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op));
585     PetscCallCeed(ceed, CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
586     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
587     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data));
588     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
589     PetscCallCeed(ceed, CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data));
590     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
591     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
592     op_ijacobian_vol = op;
593     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol));
594   }
595 
596   // *****************************************************************************
597   // Set up CEED objects for the exterior domain (surface)
598   // *****************************************************************************
599   CeedInt       height = 1, dim_sur = dim - height, P_sur = app_ctx->degree + 1, Q_sur = P_sur + app_ctx->q_extra;
600   const CeedInt q_data_size_sur = problem->q_data_size_sur, jac_data_size_sur = problem->jac_data_size_sur;
601 
602   // -----------------------------------------------------------------------------
603   // CEED Bases
604   // -----------------------------------------------------------------------------
605 
606   DMLabel  label   = 0;
607   PetscInt face_id = 0;
608   PetscInt field   = 0;  // Still want the normal, default field
609   PetscCall(CreateBasisFromPlex(ceed, dm, label, face_id, height, field, &ceed_data->basis_q_sur));
610   PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, face_id, height, field, &ceed_data->basis_x_sur));
611   PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x_sur, ceed_data->basis_q_sur, &ceed_data->basis_xc_sur));
612 
613   // -----------------------------------------------------------------------------
614   // CEED QFunctions
615   // -----------------------------------------------------------------------------
616   // -- Create QFunction for quadrature data
617   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, problem->setup_sur.qfunction_loc, &ceed_data->qf_setup_sur));
618   if (problem->setup_sur.qfunction_context) {
619     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_setup_sur, problem->setup_sur.qfunction_context));
620   }
621   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x * dim_sur, CEED_EVAL_GRAD));
622   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT));
623   PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
624 
625   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow,
626                               problem->apply_inflow_jacobian, &ceed_data->qf_apply_inflow, &ceed_data->qf_apply_inflow_jacobian));
627   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow,
628                               problem->apply_outflow_jacobian, &ceed_data->qf_apply_outflow, &ceed_data->qf_apply_outflow_jacobian));
629   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream,
630                               problem->apply_freestream_jacobian, &ceed_data->qf_apply_freestream, &ceed_data->qf_apply_freestream_jacobian));
631 
632   // *****************************************************************************
633   // CEED Operator Apply
634   // *****************************************************************************
635   // -- Apply CEED Operator for the geometric data
636   PetscCallCeed(ceed, CeedOperatorApply(ceed_data->op_setup_vol, ceed_data->x_coord, ceed_data->q_data, CEED_REQUEST_IMMEDIATE));
637 
638   // -- Create and apply CEED Composite Operator for the entire domain
639   if (!user->phys->implicit) {  // RHS
640     CeedOperator op_rhs;
641     PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_rhs_vol, NULL, height, P_sur, Q_sur, q_data_size_sur, 0, &op_rhs,
642                                       NULL));
643     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx));
644     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
645   } else {  // IFunction
646     PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_ifunction_vol, op_ijacobian_vol, height, P_sur, Q_sur,
647                                       q_data_size_sur, jac_data_size_sur, &user->op_ifunction, op_ijacobian_vol ? &user->op_ijacobian : NULL));
648     if (user->op_ijacobian) {
649       PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ijacobian, "ijacobian time shift", &user->phys->ijacobian_time_shift_label));
650     }
651     if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc, q_data_size_sur));
652     if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SGS_DD_ModelSetup(ceed, user, ceed_data, problem));
653   }
654 
655   if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem));
656   if (app_ctx->diff_filter_monitor) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem));
657 
658   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i));
659   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol));
660   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data));
661   PetscFunctionReturn(PETSC_SUCCESS);
662 }
663