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