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