xref: /honee/src/setuplibceed.c (revision 866f9b4a9493f304f160970a0673bcccc0fa276c)
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 
14a515125bSLeila Ghaffari #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 
93*866f9b4aSJames Wright static PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height,
94*866f9b4aSJames Wright                                        CeedInt Q_sur, CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedBasis basis_q_sur,
95*866f9b4aSJames Wright                                        CeedBasis basis_x_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian, CeedOperator op_apply,
96*866f9b4aSJames Wright                                        CeedOperator op_apply_ijacobian) {
9715c18037SJames Wright   CeedVector          q_data_sur, jac_data_sur          = NULL;
98*866f9b4aSJames 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;
100*866f9b4aSJames Wright   PetscInt            dm_field = 0;
101368c645fSJames Wright 
10206f41313SJames Wright   PetscFunctionBeginUser;
103368c645fSJames Wright 
104368c645fSJames Wright   // ---- CEED Restriction
10515c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &elem_restr_q_sur));
10615c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &elem_restr_x_sur));
10715c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, q_data_size_sur, &elem_restr_qd_i_sur));
108368c645fSJames Wright   if (jac_data_size_sur > 0) {
109368c645fSJames Wright     // State-dependent data will be passed from residual to Jacobian. This will be collocated.
11015c18037SJames Wright     PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_sur, &elem_restr_jd_i_sur));
111b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL));
112368c645fSJames Wright   }
113368c645fSJames Wright 
114*866f9b4aSJames Wright   {  // Create q_data_sur vector
115*866f9b4aSJames Wright     CeedOperator op_setup_sur;
116368c645fSJames Wright 
117*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_qd_i_sur, &q_data_sur, NULL));
118*866f9b4aSJames Wright 
119b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur));
120*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, basis_x_sur, CEED_VECTOR_ACTIVE));
121*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, basis_x_sur, CEED_VECTOR_NONE));
12258e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
123368c645fSJames Wright 
124*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE));
125*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur));
126*866f9b4aSJames Wright   }
127*866f9b4aSJames Wright 
128368c645fSJames Wright   // ----- CEED Operator for Physics
129b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc));
130*866f9b4aSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
131*866f9b4aSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
13258e1cbfdSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
133*866f9b4aSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord));
134*866f9b4aSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
135b4c37c5cSJames Wright   if (elem_restr_jd_i_sur)
13658e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
137368c645fSJames Wright 
1384b96a86bSJames Wright   if (qf_apply_bc_jacobian && elem_restr_jd_i_sur) {
139b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian));
140*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
141*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
14258e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
143*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord));
14458e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
145*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
146368c645fSJames Wright   }
147368c645fSJames Wright 
148368c645fSJames Wright   // ----- Apply Sub-Operator for Physics
149*866f9b4aSJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply, op_apply_bc));
150*866f9b4aSJames Wright   if (op_apply_bc_jacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply_ijacobian, op_apply_bc_jacobian));
151368c645fSJames Wright 
152368c645fSJames Wright   // ----- Cleanup
153b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_sur));
154b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data_sur));
155b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur));
156b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur));
157b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i_sur));
158b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i_sur));
159b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc));
160b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc_jacobian));
161d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
162368c645fSJames Wright }
163368c645fSJames Wright 
164*866f9b4aSJames Wright static PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur,
1652b916ea7SJeremy L Thompson                                         PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian,
16625988f00SJames Wright                                         CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) {
16725988f00SJames Wright   PetscFunctionBeginUser;
16825988f00SJames Wright   if (apply_bc.qfunction) {
169b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc));
170b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context));
171ebfabadfSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc, 0));
172b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP));
173b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD));
174b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
175b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP));
176b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP));
177b4c37c5cSJames Wright     if (jac_data_size_sur) PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
17825988f00SJames Wright   }
17925988f00SJames Wright   if (apply_bc_jacobian.qfunction) {
180b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian));
181b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context));
182ebfabadfSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc_jacobian, 0));
183b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP));
184b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD));
185b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
186b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP));
187b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
188b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP));
18925988f00SJames Wright   }
190d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
19125988f00SJames Wright }
19225988f00SJames Wright 
193*866f9b4aSJames Wright // Utility function to create CEED Composite Operator for the entire domain
194*866f9b4aSJames Wright static PetscErrorCode AddBCSubOperators(User user, Ceed ceed, DM dm, SimpleBC bc, ProblemData problem, CeedData ceed_data, CeedOperator op_apply,
195*866f9b4aSJames Wright                                         CeedOperator op_apply_ijacobian) {
196*866f9b4aSJames Wright   CeedInt       height = 1, num_comp_q, num_comp_x;
197*866f9b4aSJames Wright   CeedInt       dim_sur, P_sur = user->app_ctx->degree + 1, Q_sur = P_sur + user->app_ctx->q_extra;
198*866f9b4aSJames Wright   const CeedInt q_data_size_sur = problem->q_data_size_sur, jac_data_size_sur = user->phys->implicit ? problem->jac_data_size_sur : 0;
199*866f9b4aSJames Wright   PetscInt      dim;
200*866f9b4aSJames Wright   DMLabel       face_sets_label;
201*866f9b4aSJames Wright   CeedBasis     basis_q_sur, basis_x_sur;
202*866f9b4aSJames Wright 
203*866f9b4aSJames Wright   PetscFunctionBeginUser;
204*866f9b4aSJames Wright   PetscCall(DMGetDimension(dm, &dim));
205*866f9b4aSJames Wright   dim_sur = dim - height;
206*866f9b4aSJames Wright   {  // Get number of components and coordinate dimension from op_apply
207*866f9b4aSJames Wright     CeedOperator       *sub_ops;
208*866f9b4aSJames Wright     CeedOperatorField   field;
209*866f9b4aSJames Wright     PetscInt            sub_op_index = 0;  // will be 0 for the volume op
210*866f9b4aSJames Wright     CeedElemRestriction elem_restr_q, elem_restr_x;
211*866f9b4aSJames Wright 
212*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(op_apply, &sub_ops));
213*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
214*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
215*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
216*866f9b4aSJames Wright 
217*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "x", &field));
218*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_x));
219*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
220*866f9b4aSJames Wright   }
221*866f9b4aSJames Wright 
222*866f9b4aSJames Wright   {  // Get bases
223*866f9b4aSJames Wright     DM dm_coord;
224*866f9b4aSJames Wright 
225*866f9b4aSJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
226*866f9b4aSJames Wright     DMLabel  label       = NULL;
227*866f9b4aSJames Wright     PetscInt label_value = 0;
228*866f9b4aSJames Wright     PetscInt field       = 0;  // Still want the normal, default field
229*866f9b4aSJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm, label, label_value, height, field, &basis_q_sur));
230*866f9b4aSJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, label_value, height, field, &basis_x_sur));
231*866f9b4aSJames Wright   }
232*866f9b4aSJames Wright 
233*866f9b4aSJames Wright   PetscCall(DMGetLabel(dm, "Face Sets", &face_sets_label));
234*866f9b4aSJames Wright 
235*866f9b4aSJames Wright   {  // -- Create QFunction for quadrature data
236*866f9b4aSJames Wright     PetscCallCeed(ceed,
237*866f9b4aSJames Wright                   CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, problem->setup_sur.qfunction_loc, &ceed_data->qf_setup_sur));
238*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_setup_sur, problem->setup_sur.qfunction_context));
239*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(ceed_data->qf_setup_sur, 0));
240*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x * dim_sur, CEED_EVAL_GRAD));
241*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT));
242*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
243*866f9b4aSJames Wright   }
244*866f9b4aSJames Wright 
245*866f9b4aSJames Wright   {  // --- Create Sub-Operator for inflow boundaries
246*866f9b4aSJames Wright     CeedQFunction qf_apply_inflow = NULL, qf_apply_inflow_jacobian = NULL;
247*866f9b4aSJames Wright 
248*866f9b4aSJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow,
249*866f9b4aSJames Wright                                 problem->apply_inflow_jacobian, &qf_apply_inflow, &qf_apply_inflow_jacobian));
250*866f9b4aSJames Wright     for (CeedInt i = 0; i < bc->num_inflow; i++) {
251*866f9b4aSJames 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,
252*866f9b4aSJames Wright                                  basis_x_sur, qf_apply_inflow, qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian));
253*866f9b4aSJames Wright     }
254*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow));
255*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow_jacobian));
256*866f9b4aSJames Wright   }
257*866f9b4aSJames Wright 
258*866f9b4aSJames Wright   {  // --- Create Sub-Operator for outflow boundaries
259*866f9b4aSJames Wright     CeedQFunction qf_apply_outflow = NULL, qf_apply_outflow_jacobian = NULL;
260*866f9b4aSJames Wright 
261*866f9b4aSJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow,
262*866f9b4aSJames Wright                                 problem->apply_outflow_jacobian, &qf_apply_outflow, &qf_apply_outflow_jacobian));
263*866f9b4aSJames Wright     for (CeedInt i = 0; i < bc->num_outflow; i++) {
264*866f9b4aSJames Wright       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
265*866f9b4aSJames Wright                                  basis_q_sur, basis_x_sur, qf_apply_outflow, qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian));
266*866f9b4aSJames Wright     }
267*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow));
268*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow_jacobian));
269*866f9b4aSJames Wright   }
270*866f9b4aSJames Wright 
271*866f9b4aSJames Wright   {  // --- Create Sub-Operator for freestream boundaries
272*866f9b4aSJames Wright     CeedQFunction qf_apply_freestream = NULL, qf_apply_freestream_jacobian = NULL;
273*866f9b4aSJames Wright 
274*866f9b4aSJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream,
275*866f9b4aSJames Wright                                 problem->apply_freestream_jacobian, &qf_apply_freestream, &qf_apply_freestream_jacobian));
276*866f9b4aSJames Wright     for (CeedInt i = 0; i < bc->num_freestream; i++) {
277*866f9b4aSJames Wright       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
278*866f9b4aSJames Wright                                  basis_q_sur, basis_x_sur, qf_apply_freestream, qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian));
279*866f9b4aSJames Wright     }
280*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream));
281*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream_jacobian));
282*866f9b4aSJames Wright   }
283*866f9b4aSJames Wright 
284*866f9b4aSJames Wright   {  // --- Create Sub-Operator for slip boundaries
285*866f9b4aSJames Wright     CeedQFunction qf_apply_slip = NULL, qf_apply_slip_jacobian = NULL;
286*866f9b4aSJames Wright 
287*866f9b4aSJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_slip,
288*866f9b4aSJames Wright                                 problem->apply_slip_jacobian, &qf_apply_slip, &qf_apply_slip_jacobian));
289*866f9b4aSJames Wright     for (CeedInt i = 0; i < bc->num_slip; i++) {
290*866f9b4aSJames 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,
291*866f9b4aSJames Wright                                  basis_x_sur, qf_apply_slip, qf_apply_slip_jacobian, op_apply, op_apply_ijacobian));
292*866f9b4aSJames Wright     }
293*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip));
294*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip_jacobian));
295*866f9b4aSJames Wright   }
296*866f9b4aSJames Wright 
297*866f9b4aSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q_sur));
298*866f9b4aSJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_sur));
299*866f9b4aSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
300*866f9b4aSJames Wright }
301*866f9b4aSJames Wright 
302991aef52SJames Wright PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData problem, SimpleBC bc) {
303a515125bSLeila Ghaffari   PetscFunctionBeginUser;
304a515125bSLeila Ghaffari   const PetscInt num_comp_q = 5;
30594a7b3d2SKenneth E. Jansen   const CeedInt  dim = problem->dim, num_comp_x = problem->dim, q_data_size_vol = problem->q_data_size_vol;
30694a7b3d2SKenneth E. Jansen   CeedInt        jac_data_size_vol = num_comp_q + 6 + 3;
30794a7b3d2SKenneth E. Jansen 
3087d8a615bSJames Wright   if (problem->apply_vol_ifunction.qfunction && problem->uses_newtonian) {
30994a7b3d2SKenneth E. Jansen     NewtonianIdealGasContext gas;
31094a7b3d2SKenneth E. Jansen     PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas));
31194a7b3d2SKenneth E. Jansen     jac_data_size_vol += (gas->idl_enable ? 1 : 0);
31294a7b3d2SKenneth E. Jansen     PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas));
31394a7b3d2SKenneth E. Jansen   }
31494a7b3d2SKenneth E. Jansen 
315752f40e3SJed Brown   CeedElemRestriction elem_restr_jd_i;
316752f40e3SJed Brown   CeedVector          jac_data;
31767263decSKenneth E. Jansen   CeedInt             num_qpts;
31815c18037SJames Wright   DMLabel             domain_label = NULL;
31915c18037SJames Wright   PetscInt            label_value = 0, height = 0, dm_field = 0;
320a515125bSLeila Ghaffari 
32167263decSKenneth E. Jansen   DM dm_coord;
32267263decSKenneth E. Jansen   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
32367263decSKenneth E. Jansen 
32415c18037SJames Wright   PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &ceed_data->basis_q));
32515c18037SJames Wright   PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &ceed_data->basis_x));
326b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(ceed_data->basis_q, &num_qpts));
327a515125bSLeila Ghaffari 
32815c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, 0, &ceed_data->elem_restr_q));
32915c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &ceed_data->elem_restr_x));
33015c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, q_data_size_vol, &ceed_data->elem_restr_qd_i));
33115c18037SJames Wright   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_vol, &elem_restr_jd_i));
332ce8bebb6SJames Wright 
333b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL));
334b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL));
335b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL));
336b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL));
337ce8bebb6SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_qd_i, &ceed_data->q_data, NULL));
338ce8bebb6SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL));
339a515125bSLeila Ghaffari 
340ce8bebb6SJames Wright   {  // -- Copy PETSc coordinate vector into CEED vector
341a515125bSLeila Ghaffari     Vec X_loc;
342cac8aa24SJed Brown     DM  cdm;
343ce8bebb6SJames Wright 
3442b916ea7SJeremy L Thompson     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
3452b916ea7SJeremy L Thompson     if (cdm) {
3462b916ea7SJeremy L Thompson       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
3472b916ea7SJeremy L Thompson     } else {
3482b916ea7SJeremy L Thompson       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
349cac8aa24SJed Brown     }
3502b916ea7SJeremy L Thompson     PetscCall(VecScale(X_loc, problem->dm_scale));
351a7dac1d5SJames Wright     PetscCall(VecCopyPetscToCeed(X_loc, ceed_data->x_coord));
352ce8bebb6SJames Wright   }
353a515125bSLeila Ghaffari 
354ce8bebb6SJames Wright   {  // -- Create quadrature data
355ce8bebb6SJames Wright     CeedQFunction qf_setup_vol;
356ce8bebb6SJames Wright     CeedOperator  op_setup_vol;
3570b0430b7SJames Wright 
358ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, problem->setup_vol.qfunction_loc, &qf_setup_vol));
359ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_setup_vol, problem->setup_vol.qfunction_context));
360ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_setup_vol, 0));
361ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_vol, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
362ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT));
363ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
364a515125bSLeila Ghaffari 
365ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup_vol, NULL, NULL, &op_setup_vol));
366ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "dx", ceed_data->elem_restr_x, ceed_data->basis_x, CEED_VECTOR_ACTIVE));
367ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE));
368ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
369ce8bebb6SJames Wright 
370ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedOperatorApply(op_setup_vol, ceed_data->x_coord, ceed_data->q_data, CEED_REQUEST_IMMEDIATE));
371ce8bebb6SJames Wright 
372ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup_vol));
373ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_vol));
374ce8bebb6SJames Wright   }
375ce8bebb6SJames Wright 
376ce8bebb6SJames Wright   {  // -- Create QFunction for ICs
377*866f9b4aSJames Wright     CeedBasis     basis_xc;
378ce8bebb6SJames Wright     CeedQFunction qf_ics;
3798f18bb8bSJames Wright     CeedOperator  op_ics;
380ce8bebb6SJames Wright 
381*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_xc));
382ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &qf_ics));
383ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ics, problem->ics.qfunction_context));
384ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ics, 0));
385ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "x", num_comp_x, CEED_EVAL_INTERP));
386ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
387ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ics, "q0", num_comp_q, CEED_EVAL_NONE));
388ce8bebb6SJames Wright 
389ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics));
390*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
391*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
39258e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
393b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label));
3948f18bb8bSJames Wright     PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx));
395a515125bSLeila Ghaffari 
396*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_xc));
397ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ics));
398ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics));
399ce8bebb6SJames Wright   }
400ce8bebb6SJames Wright 
401ce8bebb6SJames Wright   if (problem->apply_vol_rhs.qfunction) {
402ce8bebb6SJames Wright     CeedQFunction qf_rhs_vol;
403a515125bSLeila Ghaffari     CeedOperator  op;
404ce8bebb6SJames Wright 
405ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &qf_rhs_vol));
406ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_vol, problem->apply_vol_rhs.qfunction_context));
407ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_rhs_vol, 0));
408ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP));
409ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
410ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
411ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP));
412ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP));
413ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
414ce8bebb6SJames Wright 
415ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_vol, NULL, NULL, &op));
416b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
417b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
41858e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
419b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
420b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
421b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
422a515125bSLeila Ghaffari     user->op_rhs_vol = op;
423ce8bebb6SJames Wright 
424ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_vol));
425a515125bSLeila Ghaffari   }
426a515125bSLeila Ghaffari 
427ce8bebb6SJames Wright   if (problem->apply_vol_ifunction.qfunction) {
428ce8bebb6SJames Wright     CeedQFunction qf_ifunction_vol;
429a515125bSLeila Ghaffari     CeedOperator  op;
430ce8bebb6SJames Wright 
431ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc,
432ce8bebb6SJames Wright                                                     &qf_ifunction_vol));
433ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context));
434ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ifunction_vol, 0));
435ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP));
436ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
437ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP));
438ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
439ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP));
440ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP));
441ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
442ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
443ce8bebb6SJames Wright 
444ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ifunction_vol, NULL, NULL, &op));
445b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
446b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
447b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed));
44858e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
449b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
450b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
451b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
45258e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
453752f40e3SJed Brown 
454a515125bSLeila Ghaffari     user->op_ifunction_vol = op;
455ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction_vol));
456a515125bSLeila Ghaffari   }
457a515125bSLeila Ghaffari 
458f0b65372SJed Brown   CeedOperator op_ijacobian_vol = NULL;
459ce8bebb6SJames Wright   if (problem->apply_vol_ijacobian.qfunction) {
460ce8bebb6SJames Wright     CeedQFunction qf_ijacobian_vol;
461f0b65372SJed Brown     CeedOperator  op;
462ce8bebb6SJames Wright 
463ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc,
464ce8bebb6SJames Wright                                                     &qf_ijacobian_vol));
465ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context));
466ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ijacobian_vol, 0));
467ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP));
468ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD));
469ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
470ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
471ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP));
472ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
473ce8bebb6SJames Wright 
474b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op));
475b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
476b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
47758e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
47858e1cbfdSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
479b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
480b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
481f0b65372SJed Brown     op_ijacobian_vol = op;
482ce8bebb6SJames Wright 
483b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol));
484f0b65372SJed Brown   }
485f0b65372SJed Brown 
486a515125bSLeila Ghaffari   // -- Create and apply CEED Composite Operator for the entire domain
487a515125bSLeila Ghaffari   if (!user->phys->implicit) {  // RHS
488da5fe0e4SJames Wright     CeedOperator op_rhs;
489*866f9b4aSJames Wright 
490*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_rhs));
491*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_rhs, user->op_rhs_vol));
492*866f9b4aSJames Wright     PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, op_rhs, NULL));
493*866f9b4aSJames Wright 
494da5fe0e4SJames Wright     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx));
495*866f9b4aSJames Wright 
496*866f9b4aSJames Wright     // ----- Get Context Labels for Operator
497*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "solution time", &user->phys->solution_time_label));
498*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "timestep size", &user->phys->timestep_size_label));
499*866f9b4aSJames Wright 
500b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
501a78efa86SJames Wright     PetscCall(CreateKSPMass(user, problem));
502c2376cc9SJames Wright     PetscCheck(app_ctx->sgs_model_type == SGS_MODEL_NONE, user->comm, PETSC_ERR_SUP, "SGS modeling not implemented for explicit timestepping");
503a515125bSLeila Ghaffari   } else {  // IFunction
504ebfabadfSJames Wright     CeedOperator op_ijacobian = NULL;
505ebfabadfSJames Wright 
506*866f9b4aSJames Wright     // Create Composite Operaters
507*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &user->op_ifunction));
508*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(user->op_ifunction, user->op_ifunction_vol));
509*866f9b4aSJames Wright     if (op_ijacobian_vol) {
510*866f9b4aSJames Wright       PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_ijacobian));
511*866f9b4aSJames Wright       PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_ijacobian, op_ijacobian_vol));
512*866f9b4aSJames Wright     }
513*866f9b4aSJames Wright     PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, user->op_ifunction, op_ijacobian));
514*866f9b4aSJames Wright 
515*866f9b4aSJames Wright     // ----- Get Context Labels for Operator
516*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "solution time", &user->phys->solution_time_label));
517*866f9b4aSJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "timestep size", &user->phys->timestep_size_label));
518*866f9b4aSJames Wright 
519ebfabadfSJames Wright     if (op_ijacobian) {
520ebfabadfSJames Wright       PetscCall(MatCeedCreate(user->dm, user->dm, op_ijacobian, NULL, &user->mat_ijacobian));
521ebfabadfSJames Wright       PetscCall(MatCeedSetLocalVectors(user->mat_ijacobian, user->Q_dot_loc, NULL));
522ebfabadfSJames Wright       PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ijacobian, "ijacobian time shift", &user->phys->ijacobian_time_shift_label));
523ebfabadfSJames Wright       PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian));
524f0b65372SJed Brown     }
525ad494f68SJames Wright     if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SgsDDSetup(ceed, user, ceed_data, problem));
5266d0190e2SJames Wright   }
52791933550SJames Wright 
528c2376cc9SJames Wright   if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc));
52991933550SJames Wright   if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem));
530aa0b7f76SJames Wright   if (app_ctx->diff_filter_monitor && !user->diff_filter) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem));
5311c17f66aSJames Wright   if (app_ctx->sgs_train_enable) PetscCall(SGS_DD_TrainingSetup(ceed, user, ceed_data, problem));
53291933550SJames Wright 
533b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i));
534b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol));
535b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data));
536d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
537a515125bSLeila Ghaffari }
538