xref: /honee/src/setuplibceed.c (revision 149fb5361f5198e41f87e8815a7e9dbfee84a96a)
1dc936754SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Setup libCEED for Navier-Stokes example using PETSc
10a515125bSLeila Ghaffari 
11e419654dSJeremy L Thompson #include <ceed.h>
12e419654dSJeremy L Thompson #include <petscdmplex.h>
13e419654dSJeremy L Thompson 
14*149fb536SJames Wright #include <navierstokes.h>
15a515125bSLeila Ghaffari 
16a78efa86SJames Wright // @brief Create CeedOperator for unstabilized mass KSP for explicit timestepping
17a78efa86SJames Wright static PetscErrorCode CreateKSPMassOperator_Unstabilized(User user, CeedOperator *op_mass) {
18a78efa86SJames Wright   Ceed                ceed = user->ceed;
19a78efa86SJames Wright   CeedInt             num_comp_q, q_data_size;
20a78efa86SJames Wright   CeedQFunction       qf_mass;
21a78efa86SJames Wright   CeedElemRestriction elem_restr_q, elem_restr_qd_i;
22a78efa86SJames Wright   CeedBasis           basis_q;
23a78efa86SJames Wright   CeedVector          q_data;
24a78efa86SJames Wright 
25a78efa86SJames Wright   PetscFunctionBeginUser;
26a78efa86SJames Wright   {  // Get restriction and basis from the RHS function
27a78efa86SJames Wright     CeedOperator     *sub_ops;
28a78efa86SJames Wright     CeedOperatorField field;
29a78efa86SJames Wright     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
30a78efa86SJames Wright 
31a78efa86SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops));
32a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
33a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
34a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q));
35a78efa86SJames Wright 
36a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field));
37a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i));
38a78efa86SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data));
39a78efa86SJames Wright   }
40a78efa86SJames Wright 
41a78efa86SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
42a78efa86SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size));
43a78efa86SJames Wright 
44a78efa86SJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
45a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
46a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "u", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
47a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data));
48a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
49a78efa86SJames Wright 
50a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
51a78efa86SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
52a78efa86SJames Wright }
53a78efa86SJames Wright 
54a78efa86SJames Wright // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes
55991aef52SJames Wright static PetscErrorCode CreateKSPMass(User user, ProblemData problem) {
56a78efa86SJames Wright   Ceed         ceed = user->ceed;
57a78efa86SJames Wright   DM           dm   = user->dm;
58a78efa86SJames Wright   CeedOperator op_mass;
59a78efa86SJames Wright 
60a78efa86SJames Wright   PetscFunctionBeginUser;
61a78efa86SJames Wright   if (problem->create_mass_operator) PetscCall(problem->create_mass_operator(user, &op_mass));
62a78efa86SJames Wright   else PetscCall(CreateKSPMassOperator_Unstabilized(user, &op_mass));
63a78efa86SJames Wright 
64a78efa86SJames Wright   {  // -- Setup KSP for mass operator
65a78efa86SJames Wright     Mat      mat_mass;
66a78efa86SJames Wright     Vec      Zeros_loc;
67a78efa86SJames Wright     MPI_Comm comm = PetscObjectComm((PetscObject)dm);
68a78efa86SJames Wright 
69a78efa86SJames Wright     PetscCall(DMCreateLocalVector(dm, &Zeros_loc));
70a78efa86SJames Wright     PetscCall(VecZeroEntries(Zeros_loc));
71a78efa86SJames Wright     PetscCall(MatCeedCreate(dm, dm, op_mass, NULL, &mat_mass));
72a78efa86SJames Wright     PetscCall(MatCeedSetLocalVectors(mat_mass, Zeros_loc, NULL));
73a78efa86SJames Wright 
74a78efa86SJames Wright     PetscCall(KSPCreate(comm, &user->mass_ksp));
75a78efa86SJames Wright     PetscCall(KSPSetOptionsPrefix(user->mass_ksp, "mass_"));
76a78efa86SJames Wright     {  // lumped by default
77a78efa86SJames Wright       PC pc;
78a78efa86SJames Wright       PetscCall(KSPGetPC(user->mass_ksp, &pc));
79a78efa86SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
80a78efa86SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
81a78efa86SJames Wright       PetscCall(KSPSetType(user->mass_ksp, KSPPREONLY));
82a78efa86SJames Wright     }
83a78efa86SJames Wright     PetscCall(KSPSetFromOptions_WithMatCeed(user->mass_ksp, mat_mass));
84a78efa86SJames Wright     PetscCall(KSPSetFromOptions(user->mass_ksp));
85a78efa86SJames Wright     PetscCall(VecDestroy(&Zeros_loc));
86a78efa86SJames Wright     PetscCall(MatDestroy(&mat_mass));
87a78efa86SJames Wright   }
88a78efa86SJames Wright 
89a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
90a78efa86SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
91a78efa86SJames Wright }
92a78efa86SJames Wright 
93866f9b4aSJames Wright static PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height,
94866f9b4aSJames Wright                                        CeedInt Q_sur, CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedBasis basis_q_sur,
95866f9b4aSJames Wright                                        CeedBasis basis_x_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian, CeedOperator op_apply,
96866f9b4aSJames Wright                                        CeedOperator op_apply_ijacobian) {
9715c18037SJames Wright   CeedVector          q_data_sur, jac_data_sur          = NULL;
98866f9b4aSJames Wright   CeedOperator        op_apply_bc, op_apply_bc_jacobian = NULL;
9915c18037SJames Wright   CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur = NULL;
100866f9b4aSJames Wright   PetscInt            dm_field = 0;
101368c645fSJames Wright 
10206f41313SJames Wright   PetscFunctionBeginUser;
10315c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &elem_restr_q_sur));
10415c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &elem_restr_x_sur));
105368c645fSJames Wright   if (jac_data_size_sur > 0) {
106368c645fSJames Wright     // State-dependent data will be passed from residual to Jacobian. This will be collocated.
10715c18037SJames Wright     PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_sur, &elem_restr_jd_i_sur));
108b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL));
109368c645fSJames Wright   }
110368c645fSJames Wright 
111c864c5abSJames 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,
112c864c5abSJames Wright                              &q_data_sur, &q_data_size_sur));
113866f9b4aSJames Wright 
114bbc90103SJames Wright   // CEED Operator for Physics
115b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc));
116866f9b4aSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
117866f9b4aSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
11858e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
119866f9b4aSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord));
120866f9b4aSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
121b4c37c5cSJames Wright   if (elem_restr_jd_i_sur)
12258e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
123368c645fSJames Wright 
1244b96a86bSJames Wright   if (qf_apply_bc_jacobian && elem_restr_jd_i_sur) {
125b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian));
126866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
127866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
12858e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
129866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord));
13058e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
131866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
132368c645fSJames Wright   }
133368c645fSJames Wright 
134bbc90103SJames Wright   // Apply Sub-Operator for Physics
135866f9b4aSJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply, op_apply_bc));
136866f9b4aSJames Wright   if (op_apply_bc_jacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply_ijacobian, op_apply_bc_jacobian));
137368c645fSJames Wright 
138b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_sur));
139b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data_sur));
140b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur));
141b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur));
142b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i_sur));
143b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i_sur));
144b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc));
145b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc_jacobian));
146d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
147368c645fSJames Wright }
148368c645fSJames Wright 
149866f9b4aSJames Wright static PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur,
1502b916ea7SJeremy L Thompson                                         PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian,
15125988f00SJames Wright                                         CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) {
15225988f00SJames Wright   PetscFunctionBeginUser;
15325988f00SJames Wright   if (apply_bc.qfunction) {
154b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc));
155b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context));
156ebfabadfSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc, 0));
157b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP));
158b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD));
159b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
160b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP));
161b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP));
162b4c37c5cSJames Wright     if (jac_data_size_sur) PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
16325988f00SJames Wright   }
16425988f00SJames Wright   if (apply_bc_jacobian.qfunction) {
165b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian));
166b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context));
167ebfabadfSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc_jacobian, 0));
168b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP));
169b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD));
170b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
171b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP));
172b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
173b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP));
17425988f00SJames Wright   }
175d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
17625988f00SJames Wright }
17725988f00SJames Wright 
178bbc90103SJames Wright // Utility function to add boundary operators to the composite operator
179866f9b4aSJames Wright static PetscErrorCode AddBCSubOperators(User user, Ceed ceed, DM dm, SimpleBC bc, ProblemData problem, CeedData ceed_data, CeedOperator op_apply,
180866f9b4aSJames Wright                                         CeedOperator op_apply_ijacobian) {
181866f9b4aSJames Wright   CeedInt       height = 1, num_comp_q, num_comp_x;
182c864c5abSJames Wright   CeedInt       P_sur = user->app_ctx->degree + 1, Q_sur = P_sur + user->app_ctx->q_extra, dim_sur, q_data_size_sur;
183c864c5abSJames Wright   const CeedInt jac_data_size_sur = user->phys->implicit ? problem->jac_data_size_sur : 0;
184866f9b4aSJames Wright   PetscInt      dim;
185866f9b4aSJames Wright   DMLabel       face_sets_label;
186866f9b4aSJames Wright   CeedBasis     basis_q_sur, basis_x_sur;
187866f9b4aSJames Wright 
188866f9b4aSJames Wright   PetscFunctionBeginUser;
189866f9b4aSJames Wright   PetscCall(DMGetDimension(dm, &dim));
190c864c5abSJames Wright   PetscCall(QDataBoundaryGetNumComponents(dm, &q_data_size_sur));
191866f9b4aSJames Wright   dim_sur = dim - height;
192866f9b4aSJames Wright   {  // Get number of components and coordinate dimension from op_apply
193866f9b4aSJames Wright     CeedOperator       *sub_ops;
194866f9b4aSJames Wright     CeedOperatorField   field;
195866f9b4aSJames Wright     PetscInt            sub_op_index = 0;  // will be 0 for the volume op
196866f9b4aSJames Wright     CeedElemRestriction elem_restr_q, elem_restr_x;
197866f9b4aSJames Wright 
198866f9b4aSJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(op_apply, &sub_ops));
199866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
200866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
201866f9b4aSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
202866f9b4aSJames Wright 
203866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "x", &field));
204866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_x));
205866f9b4aSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
206866f9b4aSJames Wright   }
207866f9b4aSJames Wright 
208866f9b4aSJames Wright   {  // Get bases
209866f9b4aSJames Wright     DM dm_coord;
210866f9b4aSJames Wright 
211866f9b4aSJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
212866f9b4aSJames Wright     DMLabel  label       = NULL;
213866f9b4aSJames Wright     PetscInt label_value = 0;
214bbc90103SJames Wright     PetscInt field       = 0;
215866f9b4aSJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm, label, label_value, height, field, &basis_q_sur));
216866f9b4aSJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, label_value, height, field, &basis_x_sur));
217866f9b4aSJames Wright   }
218866f9b4aSJames Wright 
219866f9b4aSJames Wright   PetscCall(DMGetLabel(dm, "Face Sets", &face_sets_label));
220866f9b4aSJames Wright 
221866f9b4aSJames Wright   {  // --- Create Sub-Operator for inflow boundaries
222866f9b4aSJames Wright     CeedQFunction qf_apply_inflow = NULL, qf_apply_inflow_jacobian = NULL;
223866f9b4aSJames Wright 
224866f9b4aSJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow,
225866f9b4aSJames Wright                                 problem->apply_inflow_jacobian, &qf_apply_inflow, &qf_apply_inflow_jacobian));
226866f9b4aSJames Wright     for (CeedInt i = 0; i < bc->num_inflow; i++) {
227866f9b4aSJames 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,
228866f9b4aSJames Wright                                  basis_x_sur, qf_apply_inflow, qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian));
229866f9b4aSJames Wright     }
230866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow));
231866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow_jacobian));
232866f9b4aSJames Wright   }
233866f9b4aSJames Wright 
234866f9b4aSJames Wright   {  // --- Create Sub-Operator for outflow boundaries
235866f9b4aSJames Wright     CeedQFunction qf_apply_outflow = NULL, qf_apply_outflow_jacobian = NULL;
236866f9b4aSJames Wright 
237866f9b4aSJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow,
238866f9b4aSJames Wright                                 problem->apply_outflow_jacobian, &qf_apply_outflow, &qf_apply_outflow_jacobian));
239866f9b4aSJames Wright     for (CeedInt i = 0; i < bc->num_outflow; i++) {
240866f9b4aSJames Wright       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
241866f9b4aSJames Wright                                  basis_q_sur, basis_x_sur, qf_apply_outflow, qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian));
242866f9b4aSJames Wright     }
243866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow));
244866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow_jacobian));
245866f9b4aSJames Wright   }
246866f9b4aSJames Wright 
247866f9b4aSJames Wright   {  // --- Create Sub-Operator for freestream boundaries
248866f9b4aSJames Wright     CeedQFunction qf_apply_freestream = NULL, qf_apply_freestream_jacobian = NULL;
249866f9b4aSJames Wright 
250866f9b4aSJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream,
251866f9b4aSJames Wright                                 problem->apply_freestream_jacobian, &qf_apply_freestream, &qf_apply_freestream_jacobian));
252866f9b4aSJames Wright     for (CeedInt i = 0; i < bc->num_freestream; i++) {
253866f9b4aSJames Wright       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
254866f9b4aSJames Wright                                  basis_q_sur, basis_x_sur, qf_apply_freestream, qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian));
255866f9b4aSJames Wright     }
256866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream));
257866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream_jacobian));
258866f9b4aSJames Wright   }
259866f9b4aSJames Wright 
260866f9b4aSJames Wright   {  // --- Create Sub-Operator for slip boundaries
261866f9b4aSJames Wright     CeedQFunction qf_apply_slip = NULL, qf_apply_slip_jacobian = NULL;
262866f9b4aSJames Wright 
263866f9b4aSJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_slip,
264866f9b4aSJames Wright                                 problem->apply_slip_jacobian, &qf_apply_slip, &qf_apply_slip_jacobian));
265866f9b4aSJames Wright     for (CeedInt i = 0; i < bc->num_slip; i++) {
266866f9b4aSJames 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,
267866f9b4aSJames Wright                                  basis_x_sur, qf_apply_slip, qf_apply_slip_jacobian, op_apply, op_apply_ijacobian));
268866f9b4aSJames Wright     }
269866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip));
270866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip_jacobian));
271866f9b4aSJames Wright   }
272866f9b4aSJames Wright 
273866f9b4aSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q_sur));
274866f9b4aSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_sur));
275866f9b4aSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
276866f9b4aSJames Wright }
277866f9b4aSJames Wright 
278991aef52SJames Wright PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData problem, SimpleBC bc) {
279a515125bSLeila Ghaffari   const PetscInt      num_comp_q = 5;
280c864c5abSJames Wright   const CeedInt       dim = problem->dim, num_comp_x = problem->dim;
28194a7b3d2SKenneth E. Jansen   CeedInt             jac_data_size_vol = num_comp_q + 6 + 3;
282bbc90103SJames Wright   CeedElemRestriction elem_restr_jd_i;
283bbc90103SJames Wright   CeedVector          jac_data;
284bbc90103SJames Wright   CeedOperator        op_ifunction_vol = NULL, op_rhs_vol = NULL, op_ijacobian_vol = NULL;
285bbc90103SJames Wright 
286bbc90103SJames Wright   PetscFunctionBeginUser;
28794a7b3d2SKenneth E. Jansen 
2887d8a615bSJames Wright   if (problem->apply_vol_ifunction.qfunction && problem->uses_newtonian) {
28994a7b3d2SKenneth E. Jansen     NewtonianIdealGasContext gas;
29094a7b3d2SKenneth E. Jansen     PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas));
29194a7b3d2SKenneth E. Jansen     jac_data_size_vol += (gas->idl_enable ? 1 : 0);
29294a7b3d2SKenneth E. Jansen     PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas));
29394a7b3d2SKenneth E. Jansen   }
29494a7b3d2SKenneth E. Jansen 
295bbc90103SJames Wright   {  // Create bases and element restrictions
29615c18037SJames Wright     DMLabel  domain_label = NULL;
29715c18037SJames Wright     PetscInt label_value = 0, height = 0, dm_field = 0;
29867263decSKenneth E. Jansen     DM       dm_coord;
29967263decSKenneth E. Jansen 
300bbc90103SJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
30115c18037SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &ceed_data->basis_q));
30215c18037SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &ceed_data->basis_x));
303a515125bSLeila Ghaffari 
30415c18037SJames Wright     PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, 0, &ceed_data->elem_restr_q));
30515c18037SJames Wright     PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &ceed_data->elem_restr_x));
30615c18037SJames Wright     PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_vol, &elem_restr_jd_i));
307ce8bebb6SJames Wright 
308b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL));
309b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL));
310b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL));
311b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL));
312ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL));
313a515125bSLeila Ghaffari 
314ce8bebb6SJames Wright     {  // -- Copy PETSc coordinate vector into CEED vector
315a515125bSLeila Ghaffari       Vec X_loc;
316cac8aa24SJed Brown       DM  cdm;
317ce8bebb6SJames Wright 
3182b916ea7SJeremy L Thompson       PetscCall(DMGetCellCoordinateDM(dm, &cdm));
3192b916ea7SJeremy L Thompson       if (cdm) {
3202b916ea7SJeremy L Thompson         PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
3212b916ea7SJeremy L Thompson       } else {
3222b916ea7SJeremy L Thompson         PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
323cac8aa24SJed Brown       }
3242b916ea7SJeremy L Thompson       PetscCall(VecScale(X_loc, problem->dm_scale));
325a7dac1d5SJames Wright       PetscCall(VecCopyPetscToCeed(X_loc, ceed_data->x_coord));
326ce8bebb6SJames Wright     }
327a515125bSLeila Ghaffari 
328c864c5abSJames Wright     PetscCall(QDataGet(ceed, dm, domain_label, label_value, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord,
329c864c5abSJames Wright                        &ceed_data->elem_restr_qd_i, &ceed_data->q_data, &problem->q_data_size_vol));
330ce8bebb6SJames Wright   }
331ce8bebb6SJames Wright 
332ce8bebb6SJames Wright   {  // -- Create QFunction for ICs
333866f9b4aSJames Wright     CeedBasis     basis_xc;
334ce8bebb6SJames Wright     CeedQFunction qf_ics;
3358f18bb8bSJames Wright     CeedOperator  op_ics;
336ce8bebb6SJames Wright 
337866f9b4aSJames Wright     PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_xc));
338ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &qf_ics));
339ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ics, problem->ics.qfunction_context));
340ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ics, 0));
341ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "x", num_comp_x, CEED_EVAL_INTERP));
342ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
343ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ics, "q0", num_comp_q, CEED_EVAL_NONE));
344ce8bebb6SJames Wright 
345ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics));
346866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
347866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
34858e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
349b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label));
3508f18bb8bSJames Wright     PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx));
351a515125bSLeila Ghaffari 
352866f9b4aSJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_xc));
353ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ics));
354ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics));
355ce8bebb6SJames Wright   }
356ce8bebb6SJames Wright 
357ce8bebb6SJames Wright   if (problem->apply_vol_rhs.qfunction) {
358ce8bebb6SJames Wright     CeedQFunction qf_rhs_vol;
359ce8bebb6SJames Wright 
360ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &qf_rhs_vol));
361ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_vol, problem->apply_vol_rhs.qfunction_context));
362ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_rhs_vol, 0));
363ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP));
364ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
365c864c5abSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "qdata", problem->q_data_size_vol, CEED_EVAL_NONE));
366ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP));
367ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP));
368ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
369ce8bebb6SJames Wright 
370bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_vol, NULL, NULL, &op_rhs_vol));
371bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
372bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
373bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
374bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
375bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
376bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
377ce8bebb6SJames Wright 
378ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_vol));
379a515125bSLeila Ghaffari   }
380a515125bSLeila Ghaffari 
381ce8bebb6SJames Wright   if (problem->apply_vol_ifunction.qfunction) {
382ce8bebb6SJames Wright     CeedQFunction qf_ifunction_vol;
383ce8bebb6SJames Wright 
384ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc,
385ce8bebb6SJames Wright                                                     &qf_ifunction_vol));
386ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context));
387ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ifunction_vol, 0));
388ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP));
389ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
390ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP));
391c864c5abSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "qdata", problem->q_data_size_vol, CEED_EVAL_NONE));
392ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP));
393ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP));
394ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
395ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
396ce8bebb6SJames Wright 
397bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ifunction_vol, NULL, NULL, &op_ifunction_vol));
398bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
399bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
400bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed));
401bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
402bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
403bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
404bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
405bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
406752f40e3SJed Brown 
407ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction_vol));
408a515125bSLeila Ghaffari   }
409a515125bSLeila Ghaffari 
410ce8bebb6SJames Wright   if (problem->apply_vol_ijacobian.qfunction) {
411ce8bebb6SJames Wright     CeedQFunction qf_ijacobian_vol;
412ce8bebb6SJames Wright 
413ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc,
414ce8bebb6SJames Wright                                                     &qf_ijacobian_vol));
415ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context));
416ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ijacobian_vol, 0));
417ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP));
418ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD));
419c864c5abSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", problem->q_data_size_vol, CEED_EVAL_NONE));
420ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
421ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP));
422ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
423ce8bebb6SJames Wright 
424bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op_ijacobian_vol));
425bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
426bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
427bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
428bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
429bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
430bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
431ce8bebb6SJames Wright 
432b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol));
433f0b65372SJed Brown   }
434f0b65372SJed Brown 
435a515125bSLeila Ghaffari   // -- Create and apply CEED Composite Operator for the entire domain
436a515125bSLeila Ghaffari   if (!user->phys->implicit) {  // RHS
437da5fe0e4SJames Wright     CeedOperator op_rhs;
438866f9b4aSJames Wright 
439866f9b4aSJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_rhs));
440bbc90103SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_rhs, op_rhs_vol));
441866f9b4aSJames Wright     PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, op_rhs, NULL));
442866f9b4aSJames Wright 
443da5fe0e4SJames Wright     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx));
444866f9b4aSJames Wright 
445866f9b4aSJames Wright     // ----- Get Context Labels for Operator
446866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "solution time", &user->phys->solution_time_label));
447866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "timestep size", &user->phys->timestep_size_label));
448866f9b4aSJames Wright 
449b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
450a78efa86SJames Wright     PetscCall(CreateKSPMass(user, problem));
451c2376cc9SJames Wright     PetscCheck(app_ctx->sgs_model_type == SGS_MODEL_NONE, user->comm, PETSC_ERR_SUP, "SGS modeling not implemented for explicit timestepping");
452a515125bSLeila Ghaffari   } else {  // IFunction
453ebfabadfSJames Wright     CeedOperator op_ijacobian = NULL;
454ebfabadfSJames Wright 
455866f9b4aSJames Wright     // Create Composite Operaters
456866f9b4aSJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &user->op_ifunction));
457bbc90103SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(user->op_ifunction, op_ifunction_vol));
458866f9b4aSJames Wright     if (op_ijacobian_vol) {
459866f9b4aSJames Wright       PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_ijacobian));
460866f9b4aSJames Wright       PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_ijacobian, op_ijacobian_vol));
461866f9b4aSJames Wright     }
462866f9b4aSJames Wright     PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, user->op_ifunction, op_ijacobian));
463866f9b4aSJames Wright 
464866f9b4aSJames Wright     // ----- Get Context Labels for Operator
465866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "solution time", &user->phys->solution_time_label));
466866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "timestep size", &user->phys->timestep_size_label));
467866f9b4aSJames Wright 
468ebfabadfSJames Wright     if (op_ijacobian) {
469ebfabadfSJames Wright       PetscCall(MatCeedCreate(user->dm, user->dm, op_ijacobian, NULL, &user->mat_ijacobian));
470ebfabadfSJames Wright       PetscCall(MatCeedSetLocalVectors(user->mat_ijacobian, user->Q_dot_loc, NULL));
471ebfabadfSJames Wright       PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian));
472f0b65372SJed Brown     }
473ad494f68SJames Wright     if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SgsDDSetup(ceed, user, ceed_data, problem));
4746d0190e2SJames Wright   }
47591933550SJames Wright 
476c2376cc9SJames Wright   if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc));
47791933550SJames Wright   if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem));
478aa0b7f76SJames Wright   if (app_ctx->diff_filter_monitor && !user->diff_filter) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem));
4791c17f66aSJames Wright   if (app_ctx->sgs_train_enable) PetscCall(SGS_DD_TrainingSetup(ceed, user, ceed_data, problem));
48091933550SJames Wright 
481bbc90103SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data));
482b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i));
483b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol));
484bbc90103SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ifunction_vol));
485bbc90103SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_vol));
486d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
487a515125bSLeila Ghaffari }
488