xref: /honee/src/setuplibceed.c (revision 2b916ea7fa53b5c2584160b9274b1b14ca18ff4f)
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 "../navierstokes.h"
12 
13 // Utility function - essential BC dofs are encoded in closure indices as -(i+1).
14 PetscInt Involute(PetscInt i) { return i >= 0 ? i : -(i + 1); }
15 
16 // Utility function to create local CEED restriction
17 PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt value, CeedElemRestriction *elem_restr) {
18   PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets;
19 
20   PetscFunctionBeginUser;
21   PetscCall(DMPlexGetLocalOffsets(dm, domain_label, value, height, 0, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets));
22 
23   CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, elem_restr_offsets, elem_restr);
24   PetscCall(PetscFree(elem_restr_offsets));
25 
26   PetscFunctionReturn(0);
27 }
28 
29 // Utility function to get Ceed Restriction for each domain
30 PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, PetscInt value, CeedInt Q, CeedInt q_data_size,
31                                        CeedElemRestriction *elem_restr_q, CeedElemRestriction *elem_restr_x, CeedElemRestriction *elem_restr_qd_i) {
32   DM                  dm_coord;
33   CeedInt             dim, loc_num_elem;
34   CeedInt             Q_dim;
35   CeedElemRestriction elem_restr_tmp;
36   PetscFunctionBeginUser;
37 
38   PetscCall(DMGetDimension(dm, &dim));
39   dim -= height;
40   Q_dim = CeedIntPow(Q, dim);
41   PetscCall(CreateRestrictionFromPlex(ceed, dm, height, domain_label, value, &elem_restr_tmp));
42   if (elem_restr_q) *elem_restr_q = elem_restr_tmp;
43   if (elem_restr_x) {
44     PetscCall(DMGetCellCoordinateDM(dm, &dm_coord));
45     if (!dm_coord) {
46       PetscCall(DMGetCoordinateDM(dm, &dm_coord));
47     }
48     PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL));
49     PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, height, domain_label, value, elem_restr_x));
50   }
51   if (elem_restr_qd_i) {
52     CeedElemRestrictionGetNumElements(elem_restr_tmp, &loc_num_elem);
53     CeedElemRestrictionCreateStrided(ceed, loc_num_elem, Q_dim, q_data_size, q_data_size * loc_num_elem * Q_dim, CEED_STRIDES_BACKEND,
54                                      elem_restr_qd_i);
55   }
56   if (!elem_restr_q) CeedElemRestrictionDestroy(&elem_restr_tmp);
57   PetscFunctionReturn(0);
58 }
59 
60 PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt value, CeedInt height, CeedInt Q_sur,
61                                 CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian,
62                                 CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) {
63   CeedVector          q_data_sur, jac_data_sur;
64   CeedOperator        op_setup_sur, op_apply_bc, op_apply_bc_jacobian = NULL;
65   CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur;
66   CeedInt             num_qpts_sur;
67   PetscFunctionBeginUser;
68 
69   // --- Get number of quadrature points for the boundaries
70   CeedBasisGetNumQuadraturePoints(ceed_data->basis_q_sur, &num_qpts_sur);
71 
72   // ---- CEED Restriction
73   PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, value, Q_sur, q_data_size_sur, &elem_restr_q_sur, &elem_restr_x_sur,
74                                     &elem_restr_qd_i_sur));
75   if (jac_data_size_sur > 0) {
76     // State-dependent data will be passed from residual to Jacobian. This will be collocated.
77     PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, value, Q_sur, jac_data_size_sur, NULL, NULL, &elem_restr_jd_i_sur));
78     CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL);
79   } else {
80     elem_restr_jd_i_sur = NULL;
81     jac_data_sur        = NULL;
82   }
83 
84   // ---- CEED Vector
85   PetscInt loc_num_elem_sur;
86   CeedElemRestrictionGetNumElements(elem_restr_q_sur, &loc_num_elem_sur);
87   CeedVectorCreate(ceed, q_data_size_sur * loc_num_elem_sur * num_qpts_sur, &q_data_sur);
88 
89   // ---- CEED Operator
90   // ----- CEED Operator for Setup (geometric factors)
91   CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur);
92   CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, ceed_data->basis_x_sur, CEED_VECTOR_ACTIVE);
93   CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x_sur, CEED_VECTOR_NONE);
94   CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
95 
96   // ----- CEED Operator for Physics
97   CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc);
98   CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
99   CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
100   CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur);
101   CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord);
102   CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
103   if (elem_restr_jd_i_sur) CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur);
104 
105   if (qf_apply_bc_jacobian) {
106     CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian);
107     CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
108     CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
109     CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur);
110     CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord);
111     CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur);
112     CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
113   }
114 
115   // ----- Apply CEED operator for Setup
116   CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE);
117 
118   // ----- Apply Sub-Operator for Physics
119   CeedCompositeOperatorAddSub(*op_apply, op_apply_bc);
120   if (op_apply_bc_jacobian) CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_bc_jacobian);
121 
122   // ----- Cleanup
123   CeedVectorDestroy(&q_data_sur);
124   CeedVectorDestroy(&jac_data_sur);
125   CeedElemRestrictionDestroy(&elem_restr_q_sur);
126   CeedElemRestrictionDestroy(&elem_restr_x_sur);
127   CeedElemRestrictionDestroy(&elem_restr_qd_i_sur);
128   CeedElemRestrictionDestroy(&elem_restr_jd_i_sur);
129   CeedOperatorDestroy(&op_setup_sur);
130   CeedOperatorDestroy(&op_apply_bc);
131   CeedOperatorDestroy(&op_apply_bc_jacobian);
132 
133   PetscFunctionReturn(0);
134 }
135 
136 // Utility function to create CEED Composite Operator for the entire domain
137 PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol,
138                                        CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur,
139                                        CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) {
140   DMLabel domain_label;
141   PetscFunctionBeginUser;
142 
143   // Create Composite Operaters
144   CeedCompositeOperatorCreate(ceed, op_apply);
145   if (op_apply_ijacobian) CeedCompositeOperatorCreate(ceed, op_apply_ijacobian);
146 
147   // --Apply Sub-Operator for the volume
148   CeedCompositeOperatorAddSub(*op_apply, op_apply_vol);
149   if (op_apply_ijacobian) CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_ijacobian_vol);
150 
151   // -- Create Sub-Operator for in/outflow BCs
152   if (phys->has_neumann || 1) {
153     // --- Setup
154     PetscCall(DMGetLabel(dm, "Face Sets", &domain_label));
155 
156     // --- Create Sub-Operator for inflow boundaries
157     for (CeedInt i = 0; i < bc->num_inflow; i++) {
158       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->inflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
159                                  ceed_data->qf_apply_inflow, ceed_data->qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian));
160     }
161 
162     // --- Create Sub-Operator for outflow boundaries
163     for (CeedInt i = 0; i < bc->num_outflow; i++) {
164       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
165                                  ceed_data->qf_apply_outflow, ceed_data->qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian));
166     }
167 
168     // --- Create Sub-Operator for freestream boundaries
169     for (CeedInt i = 0; i < bc->num_freestream; i++) {
170       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
171                                  ceed_data->qf_apply_freestream, ceed_data->qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian));
172     }
173   }
174 
175   // ----- Get Context Labels for Operator
176   CeedOperatorContextGetFieldLabel(*op_apply, "solution time", &phys->solution_time_label);
177   CeedOperatorContextGetFieldLabel(*op_apply, "timestep size", &phys->timestep_size_label);
178 
179   PetscFunctionReturn(0);
180 }
181 
182 PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur,
183                                  PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian,
184                                  CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) {
185   PetscFunctionBeginUser;
186 
187   if (apply_bc.qfunction) {
188     // *INDENT-OFF*
189     CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc);
190     CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context);
191     CeedQFunctionContextDestroy(&apply_bc.qfunction_context);
192     CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP);
193     CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD);
194     CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
195     CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP);
196     CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP);
197     // *INDENT-ON*
198     if (jac_data_size_sur) CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE);
199   }
200   if (apply_bc_jacobian.qfunction) {
201     // *INDENT-OFF*
202     CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian);
203     CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context);
204     CeedQFunctionContextDestroy(&apply_bc_jacobian.qfunction_context);
205     CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP);
206     CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD);
207     CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
208     CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP);
209     CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE);
210     CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP);
211     // *INDENT-ON*
212   }
213   PetscFunctionReturn(0);
214 }
215 
216 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc) {
217   PetscFunctionBeginUser;
218 
219   // *****************************************************************************
220   // Set up CEED objects for the interior domain (volume)
221   // *****************************************************************************
222   const PetscInt num_comp_q = 5;
223   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,
224                 P = app_ctx->degree + 1, Q = P + app_ctx->q_extra;
225   CeedElemRestriction elem_restr_jd_i;
226   CeedVector          jac_data;
227 
228   // -----------------------------------------------------------------------------
229   // CEED Bases
230   // -----------------------------------------------------------------------------
231   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_q, P, Q, CEED_GAUSS, &ceed_data->basis_q);
232   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS, &ceed_data->basis_x);
233   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, P, CEED_GAUSS_LOBATTO, &ceed_data->basis_xc);
234 
235   // -----------------------------------------------------------------------------
236   // CEED Restrictions
237   // -----------------------------------------------------------------------------
238   // -- Create restriction
239   PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, q_data_size_vol, &ceed_data->elem_restr_q, &ceed_data->elem_restr_x,
240                                     &ceed_data->elem_restr_qd_i));
241 
242   PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, jac_data_size_vol, NULL, NULL, &elem_restr_jd_i));
243   // -- Create E vectors
244   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL);
245   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL);
246   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL);
247 
248   // -----------------------------------------------------------------------------
249   // CEED QFunctions
250   // -----------------------------------------------------------------------------
251   // -- Create QFunction for quadrature data
252   CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, problem->setup_vol.qfunction_loc, &ceed_data->qf_setup_vol);
253   if (problem->setup_vol.qfunction_context) {
254     CeedQFunctionSetContext(ceed_data->qf_setup_vol, problem->setup_vol.qfunction_context);
255     CeedQFunctionContextDestroy(&problem->setup_vol.qfunction_context);
256   }
257   CeedQFunctionAddInput(ceed_data->qf_setup_vol, "dx", num_comp_x * dim, CEED_EVAL_GRAD);
258   CeedQFunctionAddInput(ceed_data->qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT);
259   CeedQFunctionAddOutput(ceed_data->qf_setup_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE);
260 
261   // -- Create QFunction for ICs
262   CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &ceed_data->qf_ics);
263   CeedQFunctionSetContext(ceed_data->qf_ics, problem->ics.qfunction_context);
264   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
265   CeedQFunctionAddInput(ceed_data->qf_ics, "x", num_comp_x, CEED_EVAL_INTERP);
266   CeedQFunctionAddInput(ceed_data->qf_ics, "qdata", q_data_size_vol, CEED_EVAL_NONE);
267   CeedQFunctionAddOutput(ceed_data->qf_ics, "q0", num_comp_q, CEED_EVAL_NONE);
268 
269   // -- Create QFunction for RHS
270   if (problem->apply_vol_rhs.qfunction) {
271     CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &ceed_data->qf_rhs_vol);
272     CeedQFunctionSetContext(ceed_data->qf_rhs_vol, problem->apply_vol_rhs.qfunction_context);
273     CeedQFunctionContextDestroy(&problem->apply_vol_rhs.qfunction_context);
274     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP);
275     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD);
276     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE);
277     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP);
278     CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP);
279     CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD);
280   }
281 
282   // -- Create QFunction for IFunction
283   if (problem->apply_vol_ifunction.qfunction) {
284     CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc,
285                                 &ceed_data->qf_ifunction_vol);
286     CeedQFunctionSetContext(ceed_data->qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context);
287     CeedQFunctionContextDestroy(&problem->apply_vol_ifunction.qfunction_context);
288     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP);
289     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD);
290     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP);
291     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE);
292     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP);
293     CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP);
294     CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD);
295     CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE);
296   }
297 
298   CeedQFunction qf_ijacobian_vol = NULL;
299   if (problem->apply_vol_ijacobian.qfunction) {
300     CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc, &qf_ijacobian_vol);
301     CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context);
302     CeedQFunctionContextDestroy(&problem->apply_vol_ijacobian.qfunction_context);
303     CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP);
304     CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD);
305     CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE);
306     CeedQFunctionAddInput(qf_ijacobian_vol, "x", num_comp_x, CEED_EVAL_INTERP);
307     CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE);
308     CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP);
309     CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD);
310   }
311 
312   // ---------------------------------------------------------------------------
313   // Element coordinates
314   // ---------------------------------------------------------------------------
315   // -- Create CEED vector
316   CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL);
317 
318   // -- Copy PETSc vector in CEED vector
319   Vec                X_loc;
320   const PetscScalar *X_loc_array;
321   {
322     DM cdm;
323     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
324     if (cdm) {
325       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
326     } else {
327       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
328     }
329   }
330   PetscCall(VecScale(X_loc, problem->dm_scale));
331   PetscCall(VecGetArrayRead(X_loc, &X_loc_array));
332   CeedVectorSetArray(ceed_data->x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, (PetscScalar *)X_loc_array);
333   PetscCall(VecRestoreArrayRead(X_loc, &X_loc_array));
334 
335   // -----------------------------------------------------------------------------
336   // CEED vectors
337   // -----------------------------------------------------------------------------
338   // -- Create CEED vector for geometric data
339   CeedInt  num_qpts_vol;
340   PetscInt loc_num_elem_vol;
341   CeedBasisGetNumQuadraturePoints(ceed_data->basis_q, &num_qpts_vol);
342   CeedElemRestrictionGetNumElements(ceed_data->elem_restr_q, &loc_num_elem_vol);
343   CeedVectorCreate(ceed, q_data_size_vol * loc_num_elem_vol * num_qpts_vol, &ceed_data->q_data);
344 
345   CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL);
346   // -----------------------------------------------------------------------------
347   // CEED Operators
348   // -----------------------------------------------------------------------------
349   // -- Create CEED operator for quadrature data
350   CeedOperatorCreate(ceed, ceed_data->qf_setup_vol, NULL, NULL, &ceed_data->op_setup_vol);
351   CeedOperatorSetField(ceed_data->op_setup_vol, "dx", ceed_data->elem_restr_x, ceed_data->basis_x, CEED_VECTOR_ACTIVE);
352   CeedOperatorSetField(ceed_data->op_setup_vol, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE);
353   CeedOperatorSetField(ceed_data->op_setup_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
354 
355   // -- Create CEED operator for ICs
356   CeedOperatorCreate(ceed, ceed_data->qf_ics, NULL, NULL, &ceed_data->op_ics);
357   CeedOperatorSetField(ceed_data->op_ics, "x", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE);
358   CeedOperatorSetField(ceed_data->op_ics, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
359   CeedOperatorSetField(ceed_data->op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
360   CeedOperatorContextGetFieldLabel(ceed_data->op_ics, "evaluation time", &user->phys->ics_time_label);
361 
362   // Create CEED operator for RHS
363   if (ceed_data->qf_rhs_vol) {
364     CeedOperator op;
365     CeedOperatorCreate(ceed, ceed_data->qf_rhs_vol, NULL, NULL, &op);
366     CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
367     CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
368     CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
369     CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
370     CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
371     CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
372     user->op_rhs_vol = op;
373   }
374 
375   // -- CEED operator for IFunction
376   if (ceed_data->qf_ifunction_vol) {
377     CeedOperator op;
378     CeedOperatorCreate(ceed, ceed_data->qf_ifunction_vol, NULL, NULL, &op);
379     CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
380     CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
381     CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed);
382     CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
383     CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
384     CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
385     CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
386     CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data);
387 
388     user->op_ifunction_vol = op;
389   }
390 
391   CeedOperator op_ijacobian_vol = NULL;
392   if (qf_ijacobian_vol) {
393     CeedOperator op;
394     CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op);
395     CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
396     CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
397     CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
398     CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
399     CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data);
400     CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
401     CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
402     op_ijacobian_vol = op;
403     CeedQFunctionDestroy(&qf_ijacobian_vol);
404   }
405 
406   // *****************************************************************************
407   // Set up CEED objects for the exterior domain (surface)
408   // *****************************************************************************
409   CeedInt       height = 1, dim_sur = dim - height, P_sur = app_ctx->degree + 1, Q_sur = P_sur + app_ctx->q_extra;
410   const CeedInt q_data_size_sur = problem->q_data_size_sur, jac_data_size_sur = problem->jac_data_size_sur;
411 
412   // -----------------------------------------------------------------------------
413   // CEED Bases
414   // -----------------------------------------------------------------------------
415   CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_q, P_sur, Q_sur, CEED_GAUSS, &ceed_data->basis_q_sur);
416   CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_x, 2, Q_sur, CEED_GAUSS, &ceed_data->basis_x_sur);
417   CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_x, 2, P_sur, CEED_GAUSS_LOBATTO, &ceed_data->basis_xc_sur);
418 
419   // -----------------------------------------------------------------------------
420   // CEED QFunctions
421   // -----------------------------------------------------------------------------
422   // -- Create QFunction for quadrature data
423   CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, problem->setup_sur.qfunction_loc, &ceed_data->qf_setup_sur);
424   if (problem->setup_sur.qfunction_context) {
425     CeedQFunctionSetContext(ceed_data->qf_setup_sur, problem->setup_sur.qfunction_context);
426     CeedQFunctionContextDestroy(&problem->setup_sur.qfunction_context);
427   }
428   CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x * dim_sur, CEED_EVAL_GRAD);
429   CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT);
430   CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
431 
432   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow,
433                               problem->apply_inflow_jacobian, &ceed_data->qf_apply_inflow, &ceed_data->qf_apply_inflow_jacobian));
434 
435   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow,
436                               problem->apply_outflow_jacobian, &ceed_data->qf_apply_outflow, &ceed_data->qf_apply_outflow_jacobian));
437 
438   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream,
439                               problem->apply_freestream_jacobian, &ceed_data->qf_apply_freestream, &ceed_data->qf_apply_freestream_jacobian));
440 
441   // *****************************************************************************
442   // CEED Operator Apply
443   // *****************************************************************************
444   // -- Apply CEED Operator for the geometric data
445   CeedOperatorApply(ceed_data->op_setup_vol, ceed_data->x_coord, ceed_data->q_data, CEED_REQUEST_IMMEDIATE);
446 
447   // -- Create and apply CEED Composite Operator for the entire domain
448   if (!user->phys->implicit) {  // RHS
449     PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_rhs_vol, NULL, height, P_sur, Q_sur, q_data_size_sur, 0,
450                                       &user->op_rhs, NULL));
451   } else {  // IFunction
452     PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_ifunction_vol, op_ijacobian_vol, height, P_sur, Q_sur,
453                                       q_data_size_sur, jac_data_size_sur, &user->op_ifunction, op_ijacobian_vol ? &user->op_ijacobian : NULL));
454     if (user->op_ijacobian) {
455       CeedOperatorContextGetFieldLabel(user->op_ijacobian, "ijacobian time shift", &user->phys->ijacobian_time_shift_label);
456     }
457     if (problem->use_dirichlet_ceed) {
458       PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, app_ctx, problem, bc, Q_sur, q_data_size_sur));
459     }
460   }
461   CeedElemRestrictionDestroy(&elem_restr_jd_i);
462   CeedOperatorDestroy(&op_ijacobian_vol);
463   CeedVectorDestroy(&jac_data);
464   PetscFunctionReturn(0);
465 }
466