xref: /honee/src/setuplibceed.c (revision b30619f6337e4aeddfabbcf7d00859553914ef7b)
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 PetscErrorCode SetupLibceed(Ceed ceed, DM dm, Honee honee, AppCtx app_ctx, ProblemData problem) {
92   const PetscInt      num_comp_q = 5;
93   PetscInt            dim;
94   CeedInt             num_comps_jac_data = problem->num_comps_jac_data, num_comp_x, q_data_size_vol;
95   CeedElemRestriction elem_restr_jd_i    = NULL, elem_restr_qd;
96   CeedVector          jac_data           = NULL, q_data;
97   CeedOperator        op_ifunction_vol = NULL, op_rhs_vol = NULL, op_ijacobian_vol = NULL;
98 
99   PetscFunctionBeginUser;
100   PetscCall(DMGetDimension(dm, &dim));
101   num_comp_x = dim;
102 
103   CeedElemRestriction elem_restr_diff_flux = NULL;
104   CeedVector          div_diff_flux_ceed   = NULL;
105   CeedBasis           basis_diff_flux      = NULL;
106   CeedEvalMode        eval_mode_diff_flux  = -1;
107   {  // Create bases and element restrictions
108     DMLabel  domain_label = NULL;
109     PetscInt label_value = 0, height = 0, dm_field = 0;
110     DM       dm_coord;
111 
112     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
113     PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &honee->basis_q));
114     PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &honee->basis_x));
115 
116     PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, 0, &honee->elem_restr_q));
117     PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &honee->elem_restr_x));
118     if (num_comps_jac_data) {
119       PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, num_comps_jac_data, &elem_restr_jd_i));
120       PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL));
121     }
122 
123     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_q, &honee->q_ceed, NULL));
124     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_q, &honee->q_dot_ceed, NULL));
125     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_q, &honee->g_ceed, NULL));
126     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_x, &honee->x_coord, NULL));
127 
128     {  // -- Copy PETSc coordinate vector into CEED vector
129       Vec X_loc;
130       DM  cdm;
131 
132       PetscCall(DMGetCellCoordinateDM(dm, &cdm));
133       if (cdm) {
134         PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
135       } else {
136         PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
137       }
138       PetscCall(VecScale(X_loc, honee->units->meter));
139       PetscCall(VecCopyPetscToCeed(X_loc, honee->x_coord));
140     }
141 
142     PetscCall(QDataGet(ceed, dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data,
143                        &q_data_size_vol));
144   }
145 
146   if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) {
147     PetscCheck(honee->diff_flux_proj, honee->comm, PETSC_ERR_ARG_WRONGSTATE,
148                "Divergence of diffusive flux projection requested but object not created");
149     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(honee->diff_flux_proj, &elem_restr_diff_flux, &basis_diff_flux, &div_diff_flux_ceed,
150                                                         &eval_mode_diff_flux));
151   }
152 
153   {  // -- Create QFunction for ICs
154     CeedBasis     basis_xc;
155     CeedQFunction qf_ics;
156     CeedOperator  op_ics;
157 
158     PetscCallCeed(ceed, CeedBasisCreateProjection(honee->basis_x, honee->basis_q, &basis_xc));
159     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qf_func_ptr, problem->ics.qf_loc, &qf_ics));
160     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ics, problem->ics.qfctx));
161     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ics, 0));
162     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "x", num_comp_x, CEED_EVAL_INTERP));
163     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
164     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ics, "q0", num_comp_q, CEED_EVAL_NONE));
165 
166     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics));
167     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", honee->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
168     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", honee->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
169     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", honee->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
170     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &honee->phys->ics_time_label));
171     PetscCall(OperatorApplyContextCreate(NULL, dm, honee->ceed, op_ics, honee->x_coord, NULL, NULL, honee->Q_loc, &honee->op_ics_ctx));
172 
173     PetscCallCeed(ceed, CeedBasisDestroy(&basis_xc));
174     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ics));
175     PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics));
176   }
177 
178   if (problem->apply_vol_rhs.qf_func_ptr) {
179     CeedQFunction qf_rhs_vol;
180 
181     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qf_func_ptr, problem->apply_vol_rhs.qf_loc, &qf_rhs_vol));
182     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_vol, problem->apply_vol_rhs.qfctx));
183     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_rhs_vol, 0));
184     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP));
185     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
186     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
187     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP));
188     if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE)
189       PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "div F_diff", honee->diff_flux_proj->num_diff_flux_comps, eval_mode_diff_flux));
190     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP));
191     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
192 
193     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_vol, NULL, NULL, &op_rhs_vol));
194     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
195     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
196     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
197     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "x", honee->elem_restr_x, honee->basis_x, honee->x_coord));
198     if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE)
199       PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "div F_diff", elem_restr_diff_flux, basis_diff_flux, div_diff_flux_ceed));
200     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
201     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
202 
203     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_vol));
204   }
205 
206   if (problem->apply_vol_ifunction.qf_func_ptr) {
207     CeedQFunction qf_ifunction_vol;
208 
209     PetscCallCeed(
210         ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qf_func_ptr, problem->apply_vol_ifunction.qf_loc, &qf_ifunction_vol));
211     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ifunction_vol, problem->apply_vol_ifunction.qfctx));
212     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ifunction_vol, 0));
213     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP));
214     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
215     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP));
216     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
217     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP));
218     if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE)
219       PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "div F_diff", honee->diff_flux_proj->num_diff_flux_comps, eval_mode_diff_flux));
220     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP));
221     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
222     if (num_comps_jac_data) PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "jac_data", num_comps_jac_data, CEED_EVAL_NONE));
223 
224     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ifunction_vol, NULL, NULL, &op_ifunction_vol));
225     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
226     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
227     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q dot", honee->elem_restr_q, honee->basis_q, honee->q_dot_ceed));
228     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
229     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "x", honee->elem_restr_x, honee->basis_x, honee->x_coord));
230     if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE)
231       PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "div F_diff", elem_restr_diff_flux, basis_diff_flux, div_diff_flux_ceed));
232     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
233     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
234     if (num_comps_jac_data) PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
235 
236     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction_vol));
237   }
238 
239   if (problem->apply_vol_ijacobian.qf_func_ptr) {
240     CeedQFunction qf_ijacobian_vol;
241 
242     PetscCallCeed(
243         ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qf_func_ptr, problem->apply_vol_ijacobian.qf_loc, &qf_ijacobian_vol));
244     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfctx));
245     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ijacobian_vol, 0));
246     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP));
247     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD));
248     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
249     if (num_comps_jac_data) PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", num_comps_jac_data, CEED_EVAL_NONE));
250     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP));
251     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
252 
253     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op_ijacobian_vol));
254     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "dq", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
255     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_dq", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
256     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
257     if (num_comps_jac_data) PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
258     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
259     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
260 
261     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol));
262   }
263 
264   // -- Create and apply CEED Composite Operator for the entire domain
265   if (!honee->phys->implicit) {  // RHS
266     CeedOperator op_rhs;
267 
268     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_rhs));
269     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_rhs, op_rhs_vol));
270     for (PetscCount b = 0; b < problem->num_bc_defs; b++) PetscCall(BCDefinitionAddOperators(problem->bc_defs[b], op_rhs, NULL));
271 
272     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, honee->q_ceed, honee->g_ceed, honee->Q_loc, NULL, &honee->op_rhs_ctx));
273 
274     // ----- Get Context Labels for Operator
275     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "solution time", &honee->phys->solution_time_label));
276     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "timestep size", &honee->phys->timestep_size_label));
277 
278     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
279     PetscCall(CreateKSPMass(honee, problem));
280     PetscCheck(app_ctx->sgs_model_type == SGS_MODEL_NONE, honee->comm, PETSC_ERR_SUP, "SGS modeling not implemented for explicit timestepping");
281   } else {  // IFunction
282     CeedOperator op_ijacobian = NULL;
283 
284     // Create Composite Operaters
285     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &honee->op_ifunction));
286     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(honee->op_ifunction, op_ifunction_vol));
287     if (op_ijacobian_vol) {
288       PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_ijacobian));
289       PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_ijacobian, op_ijacobian_vol));
290     }
291     for (PetscCount b = 0; b < problem->num_bc_defs; b++) PetscCall(BCDefinitionAddOperators(problem->bc_defs[b], honee->op_ifunction, op_ijacobian));
292 
293     // ----- Get Context Labels for Operator
294     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(honee->op_ifunction, "solution time", &honee->phys->solution_time_label));
295     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(honee->op_ifunction, "timestep size", &honee->phys->timestep_size_label));
296 
297     if (op_ijacobian) {
298       PetscCall(MatCreateCeed(honee->dm, honee->dm, op_ijacobian, NULL, &honee->mat_ijacobian));
299       PetscCall(MatCeedSetLocalVectors(honee->mat_ijacobian, honee->Q_dot_loc, NULL));
300       PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian));
301     }
302     if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SgsDDSetup(ceed, honee, problem));
303   }
304 
305   if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, dm, honee, problem));
306   if (app_ctx->sgs_train_enable) PetscCall(SGS_DD_TrainingSetup(ceed, honee, problem));
307   if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionSetup(honee, honee->diff_flux_proj));
308 
309   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
310   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data));
311   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
312   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i));
313   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux));
314   PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
315   PetscCallCeed(ceed, CeedVectorDestroy(&div_diff_flux_ceed));
316   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol));
317   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ifunction_vol));
318   PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_vol));
319   PetscFunctionReturn(PETSC_SUCCESS);
320 }
321