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