xref: /libCEED/examples/fluids/src/setuplibceed.c (revision e9c36be085200ad84b427d3f41921f3626bb31f1)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, 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 
31b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(user->op_rhs_ctx->op, &sub_ops));
32b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
33b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
34b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetBasis(field, &basis_q));
35b4e9a8f8SJames Wright 
36b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &field));
37b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_qd_i));
38b4e9a8f8SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetVector(field, &q_data));
39b4e9a8f8SJames Wright   }
40b4e9a8f8SJames Wright 
41b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
42b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd_i, &q_data_size));
43b4e9a8f8SJames Wright 
44b4e9a8f8SJames Wright   PetscCall(CreateMassQFunction(ceed, num_comp_q, q_data_size, &qf_mass));
45b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
46b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "u", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
47b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data));
48b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
49b4e9a8f8SJames Wright 
50b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
51b4e9a8f8SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
52b4e9a8f8SJames Wright }
53b4e9a8f8SJames Wright 
54b4e9a8f8SJames Wright // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes
55731c13d7SJames Wright static PetscErrorCode CreateKSPMass(User user, ProblemData problem) {
56b4e9a8f8SJames Wright   Ceed         ceed = user->ceed;
57b4e9a8f8SJames Wright   DM           dm   = user->dm;
58b4e9a8f8SJames Wright   CeedOperator op_mass;
59b4e9a8f8SJames Wright 
60b4e9a8f8SJames Wright   PetscFunctionBeginUser;
61b4e9a8f8SJames Wright   if (problem->create_mass_operator) PetscCall(problem->create_mass_operator(user, &op_mass));
62b4e9a8f8SJames Wright   else PetscCall(CreateKSPMassOperator_Unstabilized(user, &op_mass));
63b4e9a8f8SJames Wright 
64b4e9a8f8SJames Wright   {  // -- Setup KSP for mass operator
65b4e9a8f8SJames Wright     Mat      mat_mass;
66b4e9a8f8SJames Wright     Vec      Zeros_loc;
67b4e9a8f8SJames Wright     MPI_Comm comm = PetscObjectComm((PetscObject)dm);
68b4e9a8f8SJames Wright 
69b4e9a8f8SJames Wright     PetscCall(DMCreateLocalVector(dm, &Zeros_loc));
70b4e9a8f8SJames Wright     PetscCall(VecZeroEntries(Zeros_loc));
71b4e9a8f8SJames Wright     PetscCall(MatCeedCreate(dm, dm, op_mass, NULL, &mat_mass));
72b4e9a8f8SJames Wright     PetscCall(MatCeedSetLocalVectors(mat_mass, Zeros_loc, NULL));
73b4e9a8f8SJames Wright 
74b4e9a8f8SJames Wright     PetscCall(KSPCreate(comm, &user->mass_ksp));
75b4e9a8f8SJames Wright     PetscCall(KSPSetOptionsPrefix(user->mass_ksp, "mass_"));
76b4e9a8f8SJames Wright     {  // lumped by default
77b4e9a8f8SJames Wright       PC pc;
78b4e9a8f8SJames Wright       PetscCall(KSPGetPC(user->mass_ksp, &pc));
79b4e9a8f8SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
80b4e9a8f8SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
81b4e9a8f8SJames Wright       PetscCall(KSPSetType(user->mass_ksp, KSPPREONLY));
82b4e9a8f8SJames Wright     }
83b4e9a8f8SJames Wright     PetscCall(KSPSetFromOptions_WithMatCeed(user->mass_ksp, mat_mass));
84b4e9a8f8SJames Wright     PetscCall(KSPSetFromOptions(user->mass_ksp));
85b4e9a8f8SJames Wright     PetscCall(VecDestroy(&Zeros_loc));
86b4e9a8f8SJames Wright     PetscCall(MatDestroy(&mat_mass));
87b4e9a8f8SJames Wright   }
88b4e9a8f8SJames Wright 
89b4e9a8f8SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
90b4e9a8f8SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
91b4e9a8f8SJames Wright }
92b4e9a8f8SJames Wright 
93*e9c36be0SJames Wright static PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height,
94*e9c36be0SJames Wright                                        CeedInt Q_sur, CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedBasis basis_q_sur,
95*e9c36be0SJames Wright                                        CeedBasis basis_x_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian, CeedOperator op_apply,
96*e9c36be0SJames Wright                                        CeedOperator op_apply_ijacobian) {
97bb85d312SJames Wright   CeedVector          q_data_sur, jac_data_sur          = NULL;
98*e9c36be0SJames Wright   CeedOperator        op_apply_bc, op_apply_bc_jacobian = NULL;
99bb85d312SJames Wright   CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur = NULL;
100*e9c36be0SJames Wright   PetscInt            dm_field = 0;
1019eb9c072SJames Wright 
102f17d818dSJames Wright   PetscFunctionBeginUser;
1039eb9c072SJames Wright 
1049eb9c072SJames Wright   // ---- CEED Restriction
105bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &elem_restr_q_sur));
106bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &elem_restr_x_sur));
107bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, q_data_size_sur, &elem_restr_qd_i_sur));
1089eb9c072SJames Wright   if (jac_data_size_sur > 0) {
1099eb9c072SJames Wright     // State-dependent data will be passed from residual to Jacobian. This will be collocated.
110bb85d312SJames Wright     PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_sur, &elem_restr_jd_i_sur));
111a424bcd0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL));
1129eb9c072SJames Wright   }
1139eb9c072SJames Wright 
114*e9c36be0SJames Wright   {  // Create q_data_sur vector
115*e9c36be0SJames Wright     CeedOperator op_setup_sur;
1169eb9c072SJames Wright 
117*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_qd_i_sur, &q_data_sur, NULL));
118*e9c36be0SJames Wright 
119a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur));
120*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, basis_x_sur, CEED_VECTOR_ACTIVE));
121*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, basis_x_sur, CEED_VECTOR_NONE));
122356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
1239eb9c072SJames Wright 
124*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE));
125*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur));
126*e9c36be0SJames Wright   }
127*e9c36be0SJames Wright 
1289eb9c072SJames Wright   // ----- CEED Operator for Physics
129a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc));
130*e9c36be0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
131*e9c36be0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
132356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
133*e9c36be0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord));
134*e9c36be0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
135a424bcd0SJames Wright   if (elem_restr_jd_i_sur)
136356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
1379eb9c072SJames Wright 
138f21e6b1cSJames Wright   if (qf_apply_bc_jacobian && elem_restr_jd_i_sur) {
139a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian));
140*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
141*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
142356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
143*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord));
144356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
145*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
1469eb9c072SJames Wright   }
1479eb9c072SJames Wright 
1489eb9c072SJames Wright   // ----- Apply Sub-Operator for Physics
149*e9c36be0SJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply, op_apply_bc));
150*e9c36be0SJames Wright   if (op_apply_bc_jacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply_ijacobian, op_apply_bc_jacobian));
1519eb9c072SJames Wright 
1529eb9c072SJames Wright   // ----- Cleanup
153a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_sur));
154a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data_sur));
155a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur));
156a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur));
157a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i_sur));
158a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i_sur));
159a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc));
160a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc_jacobian));
161ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1629eb9c072SJames Wright }
1639eb9c072SJames Wright 
164*e9c36be0SJames Wright static PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur,
1652b730f8bSJeremy L Thompson                                         PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian,
1665b881d11SJames Wright                                         CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) {
1675b881d11SJames Wright   PetscFunctionBeginUser;
1685b881d11SJames Wright   if (apply_bc.qfunction) {
169a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc));
170a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context));
17191c97f41SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc, 0));
172a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP));
173a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD));
174a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
175a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP));
176a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP));
177a424bcd0SJames Wright     if (jac_data_size_sur) PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
1785b881d11SJames Wright   }
1795b881d11SJames Wright   if (apply_bc_jacobian.qfunction) {
180a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian));
181a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context));
18291c97f41SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc_jacobian, 0));
183a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP));
184a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD));
185a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
186a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP));
187a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
188a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP));
1895b881d11SJames Wright   }
190ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1915b881d11SJames Wright }
1925b881d11SJames Wright 
193*e9c36be0SJames Wright // Utility function to create CEED Composite Operator for the entire domain
194*e9c36be0SJames Wright static PetscErrorCode AddBCSubOperators(User user, Ceed ceed, DM dm, SimpleBC bc, ProblemData problem, CeedData ceed_data, CeedOperator op_apply,
195*e9c36be0SJames Wright                                         CeedOperator op_apply_ijacobian) {
196*e9c36be0SJames Wright   CeedInt       height = 1, num_comp_q, num_comp_x;
197*e9c36be0SJames Wright   CeedInt       dim_sur, P_sur = user->app_ctx->degree + 1, Q_sur = P_sur + user->app_ctx->q_extra;
198*e9c36be0SJames 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*e9c36be0SJames Wright   PetscInt      dim;
200*e9c36be0SJames Wright   DMLabel       face_sets_label;
201*e9c36be0SJames Wright   CeedBasis     basis_q_sur, basis_x_sur;
202*e9c36be0SJames Wright 
203*e9c36be0SJames Wright   PetscFunctionBeginUser;
204*e9c36be0SJames Wright   PetscCall(DMGetDimension(dm, &dim));
205*e9c36be0SJames Wright   dim_sur = dim - height;
206*e9c36be0SJames Wright   {  // Get number of components and coordinate dimension from op_apply
207*e9c36be0SJames Wright     CeedOperator       *sub_ops;
208*e9c36be0SJames Wright     CeedOperatorField   field;
209*e9c36be0SJames Wright     PetscInt            sub_op_index = 0;  // will be 0 for the volume op
210*e9c36be0SJames Wright     CeedElemRestriction elem_restr_q, elem_restr_x;
211*e9c36be0SJames Wright 
212*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(op_apply, &sub_ops));
213*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
214*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
215*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
216*e9c36be0SJames Wright 
217*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "x", &field));
218*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_x));
219*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
220*e9c36be0SJames Wright   }
221*e9c36be0SJames Wright 
222*e9c36be0SJames Wright   {  // Get bases
223*e9c36be0SJames Wright     DM dm_coord;
224*e9c36be0SJames Wright 
225*e9c36be0SJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
226*e9c36be0SJames Wright     DMLabel  label       = NULL;
227*e9c36be0SJames Wright     PetscInt label_value = 0;
228*e9c36be0SJames Wright     PetscInt field       = 0;  // Still want the normal, default field
229*e9c36be0SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm, label, label_value, height, field, &basis_q_sur));
230*e9c36be0SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, label_value, height, field, &basis_x_sur));
231*e9c36be0SJames Wright   }
232*e9c36be0SJames Wright 
233*e9c36be0SJames Wright   PetscCall(DMGetLabel(dm, "Face Sets", &face_sets_label));
234*e9c36be0SJames Wright 
235*e9c36be0SJames Wright   {  // -- Create QFunction for quadrature data
236*e9c36be0SJames Wright     PetscCallCeed(ceed,
237*e9c36be0SJames Wright                   CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, problem->setup_sur.qfunction_loc, &ceed_data->qf_setup_sur));
238*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_setup_sur, problem->setup_sur.qfunction_context));
239*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(ceed_data->qf_setup_sur, 0));
240*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x * dim_sur, CEED_EVAL_GRAD));
241*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT));
242*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
243*e9c36be0SJames Wright   }
244*e9c36be0SJames Wright 
245*e9c36be0SJames Wright   {  // --- Create Sub-Operator for inflow boundaries
246*e9c36be0SJames Wright     CeedQFunction qf_apply_inflow = NULL, qf_apply_inflow_jacobian = NULL;
247*e9c36be0SJames Wright 
248*e9c36be0SJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow,
249*e9c36be0SJames Wright                                 problem->apply_inflow_jacobian, &qf_apply_inflow, &qf_apply_inflow_jacobian));
250*e9c36be0SJames Wright     for (CeedInt i = 0; i < bc->num_inflow; i++) {
251*e9c36be0SJames 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*e9c36be0SJames Wright                                  basis_x_sur, qf_apply_inflow, qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian));
253*e9c36be0SJames Wright     }
254*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow));
255*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow_jacobian));
256*e9c36be0SJames Wright   }
257*e9c36be0SJames Wright 
258*e9c36be0SJames Wright   {  // --- Create Sub-Operator for outflow boundaries
259*e9c36be0SJames Wright     CeedQFunction qf_apply_outflow = NULL, qf_apply_outflow_jacobian = NULL;
260*e9c36be0SJames Wright 
261*e9c36be0SJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow,
262*e9c36be0SJames Wright                                 problem->apply_outflow_jacobian, &qf_apply_outflow, &qf_apply_outflow_jacobian));
263*e9c36be0SJames Wright     for (CeedInt i = 0; i < bc->num_outflow; i++) {
264*e9c36be0SJames 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*e9c36be0SJames Wright                                  basis_q_sur, basis_x_sur, qf_apply_outflow, qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian));
266*e9c36be0SJames Wright     }
267*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow));
268*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow_jacobian));
269*e9c36be0SJames Wright   }
270*e9c36be0SJames Wright 
271*e9c36be0SJames Wright   {  // --- Create Sub-Operator for freestream boundaries
272*e9c36be0SJames Wright     CeedQFunction qf_apply_freestream = NULL, qf_apply_freestream_jacobian = NULL;
273*e9c36be0SJames Wright 
274*e9c36be0SJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream,
275*e9c36be0SJames Wright                                 problem->apply_freestream_jacobian, &qf_apply_freestream, &qf_apply_freestream_jacobian));
276*e9c36be0SJames Wright     for (CeedInt i = 0; i < bc->num_freestream; i++) {
277*e9c36be0SJames 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*e9c36be0SJames Wright                                  basis_q_sur, basis_x_sur, qf_apply_freestream, qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian));
279*e9c36be0SJames Wright     }
280*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream));
281*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream_jacobian));
282*e9c36be0SJames Wright   }
283*e9c36be0SJames Wright 
284*e9c36be0SJames Wright   {  // --- Create Sub-Operator for slip boundaries
285*e9c36be0SJames Wright     CeedQFunction qf_apply_slip = NULL, qf_apply_slip_jacobian = NULL;
286*e9c36be0SJames Wright 
287*e9c36be0SJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_slip,
288*e9c36be0SJames Wright                                 problem->apply_slip_jacobian, &qf_apply_slip, &qf_apply_slip_jacobian));
289*e9c36be0SJames Wright     for (CeedInt i = 0; i < bc->num_slip; i++) {
290*e9c36be0SJames 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*e9c36be0SJames Wright                                  basis_x_sur, qf_apply_slip, qf_apply_slip_jacobian, op_apply, op_apply_ijacobian));
292*e9c36be0SJames Wright     }
293*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip));
294*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip_jacobian));
295*e9c36be0SJames Wright   }
296*e9c36be0SJames Wright 
297*e9c36be0SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q_sur));
298*e9c36be0SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_sur));
299*e9c36be0SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
300*e9c36be0SJames Wright }
301*e9c36be0SJames Wright 
302731c13d7SJames Wright PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData problem, SimpleBC bc) {
30377841947SLeila Ghaffari   PetscFunctionBeginUser;
30477841947SLeila Ghaffari   const PetscInt num_comp_q = 5;
3051d2a9659SKenneth E. Jansen   const CeedInt  dim = problem->dim, num_comp_x = problem->dim, q_data_size_vol = problem->q_data_size_vol;
3061d2a9659SKenneth E. Jansen   CeedInt        jac_data_size_vol = num_comp_q + 6 + 3;
3071d2a9659SKenneth E. Jansen 
308752a08a7SJames Wright   if (problem->apply_vol_ifunction.qfunction && problem->uses_newtonian) {
3091d2a9659SKenneth E. Jansen     NewtonianIdealGasContext gas;
3101d2a9659SKenneth E. Jansen     PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas));
3111d2a9659SKenneth E. Jansen     jac_data_size_vol += (gas->idl_enable ? 1 : 0);
3121d2a9659SKenneth E. Jansen     PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas));
3131d2a9659SKenneth E. Jansen   }
3141d2a9659SKenneth E. Jansen 
315a3ae0734SJed Brown   CeedElemRestriction elem_restr_jd_i;
316a3ae0734SJed Brown   CeedVector          jac_data;
3170814d5a7SKenneth E. Jansen   CeedInt             num_qpts;
318bb85d312SJames Wright   DMLabel             domain_label = NULL;
319bb85d312SJames Wright   PetscInt            label_value = 0, height = 0, dm_field = 0;
32077841947SLeila Ghaffari 
3210814d5a7SKenneth E. Jansen   DM dm_coord;
3220814d5a7SKenneth E. Jansen   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
3230814d5a7SKenneth E. Jansen 
324bb85d312SJames Wright   PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &ceed_data->basis_q));
325bb85d312SJames Wright   PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &ceed_data->basis_x));
326a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(ceed_data->basis_q, &num_qpts));
32777841947SLeila Ghaffari 
328bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, 0, &ceed_data->elem_restr_q));
329bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &ceed_data->elem_restr_x));
330bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, q_data_size_vol, &ceed_data->elem_restr_qd_i));
331bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_vol, &elem_restr_jd_i));
33264f98e98SJames Wright 
333a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL));
334a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL));
335a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL));
336a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL));
33764f98e98SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_qd_i, &ceed_data->q_data, NULL));
33864f98e98SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL));
33977841947SLeila Ghaffari 
34064f98e98SJames Wright   {  // -- Copy PETSc coordinate vector into CEED vector
34177841947SLeila Ghaffari     Vec X_loc;
3423796c488SJed Brown     DM  cdm;
34364f98e98SJames Wright 
3442b730f8bSJeremy L Thompson     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
3452b730f8bSJeremy L Thompson     if (cdm) {
3462b730f8bSJeremy L Thompson       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
3472b730f8bSJeremy L Thompson     } else {
3482b730f8bSJeremy L Thompson       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
3493796c488SJed Brown     }
3502b730f8bSJeremy L Thompson     PetscCall(VecScale(X_loc, problem->dm_scale));
351d0593705SJames Wright     PetscCall(VecCopyPetscToCeed(X_loc, ceed_data->x_coord));
35264f98e98SJames Wright   }
35377841947SLeila Ghaffari 
35464f98e98SJames Wright   {  // -- Create quadrature data
35564f98e98SJames Wright     CeedQFunction qf_setup_vol;
35664f98e98SJames Wright     CeedOperator  op_setup_vol;
357d52d2babSJames Wright 
35864f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, problem->setup_vol.qfunction_loc, &qf_setup_vol));
35964f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_setup_vol, problem->setup_vol.qfunction_context));
36064f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_setup_vol, 0));
36164f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_vol, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
36264f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT));
36364f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
36477841947SLeila Ghaffari 
36564f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup_vol, NULL, NULL, &op_setup_vol));
36664f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "dx", ceed_data->elem_restr_x, ceed_data->basis_x, CEED_VECTOR_ACTIVE));
36764f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE));
36864f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
36964f98e98SJames Wright 
37064f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorApply(op_setup_vol, ceed_data->x_coord, ceed_data->q_data, CEED_REQUEST_IMMEDIATE));
37164f98e98SJames Wright 
37264f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup_vol));
37364f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_vol));
37464f98e98SJames Wright   }
37564f98e98SJames Wright 
37664f98e98SJames Wright   {  // -- Create QFunction for ICs
377*e9c36be0SJames Wright     CeedBasis     basis_xc;
37864f98e98SJames Wright     CeedQFunction qf_ics;
3795263e9c6SJames Wright     CeedOperator  op_ics;
38064f98e98SJames Wright 
381*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_xc));
38264f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &qf_ics));
38364f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ics, problem->ics.qfunction_context));
38464f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ics, 0));
38564f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "x", num_comp_x, CEED_EVAL_INTERP));
38664f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
38764f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ics, "q0", num_comp_q, CEED_EVAL_NONE));
38864f98e98SJames Wright 
38964f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics));
390*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
391*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
392356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
393a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label));
3945263e9c6SJames Wright     PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx));
39577841947SLeila Ghaffari 
396*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_xc));
39764f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ics));
39864f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics));
39964f98e98SJames Wright   }
40064f98e98SJames Wright 
40164f98e98SJames Wright   if (problem->apply_vol_rhs.qfunction) {
40264f98e98SJames Wright     CeedQFunction qf_rhs_vol;
40377841947SLeila Ghaffari     CeedOperator  op;
40464f98e98SJames Wright 
40564f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &qf_rhs_vol));
40664f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_vol, problem->apply_vol_rhs.qfunction_context));
40764f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_rhs_vol, 0));
40864f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP));
40964f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
41064f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
41164f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP));
41264f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP));
41364f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
41464f98e98SJames Wright 
41564f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_vol, NULL, NULL, &op));
416a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
417a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
418356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
419a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
420a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
421a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
42277841947SLeila Ghaffari     user->op_rhs_vol = op;
42364f98e98SJames Wright 
42464f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_vol));
42577841947SLeila Ghaffari   }
42677841947SLeila Ghaffari 
42764f98e98SJames Wright   if (problem->apply_vol_ifunction.qfunction) {
42864f98e98SJames Wright     CeedQFunction qf_ifunction_vol;
42977841947SLeila Ghaffari     CeedOperator  op;
43064f98e98SJames Wright 
43164f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc,
43264f98e98SJames Wright                                                     &qf_ifunction_vol));
43364f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context));
43464f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ifunction_vol, 0));
43564f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP));
43664f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
43764f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP));
43864f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
43964f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP));
44064f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP));
44164f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
44264f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
44364f98e98SJames Wright 
44464f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ifunction_vol, NULL, NULL, &op));
445a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
446a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
447a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed));
448356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
449a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
450a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
451a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
452356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
453a3ae0734SJed Brown 
45477841947SLeila Ghaffari     user->op_ifunction_vol = op;
45564f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction_vol));
45677841947SLeila Ghaffari   }
45777841947SLeila Ghaffari 
458e334ad8fSJed Brown   CeedOperator op_ijacobian_vol = NULL;
45964f98e98SJames Wright   if (problem->apply_vol_ijacobian.qfunction) {
46064f98e98SJames Wright     CeedQFunction qf_ijacobian_vol;
461e334ad8fSJed Brown     CeedOperator  op;
46264f98e98SJames Wright 
46364f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc,
46464f98e98SJames Wright                                                     &qf_ijacobian_vol));
46564f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context));
46664f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ijacobian_vol, 0));
46764f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP));
46864f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD));
46964f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
47064f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
47164f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP));
47264f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
47364f98e98SJames Wright 
474a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op));
475a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
476a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
477356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
478356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
479a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
480a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
481e334ad8fSJed Brown     op_ijacobian_vol = op;
48264f98e98SJames Wright 
483a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol));
484e334ad8fSJed Brown   }
485e334ad8fSJed Brown 
48677841947SLeila Ghaffari   // -- Create and apply CEED Composite Operator for the entire domain
48777841947SLeila Ghaffari   if (!user->phys->implicit) {  // RHS
4889ad5e8e4SJames Wright     CeedOperator op_rhs;
489*e9c36be0SJames Wright 
490*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_rhs));
491*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_rhs, user->op_rhs_vol));
492*e9c36be0SJames Wright     PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, op_rhs, NULL));
493*e9c36be0SJames Wright 
4949ad5e8e4SJames Wright     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx));
495*e9c36be0SJames Wright 
496*e9c36be0SJames Wright     // ----- Get Context Labels for Operator
497*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "solution time", &user->phys->solution_time_label));
498*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "timestep size", &user->phys->timestep_size_label));
499*e9c36be0SJames Wright 
500a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
501b4e9a8f8SJames Wright     PetscCall(CreateKSPMass(user, problem));
502577c3361SJames Wright     PetscCheck(app_ctx->sgs_model_type == SGS_MODEL_NONE, user->comm, PETSC_ERR_SUP, "SGS modeling not implemented for explicit timestepping");
50377841947SLeila Ghaffari   } else {  // IFunction
50491c97f41SJames Wright     CeedOperator op_ijacobian = NULL;
50591c97f41SJames Wright 
506*e9c36be0SJames Wright     // Create Composite Operaters
507*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &user->op_ifunction));
508*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(user->op_ifunction, user->op_ifunction_vol));
509*e9c36be0SJames Wright     if (op_ijacobian_vol) {
510*e9c36be0SJames Wright       PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_ijacobian));
511*e9c36be0SJames Wright       PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_ijacobian, op_ijacobian_vol));
512*e9c36be0SJames Wright     }
513*e9c36be0SJames Wright     PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, user->op_ifunction, op_ijacobian));
514*e9c36be0SJames Wright 
515*e9c36be0SJames Wright     // ----- Get Context Labels for Operator
516*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "solution time", &user->phys->solution_time_label));
517*e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "timestep size", &user->phys->timestep_size_label));
518*e9c36be0SJames Wright 
51991c97f41SJames Wright     if (op_ijacobian) {
52091c97f41SJames Wright       PetscCall(MatCeedCreate(user->dm, user->dm, op_ijacobian, NULL, &user->mat_ijacobian));
52191c97f41SJames Wright       PetscCall(MatCeedSetLocalVectors(user->mat_ijacobian, user->Q_dot_loc, NULL));
52291c97f41SJames Wright       PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ijacobian, "ijacobian time shift", &user->phys->ijacobian_time_shift_label));
52391c97f41SJames Wright       PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian));
524e334ad8fSJed Brown     }
5258f5ab23bSJames Wright     if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SgsDDSetup(ceed, user, ceed_data, problem));
526dada6cc0SJames Wright   }
527f5452247SJames Wright 
528577c3361SJames Wright   if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc));
529f5452247SJames Wright   if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem));
530ee3b213aSJames Wright   if (app_ctx->diff_filter_monitor && !user->diff_filter) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem));
53128134cfaSJames Wright   if (app_ctx->sgs_train_enable) PetscCall(SGS_DD_TrainingSetup(ceed, user, ceed_data, problem));
532f5452247SJames Wright 
533a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i));
534a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol));
535a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data));
536ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
53777841947SLeila Ghaffari }
538