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