xref: /libCEED/examples/fluids/src/setuplibceed.c (revision a5b0ec6f41c4198e9c6bd74b9679c3c63ae70b36)
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 
93e9c36be0SJames Wright static PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height,
94e9c36be0SJames Wright                                        CeedInt Q_sur, CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedBasis basis_q_sur,
95e9c36be0SJames Wright                                        CeedBasis basis_x_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian, CeedOperator op_apply,
96e9c36be0SJames Wright                                        CeedOperator op_apply_ijacobian) {
97bb85d312SJames Wright   CeedVector          q_data_sur, jac_data_sur          = NULL;
98e9c36be0SJames 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;
100e9c36be0SJames Wright   PetscInt            dm_field = 0;
1019eb9c072SJames Wright 
102f17d818dSJames Wright   PetscFunctionBeginUser;
103bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, dm_field, &elem_restr_q_sur));
104bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &elem_restr_x_sur));
105bb85d312SJames Wright   PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, q_data_size_sur, &elem_restr_qd_i_sur));
1069eb9c072SJames Wright   if (jac_data_size_sur > 0) {
1079eb9c072SJames Wright     // State-dependent data will be passed from residual to Jacobian. This will be collocated.
108bb85d312SJames Wright     PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_sur, &elem_restr_jd_i_sur));
109a424bcd0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL));
1109eb9c072SJames Wright   }
1119eb9c072SJames Wright 
112e9c36be0SJames Wright   {  // Create q_data_sur vector
113e9c36be0SJames Wright     CeedOperator op_setup_sur;
1149eb9c072SJames Wright 
115e9c36be0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_qd_i_sur, &q_data_sur, NULL));
116e9c36be0SJames Wright 
117a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur));
118e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, basis_x_sur, CEED_VECTOR_ACTIVE));
119e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, basis_x_sur, CEED_VECTOR_NONE));
120356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
1219eb9c072SJames Wright 
122e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE));
123e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur));
124e9c36be0SJames Wright   }
125e9c36be0SJames Wright 
126*a5b0ec6fSJames Wright   // CEED Operator for Physics
127a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc));
128e9c36be0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
129e9c36be0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
130356036faSJeremy L Thompson   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
131e9c36be0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord));
132e9c36be0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
133a424bcd0SJames Wright   if (elem_restr_jd_i_sur)
134356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
1359eb9c072SJames Wright 
136f21e6b1cSJames Wright   if (qf_apply_bc_jacobian && elem_restr_jd_i_sur) {
137a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian));
138e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
139e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
140356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_NONE, q_data_sur));
141e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, basis_x_sur, ceed_data->x_coord));
142356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_NONE, jac_data_sur));
143e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, basis_q_sur, CEED_VECTOR_ACTIVE));
1449eb9c072SJames Wright   }
1459eb9c072SJames Wright 
146*a5b0ec6fSJames Wright   // Apply Sub-Operator for Physics
147e9c36be0SJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply, op_apply_bc));
148e9c36be0SJames Wright   if (op_apply_bc_jacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_apply_ijacobian, op_apply_bc_jacobian));
1499eb9c072SJames Wright 
150a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_sur));
151a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data_sur));
152a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur));
153a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur));
154a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i_sur));
155a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i_sur));
156a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc));
157a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc_jacobian));
158ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1599eb9c072SJames Wright }
1609eb9c072SJames Wright 
161e9c36be0SJames Wright static PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur,
1622b730f8bSJeremy L Thompson                                         PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian,
1635b881d11SJames Wright                                         CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) {
1645b881d11SJames Wright   PetscFunctionBeginUser;
1655b881d11SJames Wright   if (apply_bc.qfunction) {
166a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc));
167a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context));
16891c97f41SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc, 0));
169a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP));
170a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD));
171a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
172a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP));
173a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP));
174a424bcd0SJames Wright     if (jac_data_size_sur) PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
1755b881d11SJames Wright   }
1765b881d11SJames Wright   if (apply_bc_jacobian.qfunction) {
177a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian));
178a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context));
17991c97f41SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(*qf_apply_bc_jacobian, 0));
180a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP));
181a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD));
182a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
183a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP));
184a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
185a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP));
1865b881d11SJames Wright   }
187ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1885b881d11SJames Wright }
1895b881d11SJames Wright 
190*a5b0ec6fSJames Wright // Utility function to add boundary operators to the composite operator
191e9c36be0SJames Wright static PetscErrorCode AddBCSubOperators(User user, Ceed ceed, DM dm, SimpleBC bc, ProblemData problem, CeedData ceed_data, CeedOperator op_apply,
192e9c36be0SJames Wright                                         CeedOperator op_apply_ijacobian) {
193e9c36be0SJames Wright   CeedInt       height = 1, num_comp_q, num_comp_x;
194e9c36be0SJames Wright   CeedInt       dim_sur, P_sur = user->app_ctx->degree + 1, Q_sur = P_sur + user->app_ctx->q_extra;
195e9c36be0SJames 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;
196e9c36be0SJames Wright   PetscInt      dim;
197e9c36be0SJames Wright   DMLabel       face_sets_label;
198e9c36be0SJames Wright   CeedBasis     basis_q_sur, basis_x_sur;
199e9c36be0SJames Wright 
200e9c36be0SJames Wright   PetscFunctionBeginUser;
201e9c36be0SJames Wright   PetscCall(DMGetDimension(dm, &dim));
202e9c36be0SJames Wright   dim_sur = dim - height;
203e9c36be0SJames Wright   {  // Get number of components and coordinate dimension from op_apply
204e9c36be0SJames Wright     CeedOperator       *sub_ops;
205e9c36be0SJames Wright     CeedOperatorField   field;
206e9c36be0SJames Wright     PetscInt            sub_op_index = 0;  // will be 0 for the volume op
207e9c36be0SJames Wright     CeedElemRestriction elem_restr_q, elem_restr_x;
208e9c36be0SJames Wright 
209e9c36be0SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(op_apply, &sub_ops));
210e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &field));
211e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_q));
212e9c36be0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
213e9c36be0SJames Wright 
214e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "x", &field));
215e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetElemRestriction(field, &elem_restr_x));
216e9c36be0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_x, &num_comp_x));
217e9c36be0SJames Wright   }
218e9c36be0SJames Wright 
219e9c36be0SJames Wright   {  // Get bases
220e9c36be0SJames Wright     DM dm_coord;
221e9c36be0SJames Wright 
222e9c36be0SJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
223e9c36be0SJames Wright     DMLabel  label       = NULL;
224e9c36be0SJames Wright     PetscInt label_value = 0;
225*a5b0ec6fSJames Wright     PetscInt field       = 0;
226e9c36be0SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm, label, label_value, height, field, &basis_q_sur));
227e9c36be0SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, label_value, height, field, &basis_x_sur));
228e9c36be0SJames Wright   }
229e9c36be0SJames Wright 
230e9c36be0SJames Wright   PetscCall(DMGetLabel(dm, "Face Sets", &face_sets_label));
231e9c36be0SJames Wright 
232e9c36be0SJames Wright   {  // -- Create QFunction for quadrature data
233e9c36be0SJames Wright     PetscCallCeed(ceed,
234e9c36be0SJames Wright                   CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, problem->setup_sur.qfunction_loc, &ceed_data->qf_setup_sur));
235e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_setup_sur, problem->setup_sur.qfunction_context));
236e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(ceed_data->qf_setup_sur, 0));
237e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x * dim_sur, CEED_EVAL_GRAD));
238e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT));
239e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
240e9c36be0SJames Wright   }
241e9c36be0SJames Wright 
242e9c36be0SJames Wright   {  // --- Create Sub-Operator for inflow boundaries
243e9c36be0SJames Wright     CeedQFunction qf_apply_inflow = NULL, qf_apply_inflow_jacobian = NULL;
244e9c36be0SJames Wright 
245e9c36be0SJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow,
246e9c36be0SJames Wright                                 problem->apply_inflow_jacobian, &qf_apply_inflow, &qf_apply_inflow_jacobian));
247e9c36be0SJames Wright     for (CeedInt i = 0; i < bc->num_inflow; i++) {
248e9c36be0SJames 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,
249e9c36be0SJames Wright                                  basis_x_sur, qf_apply_inflow, qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian));
250e9c36be0SJames Wright     }
251e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow));
252e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_inflow_jacobian));
253e9c36be0SJames Wright   }
254e9c36be0SJames Wright 
255e9c36be0SJames Wright   {  // --- Create Sub-Operator for outflow boundaries
256e9c36be0SJames Wright     CeedQFunction qf_apply_outflow = NULL, qf_apply_outflow_jacobian = NULL;
257e9c36be0SJames Wright 
258e9c36be0SJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow,
259e9c36be0SJames Wright                                 problem->apply_outflow_jacobian, &qf_apply_outflow, &qf_apply_outflow_jacobian));
260e9c36be0SJames Wright     for (CeedInt i = 0; i < bc->num_outflow; i++) {
261e9c36be0SJames Wright       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
262e9c36be0SJames Wright                                  basis_q_sur, basis_x_sur, qf_apply_outflow, qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian));
263e9c36be0SJames Wright     }
264e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow));
265e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_outflow_jacobian));
266e9c36be0SJames Wright   }
267e9c36be0SJames Wright 
268e9c36be0SJames Wright   {  // --- Create Sub-Operator for freestream boundaries
269e9c36be0SJames Wright     CeedQFunction qf_apply_freestream = NULL, qf_apply_freestream_jacobian = NULL;
270e9c36be0SJames Wright 
271e9c36be0SJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream,
272e9c36be0SJames Wright                                 problem->apply_freestream_jacobian, &qf_apply_freestream, &qf_apply_freestream_jacobian));
273e9c36be0SJames Wright     for (CeedInt i = 0; i < bc->num_freestream; i++) {
274e9c36be0SJames Wright       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, face_sets_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
275e9c36be0SJames Wright                                  basis_q_sur, basis_x_sur, qf_apply_freestream, qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian));
276e9c36be0SJames Wright     }
277e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream));
278e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_freestream_jacobian));
279e9c36be0SJames Wright   }
280e9c36be0SJames Wright 
281e9c36be0SJames Wright   {  // --- Create Sub-Operator for slip boundaries
282e9c36be0SJames Wright     CeedQFunction qf_apply_slip = NULL, qf_apply_slip_jacobian = NULL;
283e9c36be0SJames Wright 
284e9c36be0SJames Wright     PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_slip,
285e9c36be0SJames Wright                                 problem->apply_slip_jacobian, &qf_apply_slip, &qf_apply_slip_jacobian));
286e9c36be0SJames Wright     for (CeedInt i = 0; i < bc->num_slip; i++) {
287e9c36be0SJames 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,
288e9c36be0SJames Wright                                  basis_x_sur, qf_apply_slip, qf_apply_slip_jacobian, op_apply, op_apply_ijacobian));
289e9c36be0SJames Wright     }
290e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip));
291e9c36be0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_apply_slip_jacobian));
292e9c36be0SJames Wright   }
293e9c36be0SJames Wright 
294e9c36be0SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q_sur));
295e9c36be0SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_x_sur));
296e9c36be0SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
297e9c36be0SJames Wright }
298e9c36be0SJames Wright 
299731c13d7SJames Wright PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData problem, SimpleBC bc) {
30077841947SLeila Ghaffari   const PetscInt      num_comp_q = 5;
3011d2a9659SKenneth E. Jansen   const CeedInt       dim = problem->dim, num_comp_x = problem->dim, q_data_size_vol = problem->q_data_size_vol;
3021d2a9659SKenneth E. Jansen   CeedInt             jac_data_size_vol = num_comp_q + 6 + 3;
303*a5b0ec6fSJames Wright   CeedElemRestriction elem_restr_jd_i;
304*a5b0ec6fSJames Wright   CeedVector          jac_data;
305*a5b0ec6fSJames Wright   CeedOperator        op_ifunction_vol = NULL, op_rhs_vol = NULL, op_ijacobian_vol = NULL;
306*a5b0ec6fSJames Wright 
307*a5b0ec6fSJames Wright   PetscFunctionBeginUser;
3081d2a9659SKenneth E. Jansen 
309752a08a7SJames Wright   if (problem->apply_vol_ifunction.qfunction && problem->uses_newtonian) {
3101d2a9659SKenneth E. Jansen     NewtonianIdealGasContext gas;
3111d2a9659SKenneth E. Jansen     PetscCallCeed(ceed, CeedQFunctionContextGetDataRead(problem->apply_vol_ifunction.qfunction_context, CEED_MEM_HOST, &gas));
3121d2a9659SKenneth E. Jansen     jac_data_size_vol += (gas->idl_enable ? 1 : 0);
3131d2a9659SKenneth E. Jansen     PetscCallCeed(ceed, CeedQFunctionContextRestoreDataRead(problem->apply_vol_ifunction.qfunction_context, &gas));
3141d2a9659SKenneth E. Jansen   }
3151d2a9659SKenneth E. Jansen 
316*a5b0ec6fSJames Wright   {  // Create bases and element restrictions
317bb85d312SJames Wright     DMLabel  domain_label = NULL;
318bb85d312SJames Wright     PetscInt label_value = 0, height = 0, dm_field = 0;
3190814d5a7SKenneth E. Jansen     DM       dm_coord;
3200814d5a7SKenneth E. Jansen 
321*a5b0ec6fSJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
322bb85d312SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &ceed_data->basis_q));
323bb85d312SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &ceed_data->basis_x));
32477841947SLeila Ghaffari 
325bb85d312SJames Wright     PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, 0, &ceed_data->elem_restr_q));
326bb85d312SJames Wright     PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &ceed_data->elem_restr_x));
327bb85d312SJames Wright     PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, q_data_size_vol, &ceed_data->elem_restr_qd_i));
328bb85d312SJames Wright     PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, jac_data_size_vol, &elem_restr_jd_i));
32964f98e98SJames Wright 
330a424bcd0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL));
331a424bcd0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL));
332a424bcd0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL));
333a424bcd0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL));
33464f98e98SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_qd_i, &ceed_data->q_data, NULL));
33564f98e98SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL));
336*a5b0ec6fSJames Wright   }
33777841947SLeila Ghaffari 
33864f98e98SJames Wright   {  // -- Copy PETSc coordinate vector into CEED vector
33977841947SLeila Ghaffari     Vec X_loc;
3403796c488SJed Brown     DM  cdm;
34164f98e98SJames Wright 
3422b730f8bSJeremy L Thompson     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
3432b730f8bSJeremy L Thompson     if (cdm) {
3442b730f8bSJeremy L Thompson       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
3452b730f8bSJeremy L Thompson     } else {
3462b730f8bSJeremy L Thompson       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
3473796c488SJed Brown     }
3482b730f8bSJeremy L Thompson     PetscCall(VecScale(X_loc, problem->dm_scale));
349d0593705SJames Wright     PetscCall(VecCopyPetscToCeed(X_loc, ceed_data->x_coord));
35064f98e98SJames Wright   }
35177841947SLeila Ghaffari 
35264f98e98SJames Wright   {  // -- Create quadrature data
35364f98e98SJames Wright     CeedQFunction qf_setup_vol;
35464f98e98SJames Wright     CeedOperator  op_setup_vol;
355d52d2babSJames Wright 
35664f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, problem->setup_vol.qfunction_loc, &qf_setup_vol));
35764f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_setup_vol, problem->setup_vol.qfunction_context));
35864f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_setup_vol, 0));
35964f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_vol, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
36064f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT));
36164f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_setup_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
36277841947SLeila Ghaffari 
36364f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_setup_vol, NULL, NULL, &op_setup_vol));
36464f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "dx", ceed_data->elem_restr_x, ceed_data->basis_x, CEED_VECTOR_ACTIVE));
36564f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE));
36664f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_setup_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
36764f98e98SJames Wright 
36864f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorApply(op_setup_vol, ceed_data->x_coord, ceed_data->q_data, CEED_REQUEST_IMMEDIATE));
36964f98e98SJames Wright 
37064f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_setup_vol));
37164f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_vol));
37264f98e98SJames Wright   }
37364f98e98SJames Wright 
37464f98e98SJames Wright   {  // -- Create QFunction for ICs
375e9c36be0SJames Wright     CeedBasis     basis_xc;
37664f98e98SJames Wright     CeedQFunction qf_ics;
3775263e9c6SJames Wright     CeedOperator  op_ics;
37864f98e98SJames Wright 
379e9c36be0SJames Wright     PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &basis_xc));
38064f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &qf_ics));
38164f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ics, problem->ics.qfunction_context));
38264f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ics, 0));
38364f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "x", num_comp_x, CEED_EVAL_INTERP));
38464f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
38564f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ics, "q0", num_comp_q, CEED_EVAL_NONE));
38664f98e98SJames Wright 
38764f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics));
388e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
389e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
390356036faSJeremy L Thompson     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
391a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label));
3925263e9c6SJames Wright     PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx));
39377841947SLeila Ghaffari 
394e9c36be0SJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_xc));
39564f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ics));
39664f98e98SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics));
39764f98e98SJames Wright   }
39864f98e98SJames Wright 
39964f98e98SJames Wright   if (problem->apply_vol_rhs.qfunction) {
40064f98e98SJames Wright     CeedQFunction qf_rhs_vol;
40164f98e98SJames Wright 
40264f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &qf_rhs_vol));
40364f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_vol, problem->apply_vol_rhs.qfunction_context));
40464f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_rhs_vol, 0));
40564f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP));
40664f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
40764f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
40864f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP));
40964f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP));
41064f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
41164f98e98SJames Wright 
412*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_vol, NULL, NULL, &op_rhs_vol));
413*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
414*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
415*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
416*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
417*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
418*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
41964f98e98SJames Wright 
42064f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_vol));
42177841947SLeila Ghaffari   }
42277841947SLeila Ghaffari 
42364f98e98SJames Wright   if (problem->apply_vol_ifunction.qfunction) {
42464f98e98SJames Wright     CeedQFunction qf_ifunction_vol;
42564f98e98SJames Wright 
42664f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc,
42764f98e98SJames Wright                                                     &qf_ifunction_vol));
42864f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context));
42964f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ifunction_vol, 0));
43064f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP));
43164f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
43264f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP));
43364f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
43464f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP));
43564f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP));
43664f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
43764f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
43864f98e98SJames Wright 
439*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ifunction_vol, NULL, NULL, &op_ifunction_vol));
440*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
441*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
442*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed));
443*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
444*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
445*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
446*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
447*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
448a3ae0734SJed Brown 
44964f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction_vol));
45077841947SLeila Ghaffari   }
45177841947SLeila Ghaffari 
45264f98e98SJames Wright   if (problem->apply_vol_ijacobian.qfunction) {
45364f98e98SJames Wright     CeedQFunction qf_ijacobian_vol;
45464f98e98SJames Wright 
45564f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc,
45664f98e98SJames Wright                                                     &qf_ijacobian_vol));
45764f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context));
45864f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ijacobian_vol, 0));
45964f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP));
46064f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD));
46164f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
46264f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
46364f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP));
46464f98e98SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
46564f98e98SJames Wright 
466*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op_ijacobian_vol));
467*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
468*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
469*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data));
470*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
471*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
472*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
47364f98e98SJames Wright 
474a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol));
475e334ad8fSJed Brown   }
476e334ad8fSJed Brown 
47777841947SLeila Ghaffari   // -- Create and apply CEED Composite Operator for the entire domain
47877841947SLeila Ghaffari   if (!user->phys->implicit) {  // RHS
4799ad5e8e4SJames Wright     CeedOperator op_rhs;
480e9c36be0SJames Wright 
481e9c36be0SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_rhs));
482*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_rhs, op_rhs_vol));
483e9c36be0SJames Wright     PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, op_rhs, NULL));
484e9c36be0SJames Wright 
4859ad5e8e4SJames Wright     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx));
486e9c36be0SJames Wright 
487e9c36be0SJames Wright     // ----- Get Context Labels for Operator
488e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "solution time", &user->phys->solution_time_label));
489e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "timestep size", &user->phys->timestep_size_label));
490e9c36be0SJames Wright 
491a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
492b4e9a8f8SJames Wright     PetscCall(CreateKSPMass(user, problem));
493577c3361SJames Wright     PetscCheck(app_ctx->sgs_model_type == SGS_MODEL_NONE, user->comm, PETSC_ERR_SUP, "SGS modeling not implemented for explicit timestepping");
49477841947SLeila Ghaffari   } else {  // IFunction
49591c97f41SJames Wright     CeedOperator op_ijacobian = NULL;
49691c97f41SJames Wright 
497e9c36be0SJames Wright     // Create Composite Operaters
498e9c36be0SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &user->op_ifunction));
499*a5b0ec6fSJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(user->op_ifunction, op_ifunction_vol));
500e9c36be0SJames Wright     if (op_ijacobian_vol) {
501e9c36be0SJames Wright       PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_ijacobian));
502e9c36be0SJames Wright       PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_ijacobian, op_ijacobian_vol));
503e9c36be0SJames Wright     }
504e9c36be0SJames Wright     PetscCall(AddBCSubOperators(user, ceed, dm, bc, problem, ceed_data, user->op_ifunction, op_ijacobian));
505e9c36be0SJames Wright 
506e9c36be0SJames Wright     // ----- Get Context Labels for Operator
507e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "solution time", &user->phys->solution_time_label));
508e9c36be0SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ifunction, "timestep size", &user->phys->timestep_size_label));
509e9c36be0SJames Wright 
51091c97f41SJames Wright     if (op_ijacobian) {
51191c97f41SJames Wright       PetscCall(MatCeedCreate(user->dm, user->dm, op_ijacobian, NULL, &user->mat_ijacobian));
51291c97f41SJames Wright       PetscCall(MatCeedSetLocalVectors(user->mat_ijacobian, user->Q_dot_loc, NULL));
51391c97f41SJames Wright       PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ijacobian, "ijacobian time shift", &user->phys->ijacobian_time_shift_label));
51491c97f41SJames Wright       PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian));
515e334ad8fSJed Brown     }
5168f5ab23bSJames Wright     if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SgsDDSetup(ceed, user, ceed_data, problem));
517dada6cc0SJames Wright   }
518f5452247SJames Wright 
519577c3361SJames Wright   if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc));
520f5452247SJames Wright   if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem));
521ee3b213aSJames Wright   if (app_ctx->diff_filter_monitor && !user->diff_filter) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem));
52228134cfaSJames Wright   if (app_ctx->sgs_train_enable) PetscCall(SGS_DD_TrainingSetup(ceed, user, ceed_data, problem));
523f5452247SJames Wright 
524*a5b0ec6fSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data));
525a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i));
526a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol));
527*a5b0ec6fSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ifunction_vol));
528*a5b0ec6fSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_vol));
529ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
53077841947SLeila Ghaffari }
531