xref: /honee/src/setuplibceed.c (revision 1abc283710a9eda9e8040cb5f0ed97f1cda12b56)
1ae2b091fSJames Wright // SPDX-FileCopyrightText: Copyright (c) 2017-2024, HONEE contributors.
2ae2b091fSJames Wright // SPDX-License-Identifier: Apache-2.0 OR BSD-2-Clause
3a515125bSLeila Ghaffari 
4a515125bSLeila Ghaffari /// @file
5ea615d4cSJames Wright /// Setup libCEED Operators for HONEE
6a515125bSLeila Ghaffari 
7e419654dSJeremy L Thompson #include <ceed.h>
8e419654dSJeremy L Thompson #include <petscdmplex.h>
9e419654dSJeremy L Thompson 
10149fb536SJames Wright #include <navierstokes.h>
11a515125bSLeila Ghaffari 
12a78efa86SJames Wright // @brief Create CeedOperator for unstabilized mass KSP for explicit timestepping
130c373b74SJames Wright static PetscErrorCode CreateKSPMassOperator_Unstabilized(Honee honee, CeedOperator *op_mass) {
140c373b74SJames Wright   Ceed                ceed = honee->ceed;
15a78efa86SJames Wright   CeedInt             num_comp_q, q_data_size;
16a78efa86SJames Wright   CeedQFunction       qf_mass;
1701e19bfaSJames Wright   CeedElemRestriction elem_restr_q, elem_restr_qd;
18a78efa86SJames Wright   CeedBasis           basis_q;
19a78efa86SJames Wright   CeedVector          q_data;
20a78efa86SJames Wright 
21a78efa86SJames Wright   PetscFunctionBeginUser;
22a78efa86SJames Wright   {  // Get restriction and basis from the RHS function
23a78efa86SJames Wright     CeedOperator     *sub_ops;
2401e19bfaSJames Wright     CeedOperatorField op_field;
25a78efa86SJames Wright     PetscInt          sub_op_index = 0;  // will be 0 for the volume op
26a78efa86SJames Wright 
270c373b74SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorGetSubList(honee->op_rhs_ctx->op, &sub_ops));
2801e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "q", &op_field));
2901e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_q, &basis_q, NULL));
3001e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorGetFieldByName(sub_ops[sub_op_index], "qdata", &op_field));
3101e19bfaSJames Wright     PetscCallCeed(ceed, CeedOperatorFieldGetData(op_field, NULL, &elem_restr_qd, NULL, &q_data));
32a78efa86SJames Wright   }
33a78efa86SJames Wright 
34a78efa86SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_q, &num_comp_q));
3501e19bfaSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumComponents(elem_restr_qd, &q_data_size));
36a78efa86SJames Wright 
3764dd23feSJames Wright   PetscCall(HoneeMassQFunctionCreate(ceed, num_comp_q, q_data_size, &qf_mass));
38a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_mass, NULL, NULL, op_mass));
3971f2ed29SJames Wright   PetscCallCeed(ceed, CeedOperatorSetName(*op_mass, "RHS Mass Operator"));
40a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "u", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
4101e19bfaSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
42a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(*op_mass, "v", elem_restr_q, basis_q, CEED_VECTOR_ACTIVE));
43a78efa86SJames Wright 
44fff85bd3SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
45fff85bd3SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_q));
46fff85bd3SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q));
4701e19bfaSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
48a78efa86SJames Wright   PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_mass));
49a78efa86SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
50a78efa86SJames Wright }
51a78efa86SJames Wright 
52a78efa86SJames Wright // @brief Create KSP to solve the inverse mass operator for explicit time stepping schemes
530c373b74SJames Wright static PetscErrorCode CreateKSPMass(Honee honee, ProblemData problem) {
540c373b74SJames Wright   Ceed         ceed = honee->ceed;
550c373b74SJames Wright   DM           dm   = honee->dm;
56a78efa86SJames Wright   CeedOperator op_mass;
57a78efa86SJames Wright 
58a78efa86SJames Wright   PetscFunctionBeginUser;
590c373b74SJames Wright   if (problem->create_mass_operator) PetscCall(problem->create_mass_operator(honee, &op_mass));
600c373b74SJames Wright   else PetscCall(CreateKSPMassOperator_Unstabilized(honee, &op_mass));
61a78efa86SJames Wright 
62a78efa86SJames Wright   {  // -- Setup KSP for mass operator
63a78efa86SJames Wright     Mat      mat_mass;
64a78efa86SJames Wright     Vec      Zeros_loc;
65a78efa86SJames Wright     MPI_Comm comm = PetscObjectComm((PetscObject)dm);
66a78efa86SJames Wright 
67a78efa86SJames Wright     PetscCall(DMCreateLocalVector(dm, &Zeros_loc));
68a78efa86SJames Wright     PetscCall(VecZeroEntries(Zeros_loc));
69000d2032SJeremy L Thompson     PetscCall(MatCreateCeed(dm, dm, op_mass, NULL, &mat_mass));
70a78efa86SJames Wright     PetscCall(MatCeedSetLocalVectors(mat_mass, Zeros_loc, NULL));
71a78efa86SJames Wright 
720c373b74SJames Wright     PetscCall(KSPCreate(comm, &honee->mass_ksp));
730c373b74SJames Wright     PetscCall(KSPSetOptionsPrefix(honee->mass_ksp, "mass_"));
7471f2ed29SJames Wright     PetscCall(PetscObjectSetName((PetscObject)honee->mass_ksp, "Explicit Mass"));
75a78efa86SJames Wright     {  // lumped by default
76a78efa86SJames Wright       PC pc;
770c373b74SJames Wright       PetscCall(KSPGetPC(honee->mass_ksp, &pc));
78a78efa86SJames Wright       PetscCall(PCSetType(pc, PCJACOBI));
79a78efa86SJames Wright       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
800c373b74SJames Wright       PetscCall(KSPSetType(honee->mass_ksp, KSPPREONLY));
81a78efa86SJames Wright     }
820c373b74SJames Wright     PetscCall(KSPSetFromOptions_WithMatCeed(honee->mass_ksp, mat_mass));
83a78efa86SJames Wright     PetscCall(VecDestroy(&Zeros_loc));
84a78efa86SJames Wright     PetscCall(MatDestroy(&mat_mass));
85a78efa86SJames Wright   }
86a78efa86SJames Wright 
87a78efa86SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_mass));
88a78efa86SJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
89a78efa86SJames Wright }
90a78efa86SJames Wright 
91d3c60affSJames Wright PetscErrorCode SetupLibceed(Ceed ceed, DM dm, Honee honee, AppCtx app_ctx, ProblemData problem) {
92a515125bSLeila Ghaffari   const PetscInt      num_comp_q = 5;
9366d54740SJames Wright   PetscInt            dim;
94*1abc2837SJames Wright   CeedInt             num_comps_jac_data = problem->num_comps_jac_data, num_comp_x, q_data_size_vol;
95be29160dSJames Wright   CeedElemRestriction elem_restr_jd_i    = NULL, elem_restr_qd;
96be29160dSJames Wright   CeedVector          jac_data           = NULL, q_data;
97bbc90103SJames Wright   CeedOperator        op_ifunction_vol = NULL, op_rhs_vol = NULL, op_ijacobian_vol = NULL;
98bbc90103SJames Wright 
99bbc90103SJames Wright   PetscFunctionBeginUser;
10066d54740SJames Wright   PetscCall(DMGetDimension(dm, &dim));
10166d54740SJames Wright   num_comp_x = dim;
10294a7b3d2SKenneth E. Jansen 
1038c85b835SJames Wright   CeedElemRestriction elem_restr_diff_flux = NULL;
1040880fbb6SJames Wright   CeedVector          div_diff_flux_ceed   = NULL;
1058c85b835SJames Wright   CeedBasis           basis_diff_flux      = NULL;
1068c85b835SJames Wright   CeedEvalMode        eval_mode_diff_flux  = -1;
107bbc90103SJames Wright   {  // Create bases and element restrictions
10815c18037SJames Wright     DMLabel  domain_label = NULL;
10915c18037SJames Wright     PetscInt label_value = 0, height = 0, dm_field = 0;
11067263decSKenneth E. Jansen     DM       dm_coord;
11167263decSKenneth E. Jansen 
112bbc90103SJames Wright     PetscCall(DMGetCoordinateDM(dm, &dm_coord));
113e3663b90SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm, domain_label, label_value, height, dm_field, &honee->basis_q));
114e3663b90SJames Wright     PetscCall(CreateBasisFromPlex(ceed, dm_coord, domain_label, label_value, height, dm_field, &honee->basis_x));
115a515125bSLeila Ghaffari 
116e3663b90SJames Wright     PetscCall(DMPlexCeedElemRestrictionCreate(ceed, dm, domain_label, label_value, height, 0, &honee->elem_restr_q));
117e3663b90SJames Wright     PetscCall(DMPlexCeedElemRestrictionCoordinateCreate(ceed, dm, domain_label, label_value, height, &honee->elem_restr_x));
118*1abc2837SJames Wright     if (num_comps_jac_data) {
119*1abc2837SJames Wright       PetscCall(DMPlexCeedElemRestrictionQDataCreate(ceed, dm, domain_label, label_value, height, num_comps_jac_data, &elem_restr_jd_i));
12028160fc2SJames Wright       PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL));
12128160fc2SJames Wright     }
122ce8bebb6SJames Wright 
123e3663b90SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_q, &honee->q_ceed, NULL));
124e3663b90SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_q, &honee->q_dot_ceed, NULL));
125e3663b90SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_q, &honee->g_ceed, NULL));
126e3663b90SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(honee->elem_restr_x, &honee->x_coord, NULL));
127a515125bSLeila Ghaffari 
128ce8bebb6SJames Wright     {  // -- Copy PETSc coordinate vector into CEED vector
129a515125bSLeila Ghaffari       Vec X_loc;
130cac8aa24SJed Brown       DM  cdm;
131ce8bebb6SJames Wright 
1322b916ea7SJeremy L Thompson       PetscCall(DMGetCellCoordinateDM(dm, &cdm));
1332b916ea7SJeremy L Thompson       if (cdm) {
1342b916ea7SJeremy L Thompson         PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
1352b916ea7SJeremy L Thompson       } else {
1362b916ea7SJeremy L Thompson         PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
137cac8aa24SJed Brown       }
1380c373b74SJames Wright       PetscCall(VecScale(X_loc, honee->units->meter));
139e3663b90SJames Wright       PetscCall(VecCopyPetscToCeed(X_loc, honee->x_coord));
140ce8bebb6SJames Wright     }
141a515125bSLeila Ghaffari 
142e3663b90SJames Wright     PetscCall(QDataGet(ceed, dm, domain_label, label_value, honee->elem_restr_x, honee->basis_x, honee->x_coord, &elem_restr_qd, &q_data,
143be29160dSJames Wright                        &q_data_size_vol));
144ce8bebb6SJames Wright   }
145ce8bebb6SJames Wright 
1467cea0f50SJames Wright   if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) {
1470c373b74SJames Wright     PetscCheck(honee->diff_flux_proj, honee->comm, PETSC_ERR_ARG_WRONGSTATE,
1487cea0f50SJames Wright                "Divergence of diffusive flux projection requested but object not created");
1490c373b74SJames Wright     PetscCall(DivDiffFluxProjectionGetOperatorFieldData(honee->diff_flux_proj, &elem_restr_diff_flux, &basis_diff_flux, &div_diff_flux_ceed,
1500880fbb6SJames Wright                                                         &eval_mode_diff_flux));
1517cea0f50SJames Wright   }
1528c85b835SJames Wright 
153ce8bebb6SJames Wright   {  // -- Create QFunction for ICs
154866f9b4aSJames Wright     CeedBasis     basis_xc;
155ce8bebb6SJames Wright     CeedQFunction qf_ics;
1568f18bb8bSJames Wright     CeedOperator  op_ics;
157ce8bebb6SJames Wright 
158e3663b90SJames Wright     PetscCallCeed(ceed, CeedBasisCreateProjection(honee->basis_x, honee->basis_q, &basis_xc));
159e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qf_func_ptr, problem->ics.qf_loc, &qf_ics));
160e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ics, problem->ics.qfctx));
161ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ics, 0));
162ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "x", num_comp_x, CEED_EVAL_INTERP));
163ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
164ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ics, "q0", num_comp_q, CEED_EVAL_NONE));
165ce8bebb6SJames Wright 
166ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ics, NULL, NULL, &op_ics));
167e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", honee->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
168e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", honee->elem_restr_x, basis_xc, CEED_VECTOR_ACTIVE));
169e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", honee->elem_restr_q, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE));
1700c373b74SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &honee->phys->ics_time_label));
171e3663b90SJames Wright     PetscCall(OperatorApplyContextCreate(NULL, dm, honee->ceed, op_ics, honee->x_coord, NULL, NULL, honee->Q_loc, &honee->op_ics_ctx));
172a515125bSLeila Ghaffari 
173866f9b4aSJames Wright     PetscCallCeed(ceed, CeedBasisDestroy(&basis_xc));
174ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ics));
175ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics));
176ce8bebb6SJames Wright   }
177ce8bebb6SJames Wright 
178e07531f7SJames Wright   if (problem->apply_vol_rhs.qf_func_ptr) {
179ce8bebb6SJames Wright     CeedQFunction qf_rhs_vol;
180ce8bebb6SJames Wright 
181e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qf_func_ptr, problem->apply_vol_rhs.qf_loc, &qf_rhs_vol));
182e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_rhs_vol, problem->apply_vol_rhs.qfctx));
183ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_rhs_vol, 0));
184ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP));
185ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
18689de3142SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
187ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP));
1885f952e8dSJames Wright     if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE)
1890c373b74SJames Wright       PetscCallCeed(ceed, CeedQFunctionAddInput(qf_rhs_vol, "div F_diff", honee->diff_flux_proj->num_diff_flux_comps, eval_mode_diff_flux));
190ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP));
191ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
192ce8bebb6SJames Wright 
193bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_rhs_vol, NULL, NULL, &op_rhs_vol));
194e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
195e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
196be29160dSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
197e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "x", honee->elem_restr_x, honee->basis_x, honee->x_coord));
1985f952e8dSJames Wright     if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE)
1995f952e8dSJames Wright       PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "div F_diff", elem_restr_diff_flux, basis_diff_flux, div_diff_flux_ceed));
200e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
201e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_rhs_vol, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
202ce8bebb6SJames Wright 
203ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_rhs_vol));
204a515125bSLeila Ghaffari   }
205a515125bSLeila Ghaffari 
206e07531f7SJames Wright   if (problem->apply_vol_ifunction.qf_func_ptr) {
207ce8bebb6SJames Wright     CeedQFunction qf_ifunction_vol;
208ce8bebb6SJames Wright 
209e07531f7SJames Wright     PetscCallCeed(
210e07531f7SJames Wright         ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qf_func_ptr, problem->apply_vol_ifunction.qf_loc, &qf_ifunction_vol));
211e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ifunction_vol, problem->apply_vol_ifunction.qfctx));
212ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ifunction_vol, 0));
213ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP));
214ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
215ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP));
21689de3142SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
217ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP));
2180880fbb6SJames Wright     if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE)
2190c373b74SJames Wright       PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ifunction_vol, "div F_diff", honee->diff_flux_proj->num_diff_flux_comps, eval_mode_diff_flux));
220ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP));
221ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
222*1abc2837SJames Wright     if (num_comps_jac_data) PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ifunction_vol, "jac_data", num_comps_jac_data, CEED_EVAL_NONE));
223ce8bebb6SJames Wright 
224bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ifunction_vol, NULL, NULL, &op_ifunction_vol));
225e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
226e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_q", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
227e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "q dot", honee->elem_restr_q, honee->basis_q, honee->q_dot_ceed));
228be29160dSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
229e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "x", honee->elem_restr_x, honee->basis_x, honee->x_coord));
2308c85b835SJames Wright     if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE)
2310880fbb6SJames Wright       PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "div F_diff", elem_restr_diff_flux, basis_diff_flux, div_diff_flux_ceed));
232e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
233e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
234*1abc2837SJames Wright     if (num_comps_jac_data) PetscCallCeed(ceed, CeedOperatorSetField(op_ifunction_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
235752f40e3SJed Brown 
236ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ifunction_vol));
237a515125bSLeila Ghaffari   }
238a515125bSLeila Ghaffari 
239e07531f7SJames Wright   if (problem->apply_vol_ijacobian.qf_func_ptr) {
240ce8bebb6SJames Wright     CeedQFunction qf_ijacobian_vol;
241ce8bebb6SJames Wright 
242e07531f7SJames Wright     PetscCallCeed(
243e07531f7SJames Wright         ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qf_func_ptr, problem->apply_vol_ijacobian.qf_loc, &qf_ijacobian_vol));
244e07531f7SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfctx));
245ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetUserFlopsEstimate(qf_ijacobian_vol, 0));
246ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP));
247ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD));
24889de3142SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
249*1abc2837SJames Wright     if (num_comps_jac_data) PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", num_comps_jac_data, CEED_EVAL_NONE));
250ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP));
251ce8bebb6SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
252ce8bebb6SJames Wright 
253bbc90103SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op_ijacobian_vol));
254e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "dq", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
255e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_dq", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
256be29160dSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "qdata", elem_restr_qd, CEED_BASIS_NONE, q_data));
257*1abc2837SJames Wright     if (num_comps_jac_data) PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "jac_data", elem_restr_jd_i, CEED_BASIS_NONE, jac_data));
258e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
259e3663b90SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_ijacobian_vol, "Grad_v", honee->elem_restr_q, honee->basis_q, CEED_VECTOR_ACTIVE));
260ce8bebb6SJames Wright 
261b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol));
262f0b65372SJed Brown   }
263f0b65372SJed Brown 
264a515125bSLeila Ghaffari   // -- Create and apply CEED Composite Operator for the entire domain
2650c373b74SJames Wright   if (!honee->phys->implicit) {  // RHS
266da5fe0e4SJames Wright     CeedOperator op_rhs;
267866f9b4aSJames Wright 
268866f9b4aSJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_rhs));
269bbc90103SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_rhs, op_rhs_vol));
270d6cac220SJames Wright     for (PetscCount b = 0; b < problem->num_bc_defs; b++) PetscCall(BCDefinitionAddOperators(problem->bc_defs[b], op_rhs, NULL));
271866f9b4aSJames Wright 
2720c373b74SJames Wright     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, honee->q_ceed, honee->g_ceed, honee->Q_loc, NULL, &honee->op_rhs_ctx));
273866f9b4aSJames Wright 
274866f9b4aSJames Wright     // ----- Get Context Labels for Operator
2750c373b74SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "solution time", &honee->phys->solution_time_label));
2760c373b74SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_rhs, "timestep size", &honee->phys->timestep_size_label));
277866f9b4aSJames Wright 
278b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
2790c373b74SJames Wright     PetscCall(CreateKSPMass(honee, problem));
2800c373b74SJames Wright     PetscCheck(app_ctx->sgs_model_type == SGS_MODEL_NONE, honee->comm, PETSC_ERR_SUP, "SGS modeling not implemented for explicit timestepping");
281a515125bSLeila Ghaffari   } else {  // IFunction
282ebfabadfSJames Wright     CeedOperator op_ijacobian = NULL;
283ebfabadfSJames Wright 
284866f9b4aSJames Wright     // Create Composite Operaters
2850c373b74SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &honee->op_ifunction));
2860c373b74SJames Wright     PetscCallCeed(ceed, CeedCompositeOperatorAddSub(honee->op_ifunction, op_ifunction_vol));
287866f9b4aSJames Wright     if (op_ijacobian_vol) {
288866f9b4aSJames Wright       PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, &op_ijacobian));
289866f9b4aSJames Wright       PetscCallCeed(ceed, CeedCompositeOperatorAddSub(op_ijacobian, op_ijacobian_vol));
290866f9b4aSJames Wright     }
291d6cac220SJames Wright     for (PetscCount b = 0; b < problem->num_bc_defs; b++) PetscCall(BCDefinitionAddOperators(problem->bc_defs[b], honee->op_ifunction, op_ijacobian));
292866f9b4aSJames Wright 
293866f9b4aSJames Wright     // ----- Get Context Labels for Operator
2940c373b74SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(honee->op_ifunction, "solution time", &honee->phys->solution_time_label));
2950c373b74SJames Wright     PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(honee->op_ifunction, "timestep size", &honee->phys->timestep_size_label));
296866f9b4aSJames Wright 
297ebfabadfSJames Wright     if (op_ijacobian) {
2980c373b74SJames Wright       PetscCall(MatCreateCeed(honee->dm, honee->dm, op_ijacobian, NULL, &honee->mat_ijacobian));
2990c373b74SJames Wright       PetscCall(MatCeedSetLocalVectors(honee->mat_ijacobian, honee->Q_dot_loc, NULL));
300ebfabadfSJames Wright       PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian));
301f0b65372SJed Brown     }
302e3663b90SJames Wright     if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SgsDDSetup(ceed, honee, problem));
3036d0190e2SJames Wright   }
30491933550SJames Wright 
305d3c60affSJames Wright   if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, dm, honee, problem));
306e3663b90SJames Wright   if (app_ctx->diff_filter_monitor && !honee->diff_filter) PetscCall(DifferentialFilterSetup(ceed, honee, problem));
307e3663b90SJames Wright   if (app_ctx->sgs_train_enable) PetscCall(SGS_DD_TrainingSetup(ceed, honee, problem));
308e3663b90SJames Wright   if (app_ctx->divFdiffproj_method != DIV_DIFF_FLUX_PROJ_NONE) PetscCall(DivDiffFluxProjectionSetup(honee, honee->diff_flux_proj));
30991933550SJames Wright 
310be29160dSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data));
3118c85b835SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data));
312be29160dSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd));
3138c85b835SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i));
3148c85b835SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_diff_flux));
3158c85b835SJames Wright   PetscCallCeed(ceed, CeedBasisDestroy(&basis_diff_flux));
3160880fbb6SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&div_diff_flux_ceed));
317b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol));
318bbc90103SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ifunction_vol));
319bbc90103SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs_vol));
320d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
321a515125bSLeila Ghaffari }
322