xref: /libCEED/examples/fluids/src/setuplibceed.c (revision 9ba83ac0e4b1fca39d6fa6737a318a9f0cbc172d)
1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
377841947SLeila Ghaffari //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
577841947SLeila Ghaffari //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
777841947SLeila Ghaffari 
877841947SLeila Ghaffari /// @file
977841947SLeila Ghaffari /// Setup libCEED for Navier-Stokes example using PETSc
1077841947SLeila Ghaffari 
1149aac155SJeremy L Thompson #include <ceed.h>
1249aac155SJeremy L Thompson #include <petscdmplex.h>
1349aac155SJeremy L Thompson 
1477841947SLeila Ghaffari #include "../navierstokes.h"
1577841947SLeila Ghaffari 
16b4e9a8f8SJames Wright // @brief Create CeedOperator for unstabilized mass KSP for explicit timestepping
17b4e9a8f8SJames Wright static PetscErrorCode CreateKSPMassOperator_Unstabilized(User user, CeedOperator *op_mass) {
18b4e9a8f8SJames Wright   Ceed                ceed = user->ceed;
19b4e9a8f8SJames Wright   CeedInt             num_comp_q, q_data_size;
20b4e9a8f8SJames Wright   CeedQFunction       qf_mass;
21b4e9a8f8SJames Wright   CeedElemRestriction elem_restr_q, elem_restr_qd_i;
22b4e9a8f8SJames Wright   CeedBasis           basis_q;
23b4e9a8f8SJames Wright   CeedVector          q_data;
24b4e9a8f8SJames Wright 
25b4e9a8f8SJames Wright   PetscFunctionBeginUser;
26b4e9a8f8SJames Wright   {  // Get restriction and basis from the RHS function
27b4e9a8f8SJames Wright     CeedOperator     *sub_ops;
28b4e9a8f8SJames Wright     CeedOperatorField field;
29b4e9a8f8SJames Wright     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
30b4e9a8f8SJames Wright 
31ed094490SJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(user->op_rhs_ctx->op, &sub_ops));
32b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
33681d0ea7SJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorFieldGetData(field, NULL, &elem_restr_q, &basis_q, NULL));
34b4e9a8f8SJames Wright 
35b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field));
36681d0ea7SJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorFieldGetData(field, NULL, &elem_restr_qd_i, NULL, &q_data));
37b4e9a8f8SJames Wright   }
38b4e9a8f8SJames Wright 
39b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
40b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size));
41b4e9a8f8SJames Wright 
42b4e9a8f8SJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
43b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
44b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "u", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
45b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data));
46b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
47b4e9a8f8SJames Wright 
48681d0ea7SJeremy L Thompson   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
49681d0ea7SJeremy L Thompson   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
50681d0ea7SJeremy L Thompson   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i));
51681d0ea7SJeremy L Thompson   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
52b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
53b4e9a8f8SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
54b4e9a8f8SJames Wright }
55b4e9a8f8SJames Wright 
56b4e9a8f8SJames Wright // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes
57731c13d7SJames Wright static PetscErrorCode CreateKSPMass(User user, ProblemData problem) {
58b4e9a8f8SJames Wright   Ceed         ceed = user->ceed;
59b4e9a8f8SJames Wright   DM           dm   = user->dm;
60b4e9a8f8SJames Wright   CeedOperator op_mass;
61b4e9a8f8SJames Wright 
62b4e9a8f8SJames Wright   PetscFunctionBeginUser;
63b4e9a8f8SJames Wright   if (problem->create_mass_operator) PetscCall(problem->create_mass_operator(user, &op_mass));
64b4e9a8f8SJames Wright   else PetscCall(CreateKSPMassOperator_Unstabilized(user, &op_mass));
65b4e9a8f8SJames Wright 
66b4e9a8f8SJames Wright   {  // -- Setup KSP for mass operator
67b4e9a8f8SJames Wright     Mat      mat_mass;
68b4e9a8f8SJames Wright     Vec      Zeros_loc;
69b4e9a8f8SJames Wright     MPI_Comm comm = PetscObjectComm((PetscObject)dm);
70b4e9a8f8SJames Wright 
71b4e9a8f8SJames Wright     PetscCall(DMCreateLocalVector(dm, &Zeros_loc));
72b4e9a8f8SJames Wright     PetscCall(VecZeroEntries(Zeros_loc));
73243afec9SJames Wright     PetscCall(MatCreateCeed(dm, dm, op_mass, NULL, &mat_mass));
74b4e9a8f8SJames Wright     PetscCall(MatCeedSetLocalVectors(mat_mass, Zeros_loc, NULL));
75b4e9a8f8SJames Wright 
76b4e9a8f8SJames Wright     PetscCall(KSPCreate(comm, &user->mass_ksp));
77b4e9a8f8SJames Wright     PetscCall(KSPSetOptionsPrefix(user->mass_ksp, "mass_"));
78b4e9a8f8SJames Wright     {  // lumped by default
79b4e9a8f8SJames Wright       PC pc;
80b4e9a8f8SJames Wright       PetscCall(KSPGetPC(user->mass_ksp, &pc));
81b4e9a8f8SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
82b4e9a8f8SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
83b4e9a8f8SJames Wright       PetscCall(KSPSetType(user->mass_ksp, KSPPREONLY));
84b4e9a8f8SJames Wright     }
85b4e9a8f8SJames Wright     PetscCall(KSPSetFromOptions_WithMatCeed(user->mass_ksp, mat_mass));
86b4e9a8f8SJames Wright     PetscCall(VecDestroy(&Zeros_loc));
87b4e9a8f8SJames Wright     PetscCall(MatDestroy(&mat_mass));
88b4e9a8f8SJames Wright   }
89b4e9a8f8SJames Wright 
90b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
91b4e9a8f8SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
92b4e9a8f8SJames Wright }
93b4e9a8f8SJames Wright 
94e9c36be0SJames Wright static PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height,
95e9c36be0SJames Wright                                        CeedInt Q_sur, CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedBasis basis_q_sur,
96e9c36be0SJames Wright                                        CeedBasis basis_x_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian, CeedOperator op_apply,
97e9c36be0SJames Wright                                        CeedOperator op_apply_ijacobian) {
98bb85d312SJames Wright   CeedVector          q_data_sur, jac_data_sur          = NULL;
99e9c36be0SJames Wright   CeedOperator        op_apply_bc, op_apply_bc_jacobian = NULL;
100bb85d312SJames Wright   CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur = NULL;
101e9c36be0SJames Wright   PetscInt            dm_field = 0;
1029eb9c072SJames Wright 
103f17d818dSJames Wright   PetscFunctionBeginUser;
104bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &elem_restr_q_sur));
105bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &elem_restr_x_sur));
1069eb9c072SJames Wright   if (jac_data_size_sur > 0) {
1079eb9c072SJames Wright     // State-dependent data will be passed from residual to Jacobian. This will be collocated.
108bb85d312SJames Wright     PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_sur, &elem_restr_jd_i_sur));
109a424bcd0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL));
1109eb9c072SJames Wright   }
1119eb9c072SJames Wright 
112cfb075a4SJames 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,
113cfb075a4SJames Wright                              &q_data_sur, &q_data_size_sur));
114e9c36be0SJames Wright 
115a5b0ec6fSJames Wright   // CEED Operator for Physics
116a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc));
117e9c36be0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
118e9c36be0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
119356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
120e9c36be0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord));
121e9c36be0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
1226c10af5dSJeremy L Thompson   if (elem_restr_jd_i_sur) {
123356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
1246c10af5dSJeremy L Thompson   }
1259eb9c072SJames Wright 
126f21e6b1cSJames Wright   if (qf_apply_bc_jacobian && elem_restr_jd_i_sur) {
127a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian));
128e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
129e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
130356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
131e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord));
132356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
133e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
1349eb9c072SJames Wright   }
1359eb9c072SJames Wright 
136a5b0ec6fSJames Wright   // Apply Sub-Operator for Physics
137ed094490SJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorCompositeAddSub(op_apply, op_apply_bc));
138ed094490SJeremy L Thompson   if (op_apply_bc_jacobian) PetscCallCeed(ceed, CeedOperatorCompositeAddSub(op_apply_ijacobian, op_apply_bc_jacobian));
1399eb9c072SJames Wright 
140a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_sur));
141a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data_sur));
142a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur));
143a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur));
144a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i_sur));
145a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i_sur));
146a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc));
147a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc_jacobian));
148ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1499eb9c072SJames Wright }
1509eb9c072SJames Wright 
151e9c36be0SJames Wright static PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur,
1522b730f8bSJeremy L Thompson                                         PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian,
1535b881d11SJames Wright                                         CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) {
1545b881d11SJames Wright   PetscFunctionBeginUser;
1555b881d11SJames Wright   if (apply_bc.qfunction) {
156a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc));
157a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context));
15891c97f41SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc, 0));
159a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP));
160a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD));
161a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
162a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP));
163a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP));
164a424bcd0SJames Wright     if (jac_data_size_sur) PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
1655b881d11SJames Wright   }
1665b881d11SJames Wright   if (apply_bc_jacobian.qfunction) {
167a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian));
168a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context));
16991c97f41SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc_jacobian, 0));
170a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP));
171a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD));
172a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
173a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP));
174a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
175a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP));
1765b881d11SJames Wright   }
177ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1785b881d11SJames Wright }
1795b881d11SJames Wright 
180a5b0ec6fSJames Wright // Utility function to add boundary operators to the composite operator
181e9c36be0SJames Wright static PetscErrorCode AddBCSubOperators(User user, Ceed ceed, DM dm, SimpleBC bc, ProblemData problem, CeedData ceed_data, CeedOperator op_apply,
182e9c36be0SJames Wright                                         CeedOperator op_apply_ijacobian) {
183e9c36be0SJames Wright   CeedInt       height = 1, num_comp_q, num_comp_x;
184cfb075a4SJames Wright   CeedInt       P_sur = user->app_ctx->degree + 1, Q_sur = P_sur + user->app_ctx->q_extra, dim_sur, q_data_size_sur;
185cfb075a4SJames Wright   const CeedInt jac_data_size_sur = user->phys->implicit ? problem->jac_data_size_sur : 0;
186e9c36be0SJames Wright   PetscInt      dim;
187e9c36be0SJames Wright   DMLabel       face_sets_label;
188e9c36be0SJames Wright   CeedBasis     basis_q_sur, basis_x_sur;
189e9c36be0SJames Wright 
190e9c36be0SJames Wright   PetscFunctionBeginUser;
191e9c36be0SJames Wright   PetscCall(DMGetDimension(dm, &dim));
192cfb075a4SJames Wright   PetscCall(QDataBoundaryGetNumComponents(dm, &q_data_size_sur));
193e9c36be0SJames Wright   dim_sur = dim - height;
194e9c36be0SJames Wright   {  // Get number of components and coordinate dimension from op_apply
195e9c36be0SJames Wright     CeedOperator       *sub_ops;
196e9c36be0SJames Wright     CeedOperatorField   field;
197e9c36be0SJames Wright     PetscInt            sub_op_index = 0;  // will be 0 for the volume op
198e9c36be0SJames Wright     CeedElemRestriction elem_restr_q, elem_restr_x;
199e9c36be0SJames Wright 
200ed094490SJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorCompositeGetSubList(op_apply, &sub_ops));
201e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
202e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
203e9c36be0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
204681d0ea7SJeremy L Thompson     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
205e9c36be0SJames Wright 
206e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "x", &field));
207e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_x));
208e9c36be0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
209681d0ea7SJeremy L Thompson     PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x));
210e9c36be0SJames Wright   }
211e9c36be0SJames Wright 
212e9c36be0SJames Wright   {  // Get bases
213e9c36be0SJames Wright     DM dm_coord;
214e9c36be0SJames Wright 
215e9c36be0SJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
216e9c36be0SJames Wright     DMLabel  label       = NULL;
217e9c36be0SJames Wright     PetscInt label_value = 0;
218a5b0ec6fSJames Wright     PetscInt field       = 0;
219e9c36be0SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm, label, label_value, height, field, &basis_q_sur));
220e9c36be0SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, label_value, height, field, &basis_x_sur));
221e9c36be0SJames Wright   }
222e9c36be0SJames Wright 
223e9c36be0SJames Wright   PetscCall(DMGetLabel(dm, "Face Sets", &face_sets_label));
224e9c36be0SJames Wright 
225e9c36be0SJames Wright   {  // --- Create Sub-Operator for inflow boundaries
226e9c36be0SJames Wright     CeedQFunction qf_apply_inflow = NULL, qf_apply_inflow_jacobian = NULL;
227e9c36be0SJames Wright 
228e9c36be0SJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow,
229e9c36be0SJames Wright                                 problem->apply_inflow_jacobian, &qf_apply_inflow, &qf_apply_inflow_jacobian));
230e9c36be0SJames Wright     for (CeedInt i = 0; i < bc->num_inflow; i++) {
231e9c36be0SJames 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,
232e9c36be0SJames Wright                                  basis_x_sur, qf_apply_inflow, qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian));
233e9c36be0SJames Wright     }
234e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow));
235e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow_jacobian));
236e9c36be0SJames Wright   }
237e9c36be0SJames Wright 
238e9c36be0SJames Wright   {  // --- Create Sub-Operator for outflow boundaries
239e9c36be0SJames Wright     CeedQFunction qf_apply_outflow = NULL, qf_apply_outflow_jacobian = NULL;
240e9c36be0SJames Wright 
241e9c36be0SJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow,
242e9c36be0SJames Wright                                 problem->apply_outflow_jacobian, &qf_apply_outflow, &qf_apply_outflow_jacobian));
243e9c36be0SJames Wright     for (CeedInt i = 0; i < bc->num_outflow; i++) {
244e9c36be0SJames Wright       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
245e9c36be0SJames Wright                                  basis_q_sur, basis_x_sur, qf_apply_outflow, qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian));
246e9c36be0SJames Wright     }
247e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow));
248e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow_jacobian));
249e9c36be0SJames Wright   }
250e9c36be0SJames Wright 
251e9c36be0SJames Wright   {  // --- Create Sub-Operator for freestream boundaries
252e9c36be0SJames Wright     CeedQFunction qf_apply_freestream = NULL, qf_apply_freestream_jacobian = NULL;
253e9c36be0SJames Wright 
254e9c36be0SJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream,
255e9c36be0SJames Wright                                 problem->apply_freestream_jacobian, &qf_apply_freestream, &qf_apply_freestream_jacobian));
256e9c36be0SJames Wright     for (CeedInt i = 0; i < bc->num_freestream; i++) {
257e9c36be0SJames Wright       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
258e9c36be0SJames Wright                                  basis_q_sur, basis_x_sur, qf_apply_freestream, qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian));
259e9c36be0SJames Wright     }
260e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream));
261e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream_jacobian));
262e9c36be0SJames Wright   }
263e9c36be0SJames Wright 
264e9c36be0SJames Wright   {  // --- Create Sub-Operator for slip boundaries
265e9c36be0SJames Wright     CeedQFunction qf_apply_slip = NULL, qf_apply_slip_jacobian = NULL;
266e9c36be0SJames Wright 
267e9c36be0SJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_slip,
268e9c36be0SJames Wright                                 problem->apply_slip_jacobian, &qf_apply_slip, &qf_apply_slip_jacobian));
269e9c36be0SJames Wright     for (CeedInt i = 0; i < bc->num_slip; i++) {
270e9c36be0SJames 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,
271e9c36be0SJames Wright                                  basis_x_sur, qf_apply_slip, qf_apply_slip_jacobian, op_apply, op_apply_ijacobian));
272e9c36be0SJames Wright     }
273e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip));
274e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip_jacobian));
275e9c36be0SJames Wright   }
276e9c36be0SJames Wright 
277e9c36be0SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q_sur));
278e9c36be0SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_sur));
279e9c36be0SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
280e9c36be0SJames Wright }
281e9c36be0SJames Wright 
282731c13d7SJames Wright PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData problem, SimpleBC bc) {
28377841947SLeila Ghaffari   const PetscInt      num_comp_q = 5;
284cfb075a4SJames Wright   const CeedInt       dim = problem->dim, num_comp_x = problem->dim;
2851d2a9659SKenneth E. Jansen   CeedInt             jac_data_size_vol = num_comp_q + 6 + 3;
286a5b0ec6fSJames Wright   CeedElemRestriction elem_restr_jd_i;
287a5b0ec6fSJames Wright   CeedVector          jac_data;
288a5b0ec6fSJames Wright   CeedOperator        op_ifunction_vol = NULL, op_rhs_vol = NULL, op_ijacobian_vol = NULL;
289a5b0ec6fSJames Wright 
290a5b0ec6fSJames Wright   PetscFunctionBeginUser;
2911d2a9659SKenneth E. Jansen 
292752a08a7SJames Wright   if (problem->apply_vol_ifunction.qfunction && problem->uses_newtonian) {
2931d2a9659SKenneth E. Jansen     NewtonianIdealGasContext gas;
2941d2a9659SKenneth E. Jansen     PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas));
2951d2a9659SKenneth E. Jansen     jac_data_size_vol += (gas->idl_enable ? 1 : 0);
2961d2a9659SKenneth E. Jansen     PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas));
2971d2a9659SKenneth E. Jansen   }
2981d2a9659SKenneth E. Jansen 
299a5b0ec6fSJames Wright   {  // Create bases and element restrictions
300bb85d312SJames Wright     DMLabel  domain_label = NULL;
301bb85d312SJames Wright     PetscInt label_value = 0, height = 0, dm_field = 0;
3020814d5a7SKenneth E. Jansen     DM       dm_coord;
3030814d5a7SKenneth E. Jansen 
304a5b0ec6fSJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
305bb85d312SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &ceed_data->basis_q));
306bb85d312SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &ceed_data->basis_x));
30777841947SLeila Ghaffari 
308bb85d312SJames Wright     PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, 0, &ceed_data->elem_restr_q));
309bb85d312SJames Wright     PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &ceed_data->elem_restr_x));
310bb85d312SJames Wright     PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_vol, &elem_restr_jd_i));
31164f98e98SJames Wright 
312a424bcd0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL));
313a424bcd0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL));
314a424bcd0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL));
315a424bcd0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL));
31664f98e98SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL));
31777841947SLeila Ghaffari 
31864f98e98SJames Wright     {  // -- Copy PETSc coordinate vector into CEED vector
31977841947SLeila Ghaffari       Vec X_loc;
3203796c488SJed Brown       DM  cdm;
32164f98e98SJames Wright 
3222b730f8bSJeremy L Thompson       PetscCall(DMGetCellCoordinateDM(dm, &cdm));
3232b730f8bSJeremy L Thompson       if (cdm) {
3242b730f8bSJeremy L Thompson         PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
3252b730f8bSJeremy L Thompson       } else {
3262b730f8bSJeremy L Thompson         PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
3273796c488SJed Brown       }
3282b730f8bSJeremy L Thompson       PetscCall(VecScale(X_loc, problem->dm_scale));
329d0593705SJames Wright       PetscCall(VecCopyPetscToCeed(X_loc, ceed_data->x_coord));
33064f98e98SJames Wright     }
33177841947SLeila Ghaffari 
332cfb075a4SJames Wright     PetscCall(QDataGet(ceed, dm, domain_label, label_value, ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord,
333cfb075a4SJames Wright                        &ceed_data->elem_restr_qd_i, &ceed_data->q_data, &problem->q_data_size_vol));
33464f98e98SJames Wright   }
33564f98e98SJames Wright 
33664f98e98SJames Wright   {  // -- Create QFunction for ICs
337e9c36be0SJames Wright     CeedBasis     basis_xc;
33864f98e98SJames Wright     CeedQFunction qf_ics;
3395263e9c6SJames Wright     CeedOperator  op_ics;
34064f98e98SJames Wright 
341e9c36be0SJames Wright     PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_xc));
34264f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &qf_ics));
34364f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ics, problem->ics.qfunction_context));
34464f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ics, 0));
34564f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "x", num_comp_x, CEED_EVAL_INTERP));
34664f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
34764f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ics, "q0", num_comp_q, CEED_EVAL_NONE));
34864f98e98SJames Wright 
34964f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics));
350e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
351e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
352356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
353a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label));
3545263e9c6SJames Wright     PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx));
35577841947SLeila Ghaffari 
356e9c36be0SJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_xc));
35764f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ics));
35864f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics));
35964f98e98SJames Wright   }
36064f98e98SJames Wright 
36164f98e98SJames Wright   if (problem->apply_vol_rhs.qfunction) {
36264f98e98SJames Wright     CeedQFunction qf_rhs_vol;
36364f98e98SJames Wright 
36464f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &qf_rhs_vol));
36564f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_vol, problem->apply_vol_rhs.qfunction_context));
36664f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_rhs_vol, 0));
36764f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP));
36864f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
369cfb075a4SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "qdata", problem->q_data_size_vol, CEED_EVAL_NONE));
37064f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP));
37164f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP));
37264f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
37364f98e98SJames Wright 
374a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_vol, NULL, NULL, &op_rhs_vol));
375a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
376a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
377a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
378a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
379a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
380a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
38164f98e98SJames Wright 
38264f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_vol));
38377841947SLeila Ghaffari   }
38477841947SLeila Ghaffari 
38564f98e98SJames Wright   if (problem->apply_vol_ifunction.qfunction) {
38664f98e98SJames Wright     CeedQFunction qf_ifunction_vol;
38764f98e98SJames Wright 
38864f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc,
38964f98e98SJames Wright                                                     &qf_ifunction_vol));
39064f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context));
39164f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ifunction_vol, 0));
39264f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP));
39364f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
39464f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP));
395cfb075a4SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "qdata", problem->q_data_size_vol, CEED_EVAL_NONE));
39664f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP));
39764f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP));
39864f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
39964f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
40064f98e98SJames Wright 
401a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ifunction_vol, NULL, NULL, &op_ifunction_vol));
402a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
403a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
404a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed));
405a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
406a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
407a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
408a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
409a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
410a3ae0734SJed Brown 
41164f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction_vol));
41277841947SLeila Ghaffari   }
41377841947SLeila Ghaffari 
41464f98e98SJames Wright   if (problem->apply_vol_ijacobian.qfunction) {
41564f98e98SJames Wright     CeedQFunction qf_ijacobian_vol;
41664f98e98SJames Wright 
41764f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc,
41864f98e98SJames Wright                                                     &qf_ijacobian_vol));
41964f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context));
42064f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ijacobian_vol, 0));
42164f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP));
42264f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD));
423cfb075a4SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", problem->q_data_size_vol, CEED_EVAL_NONE));
42464f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
42564f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP));
42664f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
42764f98e98SJames Wright 
428a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op_ijacobian_vol));
429a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
430a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
431a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
432a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
433a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
434a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
43564f98e98SJames Wright 
436a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol));
437e334ad8fSJed Brown   }
438e334ad8fSJed Brown 
43977841947SLeila Ghaffari   // -- Create and apply CEED Composite Operator for the entire domain
44077841947SLeila Ghaffari   if (!user->phys->implicit) {  // RHS
4419ad5e8e4SJames Wright     CeedOperator op_rhs;
442e9c36be0SJames Wright 
443ed094490SJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, &op_rhs));
444ed094490SJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorCompositeAddSub(op_rhs, op_rhs_vol));
445e9c36be0SJames Wright     PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, op_rhs, NULL));
446e9c36be0SJames Wright 
4479ad5e8e4SJames Wright     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx));
448e9c36be0SJames Wright 
449e9c36be0SJames Wright     // ----- Get Context Labels for Operator
450e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "solution time", &user->phys->solution_time_label));
451e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "timestep size", &user->phys->timestep_size_label));
452e9c36be0SJames Wright 
453a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
454b4e9a8f8SJames Wright     PetscCall(CreateKSPMass(user, problem));
45577841947SLeila Ghaffari   } else {  // IFunction
45691c97f41SJames Wright     CeedOperator op_ijacobian = NULL;
45791c97f41SJames Wright 
458e9c36be0SJames Wright     // Create Composite Operaters
459ed094490SJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, &user->op_ifunction));
460ed094490SJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorCompositeAddSub(user->op_ifunction, op_ifunction_vol));
461e9c36be0SJames Wright     if (op_ijacobian_vol) {
462ed094490SJeremy L Thompson       PetscCallCeed(ceed, CeedOperatorCreateComposite(ceed, &op_ijacobian));
463ed094490SJeremy L Thompson       PetscCallCeed(ceed, CeedOperatorCompositeAddSub(op_ijacobian, op_ijacobian_vol));
464e9c36be0SJames Wright     }
465e9c36be0SJames Wright     PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, user->op_ifunction, op_ijacobian));
466e9c36be0SJames Wright 
467e9c36be0SJames Wright     // ----- Get Context Labels for Operator
468e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "solution time", &user->phys->solution_time_label));
469e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "timestep size", &user->phys->timestep_size_label));
470e9c36be0SJames Wright 
47191c97f41SJames Wright     if (op_ijacobian) {
472243afec9SJames Wright       PetscCall(MatCreateCeed(user->dm, user->dm, op_ijacobian, NULL, &user->mat_ijacobian));
47391c97f41SJames Wright       PetscCall(MatCeedSetLocalVectors(user->mat_ijacobian, user->Q_dot_loc, NULL));
47491c97f41SJames Wright       PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian));
475e334ad8fSJed Brown     }
476dada6cc0SJames Wright   }
477f5452247SJames Wright 
478577c3361SJames Wright   if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc));
479f5452247SJames Wright   if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem));
480ee3b213aSJames Wright   if (app_ctx->diff_filter_monitor && !user->diff_filter) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem));
481f5452247SJames Wright 
482a5b0ec6fSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data));
483a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i));
484a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol));
485a5b0ec6fSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ifunction_vol));
486a5b0ec6fSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_vol));
487ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
48877841947SLeila Ghaffari }
489