xref: /libCEED/examples/fluids/src/setuplibceed.c (revision 3c4b7af6bf1a206cad198c63d890a93b5be72708)
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       ierr = GetRestrictionForDomain(ceed, dm, height, domain_label, bc->outflows[i],
198                                      Q_sur, jac_data_size_sur, NULL, NULL,
199                                      &elem_restr_jd_i_sur);
200       CHKERRQ(ierr);
201 
202       // ---- CEED Vector
203       PetscInt loc_num_elem_sur;
204       CeedElemRestrictionGetNumElements(elem_restr_q_sur, &loc_num_elem_sur);
205       CeedVectorCreate(ceed, q_data_size_sur*loc_num_elem_sur*num_qpts_sur,
206                        &q_data_sur);
207 
208       CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL);
209 
210       // ---- CEED Operator
211       // ----- CEED Operator for Setup (geometric factors)
212       CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur);
213       CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur,
214                            ceed_data->basis_x_sur, CEED_VECTOR_ACTIVE);
215       CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE,
216                            ceed_data->basis_x_sur, CEED_VECTOR_NONE);
217       CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur,
218                            CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
219 
220       // ----- CEED Operator for Physics
221       CeedOperatorCreate(ceed, ceed_data->qf_apply_outflow, NULL, NULL,
222                          &op_apply_outflow);
223       CeedOperatorSetField(op_apply_outflow, "q", elem_restr_q_sur,
224                            ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
225       CeedOperatorSetField(op_apply_outflow, "surface qdata", elem_restr_qd_i_sur,
226                            CEED_BASIS_COLLOCATED, q_data_sur);
227       CeedOperatorSetField(op_apply_outflow, "x", elem_restr_x_sur,
228                            ceed_data->basis_x_sur, ceed_data->x_coord);
229       CeedOperatorSetField(op_apply_outflow, "v", elem_restr_q_sur,
230                            ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
231       CeedOperatorSetField(op_apply_outflow, "surface jacobian data",
232                            elem_restr_jd_i_sur,
233                            CEED_BASIS_COLLOCATED, jac_data_sur);
234 
235       if (ceed_data->qf_apply_outflow_jacobian) {
236         CeedOperatorCreate(ceed, ceed_data->qf_apply_outflow_jacobian, NULL, NULL,
237                            &op_apply_outflow_jacobian);
238         CeedOperatorSetField(op_apply_outflow_jacobian, "dq", elem_restr_q_sur,
239                              ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
240         CeedOperatorSetField(op_apply_outflow_jacobian, "surface qdata",
241                              elem_restr_qd_i_sur,
242                              CEED_BASIS_COLLOCATED, q_data_sur);
243         CeedOperatorSetField(op_apply_outflow_jacobian, "surface jacobian data",
244                              elem_restr_jd_i_sur,
245                              CEED_BASIS_COLLOCATED, jac_data_sur);
246         CeedOperatorSetField(op_apply_outflow_jacobian, "v", elem_restr_q_sur,
247                              ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
248       }
249 
250       // ----- Apply CEED operator for Setup
251       CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur,
252                         CEED_REQUEST_IMMEDIATE);
253 
254       // ----- Apply Sub-Operator for Physics
255       CeedCompositeOperatorAddSub(*op_apply, op_apply_outflow);
256       if (op_apply_ijacobian)
257         CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_outflow_jacobian);
258 
259       // ----- Cleanup
260       CeedVectorDestroy(&q_data_sur);
261       CeedVectorDestroy(&jac_data_sur);
262       CeedElemRestrictionDestroy(&elem_restr_q_sur);
263       CeedElemRestrictionDestroy(&elem_restr_x_sur);
264       CeedElemRestrictionDestroy(&elem_restr_qd_i_sur);
265       CeedElemRestrictionDestroy(&elem_restr_jd_i_sur);
266       CeedOperatorDestroy(&op_setup_sur);
267       CeedOperatorDestroy(&op_apply_outflow);
268       CeedOperatorDestroy(&op_apply_outflow_jacobian);
269     }
270   }
271 
272   // ----- Get Context Labels for Operator
273   CeedOperatorContextGetFieldLabel(*op_apply, "solution time",
274                                    &phys->solution_time_label);
275   CeedOperatorContextGetFieldLabel(*op_apply, "timestep size",
276                                    &phys->timestep_size_label);
277 
278   PetscFunctionReturn(0);
279 }
280 
281 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user,
282                             AppCtx app_ctx, ProblemData *problem, SimpleBC bc) {
283   PetscErrorCode ierr;
284   PetscFunctionBeginUser;
285 
286   // *****************************************************************************
287   // Set up CEED objects for the interior domain (volume)
288   // *****************************************************************************
289   const PetscInt num_comp_q      = 5;
290   const CeedInt  dim             = problem->dim,
291                  num_comp_x      = problem->dim,
292                  q_data_size_vol = problem->q_data_size_vol,
293                  jac_data_size_vol = num_comp_q + 6 + 3,
294                  P               = app_ctx->degree + 1,
295                  Q               = P + app_ctx->q_extra;
296   CeedElemRestriction elem_restr_jd_i;
297   CeedVector jac_data;
298 
299   // -----------------------------------------------------------------------------
300   // CEED Bases
301   // -----------------------------------------------------------------------------
302   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_q, P, Q, CEED_GAUSS,
303                                   &ceed_data->basis_q);
304   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, CEED_GAUSS,
305                                   &ceed_data->basis_x);
306   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, P,
307                                   CEED_GAUSS_LOBATTO, &ceed_data->basis_xc);
308 
309   // -----------------------------------------------------------------------------
310   // CEED Restrictions
311   // -----------------------------------------------------------------------------
312   // -- Create restriction
313   ierr = GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, q_data_size_vol,
314                                  &ceed_data->elem_restr_q, &ceed_data->elem_restr_x,
315                                  &ceed_data->elem_restr_qd_i); CHKERRQ(ierr);
316 
317   ierr = GetRestrictionForDomain(ceed, dm, 0, 0, 0, Q, jac_data_size_vol,
318                                  NULL, NULL,
319                                  &elem_restr_jd_i); CHKERRQ(ierr);
320 // -- Create E vectors
321   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL);
322   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed,
323                                   NULL);
324   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL);
325 
326   // -----------------------------------------------------------------------------
327   // CEED QFunctions
328   // -----------------------------------------------------------------------------
329   // -- Create QFunction for quadrature data
330   CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction,
331                               problem->setup_vol.qfunction_loc,
332                               &ceed_data->qf_setup_vol);
333   if (problem->setup_vol.qfunction_context) {
334     CeedQFunctionSetContext(ceed_data->qf_setup_vol,
335                             problem->setup_vol.qfunction_context);
336     CeedQFunctionContextDestroy(&problem->setup_vol.qfunction_context);
337   }
338   CeedQFunctionAddInput(ceed_data->qf_setup_vol, "dx", num_comp_x*dim,
339                         CEED_EVAL_GRAD);
340   CeedQFunctionAddInput(ceed_data->qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT);
341   CeedQFunctionAddOutput(ceed_data->qf_setup_vol, "qdata", q_data_size_vol,
342                          CEED_EVAL_NONE);
343 
344   // -- Create QFunction for ICs
345   CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction,
346                               problem->ics.qfunction_loc,
347                               &ceed_data->qf_ics);
348   CeedQFunctionSetContext(ceed_data->qf_ics, problem->ics.qfunction_context);
349   CeedQFunctionContextDestroy(&problem->ics.qfunction_context);
350   CeedQFunctionAddInput(ceed_data->qf_ics, "x", num_comp_x, CEED_EVAL_INTERP);
351   CeedQFunctionAddOutput(ceed_data->qf_ics, "q0", num_comp_q, CEED_EVAL_NONE);
352 
353   // -- Create QFunction for RHS
354   if (problem->apply_vol_rhs.qfunction) {
355     CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction,
356                                 problem->apply_vol_rhs.qfunction_loc, &ceed_data->qf_rhs_vol);
357     CeedQFunctionSetContext(ceed_data->qf_rhs_vol,
358                             problem->apply_vol_rhs.qfunction_context);
359     CeedQFunctionContextDestroy(&problem->apply_vol_rhs.qfunction_context);
360     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP);
361     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "Grad_q", num_comp_q*dim,
362                           CEED_EVAL_GRAD);
363     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "qdata", q_data_size_vol,
364                           CEED_EVAL_NONE);
365     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP);
366     CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "v", num_comp_q,
367                            CEED_EVAL_INTERP);
368     CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "Grad_v", num_comp_q*dim,
369                            CEED_EVAL_GRAD);
370   }
371 
372   // -- Create QFunction for IFunction
373   if (problem->apply_vol_ifunction.qfunction) {
374     CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction,
375                                 problem->apply_vol_ifunction.qfunction_loc, &ceed_data->qf_ifunction_vol);
376     CeedQFunctionSetContext(ceed_data->qf_ifunction_vol,
377                             problem->apply_vol_ifunction.qfunction_context);
378     CeedQFunctionContextDestroy(&problem->apply_vol_ifunction.qfunction_context);
379     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q", num_comp_q,
380                           CEED_EVAL_INTERP);
381     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "Grad_q", num_comp_q*dim,
382                           CEED_EVAL_GRAD);
383     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q dot", num_comp_q,
384                           CEED_EVAL_INTERP);
385     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "qdata", q_data_size_vol,
386                           CEED_EVAL_NONE);
387     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "x", num_comp_x,
388                           CEED_EVAL_INTERP);
389     CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "v", num_comp_q,
390                            CEED_EVAL_INTERP);
391     CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "Grad_v", num_comp_q*dim,
392                            CEED_EVAL_GRAD);
393     CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "jac_data",
394                            jac_data_size_vol, CEED_EVAL_NONE);
395   }
396 
397   CeedQFunction qf_ijacobian_vol = NULL;
398   if (problem->apply_vol_ijacobian.qfunction) {
399     CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction,
400                                 problem->apply_vol_ijacobian.qfunction_loc, &qf_ijacobian_vol);
401     CeedQFunctionSetContext(qf_ijacobian_vol,
402                             problem->apply_vol_ijacobian.qfunction_context);
403     CeedQFunctionContextDestroy(&problem->apply_vol_ijacobian.qfunction_context);
404     CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q,
405                           CEED_EVAL_INTERP);
406     CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q*dim,
407                           CEED_EVAL_GRAD);
408     CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol,
409                           CEED_EVAL_NONE);
410     CeedQFunctionAddInput(qf_ijacobian_vol, "x", num_comp_x,
411                           CEED_EVAL_INTERP);
412     CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data",
413                           jac_data_size_vol, CEED_EVAL_NONE);
414     CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q,
415                            CEED_EVAL_INTERP);
416     CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q*dim,
417                            CEED_EVAL_GRAD);
418   }
419 
420   // ---------------------------------------------------------------------------
421   // Element coordinates
422   // ---------------------------------------------------------------------------
423   // -- Create CEED vector
424   CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord,
425                                   NULL);
426 
427   // -- Copy PETSc vector in CEED vector
428   Vec               X_loc;
429   const PetscScalar *X_loc_array;
430   ierr = DMGetCoordinatesLocal(dm, &X_loc); CHKERRQ(ierr);
431   ierr = VecScale(X_loc, problem->dm_scale); CHKERRQ(ierr);
432   ierr = VecGetArrayRead(X_loc, &X_loc_array); CHKERRQ(ierr);
433   CeedVectorSetArray(ceed_data->x_coord, CEED_MEM_HOST, CEED_COPY_VALUES,
434                      (PetscScalar *)X_loc_array);
435   ierr = VecRestoreArrayRead(X_loc, &X_loc_array); CHKERRQ(ierr);
436 
437   // -----------------------------------------------------------------------------
438   // CEED vectors
439   // -----------------------------------------------------------------------------
440   // -- Create CEED vector for geometric data
441   CeedInt  num_qpts_vol;
442   PetscInt loc_num_elem_vol;
443   CeedBasisGetNumQuadraturePoints(ceed_data->basis_q, &num_qpts_vol);
444   CeedElemRestrictionGetNumElements(ceed_data->elem_restr_q, &loc_num_elem_vol);
445   CeedVectorCreate(ceed, q_data_size_vol*loc_num_elem_vol*num_qpts_vol,
446                    &ceed_data->q_data);
447 
448   CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL);
449   // -----------------------------------------------------------------------------
450   // CEED Operators
451   // -----------------------------------------------------------------------------
452   // -- Create CEED operator for quadrature data
453   CeedOperatorCreate(ceed, ceed_data->qf_setup_vol, NULL, NULL,
454                      &ceed_data->op_setup_vol);
455   CeedOperatorSetField(ceed_data->op_setup_vol, "dx", ceed_data->elem_restr_x,
456                        ceed_data->basis_x, CEED_VECTOR_ACTIVE);
457   CeedOperatorSetField(ceed_data->op_setup_vol, "weight",
458                        CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE);
459   CeedOperatorSetField(ceed_data->op_setup_vol, "qdata",
460                        ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
461 
462   // -- Create CEED operator for ICs
463   CeedOperatorCreate(ceed, ceed_data->qf_ics, NULL, NULL, &ceed_data->op_ics);
464   CeedOperatorSetField(ceed_data->op_ics, "x", ceed_data->elem_restr_x,
465                        ceed_data->basis_xc, CEED_VECTOR_ACTIVE);
466   CeedOperatorSetField(ceed_data->op_ics, "q0", ceed_data->elem_restr_q,
467                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
468   CeedOperatorContextGetFieldLabel(ceed_data->op_ics, "evaluation time",
469                                    &user->phys->ics_time_label);
470 
471   // Create CEED operator for RHS
472   if (ceed_data->qf_rhs_vol) {
473     CeedOperator op;
474     CeedOperatorCreate(ceed, ceed_data->qf_rhs_vol, NULL, NULL, &op);
475     CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q,
476                          CEED_VECTOR_ACTIVE);
477     CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q,
478                          CEED_VECTOR_ACTIVE);
479     CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i,
480                          CEED_BASIS_COLLOCATED, ceed_data->q_data);
481     CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x,
482                          ceed_data->x_coord);
483     CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q,
484                          CEED_VECTOR_ACTIVE);
485     CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q,
486                          CEED_VECTOR_ACTIVE);
487     user->op_rhs_vol = op;
488   }
489 
490   // -- CEED operator for IFunction
491   if (ceed_data->qf_ifunction_vol) {
492     CeedOperator op;
493     CeedOperatorCreate(ceed, ceed_data->qf_ifunction_vol, NULL, NULL, &op);
494     CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q,
495                          CEED_VECTOR_ACTIVE);
496     CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q,
497                          CEED_VECTOR_ACTIVE);
498     CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q,
499                          user->q_dot_ceed);
500     CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i,
501                          CEED_BASIS_COLLOCATED, ceed_data->q_data);
502     CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x,
503                          ceed_data->x_coord);
504     CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q,
505                          CEED_VECTOR_ACTIVE);
506     CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q,
507                          CEED_VECTOR_ACTIVE);
508     CeedOperatorSetField(op, "jac_data", elem_restr_jd_i,
509                          CEED_BASIS_COLLOCATED, jac_data);
510 
511     user->op_ifunction_vol = op;
512   }
513 
514   CeedOperator op_ijacobian_vol = NULL;
515   if (qf_ijacobian_vol) {
516     CeedOperator op;
517     CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op);
518     CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q,
519                          CEED_VECTOR_ACTIVE);
520     CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q,
521                          CEED_VECTOR_ACTIVE);
522     CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i,
523                          CEED_BASIS_COLLOCATED, ceed_data->q_data);
524     CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x,
525                          ceed_data->x_coord);
526     CeedOperatorSetField(op, "jac_data", elem_restr_jd_i,
527                          CEED_BASIS_COLLOCATED, jac_data);
528     CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q,
529                          CEED_VECTOR_ACTIVE);
530     CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q,
531                          CEED_VECTOR_ACTIVE);
532     op_ijacobian_vol = op;
533     CeedQFunctionDestroy(&qf_ijacobian_vol);
534   }
535 
536   // *****************************************************************************
537   // Set up CEED objects for the exterior domain (surface)
538   // *****************************************************************************
539   CeedInt height  = 1,
540           dim_sur = dim - height,
541           P_sur   = app_ctx->degree + 1,
542           Q_sur   = P_sur + app_ctx->q_extra;
543   const CeedInt q_data_size_sur = problem->q_data_size_sur,
544                 jac_data_size_sur = problem->jac_data_size_sur;
545 
546   // -----------------------------------------------------------------------------
547   // CEED Bases
548   // -----------------------------------------------------------------------------
549   CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_q, P_sur, Q_sur,
550                                   CEED_GAUSS, &ceed_data->basis_q_sur);
551   CeedBasisCreateTensorH1Lagrange(ceed, dim_sur, num_comp_x, 2, Q_sur, CEED_GAUSS,
552                                   &ceed_data->basis_x_sur);
553 
554   // -----------------------------------------------------------------------------
555   // CEED QFunctions
556   // -----------------------------------------------------------------------------
557   // -- Create QFunction for quadrature data
558   CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction,
559                               problem->setup_sur.qfunction_loc,
560                               &ceed_data->qf_setup_sur);
561   if (problem->setup_sur.qfunction_context) {
562     CeedQFunctionSetContext(ceed_data->qf_setup_sur,
563                             problem->setup_sur.qfunction_context);
564     CeedQFunctionContextDestroy(&problem->setup_sur.qfunction_context);
565   }
566   CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x*dim_sur,
567                         CEED_EVAL_GRAD);
568   CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT);
569   CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata",
570                          q_data_size_sur, CEED_EVAL_NONE);
571 
572   // -- Creat QFunction for inflow boundaries
573   if (problem->apply_inflow.qfunction) {
574     CeedQFunctionCreateInterior(ceed, 1, problem->apply_inflow.qfunction,
575                                 problem->apply_inflow.qfunction_loc, &ceed_data->qf_apply_inflow);
576     CeedQFunctionSetContext(ceed_data->qf_apply_inflow,
577                             problem->apply_inflow.qfunction_context);
578     CeedQFunctionContextDestroy(&problem->apply_inflow.qfunction_context);
579     CeedQFunctionAddInput(ceed_data->qf_apply_inflow, "q", num_comp_q,
580                           CEED_EVAL_INTERP);
581     CeedQFunctionAddInput(ceed_data->qf_apply_inflow, "surface qdata",
582                           q_data_size_sur, CEED_EVAL_NONE);
583     CeedQFunctionAddInput(ceed_data->qf_apply_inflow, "x", num_comp_x,
584                           CEED_EVAL_INTERP);
585     CeedQFunctionAddOutput(ceed_data->qf_apply_inflow, "v", num_comp_q,
586                            CEED_EVAL_INTERP);
587   }
588   if (problem->apply_inflow_jacobian.qfunction) {
589     CeedQFunctionCreateInterior(ceed, 1, problem->apply_inflow_jacobian.qfunction,
590                                 problem->apply_inflow_jacobian.qfunction_loc,
591                                 &ceed_data->qf_apply_inflow_jacobian);
592     CeedQFunctionSetContext(ceed_data->qf_apply_inflow_jacobian,
593                             problem->apply_inflow_jacobian.qfunction_context);
594     CeedQFunctionContextDestroy(&problem->apply_inflow_jacobian.qfunction_context);
595     CeedQFunctionAddInput(ceed_data->qf_apply_inflow_jacobian, "dq", num_comp_q,
596                           CEED_EVAL_INTERP);
597     CeedQFunctionAddInput(ceed_data->qf_apply_inflow_jacobian, "surface qdata",
598                           q_data_size_sur, CEED_EVAL_NONE);
599     CeedQFunctionAddInput(ceed_data->qf_apply_inflow_jacobian, "x", num_comp_x,
600                           CEED_EVAL_INTERP);
601     CeedQFunctionAddOutput(ceed_data->qf_apply_inflow_jacobian, "v", num_comp_q,
602                            CEED_EVAL_INTERP);
603   }
604 
605   // -- Creat QFunction for outflow boundaries
606   if (problem->apply_outflow.qfunction) {
607     CeedQFunctionCreateInterior(ceed, 1, problem->apply_outflow.qfunction,
608                                 problem->apply_outflow.qfunction_loc, &ceed_data->qf_apply_outflow);
609     CeedQFunctionSetContext(ceed_data->qf_apply_outflow,
610                             problem->apply_outflow.qfunction_context);
611     CeedQFunctionContextDestroy(&problem->apply_outflow.qfunction_context);
612     CeedQFunctionAddInput(ceed_data->qf_apply_outflow, "q", num_comp_q,
613                           CEED_EVAL_INTERP);
614     CeedQFunctionAddInput(ceed_data->qf_apply_outflow, "surface qdata",
615                           q_data_size_sur, CEED_EVAL_NONE);
616     CeedQFunctionAddInput(ceed_data->qf_apply_outflow, "x", num_comp_x,
617                           CEED_EVAL_INTERP);
618     CeedQFunctionAddOutput(ceed_data->qf_apply_outflow, "v", num_comp_q,
619                            CEED_EVAL_INTERP);
620     CeedQFunctionAddOutput(ceed_data->qf_apply_outflow, "surface jacobian data",
621                            jac_data_size_sur,
622                            CEED_EVAL_NONE);
623   }
624   if (problem->apply_outflow_jacobian.qfunction) {
625     CeedQFunctionCreateInterior(ceed, 1, problem->apply_outflow_jacobian.qfunction,
626                                 problem->apply_outflow_jacobian.qfunction_loc,
627                                 &ceed_data->qf_apply_outflow_jacobian);
628     CeedQFunctionSetContext(ceed_data->qf_apply_outflow_jacobian,
629                             problem->apply_outflow_jacobian.qfunction_context);
630     CeedQFunctionContextDestroy(&problem->apply_outflow_jacobian.qfunction_context);
631     CeedQFunctionAddInput(ceed_data->qf_apply_outflow_jacobian, "dq", num_comp_q,
632                           CEED_EVAL_INTERP);
633     CeedQFunctionAddInput(ceed_data->qf_apply_outflow_jacobian, "surface qdata",
634                           q_data_size_sur, CEED_EVAL_NONE);
635     CeedQFunctionAddInput(ceed_data->qf_apply_outflow_jacobian,
636                           "surface jacobian data",
637                           jac_data_size_sur, CEED_EVAL_NONE);
638     CeedQFunctionAddOutput(ceed_data->qf_apply_outflow_jacobian, "v", num_comp_q,
639                            CEED_EVAL_INTERP);
640   }
641 
642   // *****************************************************************************
643   // CEED Operator Apply
644   // *****************************************************************************
645   // -- Apply CEED Operator for the geometric data
646   CeedOperatorApply(ceed_data->op_setup_vol, ceed_data->x_coord,
647                     ceed_data->q_data, CEED_REQUEST_IMMEDIATE);
648 
649   // -- Create and apply CEED Composite Operator for the entire domain
650   if (!user->phys->implicit) { // RHS
651     ierr = CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys,
652                                    user->op_rhs_vol, NULL, height, P_sur, Q_sur,
653                                    q_data_size_sur, 0,
654                                    &user->op_rhs, NULL); CHKERRQ(ierr);
655   } else { // IFunction
656     ierr = CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys,
657                                    user->op_ifunction_vol, op_ijacobian_vol,
658                                    height, P_sur, Q_sur,
659                                    q_data_size_sur, jac_data_size_sur,
660                                    &user->op_ifunction,
661                                    op_ijacobian_vol ? &user->op_ijacobian : NULL); CHKERRQ(ierr);
662     if (user->op_ijacobian) {
663       CeedOperatorContextGetFieldLabel(user->op_ijacobian, "ijacobian time shift",
664                                        &user->phys->ijacobian_time_shift_label);
665     }
666 
667   }
668   CeedElemRestrictionDestroy(&elem_restr_jd_i);
669   CeedVectorDestroy(&jac_data);
670   PetscFunctionReturn(0);
671 }
672