xref: /honee/src/setuplibceed.c (revision b4c37c5c26bd208d7f474cbd9aa0850e82e1a854)
1727da7e7SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2727da7e7SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3a515125bSLeila Ghaffari //
4727da7e7SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5a515125bSLeila Ghaffari //
6727da7e7SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7a515125bSLeila Ghaffari 
8a515125bSLeila Ghaffari /// @file
9a515125bSLeila Ghaffari /// Setup libCEED for Navier-Stokes example using PETSc
10a515125bSLeila Ghaffari 
11e419654dSJeremy L Thompson #include <ceed.h>
12e419654dSJeremy L Thompson #include <petscdmplex.h>
13e419654dSJeremy L Thompson 
1467263decSKenneth E. Jansen #include <petscds.h>
15a515125bSLeila Ghaffari #include "../navierstokes.h"
16a515125bSLeila Ghaffari 
17a515125bSLeila Ghaffari // Utility function to create local CEED restriction
1836c6cbc8SJames Wright PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt label_value, PetscInt dm_field,
1936c6cbc8SJames Wright                                          CeedElemRestriction *elem_restr) {
20defe8520SJames Wright   PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets_petsc;
21defe8520SJames Wright   CeedInt *elem_restr_offsets_ceed;
226b1ccf21SJeremy L Thompson 
23a515125bSLeila Ghaffari   PetscFunctionBeginUser;
24defe8520SJames Wright   PetscCall(
25defe8520SJames Wright       DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets_petsc));
26a515125bSLeila Ghaffari 
27defe8520SJames Wright   PetscCall(IntArrayP2C(num_elem * elem_size, &elem_restr_offsets_petsc, &elem_restr_offsets_ceed));
28*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES,
29*b4c37c5cSJames Wright                                                 elem_restr_offsets_ceed, elem_restr));
30defe8520SJames Wright   PetscCall(PetscFree(elem_restr_offsets_ceed));
31a515125bSLeila Ghaffari 
32d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
33a515125bSLeila Ghaffari }
34a515125bSLeila Ghaffari 
35a515125bSLeila Ghaffari // Utility function to get Ceed Restriction for each domain
3636c6cbc8SJames Wright PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, PetscInt label_value, PetscInt dm_field, CeedInt Q,
3736c6cbc8SJames Wright                                        CeedInt q_data_size, CeedElemRestriction *elem_restr_q, CeedElemRestriction *elem_restr_x,
3836c6cbc8SJames Wright                                        CeedElemRestriction *elem_restr_qd_i) {
39a515125bSLeila Ghaffari   DM                  dm_coord;
4067263decSKenneth E. Jansen   CeedInt             loc_num_elem;
41defe8520SJames Wright   PetscInt            dim;
42395c96d3SJed Brown   CeedElemRestriction elem_restr_tmp;
43a515125bSLeila Ghaffari   PetscFunctionBeginUser;
44a515125bSLeila Ghaffari 
452b916ea7SJeremy L Thompson   PetscCall(DMGetDimension(dm, &dim));
46a515125bSLeila Ghaffari   dim -= height;
4736c6cbc8SJames Wright   PetscCall(CreateRestrictionFromPlex(ceed, dm, height, domain_label, label_value, dm_field, &elem_restr_tmp));
48395c96d3SJed Brown   if (elem_restr_q) *elem_restr_q = elem_restr_tmp;
49395c96d3SJed Brown   if (elem_restr_x) {
502b916ea7SJeremy L Thompson     PetscCall(DMGetCellCoordinateDM(dm, &dm_coord));
51cac8aa24SJed Brown     if (!dm_coord) {
522b916ea7SJeremy L Thompson       PetscCall(DMGetCoordinateDM(dm, &dm_coord));
53cac8aa24SJed Brown     }
542b916ea7SJeremy L Thompson     PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL));
5536c6cbc8SJames Wright     PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, height, domain_label, label_value, dm_field, elem_restr_x));
56395c96d3SJed Brown   }
57395c96d3SJed Brown   if (elem_restr_qd_i) {
58*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_tmp, &loc_num_elem));
59*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, loc_num_elem, Q, q_data_size, q_data_size * loc_num_elem * Q, CEED_STRIDES_BACKEND,
60*b4c37c5cSJames Wright                                                          elem_restr_qd_i));
61395c96d3SJed Brown   }
62*b4c37c5cSJames Wright   if (!elem_restr_q) PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_tmp));
63d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
64a515125bSLeila Ghaffari }
65a515125bSLeila Ghaffari 
6636c6cbc8SJames Wright PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height, CeedInt Q_sur,
672b916ea7SJeremy L Thompson                                 CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian,
68368c645fSJames Wright                                 CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) {
69368c645fSJames Wright   CeedVector          q_data_sur, jac_data_sur;
702b916ea7SJeremy L Thompson   CeedOperator        op_setup_sur, op_apply_bc, op_apply_bc_jacobian = NULL;
712b916ea7SJeremy L Thompson   CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur;
72368c645fSJames Wright   CeedInt             num_qpts_sur;
73368c645fSJames Wright   PetscFunctionBeginUser;
74368c645fSJames Wright 
75368c645fSJames Wright   // --- Get number of quadrature points for the boundaries
76*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(ceed_data->basis_q_sur, &num_qpts_sur));
77368c645fSJames Wright 
78368c645fSJames Wright   // ---- CEED Restriction
7967263decSKenneth E. Jansen   PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, label_value, 0, num_qpts_sur, q_data_size_sur, &elem_restr_q_sur,
8067263decSKenneth E. Jansen                                     &elem_restr_x_sur, &elem_restr_qd_i_sur));
81368c645fSJames Wright   if (jac_data_size_sur > 0) {
82368c645fSJames Wright     // State-dependent data will be passed from residual to Jacobian. This will be collocated.
8367263decSKenneth E. Jansen     PetscCall(
8467263decSKenneth E. Jansen         GetRestrictionForDomain(ceed, dm, height, domain_label, label_value, 0, num_qpts_sur, jac_data_size_sur, NULL, NULL, &elem_restr_jd_i_sur));
85*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL));
86368c645fSJames Wright   } else {
87368c645fSJames Wright     elem_restr_jd_i_sur = NULL;
88368c645fSJames Wright     jac_data_sur        = NULL;
89368c645fSJames Wright   }
90368c645fSJames Wright 
91368c645fSJames Wright   // ---- CEED Vector
92defe8520SJames Wright   CeedInt loc_num_elem_sur;
93*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_q_sur, &loc_num_elem_sur));
94*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorCreate(ceed, q_data_size_sur * loc_num_elem_sur * num_qpts_sur, &q_data_sur));
95368c645fSJames Wright 
96368c645fSJames Wright   // ---- CEED Operator
97368c645fSJames Wright   // ----- CEED Operator for Setup (geometric factors)
98*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur));
99*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, ceed_data->basis_x_sur, CEED_VECTOR_ACTIVE));
100*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x_sur, CEED_VECTOR_NONE));
101*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
102368c645fSJames Wright 
103368c645fSJames Wright   // ----- CEED Operator for Physics
104*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc));
105*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
106*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
107*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur));
108*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord));
109*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
110*b4c37c5cSJames Wright   if (elem_restr_jd_i_sur)
111*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur));
112368c645fSJames Wright 
113368c645fSJames Wright   if (qf_apply_bc_jacobian) {
114*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian));
115*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
116*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
117*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur));
118*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord));
119*b4c37c5cSJames Wright     PetscCallCeed(ceed,
120*b4c37c5cSJames Wright                   CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur));
121*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
122368c645fSJames Wright   }
123368c645fSJames Wright 
124368c645fSJames Wright   // ----- Apply CEED operator for Setup
125*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE));
126368c645fSJames Wright 
127368c645fSJames Wright   // ----- Apply Sub-Operator for Physics
128*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply, op_apply_bc));
129*b4c37c5cSJames Wright   if (op_apply_bc_jacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_bc_jacobian));
130368c645fSJames Wright 
131368c645fSJames Wright   // ----- Cleanup
132*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_sur));
133*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data_sur));
134*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur));
135*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur));
136*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i_sur));
137*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i_sur));
138*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur));
139*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc));
140*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc_jacobian));
141368c645fSJames Wright 
142d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
143368c645fSJames Wright }
144368c645fSJames Wright 
145a515125bSLeila Ghaffari // Utility function to create CEED Composite Operator for the entire domain
1462b916ea7SJeremy L Thompson PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol,
1472b916ea7SJeremy L Thompson                                        CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur,
1482b916ea7SJeremy L Thompson                                        CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) {
149a515125bSLeila Ghaffari   DMLabel domain_label;
150a515125bSLeila Ghaffari   PetscFunctionBeginUser;
151a515125bSLeila Ghaffari 
152a515125bSLeila Ghaffari   // Create Composite Operaters
153*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_apply));
154*b4c37c5cSJames Wright   if (op_apply_ijacobian) PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_apply_ijacobian));
155a515125bSLeila Ghaffari 
156a515125bSLeila Ghaffari   // --Apply Sub-Operator for the volume
157*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply, op_apply_vol));
158*b4c37c5cSJames Wright   if (op_apply_ijacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_ijacobian_vol));
159a515125bSLeila Ghaffari 
160a515125bSLeila Ghaffari   // -- Create Sub-Operator for in/outflow BCs
161bb8a0c61SJames Wright   if (phys->has_neumann || 1) {
162a515125bSLeila Ghaffari     // --- Setup
1632b916ea7SJeremy L Thompson     PetscCall(DMGetLabel(dm, "Face Sets", &domain_label));
164a515125bSLeila Ghaffari 
165002797a3SLeila Ghaffari     // --- Create Sub-Operator for inflow boundaries
166002797a3SLeila Ghaffari     for (CeedInt i = 0; i < bc->num_inflow; i++) {
1672b916ea7SJeremy L Thompson       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->inflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
1682b916ea7SJeremy L Thompson                                  ceed_data->qf_apply_inflow, ceed_data->qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian));
169a515125bSLeila Ghaffari     }
170002797a3SLeila Ghaffari     // --- Create Sub-Operator for outflow boundaries
171002797a3SLeila Ghaffari     for (CeedInt i = 0; i < bc->num_outflow; i++) {
1722b916ea7SJeremy L Thompson       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
1732b916ea7SJeremy L Thompson                                  ceed_data->qf_apply_outflow, ceed_data->qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian));
1742556a851SJed Brown     }
175df55ba5fSJames Wright     // --- Create Sub-Operator for freestream boundaries
176df55ba5fSJames Wright     for (CeedInt i = 0; i < bc->num_freestream; i++) {
1772b916ea7SJeremy L Thompson       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
1782b916ea7SJeremy L Thompson                                  ceed_data->qf_apply_freestream, ceed_data->qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian));
179002797a3SLeila Ghaffari     }
180002797a3SLeila Ghaffari   }
18192ada588SJames Wright 
18292ada588SJames Wright   // ----- Get Context Labels for Operator
183*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(*op_apply, "solution time", &phys->solution_time_label));
184*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(*op_apply, "timestep size", &phys->timestep_size_label));
18592ada588SJames Wright 
186d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
187a515125bSLeila Ghaffari }
188a515125bSLeila Ghaffari 
1892b916ea7SJeremy L Thompson PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur,
1902b916ea7SJeremy L Thompson                                  PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian,
19125988f00SJames Wright                                  CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) {
19225988f00SJames Wright   PetscFunctionBeginUser;
19325988f00SJames Wright 
19425988f00SJames Wright   if (apply_bc.qfunction) {
195*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc));
196*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context));
197*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP));
198*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD));
199*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
200*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP));
201*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP));
202*b4c37c5cSJames Wright     if (jac_data_size_sur) PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
20325988f00SJames Wright   }
20425988f00SJames Wright   if (apply_bc_jacobian.qfunction) {
205*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian));
206*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context));
207*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP));
208*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD));
209*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
210*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP));
211*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
212*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP));
21325988f00SJames Wright   }
214d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
21525988f00SJames Wright }
21625988f00SJames Wright 
21767263decSKenneth E. Jansen // -----------------------------------------------------------------------------
21867263decSKenneth E. Jansen // Convert DM field to DS field
21967263decSKenneth E. Jansen // -----------------------------------------------------------------------------
22067263decSKenneth E. Jansen PetscErrorCode DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, PetscInt *ds_field) {
22167263decSKenneth E. Jansen   PetscDS         ds;
22267263decSKenneth E. Jansen   IS              field_is;
22367263decSKenneth E. Jansen   const PetscInt *fields;
22467263decSKenneth E. Jansen   PetscInt        num_fields;
22567263decSKenneth E. Jansen 
22667263decSKenneth E. Jansen   PetscFunctionBeginUser;
22767263decSKenneth E. Jansen 
22867263decSKenneth E. Jansen   // Translate dm_field to ds_field
22967263decSKenneth E. Jansen   PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds, NULL));
23067263decSKenneth E. Jansen   PetscCall(ISGetIndices(field_is, &fields));
23167263decSKenneth E. Jansen   PetscCall(ISGetSize(field_is, &num_fields));
23267263decSKenneth E. Jansen   for (PetscInt i = 0; i < num_fields; i++) {
23367263decSKenneth E. Jansen     if (dm_field == fields[i]) {
23467263decSKenneth E. Jansen       *ds_field = i;
23567263decSKenneth E. Jansen       break;
23667263decSKenneth E. Jansen     }
23767263decSKenneth E. Jansen   }
23867263decSKenneth E. Jansen   PetscCall(ISRestoreIndices(field_is, &fields));
23967263decSKenneth E. Jansen 
24067263decSKenneth E. Jansen   if (*ds_field == -1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field);
24167263decSKenneth E. Jansen 
24267263decSKenneth E. Jansen   PetscFunctionReturn(PETSC_SUCCESS);
24367263decSKenneth E. Jansen }
24467263decSKenneth E. Jansen 
24567263decSKenneth E. Jansen // -----------------------------------------------------------------------------
24667263decSKenneth E. Jansen // Utility function - convert from DMPolytopeType to CeedElemTopology
24767263decSKenneth E. Jansen // -----------------------------------------------------------------------------
24867263decSKenneth E. Jansen CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) {
24967263decSKenneth E. Jansen   switch (cell_type) {
25067263decSKenneth E. Jansen     case DM_POLYTOPE_TRIANGLE:
25167263decSKenneth E. Jansen       return CEED_TOPOLOGY_TRIANGLE;
25267263decSKenneth E. Jansen     case DM_POLYTOPE_QUADRILATERAL:
25367263decSKenneth E. Jansen       return CEED_TOPOLOGY_QUAD;
25467263decSKenneth E. Jansen     case DM_POLYTOPE_TETRAHEDRON:
25567263decSKenneth E. Jansen       return CEED_TOPOLOGY_TET;
25667263decSKenneth E. Jansen     case DM_POLYTOPE_HEXAHEDRON:
25767263decSKenneth E. Jansen       return CEED_TOPOLOGY_HEX;
25867263decSKenneth E. Jansen     default:
25967263decSKenneth E. Jansen       return 0;
26067263decSKenneth E. Jansen   }
26167263decSKenneth E. Jansen }
26267263decSKenneth E. Jansen 
26367263decSKenneth E. Jansen // -----------------------------------------------------------------------------
26467263decSKenneth E. Jansen // Create libCEED Basis from PetscTabulation
26567263decSKenneth E. Jansen // -----------------------------------------------------------------------------
26667263decSKenneth E. Jansen PetscErrorCode BasisCreateFromTabulation(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt face, PetscFE fe,
26767263decSKenneth E. Jansen                                          PetscTabulation basis_tabulation, PetscQuadrature quadrature, CeedBasis *basis) {
26867263decSKenneth E. Jansen   PetscInt           first_point;
26967263decSKenneth E. Jansen   PetscInt           ids[1] = {label_value};
27067263decSKenneth E. Jansen   DMLabel            depth_label;
27167263decSKenneth E. Jansen   DMPolytopeType     cell_type;
27267263decSKenneth E. Jansen   CeedElemTopology   elem_topo;
27367263decSKenneth E. Jansen   PetscScalar       *q_points, *interp, *grad;
27467263decSKenneth E. Jansen   const PetscScalar *q_weights;
27567263decSKenneth E. Jansen   PetscDualSpace     dual_space;
27667263decSKenneth E. Jansen   PetscInt           num_dual_basis_vectors;
27767263decSKenneth E. Jansen   PetscInt           dim, num_comp, P, Q;
27867263decSKenneth E. Jansen 
27967263decSKenneth E. Jansen   PetscFunctionBeginUser;
28067263decSKenneth E. Jansen 
28167263decSKenneth E. Jansen   // General basis information
28267263decSKenneth E. Jansen   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
28367263decSKenneth E. Jansen   PetscCall(PetscFEGetNumComponents(fe, &num_comp));
28467263decSKenneth E. Jansen   PetscCall(PetscFEGetDualSpace(fe, &dual_space));
28567263decSKenneth E. Jansen   PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
28667263decSKenneth E. Jansen   P = num_dual_basis_vectors / num_comp;
28767263decSKenneth E. Jansen 
28867263decSKenneth E. Jansen   // Use depth label if no domain label present
28967263decSKenneth E. Jansen   if (!domain_label) {
29067263decSKenneth E. Jansen     PetscInt depth;
29167263decSKenneth E. Jansen 
29267263decSKenneth E. Jansen     PetscCall(DMPlexGetDepth(dm, &depth));
29367263decSKenneth E. Jansen     PetscCall(DMPlexGetDepthLabel(dm, &depth_label));
29467263decSKenneth E. Jansen     ids[0] = depth - height;
29567263decSKenneth E. Jansen   }
29667263decSKenneth E. Jansen 
29767263decSKenneth E. Jansen   // Get cell interp, grad, and quadrature data
29867263decSKenneth E. Jansen   PetscCall(DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label, 1, ids, height, &first_point, NULL));
29967263decSKenneth E. Jansen   PetscCall(DMPlexGetCellType(dm, first_point, &cell_type));
30067263decSKenneth E. Jansen   elem_topo = ElemTopologyP2C(cell_type);
30167263decSKenneth E. Jansen   if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMPlex topology not supported");
30267263decSKenneth E. Jansen   {
30367263decSKenneth E. Jansen     size_t             q_points_size;
30467263decSKenneth E. Jansen     const PetscScalar *q_points_petsc;
30567263decSKenneth E. Jansen     PetscInt           q_dim;
30667263decSKenneth E. Jansen 
30767263decSKenneth E. Jansen     PetscCall(PetscQuadratureGetData(quadrature, &q_dim, NULL, &Q, &q_points_petsc, &q_weights));
30867263decSKenneth E. Jansen     q_points_size = Q * dim * sizeof(CeedScalar);
30967263decSKenneth E. Jansen     PetscCall(PetscCalloc(q_points_size, &q_points));
31067263decSKenneth E. Jansen     for (PetscInt q = 0; q < Q; q++) {
31167263decSKenneth E. Jansen       for (PetscInt d = 0; d < q_dim; d++) q_points[q * dim + d] = q_points_petsc[q * q_dim + d];
31267263decSKenneth E. Jansen     }
31367263decSKenneth E. Jansen   }
31467263decSKenneth E. Jansen 
31567263decSKenneth E. Jansen   {  // Convert to libCEED orientation
31667263decSKenneth E. Jansen     PetscBool       is_simplex  = PETSC_FALSE;
31767263decSKenneth E. Jansen     IS              permutation = NULL;
31867263decSKenneth E. Jansen     const PetscInt *permutation_indices;
31967263decSKenneth E. Jansen 
32067263decSKenneth E. Jansen     PetscCall(DMPlexIsSimplex(dm, &is_simplex));
32167263decSKenneth E. Jansen     if (!is_simplex) {
32267263decSKenneth E. Jansen       PetscSection section;
32367263decSKenneth E. Jansen 
32467263decSKenneth E. Jansen       // -- Get permutation
32567263decSKenneth E. Jansen       PetscCall(DMGetLocalSection(dm, &section));
32667263decSKenneth E. Jansen       PetscCall(PetscSectionGetClosurePermutation(section, (PetscObject)dm, dim, num_comp * P, &permutation));
32767263decSKenneth E. Jansen       PetscCall(ISGetIndices(permutation, &permutation_indices));
32867263decSKenneth E. Jansen     }
32967263decSKenneth E. Jansen 
33067263decSKenneth E. Jansen     // -- Copy interp, grad matrices
33167263decSKenneth E. Jansen     PetscCall(PetscCalloc(P * Q * sizeof(CeedScalar), &interp));
33267263decSKenneth E. Jansen     PetscCall(PetscCalloc(P * Q * dim * sizeof(CeedScalar), &grad));
33367263decSKenneth E. Jansen     const CeedInt c = 0;
33467263decSKenneth E. Jansen     for (CeedInt q = 0; q < Q; q++) {
33567263decSKenneth E. Jansen       for (CeedInt p_ceed = 0; p_ceed < P; p_ceed++) {
33667263decSKenneth E. Jansen         CeedInt p_petsc = is_simplex ? (p_ceed * num_comp) : permutation_indices[p_ceed * num_comp];
33767263decSKenneth E. Jansen 
33867263decSKenneth E. Jansen         interp[q * P + p_ceed] = basis_tabulation->T[0][((face * Q + q) * P * num_comp + p_petsc) * num_comp + c];
33967263decSKenneth E. Jansen         for (CeedInt d = 0; d < dim; d++) {
34067263decSKenneth E. Jansen           grad[(d * Q + q) * P + p_ceed] = basis_tabulation->T[1][(((face * Q + q) * P * num_comp + p_petsc) * num_comp + c) * dim + d];
34167263decSKenneth E. Jansen         }
34267263decSKenneth E. Jansen       }
34367263decSKenneth E. Jansen     }
34467263decSKenneth E. Jansen 
34567263decSKenneth E. Jansen     // -- Cleanup
34667263decSKenneth E. Jansen     if (permutation) PetscCall(ISRestoreIndices(permutation, &permutation_indices));
34767263decSKenneth E. Jansen     PetscCall(ISDestroy(&permutation));
34867263decSKenneth E. Jansen   }
34967263decSKenneth E. Jansen 
35067263decSKenneth E. Jansen   // Finally, create libCEED basis
351*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points, q_weights, basis));
35267263decSKenneth E. Jansen   PetscCall(PetscFree(q_points));
35367263decSKenneth E. Jansen   PetscCall(PetscFree(interp));
35467263decSKenneth E. Jansen   PetscCall(PetscFree(grad));
35567263decSKenneth E. Jansen 
35667263decSKenneth E. Jansen   PetscFunctionReturn(PETSC_SUCCESS);
35767263decSKenneth E. Jansen }
35867263decSKenneth E. Jansen 
35967263decSKenneth E. Jansen // -----------------------------------------------------------------------------
36067263decSKenneth E. Jansen // Get CEED Basis from DMPlex
36167263decSKenneth E. Jansen // -----------------------------------------------------------------------------
36267263decSKenneth E. Jansen PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis) {
36367263decSKenneth E. Jansen   PetscDS         ds;
36467263decSKenneth E. Jansen   PetscFE         fe;
36567263decSKenneth E. Jansen   PetscQuadrature quadrature;
36667263decSKenneth E. Jansen   PetscBool       is_simplex = PETSC_TRUE;
36767263decSKenneth E. Jansen   PetscInt        ds_field   = -1;
36867263decSKenneth E. Jansen 
36967263decSKenneth E. Jansen   PetscFunctionBeginUser;
37067263decSKenneth E. Jansen 
37167263decSKenneth E. Jansen   // Get element information
37267263decSKenneth E. Jansen   PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL));
37367263decSKenneth E. Jansen   PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field));
37467263decSKenneth E. Jansen   PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe));
37567263decSKenneth E. Jansen   PetscCall(PetscFEGetHeightSubspace(fe, height, &fe));
37667263decSKenneth E. Jansen   PetscCall(PetscFEGetQuadrature(fe, &quadrature));
37767263decSKenneth E. Jansen 
37867263decSKenneth E. Jansen   // Check if simplex or tensor-product mesh
37967263decSKenneth E. Jansen   PetscCall(DMPlexIsSimplex(dm, &is_simplex));
38067263decSKenneth E. Jansen 
38167263decSKenneth E. Jansen   // Build libCEED basis
38267263decSKenneth E. Jansen   if (is_simplex) {
38367263decSKenneth E. Jansen     PetscTabulation basis_tabulation;
38467263decSKenneth E. Jansen     PetscInt        num_derivatives = 1, face = 0;
38567263decSKenneth E. Jansen 
38667263decSKenneth E. Jansen     PetscCall(PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation));
38767263decSKenneth E. Jansen     PetscCall(BasisCreateFromTabulation(ceed, dm, domain_label, label_value, height, face, fe, basis_tabulation, quadrature, basis));
38867263decSKenneth E. Jansen   } else {
38967263decSKenneth E. Jansen     PetscDualSpace dual_space;
39067263decSKenneth E. Jansen     PetscInt       num_dual_basis_vectors;
39167263decSKenneth E. Jansen     PetscInt       dim, num_comp, P, Q;
39267263decSKenneth E. Jansen 
39367263decSKenneth E. Jansen     PetscCall(PetscFEGetSpatialDimension(fe, &dim));
39467263decSKenneth E. Jansen     PetscCall(PetscFEGetNumComponents(fe, &num_comp));
39567263decSKenneth E. Jansen     PetscCall(PetscFEGetDualSpace(fe, &dual_space));
39667263decSKenneth E. Jansen     PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
39767263decSKenneth E. Jansen     P = num_dual_basis_vectors / num_comp;
39867263decSKenneth E. Jansen     PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL));
39967263decSKenneth E. Jansen 
40067263decSKenneth E. Jansen     CeedInt P_1d = (CeedInt)round(pow(P, 1.0 / dim));
40167263decSKenneth E. Jansen     CeedInt Q_1d = (CeedInt)round(pow(Q, 1.0 / dim));
40267263decSKenneth E. Jansen 
403*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d, CEED_GAUSS, basis));
40467263decSKenneth E. Jansen   }
40567263decSKenneth E. Jansen 
40667263decSKenneth E. Jansen   PetscFunctionReturn(PETSC_SUCCESS);
40767263decSKenneth E. Jansen }
40867263decSKenneth E. Jansen 
4092b916ea7SJeremy L Thompson PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc) {
410a515125bSLeila Ghaffari   PetscFunctionBeginUser;
411a515125bSLeila Ghaffari 
412a515125bSLeila Ghaffari   // *****************************************************************************
413a515125bSLeila Ghaffari   // Set up CEED objects for the interior domain (volume)
414a515125bSLeila Ghaffari   // *****************************************************************************
415a515125bSLeila Ghaffari   const PetscInt num_comp_q = 5;
41667263decSKenneth E. Jansen   const CeedInt  dim = problem->dim, num_comp_x = problem->dim, q_data_size_vol = problem->q_data_size_vol, jac_data_size_vol = num_comp_q + 6 + 3;
417752f40e3SJed Brown   CeedElemRestriction elem_restr_jd_i;
418752f40e3SJed Brown   CeedVector          jac_data;
41967263decSKenneth E. Jansen   CeedInt             num_qpts;
420a515125bSLeila Ghaffari 
421a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
422a515125bSLeila Ghaffari   // CEED Bases
423a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
42467263decSKenneth E. Jansen   DM dm_coord;
42567263decSKenneth E. Jansen   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
42667263decSKenneth E. Jansen 
42767263decSKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, &ceed_data->basis_q));
42867263decSKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, dm_coord, 0, 0, 0, 0, &ceed_data->basis_x));
429*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &ceed_data->basis_xc));
430*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(ceed_data->basis_q, &num_qpts));
431a515125bSLeila Ghaffari 
432a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
433a515125bSLeila Ghaffari   // CEED Restrictions
434a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
435a515125bSLeila Ghaffari   // -- Create restriction
43667263decSKenneth E. Jansen   PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, 0, num_qpts, q_data_size_vol, &ceed_data->elem_restr_q, &ceed_data->elem_restr_x,
4372b916ea7SJeremy L Thompson                                     &ceed_data->elem_restr_qd_i));
438752f40e3SJed Brown 
43967263decSKenneth E. Jansen   PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, 0, num_qpts, jac_data_size_vol, NULL, NULL, &elem_restr_jd_i));
440a515125bSLeila Ghaffari   // -- Create E vectors
441*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL));
442*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL));
443*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL));
444a515125bSLeila Ghaffari 
445a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
446a515125bSLeila Ghaffari   // CEED QFunctions
447a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
448a515125bSLeila Ghaffari   // -- Create QFunction for quadrature data
449*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, problem->setup_vol.qfunction_loc, &ceed_data->qf_setup_vol));
45015a3537eSJed Brown   if (problem->setup_vol.qfunction_context) {
451*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_setup_vol, problem->setup_vol.qfunction_context));
45215a3537eSJed Brown   }
453*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_vol, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
454*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT));
455*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_setup_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
456a515125bSLeila Ghaffari 
457a515125bSLeila Ghaffari   // -- Create QFunction for ICs
458*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &ceed_data->qf_ics));
459*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_ics, problem->ics.qfunction_context));
460*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ics, "x", num_comp_x, CEED_EVAL_INTERP));
461*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
462*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_ics, "q0", num_comp_q, CEED_EVAL_NONE));
463a515125bSLeila Ghaffari 
464a515125bSLeila Ghaffari   // -- Create QFunction for RHS
4659785fe93SJed Brown   if (problem->apply_vol_rhs.qfunction) {
466*b4c37c5cSJames Wright     PetscCallCeed(
467*b4c37c5cSJames Wright         ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &ceed_data->qf_rhs_vol));
468*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_rhs_vol, problem->apply_vol_rhs.qfunction_context));
469*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP));
470*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
471*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
472*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP));
473*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP));
474*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
475a515125bSLeila Ghaffari   }
476a515125bSLeila Ghaffari 
477a515125bSLeila Ghaffari   // -- Create QFunction for IFunction
4789785fe93SJed Brown   if (problem->apply_vol_ifunction.qfunction) {
479*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc,
480*b4c37c5cSJames Wright                                                     &ceed_data->qf_ifunction_vol));
481*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context));
482*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP));
483*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
484*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP));
485*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
486*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP));
487*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP));
488*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
489*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
490a515125bSLeila Ghaffari   }
491a515125bSLeila Ghaffari 
492f0b65372SJed Brown   CeedQFunction qf_ijacobian_vol = NULL;
493f0b65372SJed Brown   if (problem->apply_vol_ijacobian.qfunction) {
494*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc,
495*b4c37c5cSJames Wright                                                     &qf_ijacobian_vol));
496*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context));
497*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP));
498*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD));
499*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
500*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "x", num_comp_x, CEED_EVAL_INTERP));
501*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
502*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP));
503*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
504f0b65372SJed Brown   }
505f0b65372SJed Brown 
506a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
507a515125bSLeila Ghaffari   // Element coordinates
508a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
509a515125bSLeila Ghaffari   // -- Create CEED vector
510*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL));
511a515125bSLeila Ghaffari 
512a515125bSLeila Ghaffari   // -- Copy PETSc vector in CEED vector
513a515125bSLeila Ghaffari   Vec X_loc;
514cac8aa24SJed Brown   {
515cac8aa24SJed Brown     DM cdm;
5162b916ea7SJeremy L Thompson     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
5172b916ea7SJeremy L Thompson     if (cdm) {
5182b916ea7SJeremy L Thompson       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
5192b916ea7SJeremy L Thompson     } else {
5202b916ea7SJeremy L Thompson       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
521cac8aa24SJed Brown     }
5222b916ea7SJeremy L Thompson   }
5232b916ea7SJeremy L Thompson   PetscCall(VecScale(X_loc, problem->dm_scale));
524502d3feeSJames Wright   PetscCall(VecCopyP2C(X_loc, ceed_data->x_coord));
525a515125bSLeila Ghaffari 
526a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
527a515125bSLeila Ghaffari   // CEED vectors
528a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
529a515125bSLeila Ghaffari   // -- Create CEED vector for geometric data
530*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_qd_i, &ceed_data->q_data, NULL));
531*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL));
5320b0430b7SJames Wright 
533a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
534a515125bSLeila Ghaffari   // CEED Operators
535a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
536a515125bSLeila Ghaffari   // -- Create CEED operator for quadrature data
537*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_vol, NULL, NULL, &ceed_data->op_setup_vol));
538*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(ceed_data->op_setup_vol, "dx", ceed_data->elem_restr_x, ceed_data->basis_x, CEED_VECTOR_ACTIVE));
539*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(ceed_data->op_setup_vol, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE));
540*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(ceed_data->op_setup_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
541a515125bSLeila Ghaffari 
542a515125bSLeila Ghaffari   // -- Create CEED operator for ICs
5438f18bb8bSJames Wright   CeedOperator op_ics;
544*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_ics, NULL, NULL, &op_ics));
545*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE));
546*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE));
547*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
548*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label));
5498f18bb8bSJames Wright   PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx));
550*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics));
551a515125bSLeila Ghaffari 
552a515125bSLeila Ghaffari   // Create CEED operator for RHS
553a515125bSLeila Ghaffari   if (ceed_data->qf_rhs_vol) {
554a515125bSLeila Ghaffari     CeedOperator op;
555*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_rhs_vol, NULL, NULL, &op));
556*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
557*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
558*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data));
559*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
560*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
561*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
562a515125bSLeila Ghaffari     user->op_rhs_vol = op;
563a515125bSLeila Ghaffari   }
564a515125bSLeila Ghaffari 
565a515125bSLeila Ghaffari   // -- CEED operator for IFunction
566a515125bSLeila Ghaffari   if (ceed_data->qf_ifunction_vol) {
567a515125bSLeila Ghaffari     CeedOperator op;
568*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_ifunction_vol, NULL, NULL, &op));
569*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
570*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
571*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed));
572*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data));
573*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
574*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
575*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
576*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data));
577752f40e3SJed Brown 
578a515125bSLeila Ghaffari     user->op_ifunction_vol = op;
579a515125bSLeila Ghaffari   }
580a515125bSLeila Ghaffari 
581f0b65372SJed Brown   CeedOperator op_ijacobian_vol = NULL;
582f0b65372SJed Brown   if (qf_ijacobian_vol) {
583f0b65372SJed Brown     CeedOperator op;
584*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op));
585*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
586*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
587*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data));
588*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
589*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data));
590*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
591*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
592f0b65372SJed Brown     op_ijacobian_vol = op;
593*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol));
594f0b65372SJed Brown   }
595f0b65372SJed Brown 
596a515125bSLeila Ghaffari   // *****************************************************************************
597a515125bSLeila Ghaffari   // Set up CEED objects for the exterior domain (surface)
598a515125bSLeila Ghaffari   // *****************************************************************************
5992b916ea7SJeremy L Thompson   CeedInt       height = 1, dim_sur = dim - height, P_sur = app_ctx->degree + 1, Q_sur = P_sur + app_ctx->q_extra;
6002b916ea7SJeremy L Thompson   const CeedInt q_data_size_sur = problem->q_data_size_sur, jac_data_size_sur = problem->jac_data_size_sur;
601a515125bSLeila Ghaffari 
602a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
603a515125bSLeila Ghaffari   // CEED Bases
604a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
60567263decSKenneth E. Jansen 
60667263decSKenneth E. Jansen   DMLabel  label   = 0;
60767263decSKenneth E. Jansen   PetscInt face_id = 0;
60867263decSKenneth E. Jansen   PetscInt field   = 0;  // Still want the normal, default field
60967263decSKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, dm, label, face_id, height, field, &ceed_data->basis_q_sur));
61067263decSKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, face_id, height, field, &ceed_data->basis_x_sur));
611*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x_sur, ceed_data->basis_q_sur, &ceed_data->basis_xc_sur));
612a515125bSLeila Ghaffari 
613a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
614a515125bSLeila Ghaffari   // CEED QFunctions
615a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
616a515125bSLeila Ghaffari   // -- Create QFunction for quadrature data
617*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, problem->setup_sur.qfunction_loc, &ceed_data->qf_setup_sur));
61815a3537eSJed Brown   if (problem->setup_sur.qfunction_context) {
619*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_setup_sur, problem->setup_sur.qfunction_context));
62015a3537eSJed Brown   }
621*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x * dim_sur, CEED_EVAL_GRAD));
622*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT));
623*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
624a515125bSLeila Ghaffari 
6252b916ea7SJeremy L Thompson   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow,
6262b916ea7SJeremy L Thompson                               problem->apply_inflow_jacobian, &ceed_data->qf_apply_inflow, &ceed_data->qf_apply_inflow_jacobian));
6272b916ea7SJeremy L Thompson   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow,
6282b916ea7SJeremy L Thompson                               problem->apply_outflow_jacobian, &ceed_data->qf_apply_outflow, &ceed_data->qf_apply_outflow_jacobian));
6292b916ea7SJeremy L Thompson   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream,
6302b916ea7SJeremy L Thompson                               problem->apply_freestream_jacobian, &ceed_data->qf_apply_freestream, &ceed_data->qf_apply_freestream_jacobian));
631a515125bSLeila Ghaffari 
632a515125bSLeila Ghaffari   // *****************************************************************************
633a515125bSLeila Ghaffari   // CEED Operator Apply
634a515125bSLeila Ghaffari   // *****************************************************************************
635a515125bSLeila Ghaffari   // -- Apply CEED Operator for the geometric data
636*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorApply(ceed_data->op_setup_vol, ceed_data->x_coord, ceed_data->q_data, CEED_REQUEST_IMMEDIATE));
637a515125bSLeila Ghaffari 
638a515125bSLeila Ghaffari   // -- Create and apply CEED Composite Operator for the entire domain
639a515125bSLeila Ghaffari   if (!user->phys->implicit) {  // RHS
640da5fe0e4SJames Wright     CeedOperator op_rhs;
641da5fe0e4SJames 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,
642da5fe0e4SJames Wright                                       NULL));
643da5fe0e4SJames Wright     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx));
644*b4c37c5cSJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
645a515125bSLeila Ghaffari   } else {  // IFunction
6462b916ea7SJeremy L Thompson     PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_ifunction_vol, op_ijacobian_vol, height, P_sur, Q_sur,
6472b916ea7SJeremy L Thompson                                       q_data_size_sur, jac_data_size_sur, &user->op_ifunction, op_ijacobian_vol ? &user->op_ijacobian : NULL));
648f0b65372SJed Brown     if (user->op_ijacobian) {
649*b4c37c5cSJames Wright       PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ijacobian, "ijacobian time shift", &user->phys->ijacobian_time_shift_label));
650f0b65372SJed Brown     }
65167263decSKenneth E. Jansen     if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc, q_data_size_sur));
65291933550SJames Wright     if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SGS_DD_ModelSetup(ceed, user, ceed_data, problem));
6536d0190e2SJames Wright   }
65491933550SJames Wright 
65591933550SJames Wright   if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem));
65688b07121SJames Wright   if (app_ctx->diff_filter_monitor) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem));
65791933550SJames Wright 
658*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i));
659*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol));
660*b4c37c5cSJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data));
661d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
662a515125bSLeila Ghaffari }
663