xref: /libCEED/examples/fluids/src/setuplibceed.c (revision e9c36be085200ad84b427d3f41921f3626bb31f1)
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(KSPSetFromOptions(user->mass_ksp));
85     PetscCall(VecDestroy(&Zeros_loc));
86     PetscCall(MatDestroy(&mat_mass));
87   }
88 
89   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
90   PetscFunctionReturn(PETSC_SUCCESS);
91 }
92 
93 static PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height,
94                                        CeedInt Q_sur, CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedBasis basis_q_sur,
95                                        CeedBasis basis_x_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian, CeedOperator op_apply,
96                                        CeedOperator op_apply_ijacobian) {
97   CeedVector          q_data_sur, jac_data_sur          = NULL;
98   CeedOperator        op_apply_bc, op_apply_bc_jacobian = NULL;
99   CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur = NULL;
100   PetscInt            dm_field = 0;
101 
102   PetscFunctionBeginUser;
103 
104   // ---- CEED Restriction
105   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &elem_restr_q_sur));
106   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &elem_restr_x_sur));
107   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, q_data_size_sur, &elem_restr_qd_i_sur));
108   if (jac_data_size_sur > 0) {
109     // State-dependent data will be passed from residual to Jacobian. This will be collocated.
110     PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_sur, &elem_restr_jd_i_sur));
111     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL));
112   }
113 
114   {  // Create q_data_sur vector
115     CeedOperator op_setup_sur;
116 
117     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_qd_i_sur, &q_data_sur, NULL));
118 
119     PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur));
120     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, basis_x_sur, CEED_VECTOR_ACTIVE));
121     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, basis_x_sur, CEED_VECTOR_NONE));
122     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
123 
124     PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE));
125     PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur));
126   }
127 
128   // ----- CEED Operator for Physics
129   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc));
130   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
131   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
132   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
133   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord));
134   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
135   if (elem_restr_jd_i_sur)
136     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
137 
138   if (qf_apply_bc_jacobian && elem_restr_jd_i_sur) {
139     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian));
140     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
141     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
142     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
143     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord));
144     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
145     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
146   }
147 
148   // ----- Apply Sub-Operator for Physics
149   PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply, op_apply_bc));
150   if (op_apply_bc_jacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply_ijacobian, op_apply_bc_jacobian));
151 
152   // ----- Cleanup
153   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_sur));
154   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data_sur));
155   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur));
156   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur));
157   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i_sur));
158   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i_sur));
159   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc));
160   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc_jacobian));
161   PetscFunctionReturn(PETSC_SUCCESS);
162 }
163 
164 static PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur,
165                                         PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian,
166                                         CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) {
167   PetscFunctionBeginUser;
168   if (apply_bc.qfunction) {
169     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc));
170     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context));
171     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc, 0));
172     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP));
173     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD));
174     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
175     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP));
176     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP));
177     if (jac_data_size_sur) PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
178   }
179   if (apply_bc_jacobian.qfunction) {
180     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian));
181     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context));
182     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc_jacobian, 0));
183     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP));
184     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD));
185     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
186     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP));
187     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
188     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP));
189   }
190   PetscFunctionReturn(PETSC_SUCCESS);
191 }
192 
193 // Utility function to create CEED Composite Operator for the entire domain
194 static PetscErrorCode AddBCSubOperators(User user, Ceed ceed, DM dm, SimpleBC bc, ProblemData problem, CeedData ceed_data, CeedOperator op_apply,
195                                         CeedOperator op_apply_ijacobian) {
196   CeedInt       height = 1, num_comp_q, num_comp_x;
197   CeedInt       dim_sur, P_sur = user->app_ctx->degree + 1, Q_sur = P_sur + user->app_ctx->q_extra;
198   const CeedInt q_data_size_sur = problem->q_data_size_sur, jac_data_size_sur = user->phys->implicit ? problem->jac_data_size_sur : 0;
199   PetscInt      dim;
200   DMLabel       face_sets_label;
201   CeedBasis     basis_q_sur, basis_x_sur;
202 
203   PetscFunctionBeginUser;
204   PetscCall(DMGetDimension(dm, &dim));
205   dim_sur = dim - height;
206   {  // Get number of components and coordinate dimension from op_apply
207     CeedOperator       *sub_ops;
208     CeedOperatorField   field;
209     PetscInt            sub_op_index = 0;  // will be 0 for the volume op
210     CeedElemRestriction elem_restr_q, elem_restr_x;
211 
212     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(op_apply, &sub_ops));
213     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
214     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
215     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
216 
217     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "x", &field));
218     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_x));
219     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
220   }
221 
222   {  // Get bases
223     DM dm_coord;
224 
225     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
226     DMLabel  label       = NULL;
227     PetscInt label_value = 0;
228     PetscInt field       = 0;  // Still want the normal, default field
229     PetscCall(CreateBasisFromPlex(ceed, dm, label, label_value, height, field, &basis_q_sur));
230     PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, label_value, height, field, &basis_x_sur));
231   }
232 
233   PetscCall(DMGetLabel(dm, "Face Sets", &face_sets_label));
234 
235   {  // -- Create QFunction for quadrature data
236     PetscCallCeed(ceed,
237                   CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, problem->setup_sur.qfunction_loc, &ceed_data->qf_setup_sur));
238     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_setup_sur, problem->setup_sur.qfunction_context));
239     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(ceed_data->qf_setup_sur, 0));
240     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x * dim_sur, CEED_EVAL_GRAD));
241     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT));
242     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
243   }
244 
245   {  // --- Create Sub-Operator for inflow boundaries
246     CeedQFunction qf_apply_inflow = NULL, qf_apply_inflow_jacobian = NULL;
247 
248     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow,
249                                 problem->apply_inflow_jacobian, &qf_apply_inflow, &qf_apply_inflow_jacobian));
250     for (CeedInt i = 0; i < bc->num_inflow; i++) {
251       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,
252                                  basis_x_sur, qf_apply_inflow, qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian));
253     }
254     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow));
255     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow_jacobian));
256   }
257 
258   {  // --- Create Sub-Operator for outflow boundaries
259     CeedQFunction qf_apply_outflow = NULL, qf_apply_outflow_jacobian = NULL;
260 
261     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow,
262                                 problem->apply_outflow_jacobian, &qf_apply_outflow, &qf_apply_outflow_jacobian));
263     for (CeedInt i = 0; i < bc->num_outflow; i++) {
264       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
265                                  basis_q_sur, basis_x_sur, qf_apply_outflow, qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian));
266     }
267     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow));
268     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow_jacobian));
269   }
270 
271   {  // --- Create Sub-Operator for freestream boundaries
272     CeedQFunction qf_apply_freestream = NULL, qf_apply_freestream_jacobian = NULL;
273 
274     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream,
275                                 problem->apply_freestream_jacobian, &qf_apply_freestream, &qf_apply_freestream_jacobian));
276     for (CeedInt i = 0; i < bc->num_freestream; i++) {
277       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
278                                  basis_q_sur, basis_x_sur, qf_apply_freestream, qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian));
279     }
280     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream));
281     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream_jacobian));
282   }
283 
284   {  // --- Create Sub-Operator for slip boundaries
285     CeedQFunction qf_apply_slip = NULL, qf_apply_slip_jacobian = NULL;
286 
287     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_slip,
288                                 problem->apply_slip_jacobian, &qf_apply_slip, &qf_apply_slip_jacobian));
289     for (CeedInt i = 0; i < bc->num_slip; i++) {
290       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,
291                                  basis_x_sur, qf_apply_slip, qf_apply_slip_jacobian, op_apply, op_apply_ijacobian));
292     }
293     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip));
294     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip_jacobian));
295   }
296 
297   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q_sur));
298   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_sur));
299   PetscFunctionReturn(PETSC_SUCCESS);
300 }
301 
302 PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData problem, SimpleBC bc) {
303   PetscFunctionBeginUser;
304   const PetscInt num_comp_q = 5;
305   const CeedInt  dim = problem->dim, num_comp_x = problem->dim, q_data_size_vol = problem->q_data_size_vol;
306   CeedInt        jac_data_size_vol = num_comp_q + 6 + 3;
307 
308   if (problem->apply_vol_ifunction.qfunction && problem->uses_newtonian) {
309     NewtonianIdealGasContext gas;
310     PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas));
311     jac_data_size_vol += (gas->idl_enable ? 1 : 0);
312     PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas));
313   }
314 
315   CeedElemRestriction elem_restr_jd_i;
316   CeedVector          jac_data;
317   CeedInt             num_qpts;
318   DMLabel             domain_label = NULL;
319   PetscInt            label_value = 0, height = 0, dm_field = 0;
320 
321   DM dm_coord;
322   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
323 
324   PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &ceed_data->basis_q));
325   PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &ceed_data->basis_x));
326   PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(ceed_data->basis_q, &num_qpts));
327 
328   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, 0, &ceed_data->elem_restr_q));
329   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &ceed_data->elem_restr_x));
330   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, q_data_size_vol, &ceed_data->elem_restr_qd_i));
331   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_vol, &elem_restr_jd_i));
332 
333   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL));
334   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL));
335   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL));
336   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL));
337   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_qd_i, &ceed_data->q_data, NULL));
338   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL));
339 
340   {  // -- Copy PETSc coordinate vector into CEED vector
341     Vec X_loc;
342     DM  cdm;
343 
344     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
345     if (cdm) {
346       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
347     } else {
348       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
349     }
350     PetscCall(VecScale(X_loc, problem->dm_scale));
351     PetscCall(VecCopyPetscToCeed(X_loc, ceed_data->x_coord));
352   }
353 
354   {  // -- Create quadrature data
355     CeedQFunction qf_setup_vol;
356     CeedOperator  op_setup_vol;
357 
358     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, problem->setup_vol.qfunction_loc, &qf_setup_vol));
359     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_setup_vol, problem->setup_vol.qfunction_context));
360     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_setup_vol, 0));
361     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_vol, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
362     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT));
363     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
364 
365     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup_vol, NULL, NULL, &op_setup_vol));
366     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "dx", ceed_data->elem_restr_x, ceed_data->basis_x, CEED_VECTOR_ACTIVE));
367     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE));
368     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
369 
370     PetscCallCeed(ceed, CeedOperatorApply(op_setup_vol, ceed_data->x_coord, ceed_data->q_data, CEED_REQUEST_IMMEDIATE));
371 
372     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup_vol));
373     PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_vol));
374   }
375 
376   {  // -- Create QFunction for ICs
377     CeedBasis     basis_xc;
378     CeedQFunction qf_ics;
379     CeedOperator  op_ics;
380 
381     PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_xc));
382     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &qf_ics));
383     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ics, problem->ics.qfunction_context));
384     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ics, 0));
385     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "x", num_comp_x, CEED_EVAL_INTERP));
386     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
387     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ics, "q0", num_comp_q, CEED_EVAL_NONE));
388 
389     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics));
390     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
391     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
392     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
393     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label));
394     PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx));
395 
396     PetscCallCeed(ceed, CeedBasisDestroy(&basis_xc));
397     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ics));
398     PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics));
399   }
400 
401   if (problem->apply_vol_rhs.qfunction) {
402     CeedQFunction qf_rhs_vol;
403     CeedOperator  op;
404 
405     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &qf_rhs_vol));
406     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_vol, problem->apply_vol_rhs.qfunction_context));
407     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_rhs_vol, 0));
408     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP));
409     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
410     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
411     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP));
412     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP));
413     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
414 
415     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_vol, NULL, NULL, &op));
416     PetscCallCeed(ceed, CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
417     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
418     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
419     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
420     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
421     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
422     user->op_rhs_vol = op;
423 
424     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_vol));
425   }
426 
427   if (problem->apply_vol_ifunction.qfunction) {
428     CeedQFunction qf_ifunction_vol;
429     CeedOperator  op;
430 
431     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc,
432                                                     &qf_ifunction_vol));
433     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context));
434     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ifunction_vol, 0));
435     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP));
436     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
437     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP));
438     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
439     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP));
440     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP));
441     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
442     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
443 
444     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ifunction_vol, NULL, NULL, &op));
445     PetscCallCeed(ceed, CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
446     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
447     PetscCallCeed(ceed, CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed));
448     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
449     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
450     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
451     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
452     PetscCallCeed(ceed, CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
453 
454     user->op_ifunction_vol = op;
455     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction_vol));
456   }
457 
458   CeedOperator op_ijacobian_vol = NULL;
459   if (problem->apply_vol_ijacobian.qfunction) {
460     CeedQFunction qf_ijacobian_vol;
461     CeedOperator  op;
462 
463     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc,
464                                                     &qf_ijacobian_vol));
465     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context));
466     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ijacobian_vol, 0));
467     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP));
468     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD));
469     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
470     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
471     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP));
472     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
473 
474     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op));
475     PetscCallCeed(ceed, CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
476     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
477     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
478     PetscCallCeed(ceed, CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
479     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
480     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
481     op_ijacobian_vol = op;
482 
483     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol));
484   }
485 
486   // -- Create and apply CEED Composite Operator for the entire domain
487   if (!user->phys->implicit) {  // RHS
488     CeedOperator op_rhs;
489 
490     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_rhs));
491     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_rhs, user->op_rhs_vol));
492     PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, op_rhs, NULL));
493 
494     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx));
495 
496     // ----- Get Context Labels for Operator
497     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "solution time", &user->phys->solution_time_label));
498     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "timestep size", &user->phys->timestep_size_label));
499 
500     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
501     PetscCall(CreateKSPMass(user, problem));
502     PetscCheck(app_ctx->sgs_model_type == SGS_MODEL_NONE, user->comm, PETSC_ERR_SUP, "SGS modeling not implemented for explicit timestepping");
503   } else {  // IFunction
504     CeedOperator op_ijacobian = NULL;
505 
506     // Create Composite Operaters
507     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &user->op_ifunction));
508     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(user->op_ifunction, user->op_ifunction_vol));
509     if (op_ijacobian_vol) {
510       PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_ijacobian));
511       PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_ijacobian, op_ijacobian_vol));
512     }
513     PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, user->op_ifunction, op_ijacobian));
514 
515     // ----- Get Context Labels for Operator
516     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "solution time", &user->phys->solution_time_label));
517     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "timestep size", &user->phys->timestep_size_label));
518 
519     if (op_ijacobian) {
520       PetscCall(MatCeedCreate(user->dm, user->dm, op_ijacobian, NULL, &user->mat_ijacobian));
521       PetscCall(MatCeedSetLocalVectors(user->mat_ijacobian, user->Q_dot_loc, NULL));
522       PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ijacobian, "ijacobian time shift", &user->phys->ijacobian_time_shift_label));
523       PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian));
524     }
525     if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SgsDDSetup(ceed, user, ceed_data, problem));
526   }
527 
528   if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc));
529   if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem));
530   if (app_ctx->diff_filter_monitor && !user->diff_filter) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem));
531   if (app_ctx->sgs_train_enable) PetscCall(SGS_DD_TrainingSetup(ceed, user, ceed_data, problem));
532 
533   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i));
534   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol));
535   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data));
536   PetscFunctionReturn(PETSC_SUCCESS);
537 }
538