xref: /libCEED/examples/fluids/src/setuplibceed.c (revision 830b33163ef00722eb808495168d4e17b4167d37)
1 // Copyright (c) 2017-2024, 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 // @brief Create CeedOperator for unstabilized mass KSP for explicit timestepping
17 static PetscErrorCode CreateKSPMassOperator_Unstabilized(User user, CeedOperator *op_mass) {
18   Ceed                ceed = user->ceed;
19   CeedInt             num_comp_q, q_data_size;
20   CeedQFunction       qf_mass;
21   CeedElemRestriction elem_restr_q, elem_restr_qd_i;
22   CeedBasis           basis_q;
23   CeedVector          q_data;
24 
25   PetscFunctionBeginUser;
26   {  // Get restriction and basis from the RHS function
27     CeedOperator     *sub_ops;
28     CeedOperatorField field;
29     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
30 
31     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops));
32     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
33     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
34     PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q));
35 
36     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field));
37     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i));
38     PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data));
39   }
40 
41   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
42   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size));
43 
44   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
45   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
46   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "u", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
47   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data));
48   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
49 
50   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
51   PetscFunctionReturn(PETSC_SUCCESS);
52 }
53 
54 // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes
55 static PetscErrorCode CreateKSPMass(User user, ProblemData problem) {
56   Ceed         ceed = user->ceed;
57   DM           dm   = user->dm;
58   CeedOperator op_mass;
59 
60   PetscFunctionBeginUser;
61   if (problem->create_mass_operator) PetscCall(problem->create_mass_operator(user, &op_mass));
62   else PetscCall(CreateKSPMassOperator_Unstabilized(user, &op_mass));
63 
64   {  // -- Setup KSP for mass operator
65     Mat      mat_mass;
66     Vec      Zeros_loc;
67     MPI_Comm comm = PetscObjectComm((PetscObject)dm);
68 
69     PetscCall(DMCreateLocalVector(dm, &Zeros_loc));
70     PetscCall(VecZeroEntries(Zeros_loc));
71     PetscCall(MatCeedCreate(dm, dm, op_mass, NULL, &mat_mass));
72     PetscCall(MatCeedSetLocalVectors(mat_mass, Zeros_loc, NULL));
73 
74     PetscCall(KSPCreate(comm, &user->mass_ksp));
75     PetscCall(KSPSetOptionsPrefix(user->mass_ksp, "mass_"));
76     {  // lumped by default
77       PC pc;
78       PetscCall(KSPGetPC(user->mass_ksp, &pc));
79       PetscCall(PCSetType(pc, PCJACOBI));
80       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
81       PetscCall(KSPSetType(user->mass_ksp, KSPPREONLY));
82     }
83     PetscCall(KSPSetFromOptions_WithMatCeed(user->mass_ksp, mat_mass));
84     PetscCall(VecDestroy(&Zeros_loc));
85     PetscCall(MatDestroy(&mat_mass));
86   }
87 
88   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
89   PetscFunctionReturn(PETSC_SUCCESS);
90 }
91 
92 static PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height,
93                                        CeedInt Q_sur, CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedBasis basis_q_sur,
94                                        CeedBasis basis_x_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian, CeedOperator op_apply,
95                                        CeedOperator op_apply_ijacobian) {
96   CeedVector          q_data_sur, jac_data_sur          = NULL;
97   CeedOperator        op_apply_bc, op_apply_bc_jacobian = NULL;
98   CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur = NULL;
99   PetscInt            dm_field = 0;
100 
101   PetscFunctionBeginUser;
102   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &elem_restr_q_sur));
103   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &elem_restr_x_sur));
104   if (jac_data_size_sur > 0) {
105     // State-dependent data will be passed from residual to Jacobian. This will be collocated.
106     PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_sur, &elem_restr_jd_i_sur));
107     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL));
108   }
109 
110   PetscCall(QDataBoundaryGet(ceed, dm, domain_label, label_value, elem_restr_x_sur, basis_x_sur, ceed_data->x_coord, &elem_restr_qd_i_sur,
111                              &q_data_sur, &q_data_size_sur));
112 
113   // CEED Operator for Physics
114   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc));
115   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
116   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
117   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
118   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord));
119   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
120   if (elem_restr_jd_i_sur)
121     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
122 
123   if (qf_apply_bc_jacobian && elem_restr_jd_i_sur) {
124     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian));
125     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
126     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
127     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
128     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord));
129     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
130     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
131   }
132 
133   // Apply Sub-Operator for Physics
134   PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply, op_apply_bc));
135   if (op_apply_bc_jacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply_ijacobian, op_apply_bc_jacobian));
136 
137   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_sur));
138   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data_sur));
139   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur));
140   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur));
141   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i_sur));
142   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i_sur));
143   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc));
144   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc_jacobian));
145   PetscFunctionReturn(PETSC_SUCCESS);
146 }
147 
148 static PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur,
149                                         PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian,
150                                         CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) {
151   PetscFunctionBeginUser;
152   if (apply_bc.qfunction) {
153     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc));
154     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context));
155     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc, 0));
156     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP));
157     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD));
158     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
159     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP));
160     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP));
161     if (jac_data_size_sur) PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
162   }
163   if (apply_bc_jacobian.qfunction) {
164     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian));
165     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context));
166     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc_jacobian, 0));
167     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP));
168     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD));
169     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
170     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP));
171     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
172     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP));
173   }
174   PetscFunctionReturn(PETSC_SUCCESS);
175 }
176 
177 // Utility function to add boundary operators to the composite operator
178 static PetscErrorCode AddBCSubOperators(User user, Ceed ceed, DM dm, SimpleBC bc, ProblemData problem, CeedData ceed_data, CeedOperator op_apply,
179                                         CeedOperator op_apply_ijacobian) {
180   CeedInt       height = 1, num_comp_q, num_comp_x;
181   CeedInt       P_sur = user->app_ctx->degree + 1, Q_sur = P_sur + user->app_ctx->q_extra, dim_sur, q_data_size_sur;
182   const CeedInt jac_data_size_sur = user->phys->implicit ? problem->jac_data_size_sur : 0;
183   PetscInt      dim;
184   DMLabel       face_sets_label;
185   CeedBasis     basis_q_sur, basis_x_sur;
186 
187   PetscFunctionBeginUser;
188   PetscCall(DMGetDimension(dm, &dim));
189   PetscCall(QDataBoundaryGetNumComponents(dm, &q_data_size_sur));
190   dim_sur = dim - height;
191   {  // Get number of components and coordinate dimension from op_apply
192     CeedOperator       *sub_ops;
193     CeedOperatorField   field;
194     PetscInt            sub_op_index = 0;  // will be 0 for the volume op
195     CeedElemRestriction elem_restr_q, elem_restr_x;
196 
197     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(op_apply, &sub_ops));
198     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
199     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
200     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
201 
202     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "x", &field));
203     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_x));
204     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
205   }
206 
207   {  // Get bases
208     DM dm_coord;
209 
210     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
211     DMLabel  label       = NULL;
212     PetscInt label_value = 0;
213     PetscInt field       = 0;
214     PetscCall(CreateBasisFromPlex(ceed, dm, label, label_value, height, field, &basis_q_sur));
215     PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, label_value, height, field, &basis_x_sur));
216   }
217 
218   PetscCall(DMGetLabel(dm, "Face Sets", &face_sets_label));
219 
220   {  // --- Create Sub-Operator for inflow boundaries
221     CeedQFunction qf_apply_inflow = NULL, qf_apply_inflow_jacobian = NULL;
222 
223     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow,
224                                 problem->apply_inflow_jacobian, &qf_apply_inflow, &qf_apply_inflow_jacobian));
225     for (CeedInt i = 0; i < bc->num_inflow; i++) {
226       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->inflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, basis_q_sur,
227                                  basis_x_sur, qf_apply_inflow, qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian));
228     }
229     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow));
230     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow_jacobian));
231   }
232 
233   {  // --- Create Sub-Operator for outflow boundaries
234     CeedQFunction qf_apply_outflow = NULL, qf_apply_outflow_jacobian = NULL;
235 
236     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow,
237                                 problem->apply_outflow_jacobian, &qf_apply_outflow, &qf_apply_outflow_jacobian));
238     for (CeedInt i = 0; i < bc->num_outflow; i++) {
239       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
240                                  basis_q_sur, basis_x_sur, qf_apply_outflow, qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian));
241     }
242     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow));
243     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow_jacobian));
244   }
245 
246   {  // --- Create Sub-Operator for freestream boundaries
247     CeedQFunction qf_apply_freestream = NULL, qf_apply_freestream_jacobian = NULL;
248 
249     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream,
250                                 problem->apply_freestream_jacobian, &qf_apply_freestream, &qf_apply_freestream_jacobian));
251     for (CeedInt i = 0; i < bc->num_freestream; i++) {
252       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
253                                  basis_q_sur, basis_x_sur, qf_apply_freestream, qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian));
254     }
255     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream));
256     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream_jacobian));
257   }
258 
259   {  // --- Create Sub-Operator for slip boundaries
260     CeedQFunction qf_apply_slip = NULL, qf_apply_slip_jacobian = NULL;
261 
262     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_slip,
263                                 problem->apply_slip_jacobian, &qf_apply_slip, &qf_apply_slip_jacobian));
264     for (CeedInt i = 0; i < bc->num_slip; i++) {
265       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->slips[i], height, Q_sur, q_data_size_sur, jac_data_size_sur, basis_q_sur,
266                                  basis_x_sur, qf_apply_slip, qf_apply_slip_jacobian, op_apply, op_apply_ijacobian));
267     }
268     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip));
269     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip_jacobian));
270   }
271 
272   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q_sur));
273   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_sur));
274   PetscFunctionReturn(PETSC_SUCCESS);
275 }
276 
277 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData problem, SimpleBC bc) {
278   const PetscInt      num_comp_q = 5;
279   const CeedInt       dim = problem->dim, num_comp_x = problem->dim;
280   CeedInt             jac_data_size_vol = num_comp_q + 6 + 3;
281   CeedElemRestriction elem_restr_jd_i;
282   CeedVector          jac_data;
283   CeedOperator        op_ifunction_vol = NULL, op_rhs_vol = NULL, op_ijacobian_vol = NULL;
284 
285   PetscFunctionBeginUser;
286 
287   if (problem->apply_vol_ifunction.qfunction && problem->uses_newtonian) {
288     NewtonianIdealGasContext gas;
289     PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas));
290     jac_data_size_vol += (gas->idl_enable ? 1 : 0);
291     PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas));
292   }
293 
294   {  // Create bases and element restrictions
295     DMLabel  domain_label = NULL;
296     PetscInt label_value = 0, height = 0, dm_field = 0;
297     DM       dm_coord;
298 
299     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
300     PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &ceed_data->basis_q));
301     PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &ceed_data->basis_x));
302 
303     PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, 0, &ceed_data->elem_restr_q));
304     PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &ceed_data->elem_restr_x));
305     PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_vol, &elem_restr_jd_i));
306 
307     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL));
308     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL));
309     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL));
310     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL));
311     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL));
312 
313     {  // -- Copy PETSc coordinate vector into CEED vector
314       Vec X_loc;
315       DM  cdm;
316 
317       PetscCall(DMGetCellCoordinateDM(dm, &cdm));
318       if (cdm) {
319         PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
320       } else {
321         PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
322       }
323       PetscCall(VecScale(X_loc, problem->dm_scale));
324       PetscCall(VecCopyPetscToCeed(X_loc, ceed_data->x_coord));
325     }
326 
327     PetscCall(QDataGet(ceed, dm, domain_label, label_value, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord,
328                        &ceed_data->elem_restr_qd_i, &ceed_data->q_data, &problem->q_data_size_vol));
329   }
330 
331   {  // -- Create QFunction for ICs
332     CeedBasis     basis_xc;
333     CeedQFunction qf_ics;
334     CeedOperator  op_ics;
335 
336     PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_xc));
337     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &qf_ics));
338     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ics, problem->ics.qfunction_context));
339     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ics, 0));
340     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "x", num_comp_x, CEED_EVAL_INTERP));
341     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
342     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ics, "q0", num_comp_q, CEED_EVAL_NONE));
343 
344     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics));
345     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
346     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
347     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
348     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label));
349     PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx));
350 
351     PetscCallCeed(ceed, CeedBasisDestroy(&basis_xc));
352     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ics));
353     PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics));
354   }
355 
356   if (problem->apply_vol_rhs.qfunction) {
357     CeedQFunction qf_rhs_vol;
358 
359     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &qf_rhs_vol));
360     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_vol, problem->apply_vol_rhs.qfunction_context));
361     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_rhs_vol, 0));
362     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP));
363     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
364     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "qdata", problem->q_data_size_vol, CEED_EVAL_NONE));
365     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP));
366     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP));
367     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
368 
369     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_vol, NULL, NULL, &op_rhs_vol));
370     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
371     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
372     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
373     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
374     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
375     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
376 
377     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_vol));
378   }
379 
380   if (problem->apply_vol_ifunction.qfunction) {
381     CeedQFunction qf_ifunction_vol;
382 
383     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc,
384                                                     &qf_ifunction_vol));
385     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context));
386     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ifunction_vol, 0));
387     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP));
388     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
389     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP));
390     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "qdata", problem->q_data_size_vol, CEED_EVAL_NONE));
391     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP));
392     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP));
393     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
394     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
395 
396     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ifunction_vol, NULL, NULL, &op_ifunction_vol));
397     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
398     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
399     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed));
400     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
401     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
402     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
403     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
404     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
405 
406     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction_vol));
407   }
408 
409   if (problem->apply_vol_ijacobian.qfunction) {
410     CeedQFunction qf_ijacobian_vol;
411 
412     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc,
413                                                     &qf_ijacobian_vol));
414     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context));
415     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ijacobian_vol, 0));
416     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP));
417     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD));
418     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", problem->q_data_size_vol, CEED_EVAL_NONE));
419     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
420     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP));
421     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
422 
423     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op_ijacobian_vol));
424     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
425     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
426     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
427     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
428     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
429     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
430 
431     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol));
432   }
433 
434   // -- Create and apply CEED Composite Operator for the entire domain
435   if (!user->phys->implicit) {  // RHS
436     CeedOperator op_rhs;
437 
438     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_rhs));
439     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_rhs, op_rhs_vol));
440     PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, op_rhs, NULL));
441 
442     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx));
443 
444     // ----- Get Context Labels for Operator
445     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "solution time", &user->phys->solution_time_label));
446     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "timestep size", &user->phys->timestep_size_label));
447 
448     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
449     PetscCall(CreateKSPMass(user, problem));
450   } else {  // IFunction
451     CeedOperator op_ijacobian = NULL;
452 
453     // Create Composite Operaters
454     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &user->op_ifunction));
455     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(user->op_ifunction, op_ifunction_vol));
456     if (op_ijacobian_vol) {
457       PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_ijacobian));
458       PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_ijacobian, op_ijacobian_vol));
459     }
460     PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, user->op_ifunction, op_ijacobian));
461 
462     // ----- Get Context Labels for Operator
463     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "solution time", &user->phys->solution_time_label));
464     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "timestep size", &user->phys->timestep_size_label));
465 
466     if (op_ijacobian) {
467       PetscCall(MatCeedCreate(user->dm, user->dm, op_ijacobian, NULL, &user->mat_ijacobian));
468       PetscCall(MatCeedSetLocalVectors(user->mat_ijacobian, user->Q_dot_loc, NULL));
469       PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian));
470     }
471   }
472 
473   if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc));
474   if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem));
475   if (app_ctx->diff_filter_monitor && !user->diff_filter) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem));
476 
477   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data));
478   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i));
479   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol));
480   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ifunction_vol));
481   PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_vol));
482   PetscFunctionReturn(PETSC_SUCCESS);
483 }
484