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