xref: /libCEED/examples/fluids/src/setuplibceed.c (revision 64f98e986d64302820caec520c05deb68bf27a53)
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 
93431cd09aSJames Wright PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height, CeedInt Q_sur,
942b730f8bSJeremy L Thompson                                 CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian,
959eb9c072SJames Wright                                 CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) {
96bb85d312SJames Wright   CeedVector          q_data_sur, jac_data_sur = NULL;
972b730f8bSJeremy L Thompson   CeedOperator        op_setup_sur, op_apply_bc, op_apply_bc_jacobian = NULL;
98bb85d312SJames Wright   CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur = NULL;
99bb85d312SJames Wright   CeedInt             num_qpts_sur, dm_field = 0;
1009eb9c072SJames Wright 
101f17d818dSJames Wright   PetscFunctionBeginUser;
1029eb9c072SJames Wright   // --- Get number of quadrature points for the boundaries
103a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(ceed_data->basis_q_sur, &num_qpts_sur));
1049eb9c072SJames Wright 
1059eb9c072SJames Wright   // ---- CEED Restriction
106bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &elem_restr_q_sur));
107bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &elem_restr_x_sur));
108bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, q_data_size_sur, &elem_restr_qd_i_sur));
1099eb9c072SJames Wright   if (jac_data_size_sur > 0) {
1109eb9c072SJames Wright     // State-dependent data will be passed from residual to Jacobian. This will be collocated.
111bb85d312SJames Wright     PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_sur, &elem_restr_jd_i_sur));
112a424bcd0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL));
1139eb9c072SJames Wright   }
1149eb9c072SJames Wright 
1159eb9c072SJames Wright   // ---- CEED Vector
116457e73b2SJames Wright   CeedInt loc_num_elem_sur;
117a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_q_sur, &loc_num_elem_sur));
118a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorCreate(ceed, q_data_size_sur * loc_num_elem_sur * num_qpts_sur, &q_data_sur));
1199eb9c072SJames Wright 
1209eb9c072SJames Wright   // ---- CEED Operator
1219eb9c072SJames Wright   // ----- CEED Operator for Setup (geometric factors)
122a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur));
123a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, ceed_data->basis_x_sur, CEED_VECTOR_ACTIVE));
124a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x_sur, CEED_VECTOR_NONE));
125356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
1269eb9c072SJames Wright 
1279eb9c072SJames Wright   // ----- CEED Operator for Physics
128a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc));
129a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
130a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
131356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
132a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord));
133a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
134a424bcd0SJames Wright   if (elem_restr_jd_i_sur)
135356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
1369eb9c072SJames Wright 
137f21e6b1cSJames Wright   if (qf_apply_bc_jacobian && elem_restr_jd_i_sur) {
138a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian));
139a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
140a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
141356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
142a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord));
143356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
144a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
1459eb9c072SJames Wright   }
1469eb9c072SJames Wright 
1479eb9c072SJames Wright   // ----- Apply CEED operator for Setup
148a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE));
1499eb9c072SJames Wright 
1509eb9c072SJames Wright   // ----- Apply Sub-Operator for Physics
151a424bcd0SJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply, op_apply_bc));
152a424bcd0SJames Wright   if (op_apply_bc_jacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_bc_jacobian));
1539eb9c072SJames Wright 
1549eb9c072SJames Wright   // ----- Cleanup
155a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_sur));
156a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data_sur));
157a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur));
158a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur));
159a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i_sur));
160a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i_sur));
161a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur));
162a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc));
163a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc_jacobian));
164ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1659eb9c072SJames Wright }
1669eb9c072SJames Wright 
16777841947SLeila Ghaffari // Utility function to create CEED Composite Operator for the entire domain
1682b730f8bSJeremy L Thompson PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol,
1692b730f8bSJeremy L Thompson                                        CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur,
1702b730f8bSJeremy L Thompson                                        CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) {
17177841947SLeila Ghaffari   DMLabel domain_label;
17277841947SLeila Ghaffari 
173f17d818dSJames Wright   PetscFunctionBeginUser;
17477841947SLeila Ghaffari   // Create Composite Operaters
175a424bcd0SJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_apply));
176a424bcd0SJames Wright   if (op_apply_ijacobian) PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_apply_ijacobian));
17777841947SLeila Ghaffari 
17877841947SLeila Ghaffari   // --Apply Sub-Operator for the volume
179a424bcd0SJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply, op_apply_vol));
180a424bcd0SJames Wright   if (op_apply_ijacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_ijacobian_vol));
18177841947SLeila Ghaffari 
18277841947SLeila Ghaffari   // -- Create Sub-Operator for in/outflow BCs
1832b730f8bSJeremy L Thompson   PetscCall(DMGetLabel(dm, "Face Sets", &domain_label));
18477841947SLeila Ghaffari 
1852fe7aee7SLeila Ghaffari   // --- Create Sub-Operator for inflow boundaries
1862fe7aee7SLeila Ghaffari   for (CeedInt i = 0; i < bc->num_inflow; i++) {
1872b730f8bSJeremy L Thompson     PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->inflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
1882b730f8bSJeremy L Thompson                                ceed_data->qf_apply_inflow, ceed_data->qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian));
18977841947SLeila Ghaffari   }
1902fe7aee7SLeila Ghaffari   // --- Create Sub-Operator for outflow boundaries
1912fe7aee7SLeila Ghaffari   for (CeedInt i = 0; i < bc->num_outflow; i++) {
1922b730f8bSJeremy L Thompson     PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
1932b730f8bSJeremy L Thompson                                ceed_data->qf_apply_outflow, ceed_data->qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian));
19439c69132SJed Brown   }
195f4ca79c2SJames Wright   // --- Create Sub-Operator for freestream boundaries
196f4ca79c2SJames Wright   for (CeedInt i = 0; i < bc->num_freestream; i++) {
1972b730f8bSJeremy L Thompson     PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
1982b730f8bSJeremy L Thompson                                ceed_data->qf_apply_freestream, ceed_data->qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian));
1992fe7aee7SLeila Ghaffari   }
2009f844368SJames Wright   // --- Create Sub-Operator for slip boundaries
2019f844368SJames Wright   for (CeedInt i = 0; i < bc->num_slip; i++) {
2029f844368SJames Wright     PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->slips[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
2039f844368SJames Wright                                ceed_data->qf_apply_slip, ceed_data->qf_apply_slip_jacobian, op_apply, op_apply_ijacobian));
2049f844368SJames Wright   }
20511436a05SJames Wright 
20611436a05SJames Wright   // ----- Get Context Labels for Operator
207a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(*op_apply, "solution time", &phys->solution_time_label));
208a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(*op_apply, "timestep size", &phys->timestep_size_label));
209ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
21077841947SLeila Ghaffari }
21177841947SLeila Ghaffari 
2122b730f8bSJeremy L Thompson PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur,
2132b730f8bSJeremy L Thompson                                  PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian,
2145b881d11SJames Wright                                  CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) {
2155b881d11SJames Wright   PetscFunctionBeginUser;
2165b881d11SJames Wright   if (apply_bc.qfunction) {
217a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc));
218a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context));
21991c97f41SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc, 0));
220a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP));
221a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD));
222a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
223a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP));
224a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP));
225a424bcd0SJames Wright     if (jac_data_size_sur) PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
2265b881d11SJames Wright   }
2275b881d11SJames Wright   if (apply_bc_jacobian.qfunction) {
228a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian));
229a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context));
23091c97f41SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc_jacobian, 0));
231a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP));
232a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD));
233a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
234a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP));
235a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
236a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP));
2375b881d11SJames Wright   }
238ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2395b881d11SJames Wright }
2405b881d11SJames Wright 
241731c13d7SJames Wright PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData problem, SimpleBC bc) {
24277841947SLeila Ghaffari   PetscFunctionBeginUser;
24377841947SLeila Ghaffari   const PetscInt num_comp_q = 5;
2441d2a9659SKenneth E. Jansen   const CeedInt  dim = problem->dim, num_comp_x = problem->dim, q_data_size_vol = problem->q_data_size_vol;
2451d2a9659SKenneth E. Jansen   CeedInt        jac_data_size_vol = num_comp_q + 6 + 3;
2461d2a9659SKenneth E. Jansen 
247752a08a7SJames Wright   if (problem->apply_vol_ifunction.qfunction && problem->uses_newtonian) {
2481d2a9659SKenneth E. Jansen     NewtonianIdealGasContext gas;
2491d2a9659SKenneth E. Jansen     PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas));
2501d2a9659SKenneth E. Jansen     jac_data_size_vol += (gas->idl_enable ? 1 : 0);
2511d2a9659SKenneth E. Jansen     PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas));
2521d2a9659SKenneth E. Jansen   }
2531d2a9659SKenneth E. Jansen 
254a3ae0734SJed Brown   CeedElemRestriction elem_restr_jd_i;
255a3ae0734SJed Brown   CeedVector          jac_data;
2560814d5a7SKenneth E. Jansen   CeedInt             num_qpts;
257bb85d312SJames Wright   DMLabel             domain_label = NULL;
258bb85d312SJames Wright   PetscInt            label_value = 0, height = 0, dm_field = 0;
25977841947SLeila Ghaffari 
2600814d5a7SKenneth E. Jansen   DM dm_coord;
2610814d5a7SKenneth E. Jansen   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
2620814d5a7SKenneth E. Jansen 
263bb85d312SJames Wright   PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &ceed_data->basis_q));
264bb85d312SJames Wright   PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &ceed_data->basis_x));
265a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &ceed_data->basis_xc));
266a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(ceed_data->basis_q, &num_qpts));
26777841947SLeila Ghaffari 
268bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, 0, &ceed_data->elem_restr_q));
269bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &ceed_data->elem_restr_x));
270bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, q_data_size_vol, &ceed_data->elem_restr_qd_i));
271bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_vol, &elem_restr_jd_i));
272*64f98e98SJames Wright 
273a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL));
274a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL));
275a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL));
276a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL));
277*64f98e98SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_qd_i, &ceed_data->q_data, NULL));
278*64f98e98SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL));
27977841947SLeila Ghaffari 
280*64f98e98SJames Wright   {  // -- Copy PETSc coordinate vector into CEED vector
28177841947SLeila Ghaffari     Vec X_loc;
2823796c488SJed Brown     DM  cdm;
283*64f98e98SJames Wright 
2842b730f8bSJeremy L Thompson     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
2852b730f8bSJeremy L Thompson     if (cdm) {
2862b730f8bSJeremy L Thompson       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
2872b730f8bSJeremy L Thompson     } else {
2882b730f8bSJeremy L Thompson       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
2893796c488SJed Brown     }
2902b730f8bSJeremy L Thompson     PetscCall(VecScale(X_loc, problem->dm_scale));
291d0593705SJames Wright     PetscCall(VecCopyPetscToCeed(X_loc, ceed_data->x_coord));
292*64f98e98SJames Wright   }
29377841947SLeila Ghaffari 
294*64f98e98SJames Wright   {  // -- Create quadrature data
295*64f98e98SJames Wright     CeedQFunction qf_setup_vol;
296*64f98e98SJames Wright     CeedOperator  op_setup_vol;
297d52d2babSJames Wright 
298*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, problem->setup_vol.qfunction_loc, &qf_setup_vol));
299*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_setup_vol, problem->setup_vol.qfunction_context));
300*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_setup_vol, 0));
301*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_vol, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
302*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT));
303*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
30477841947SLeila Ghaffari 
305*64f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup_vol, NULL, NULL, &op_setup_vol));
306*64f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "dx", ceed_data->elem_restr_x, ceed_data->basis_x, CEED_VECTOR_ACTIVE));
307*64f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE));
308*64f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
309*64f98e98SJames Wright 
310*64f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorApply(op_setup_vol, ceed_data->x_coord, ceed_data->q_data, CEED_REQUEST_IMMEDIATE));
311*64f98e98SJames Wright 
312*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup_vol));
313*64f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_vol));
314*64f98e98SJames Wright   }
315*64f98e98SJames Wright 
316*64f98e98SJames Wright   {  // -- Create QFunction for ICs
317*64f98e98SJames Wright     CeedQFunction qf_ics;
3185263e9c6SJames Wright     CeedOperator  op_ics;
319*64f98e98SJames Wright 
320*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &qf_ics));
321*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ics, problem->ics.qfunction_context));
322*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ics, 0));
323*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "x", num_comp_x, CEED_EVAL_INTERP));
324*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
325*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ics, "q0", num_comp_q, CEED_EVAL_NONE));
326*64f98e98SJames Wright 
327*64f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics));
328a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE));
329a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE));
330356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
331a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label));
3325263e9c6SJames Wright     PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx));
33377841947SLeila Ghaffari 
334*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ics));
335*64f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics));
336*64f98e98SJames Wright   }
337*64f98e98SJames Wright 
338*64f98e98SJames Wright   if (problem->apply_vol_rhs.qfunction) {
339*64f98e98SJames Wright     CeedQFunction qf_rhs_vol;
34077841947SLeila Ghaffari     CeedOperator  op;
341*64f98e98SJames Wright 
342*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &qf_rhs_vol));
343*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_vol, problem->apply_vol_rhs.qfunction_context));
344*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_rhs_vol, 0));
345*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP));
346*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
347*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
348*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP));
349*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP));
350*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
351*64f98e98SJames Wright 
352*64f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_vol, NULL, NULL, &op));
353a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
354a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
355356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
356a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
357a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
358a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
35977841947SLeila Ghaffari     user->op_rhs_vol = op;
360*64f98e98SJames Wright 
361*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_vol));
36277841947SLeila Ghaffari   }
36377841947SLeila Ghaffari 
364*64f98e98SJames Wright   if (problem->apply_vol_ifunction.qfunction) {
365*64f98e98SJames Wright     CeedQFunction qf_ifunction_vol;
36677841947SLeila Ghaffari     CeedOperator  op;
367*64f98e98SJames Wright 
368*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc,
369*64f98e98SJames Wright                                                     &qf_ifunction_vol));
370*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context));
371*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ifunction_vol, 0));
372*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP));
373*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
374*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP));
375*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
376*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP));
377*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP));
378*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
379*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
380*64f98e98SJames Wright 
381*64f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ifunction_vol, NULL, NULL, &op));
382a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
383a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
384a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed));
385356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
386a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
387a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
388a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
389356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
390a3ae0734SJed Brown 
39177841947SLeila Ghaffari     user->op_ifunction_vol = op;
392*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction_vol));
39377841947SLeila Ghaffari   }
39477841947SLeila Ghaffari 
395e334ad8fSJed Brown   CeedOperator op_ijacobian_vol = NULL;
396*64f98e98SJames Wright   if (problem->apply_vol_ijacobian.qfunction) {
397*64f98e98SJames Wright     CeedQFunction qf_ijacobian_vol;
398e334ad8fSJed Brown     CeedOperator  op;
399*64f98e98SJames Wright 
400*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc,
401*64f98e98SJames Wright                                                     &qf_ijacobian_vol));
402*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context));
403*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ijacobian_vol, 0));
404*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP));
405*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD));
406*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
407*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
408*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP));
409*64f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
410*64f98e98SJames Wright 
411a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op));
412a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
413a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
414356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
415356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
416a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
417a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
418e334ad8fSJed Brown     op_ijacobian_vol = op;
419*64f98e98SJames Wright 
420a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol));
421e334ad8fSJed Brown   }
422e334ad8fSJed Brown 
42377841947SLeila Ghaffari   // *****************************************************************************
42477841947SLeila Ghaffari   // Set up CEED objects for the exterior domain (surface)
42577841947SLeila Ghaffari   // *****************************************************************************
426bb85d312SJames Wright   height                = 1;
427bb85d312SJames Wright   CeedInt       dim_sur = dim - height, P_sur = app_ctx->degree + 1, Q_sur = P_sur + app_ctx->q_extra;
428f21e6b1cSJames 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;
42977841947SLeila Ghaffari 
43077841947SLeila Ghaffari   // -----------------------------------------------------------------------------
43177841947SLeila Ghaffari   // CEED Bases
43277841947SLeila Ghaffari   // -----------------------------------------------------------------------------
4330814d5a7SKenneth E. Jansen 
434*64f98e98SJames Wright   {
435*64f98e98SJames Wright     DM dm_coord;
436*64f98e98SJames Wright 
437*64f98e98SJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
438*64f98e98SJames Wright     DMLabel  label   = NULL;
4390814d5a7SKenneth E. Jansen     PetscInt face_id = 0;
4400814d5a7SKenneth E. Jansen     PetscInt field   = 0;  // Still want the normal, default field
4410814d5a7SKenneth E. Jansen     PetscCall(CreateBasisFromPlex(ceed, dm, label, face_id, height, field, &ceed_data->basis_q_sur));
4420814d5a7SKenneth E. Jansen     PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, face_id, height, field, &ceed_data->basis_x_sur));
443*64f98e98SJames Wright   }
44477841947SLeila Ghaffari 
44577841947SLeila Ghaffari   // -----------------------------------------------------------------------------
44677841947SLeila Ghaffari   // CEED QFunctions
44777841947SLeila Ghaffari   // -----------------------------------------------------------------------------
44877841947SLeila Ghaffari   // -- Create QFunction for quadrature data
449a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, problem->setup_sur.qfunction_loc, &ceed_data->qf_setup_sur));
450841e4c73SJed Brown   if (problem->setup_sur.qfunction_context) {
451a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_setup_sur, problem->setup_sur.qfunction_context));
452841e4c73SJed Brown   }
45391c97f41SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(ceed_data->qf_setup_sur, 0));
454a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x * dim_sur, CEED_EVAL_GRAD));
455a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT));
456a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
45777841947SLeila Ghaffari 
4582b730f8bSJeremy L Thompson   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow,
4592b730f8bSJeremy L Thompson                               problem->apply_inflow_jacobian, &ceed_data->qf_apply_inflow, &ceed_data->qf_apply_inflow_jacobian));
4602b730f8bSJeremy L Thompson   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow,
4612b730f8bSJeremy L Thompson                               problem->apply_outflow_jacobian, &ceed_data->qf_apply_outflow, &ceed_data->qf_apply_outflow_jacobian));
4622b730f8bSJeremy L Thompson   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream,
4632b730f8bSJeremy L Thompson                               problem->apply_freestream_jacobian, &ceed_data->qf_apply_freestream, &ceed_data->qf_apply_freestream_jacobian));
4649f844368SJames Wright   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_slip,
4659f844368SJames Wright                               problem->apply_slip_jacobian, &ceed_data->qf_apply_slip, &ceed_data->qf_apply_slip_jacobian));
46677841947SLeila Ghaffari 
46777841947SLeila Ghaffari   // *****************************************************************************
46877841947SLeila Ghaffari   // CEED Operator Apply
46977841947SLeila Ghaffari   // *****************************************************************************
47077841947SLeila Ghaffari   // -- Create and apply CEED Composite Operator for the entire domain
47177841947SLeila Ghaffari   if (!user->phys->implicit) {  // RHS
4729ad5e8e4SJames Wright     CeedOperator op_rhs;
4739ad5e8e4SJames Wright     PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_rhs_vol, NULL, height, P_sur, Q_sur, q_data_size_sur, 0, &op_rhs,
4749ad5e8e4SJames Wright                                       NULL));
4759ad5e8e4SJames Wright     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx));
476a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
477b4e9a8f8SJames Wright     PetscCall(CreateKSPMass(user, problem));
478577c3361SJames Wright     PetscCheck(app_ctx->sgs_model_type == SGS_MODEL_NONE, user->comm, PETSC_ERR_SUP, "SGS modeling not implemented for explicit timestepping");
47977841947SLeila Ghaffari   } else {  // IFunction
48091c97f41SJames Wright     CeedOperator op_ijacobian = NULL;
48191c97f41SJames Wright 
4822b730f8bSJeremy L Thompson     PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_ifunction_vol, op_ijacobian_vol, height, P_sur, Q_sur,
48391c97f41SJames Wright                                       q_data_size_sur, jac_data_size_sur, &user->op_ifunction, op_ijacobian_vol ? &op_ijacobian : NULL));
48491c97f41SJames Wright     if (op_ijacobian) {
48591c97f41SJames Wright       PetscCall(MatCeedCreate(user->dm, user->dm, op_ijacobian, NULL, &user->mat_ijacobian));
48691c97f41SJames Wright       PetscCall(MatCeedSetLocalVectors(user->mat_ijacobian, user->Q_dot_loc, NULL));
48791c97f41SJames Wright       PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ijacobian, "ijacobian time shift", &user->phys->ijacobian_time_shift_label));
48891c97f41SJames Wright       PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian));
489e334ad8fSJed Brown     }
4908f5ab23bSJames Wright     if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SgsDDSetup(ceed, user, ceed_data, problem));
491dada6cc0SJames Wright   }
492f5452247SJames Wright 
493577c3361SJames Wright   if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc));
494f5452247SJames Wright   if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem));
495ee3b213aSJames Wright   if (app_ctx->diff_filter_monitor && !user->diff_filter) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem));
49628134cfaSJames Wright   if (app_ctx->sgs_train_enable) PetscCall(SGS_DD_TrainingSetup(ceed, user, ceed_data, problem));
497f5452247SJames Wright 
498a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i));
499a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol));
500a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data));
501ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
50277841947SLeila Ghaffari }
503