xref: /honee/src/setuplibceed.c (revision 8a02cd4c87f0d986fdd337778674e558be30ecff)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3a515125bSLeila Ghaffari 
4a515125bSLeila Ghaffari /// @file
5a515125bSLeila Ghaffari /// Setup libCEED for Navier-Stokes example using PETSc
6a515125bSLeila Ghaffari 
7e419654dSJeremy L Thompson #include <ceed.h>
8e419654dSJeremy L Thompson #include <petscdmplex.h>
9e419654dSJeremy L Thompson 
10149fb536SJames Wright #include <navierstokes.h>
11a515125bSLeila Ghaffari 
12a78efa86SJames Wright // @brief Create CeedOperator for unstabilized mass KSP for explicit timestepping
13a78efa86SJames Wright static PetscErrorCode CreateKSPMassOperator_Unstabilized(User user, CeedOperator *op_mass) {
14a78efa86SJames Wright   Ceed                ceed = user->ceed;
15a78efa86SJames Wright   CeedInt             num_comp_q, q_data_size;
16a78efa86SJames Wright   CeedQFunction       qf_mass;
17a78efa86SJames Wright   CeedElemRestriction elem_restr_q, elem_restr_qd_i;
18a78efa86SJames Wright   CeedBasis           basis_q;
19a78efa86SJames Wright   CeedVector          q_data;
20a78efa86SJames Wright 
21a78efa86SJames Wright   PetscFunctionBeginUser;
22a78efa86SJames Wright   {  // Get restriction and basis from the RHS function
23a78efa86SJames Wright     CeedOperator     *sub_ops;
24a78efa86SJames Wright     CeedOperatorField field;
25a78efa86SJames Wright     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
26a78efa86SJames Wright 
27a78efa86SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops));
28a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
29a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
30a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q));
31a78efa86SJames Wright 
32a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field));
33a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i));
34a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data));
35a78efa86SJames Wright   }
36a78efa86SJames Wright 
37a78efa86SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
38a78efa86SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size));
39a78efa86SJames Wright 
40a78efa86SJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
41a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
42a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "u", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
43a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data));
44a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
45a78efa86SJames Wright 
46a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
47a78efa86SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
48a78efa86SJames Wright }
49a78efa86SJames Wright 
50a78efa86SJames Wright // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes
51991aef52SJames Wright static PetscErrorCode CreateKSPMass(User user, ProblemData problem) {
52a78efa86SJames Wright   Ceed         ceed = user->ceed;
53a78efa86SJames Wright   DM           dm   = user->dm;
54a78efa86SJames Wright   CeedOperator op_mass;
55a78efa86SJames Wright 
56a78efa86SJames Wright   PetscFunctionBeginUser;
57a78efa86SJames Wright   if (problem->create_mass_operator) PetscCall(problem->create_mass_operator(user, &op_mass));
58a78efa86SJames Wright   else PetscCall(CreateKSPMassOperator_Unstabilized(user, &op_mass));
59a78efa86SJames Wright 
60a78efa86SJames Wright   {  // -- Setup KSP for mass operator
61a78efa86SJames Wright     Mat      mat_mass;
62a78efa86SJames Wright     Vec      Zeros_loc;
63a78efa86SJames Wright     MPI_Comm comm = PetscObjectComm((PetscObject)dm);
64a78efa86SJames Wright 
65a78efa86SJames Wright     PetscCall(DMCreateLocalVector(dm, &Zeros_loc));
66a78efa86SJames Wright     PetscCall(VecZeroEntries(Zeros_loc));
67000d2032SJeremy L Thompson     PetscCall(MatCreateCeed(dm, dm, op_mass, NULL, &mat_mass));
68a78efa86SJames Wright     PetscCall(MatCeedSetLocalVectors(mat_mass, Zeros_loc, NULL));
69a78efa86SJames Wright 
70a78efa86SJames Wright     PetscCall(KSPCreate(comm, &user->mass_ksp));
71a78efa86SJames Wright     PetscCall(KSPSetOptionsPrefix(user->mass_ksp, "mass_"));
72a78efa86SJames Wright     {  // lumped by default
73a78efa86SJames Wright       PC pc;
74a78efa86SJames Wright       PetscCall(KSPGetPC(user->mass_ksp, &pc));
75a78efa86SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
76a78efa86SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
77a78efa86SJames Wright       PetscCall(KSPSetType(user->mass_ksp, KSPPREONLY));
78a78efa86SJames Wright     }
79a78efa86SJames Wright     PetscCall(KSPSetFromOptions_WithMatCeed(user->mass_ksp, mat_mass));
80a78efa86SJames Wright     PetscCall(VecDestroy(&Zeros_loc));
81a78efa86SJames Wright     PetscCall(MatDestroy(&mat_mass));
82a78efa86SJames Wright   }
83a78efa86SJames Wright 
84a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
85a78efa86SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
86a78efa86SJames Wright }
87a78efa86SJames Wright 
88866f9b4aSJames Wright static PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height,
89866f9b4aSJames Wright                                        CeedInt Q_sur, CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedBasis basis_q_sur,
90866f9b4aSJames Wright                                        CeedBasis basis_x_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian, CeedOperator op_apply,
91866f9b4aSJames Wright                                        CeedOperator op_apply_ijacobian) {
9215c18037SJames Wright   CeedVector          q_data_sur, jac_data_sur          = NULL;
93866f9b4aSJames Wright   CeedOperator        op_apply_bc, op_apply_bc_jacobian = NULL;
9415c18037SJames Wright   CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur = NULL;
95866f9b4aSJames Wright   PetscInt            dm_field = 0;
96368c645fSJames Wright 
9706f41313SJames Wright   PetscFunctionBeginUser;
9815c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &elem_restr_q_sur));
9915c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &elem_restr_x_sur));
100368c645fSJames Wright   if (jac_data_size_sur > 0) {
101368c645fSJames Wright     // State-dependent data will be passed from residual to Jacobian. This will be collocated.
10215c18037SJames Wright     PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_sur, &elem_restr_jd_i_sur));
103b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL));
104368c645fSJames Wright   }
105368c645fSJames Wright 
106c864c5abSJames Wright   PetscCall(QDataBoundaryGet(ceed, dm, domain_label, label_value, elem_restr_x_sur, basis_x_sur, ceed_data->x_coord, &elem_restr_qd_i_sur,
107c864c5abSJames Wright                              &q_data_sur, &q_data_size_sur));
108866f9b4aSJames Wright 
109bbc90103SJames Wright   // CEED Operator for Physics
110b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc));
111866f9b4aSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
112866f9b4aSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
11358e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
114866f9b4aSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord));
115866f9b4aSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
116b4c37c5cSJames Wright   if (elem_restr_jd_i_sur)
11758e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
118368c645fSJames Wright 
1194b96a86bSJames Wright   if (qf_apply_bc_jacobian && elem_restr_jd_i_sur) {
120b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian));
121866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
122866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
12358e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
124866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord));
12558e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
126866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
127368c645fSJames Wright   }
128368c645fSJames Wright 
129bbc90103SJames Wright   // Apply Sub-Operator for Physics
130866f9b4aSJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply, op_apply_bc));
131866f9b4aSJames Wright   if (op_apply_bc_jacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply_ijacobian, op_apply_bc_jacobian));
132368c645fSJames Wright 
133b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_sur));
134b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data_sur));
135b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur));
136b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur));
137b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i_sur));
138b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i_sur));
139b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc));
140b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc_jacobian));
141d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
142368c645fSJames Wright }
143368c645fSJames Wright 
144866f9b4aSJames Wright static PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur,
1452b916ea7SJeremy L Thompson                                         PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian,
14625988f00SJames Wright                                         CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) {
14725988f00SJames Wright   PetscFunctionBeginUser;
148e07531f7SJames Wright   if (apply_bc.qf_func_ptr) {
149e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc.qf_func_ptr, apply_bc.qf_loc, qf_apply_bc));
150e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfctx));
151ebfabadfSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc, 0));
152b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP));
153b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD));
154b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
155b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP));
156b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP));
157b4c37c5cSJames Wright     if (jac_data_size_sur) PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
15825988f00SJames Wright   }
159e07531f7SJames Wright   if (apply_bc_jacobian.qf_func_ptr) {
160e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qf_func_ptr, apply_bc_jacobian.qf_loc, qf_apply_bc_jacobian));
161e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfctx));
162ebfabadfSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc_jacobian, 0));
163b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP));
164b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD));
165b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
166b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP));
167b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
168b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP));
16925988f00SJames Wright   }
170d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
17125988f00SJames Wright }
17225988f00SJames Wright 
173bbc90103SJames Wright // Utility function to add boundary operators to the composite operator
174866f9b4aSJames Wright static PetscErrorCode AddBCSubOperators(User user, Ceed ceed, DM dm, SimpleBC bc, ProblemData problem, CeedData ceed_data, CeedOperator op_apply,
175866f9b4aSJames Wright                                         CeedOperator op_apply_ijacobian) {
176866f9b4aSJames Wright   CeedInt       height = 1, num_comp_q, num_comp_x;
177c864c5abSJames Wright   CeedInt       P_sur = user->app_ctx->degree + 1, Q_sur = P_sur + user->app_ctx->q_extra, dim_sur, q_data_size_sur;
178c864c5abSJames Wright   const CeedInt jac_data_size_sur = user->phys->implicit ? problem->jac_data_size_sur : 0;
179866f9b4aSJames Wright   PetscInt      dim;
180866f9b4aSJames Wright   DMLabel       face_sets_label;
181866f9b4aSJames Wright   CeedBasis     basis_q_sur, basis_x_sur;
182866f9b4aSJames Wright 
183866f9b4aSJames Wright   PetscFunctionBeginUser;
184866f9b4aSJames Wright   PetscCall(DMGetDimension(dm, &dim));
185c864c5abSJames Wright   PetscCall(QDataBoundaryGetNumComponents(dm, &q_data_size_sur));
186866f9b4aSJames Wright   dim_sur = dim - height;
187866f9b4aSJames Wright   {  // Get number of components and coordinate dimension from op_apply
188866f9b4aSJames Wright     CeedOperator       *sub_ops;
189866f9b4aSJames Wright     CeedOperatorField   field;
190866f9b4aSJames Wright     PetscInt            sub_op_index = 0;  // will be 0 for the volume op
191866f9b4aSJames Wright     CeedElemRestriction elem_restr_q, elem_restr_x;
192866f9b4aSJames Wright 
193866f9b4aSJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(op_apply, &sub_ops));
194866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
195866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
196866f9b4aSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
197866f9b4aSJames Wright 
198866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "x", &field));
199866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_x));
200866f9b4aSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
201866f9b4aSJames Wright   }
202866f9b4aSJames Wright 
203866f9b4aSJames Wright   {  // Get bases
204866f9b4aSJames Wright     DM dm_coord;
205866f9b4aSJames Wright 
206866f9b4aSJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
207866f9b4aSJames Wright     DMLabel  label       = NULL;
208866f9b4aSJames Wright     PetscInt label_value = 0;
209bbc90103SJames Wright     PetscInt field       = 0;
210866f9b4aSJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm, label, label_value, height, field, &basis_q_sur));
211866f9b4aSJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, label_value, height, field, &basis_x_sur));
212866f9b4aSJames Wright   }
213866f9b4aSJames Wright 
214866f9b4aSJames Wright   PetscCall(DMGetLabel(dm, "Face Sets", &face_sets_label));
215866f9b4aSJames Wright 
216866f9b4aSJames Wright   {  // --- Create Sub-Operator for inflow boundaries
217866f9b4aSJames Wright     CeedQFunction qf_apply_inflow = NULL, qf_apply_inflow_jacobian = NULL;
218866f9b4aSJames Wright 
219866f9b4aSJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow,
220866f9b4aSJames Wright                                 problem->apply_inflow_jacobian, &qf_apply_inflow, &qf_apply_inflow_jacobian));
221866f9b4aSJames Wright     for (CeedInt i = 0; i < bc->num_inflow; i++) {
222866f9b4aSJames Wright       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,
223866f9b4aSJames Wright                                  basis_x_sur, qf_apply_inflow, qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian));
224866f9b4aSJames Wright     }
225866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow));
226866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow_jacobian));
227866f9b4aSJames Wright   }
228866f9b4aSJames Wright 
229866f9b4aSJames Wright   {  // --- Create Sub-Operator for outflow boundaries
230866f9b4aSJames Wright     CeedQFunction qf_apply_outflow = NULL, qf_apply_outflow_jacobian = NULL;
231866f9b4aSJames Wright 
232866f9b4aSJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow,
233866f9b4aSJames Wright                                 problem->apply_outflow_jacobian, &qf_apply_outflow, &qf_apply_outflow_jacobian));
234866f9b4aSJames Wright     for (CeedInt i = 0; i < bc->num_outflow; i++) {
235866f9b4aSJames Wright       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
236866f9b4aSJames Wright                                  basis_q_sur, basis_x_sur, qf_apply_outflow, qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian));
237866f9b4aSJames Wright     }
238866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow));
239866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow_jacobian));
240866f9b4aSJames Wright   }
241866f9b4aSJames Wright 
242866f9b4aSJames Wright   {  // --- Create Sub-Operator for freestream boundaries
243866f9b4aSJames Wright     CeedQFunction qf_apply_freestream = NULL, qf_apply_freestream_jacobian = NULL;
244866f9b4aSJames Wright 
245866f9b4aSJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream,
246866f9b4aSJames Wright                                 problem->apply_freestream_jacobian, &qf_apply_freestream, &qf_apply_freestream_jacobian));
247866f9b4aSJames Wright     for (CeedInt i = 0; i < bc->num_freestream; i++) {
248866f9b4aSJames Wright       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
249866f9b4aSJames Wright                                  basis_q_sur, basis_x_sur, qf_apply_freestream, qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian));
250866f9b4aSJames Wright     }
251866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream));
252866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream_jacobian));
253866f9b4aSJames Wright   }
254866f9b4aSJames Wright 
255866f9b4aSJames Wright   {  // --- Create Sub-Operator for slip boundaries
256866f9b4aSJames Wright     CeedQFunction qf_apply_slip = NULL, qf_apply_slip_jacobian = NULL;
257866f9b4aSJames Wright 
258866f9b4aSJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_slip,
259866f9b4aSJames Wright                                 problem->apply_slip_jacobian, &qf_apply_slip, &qf_apply_slip_jacobian));
260866f9b4aSJames Wright     for (CeedInt i = 0; i < bc->num_slip; i++) {
261866f9b4aSJames Wright       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,
262866f9b4aSJames Wright                                  basis_x_sur, qf_apply_slip, qf_apply_slip_jacobian, op_apply, op_apply_ijacobian));
263866f9b4aSJames Wright     }
264866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip));
265866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip_jacobian));
266866f9b4aSJames Wright   }
267866f9b4aSJames Wright 
268866f9b4aSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q_sur));
269866f9b4aSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_sur));
270866f9b4aSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
271866f9b4aSJames Wright }
272866f9b4aSJames Wright 
273991aef52SJames Wright PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData problem, SimpleBC bc) {
274a515125bSLeila Ghaffari   const PetscInt      num_comp_q = 5;
27566d54740SJames Wright   PetscInt            dim;
27689de3142SJames Wright   CeedInt             jac_data_size_vol = num_comp_q + 6 + 3, num_comp_x, q_data_size_vol;
277bbc90103SJames Wright   CeedElemRestriction elem_restr_jd_i;
278bbc90103SJames Wright   CeedVector          jac_data;
279bbc90103SJames Wright   CeedOperator        op_ifunction_vol = NULL, op_rhs_vol = NULL, op_ijacobian_vol = NULL;
280bbc90103SJames Wright 
281bbc90103SJames Wright   PetscFunctionBeginUser;
28266d54740SJames Wright   PetscCall(DMGetDimension(dm, &dim));
28366d54740SJames Wright   num_comp_x = dim;
284e07531f7SJames Wright   if (problem->apply_vol_ifunction.qf_func_ptr && problem->uses_newtonian) {
28594a7b3d2SKenneth E. Jansen     NewtonianIdealGasContext gas;
286e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfctx, CEED_MEM_HOST, &gas));
28794a7b3d2SKenneth E. Jansen     jac_data_size_vol += (gas->idl_enable ? 1 : 0);
288e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfctx, &gas));
28994a7b3d2SKenneth E. Jansen   }
29094a7b3d2SKenneth E. Jansen 
291bbc90103SJames Wright   {  // Create bases and element restrictions
29215c18037SJames Wright     DMLabel  domain_label = NULL;
29315c18037SJames Wright     PetscInt label_value = 0, height = 0, dm_field = 0;
29467263decSKenneth E. Jansen     DM       dm_coord;
29567263decSKenneth E. Jansen 
296bbc90103SJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
29715c18037SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &ceed_data->basis_q));
29815c18037SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &ceed_data->basis_x));
299a515125bSLeila Ghaffari 
30015c18037SJames Wright     PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, 0, &ceed_data->elem_restr_q));
30115c18037SJames Wright     PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &ceed_data->elem_restr_x));
30215c18037SJames Wright     PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_vol, &elem_restr_jd_i));
303ce8bebb6SJames Wright 
304b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL));
305b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL));
306b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL));
307b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL));
308ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL));
309a515125bSLeila Ghaffari 
310ce8bebb6SJames Wright     {  // -- Copy PETSc coordinate vector into CEED vector
311a515125bSLeila Ghaffari       Vec X_loc;
312cac8aa24SJed Brown       DM  cdm;
313ce8bebb6SJames Wright 
3142b916ea7SJeremy L Thompson       PetscCall(DMGetCellCoordinateDM(dm, &cdm));
3152b916ea7SJeremy L Thompson       if (cdm) {
3162b916ea7SJeremy L Thompson         PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
3172b916ea7SJeremy L Thompson       } else {
3182b916ea7SJeremy L Thompson         PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
319cac8aa24SJed Brown       }
320*8a02cd4cSJames Wright       PetscCall(VecScale(X_loc, user->units->meter));
321a7dac1d5SJames Wright       PetscCall(VecCopyPetscToCeed(X_loc, ceed_data->x_coord));
322ce8bebb6SJames Wright     }
323a515125bSLeila Ghaffari 
324c864c5abSJames Wright     PetscCall(QDataGet(ceed, dm, domain_label, label_value, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord,
32589de3142SJames Wright                        &ceed_data->elem_restr_qd_i, &ceed_data->q_data, &q_data_size_vol));
326ce8bebb6SJames Wright   }
327ce8bebb6SJames Wright 
328ce8bebb6SJames Wright   {  // -- Create QFunction for ICs
329866f9b4aSJames Wright     CeedBasis     basis_xc;
330ce8bebb6SJames Wright     CeedQFunction qf_ics;
3318f18bb8bSJames Wright     CeedOperator  op_ics;
332ce8bebb6SJames Wright 
333866f9b4aSJames Wright     PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_xc));
334e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qf_func_ptr, problem->ics.qf_loc, &qf_ics));
335e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ics, problem->ics.qfctx));
336ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ics, 0));
337ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "x", num_comp_x, CEED_EVAL_INTERP));
338ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
339ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ics, "q0", num_comp_q, CEED_EVAL_NONE));
340ce8bebb6SJames Wright 
341ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics));
342866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
343866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
34458e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
345b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label));
3468f18bb8bSJames Wright     PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx));
347a515125bSLeila Ghaffari 
348866f9b4aSJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_xc));
349ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ics));
350ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics));
351ce8bebb6SJames Wright   }
352ce8bebb6SJames Wright 
353e07531f7SJames Wright   if (problem->apply_vol_rhs.qf_func_ptr) {
354ce8bebb6SJames Wright     CeedQFunction qf_rhs_vol;
355ce8bebb6SJames Wright 
356e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qf_func_ptr, problem->apply_vol_rhs.qf_loc, &qf_rhs_vol));
357e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_vol, problem->apply_vol_rhs.qfctx));
358ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_rhs_vol, 0));
359ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP));
360ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
36189de3142SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
362ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP));
363ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP));
364ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
365ce8bebb6SJames Wright 
366bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_vol, NULL, NULL, &op_rhs_vol));
367bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
368bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
369bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
370bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
371bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
372bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
373ce8bebb6SJames Wright 
374ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_vol));
375a515125bSLeila Ghaffari   }
376a515125bSLeila Ghaffari 
377e07531f7SJames Wright   if (problem->apply_vol_ifunction.qf_func_ptr) {
378ce8bebb6SJames Wright     CeedQFunction qf_ifunction_vol;
379ce8bebb6SJames Wright 
380e07531f7SJames Wright     PetscCallCeed(
381e07531f7SJames Wright         ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qf_func_ptr, problem->apply_vol_ifunction.qf_loc, &qf_ifunction_vol));
382e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ifunction_vol, problem->apply_vol_ifunction.qfctx));
383ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ifunction_vol, 0));
384ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP));
385ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
386ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP));
38789de3142SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
388ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP));
389ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP));
390ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
391ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
392ce8bebb6SJames Wright 
393bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ifunction_vol, NULL, NULL, &op_ifunction_vol));
394bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
395bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
396bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed));
397bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
398bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
399bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
400bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
401bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
402752f40e3SJed Brown 
403ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction_vol));
404a515125bSLeila Ghaffari   }
405a515125bSLeila Ghaffari 
406e07531f7SJames Wright   if (problem->apply_vol_ijacobian.qf_func_ptr) {
407ce8bebb6SJames Wright     CeedQFunction qf_ijacobian_vol;
408ce8bebb6SJames Wright 
409e07531f7SJames Wright     PetscCallCeed(
410e07531f7SJames Wright         ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qf_func_ptr, problem->apply_vol_ijacobian.qf_loc, &qf_ijacobian_vol));
411e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfctx));
412ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ijacobian_vol, 0));
413ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP));
414ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD));
41589de3142SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
416ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
417ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP));
418ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
419ce8bebb6SJames Wright 
420bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op_ijacobian_vol));
421bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
422bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
423bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
424bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
425bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
426bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
427ce8bebb6SJames Wright 
428b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol));
429f0b65372SJed Brown   }
430f0b65372SJed Brown 
431a515125bSLeila Ghaffari   // -- Create and apply CEED Composite Operator for the entire domain
432a515125bSLeila Ghaffari   if (!user->phys->implicit) {  // RHS
433da5fe0e4SJames Wright     CeedOperator op_rhs;
434866f9b4aSJames Wright 
435866f9b4aSJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_rhs));
436bbc90103SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_rhs, op_rhs_vol));
437866f9b4aSJames Wright     PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, op_rhs, NULL));
438866f9b4aSJames Wright 
439da5fe0e4SJames Wright     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx));
440866f9b4aSJames Wright 
441866f9b4aSJames Wright     // ----- Get Context Labels for Operator
442866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "solution time", &user->phys->solution_time_label));
443866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "timestep size", &user->phys->timestep_size_label));
444866f9b4aSJames Wright 
445b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
446a78efa86SJames Wright     PetscCall(CreateKSPMass(user, problem));
447c2376cc9SJames Wright     PetscCheck(app_ctx->sgs_model_type == SGS_MODEL_NONE, user->comm, PETSC_ERR_SUP, "SGS modeling not implemented for explicit timestepping");
448a515125bSLeila Ghaffari   } else {  // IFunction
449ebfabadfSJames Wright     CeedOperator op_ijacobian = NULL;
450ebfabadfSJames Wright 
451866f9b4aSJames Wright     // Create Composite Operaters
452866f9b4aSJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &user->op_ifunction));
453bbc90103SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(user->op_ifunction, op_ifunction_vol));
454866f9b4aSJames Wright     if (op_ijacobian_vol) {
455866f9b4aSJames Wright       PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_ijacobian));
456866f9b4aSJames Wright       PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_ijacobian, op_ijacobian_vol));
457866f9b4aSJames Wright     }
458866f9b4aSJames Wright     PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, user->op_ifunction, op_ijacobian));
459866f9b4aSJames Wright 
460866f9b4aSJames Wright     // ----- Get Context Labels for Operator
461866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "solution time", &user->phys->solution_time_label));
462866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "timestep size", &user->phys->timestep_size_label));
463866f9b4aSJames Wright 
464ebfabadfSJames Wright     if (op_ijacobian) {
465000d2032SJeremy L Thompson       PetscCall(MatCreateCeed(user->dm, user->dm, op_ijacobian, NULL, &user->mat_ijacobian));
466ebfabadfSJames Wright       PetscCall(MatCeedSetLocalVectors(user->mat_ijacobian, user->Q_dot_loc, NULL));
467ebfabadfSJames Wright       PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian));
468f0b65372SJed Brown     }
469ad494f68SJames Wright     if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SgsDDSetup(ceed, user, ceed_data, problem));
4706d0190e2SJames Wright   }
47191933550SJames Wright 
472c2376cc9SJames Wright   if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc));
47391933550SJames Wright   if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem));
474aa0b7f76SJames Wright   if (app_ctx->diff_filter_monitor && !user->diff_filter) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem));
4751c17f66aSJames Wright   if (app_ctx->sgs_train_enable) PetscCall(SGS_DD_TrainingSetup(ceed, user, ceed_data, problem));
47691933550SJames Wright 
477bbc90103SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data));
478b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i));
479b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol));
480bbc90103SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ifunction_vol));
481bbc90103SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_vol));
482d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
483a515125bSLeila Ghaffari }
484