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