xref: /honee/src/setuplibceed.c (revision 67263decd0c197aa72b4d5e710d60290ad1c37bc)
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 
14*67263decSKenneth 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));
28defe8520SJames Wright   CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, elem_restr_offsets_ceed, elem_restr);
29defe8520SJames Wright   PetscCall(PetscFree(elem_restr_offsets_ceed));
30a515125bSLeila Ghaffari 
31d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
32a515125bSLeila Ghaffari }
33a515125bSLeila Ghaffari 
34a515125bSLeila Ghaffari // Utility function to get Ceed Restriction for each domain
3536c6cbc8SJames Wright PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, PetscInt label_value, PetscInt dm_field, CeedInt Q,
3636c6cbc8SJames Wright                                        CeedInt q_data_size, CeedElemRestriction *elem_restr_q, CeedElemRestriction *elem_restr_x,
3736c6cbc8SJames Wright                                        CeedElemRestriction *elem_restr_qd_i) {
38a515125bSLeila Ghaffari   DM                  dm_coord;
39*67263decSKenneth E. Jansen   CeedInt             loc_num_elem;
40defe8520SJames Wright   PetscInt            dim;
41395c96d3SJed Brown   CeedElemRestriction elem_restr_tmp;
42a515125bSLeila Ghaffari   PetscFunctionBeginUser;
43a515125bSLeila Ghaffari 
442b916ea7SJeremy L Thompson   PetscCall(DMGetDimension(dm, &dim));
45a515125bSLeila Ghaffari   dim -= height;
4636c6cbc8SJames Wright   PetscCall(CreateRestrictionFromPlex(ceed, dm, height, domain_label, label_value, dm_field, &elem_restr_tmp));
47395c96d3SJed Brown   if (elem_restr_q) *elem_restr_q = elem_restr_tmp;
48395c96d3SJed Brown   if (elem_restr_x) {
492b916ea7SJeremy L Thompson     PetscCall(DMGetCellCoordinateDM(dm, &dm_coord));
50cac8aa24SJed Brown     if (!dm_coord) {
512b916ea7SJeremy L Thompson       PetscCall(DMGetCoordinateDM(dm, &dm_coord));
52cac8aa24SJed Brown     }
532b916ea7SJeremy L Thompson     PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL));
5436c6cbc8SJames Wright     PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, height, domain_label, label_value, dm_field, elem_restr_x));
55395c96d3SJed Brown   }
56395c96d3SJed Brown   if (elem_restr_qd_i) {
57395c96d3SJed Brown     CeedElemRestrictionGetNumElements(elem_restr_tmp, &loc_num_elem);
58*67263decSKenneth E. Jansen     CeedElemRestrictionCreateStrided(ceed, loc_num_elem, Q, q_data_size, q_data_size * loc_num_elem * Q, CEED_STRIDES_BACKEND, elem_restr_qd_i);
59395c96d3SJed Brown   }
60395c96d3SJed Brown   if (!elem_restr_q) CeedElemRestrictionDestroy(&elem_restr_tmp);
61d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
62a515125bSLeila Ghaffari }
63a515125bSLeila Ghaffari 
6436c6cbc8SJames Wright PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height, CeedInt Q_sur,
652b916ea7SJeremy L Thompson                                 CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian,
66368c645fSJames Wright                                 CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) {
67368c645fSJames Wright   CeedVector          q_data_sur, jac_data_sur;
682b916ea7SJeremy L Thompson   CeedOperator        op_setup_sur, op_apply_bc, op_apply_bc_jacobian = NULL;
692b916ea7SJeremy L Thompson   CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur;
70368c645fSJames Wright   CeedInt             num_qpts_sur;
71368c645fSJames Wright   PetscFunctionBeginUser;
72368c645fSJames Wright 
73368c645fSJames Wright   // --- Get number of quadrature points for the boundaries
74368c645fSJames Wright   CeedBasisGetNumQuadraturePoints(ceed_data->basis_q_sur, &num_qpts_sur);
75368c645fSJames Wright 
76368c645fSJames Wright   // ---- CEED Restriction
77*67263decSKenneth E. Jansen   PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, label_value, 0, num_qpts_sur, q_data_size_sur, &elem_restr_q_sur,
78*67263decSKenneth E. Jansen                                     &elem_restr_x_sur, &elem_restr_qd_i_sur));
79368c645fSJames Wright   if (jac_data_size_sur > 0) {
80368c645fSJames Wright     // State-dependent data will be passed from residual to Jacobian. This will be collocated.
81*67263decSKenneth E. Jansen     PetscCall(
82*67263decSKenneth 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));
83368c645fSJames Wright     CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL);
84368c645fSJames Wright   } else {
85368c645fSJames Wright     elem_restr_jd_i_sur = NULL;
86368c645fSJames Wright     jac_data_sur        = NULL;
87368c645fSJames Wright   }
88368c645fSJames Wright 
89368c645fSJames Wright   // ---- CEED Vector
90defe8520SJames Wright   CeedInt loc_num_elem_sur;
91368c645fSJames Wright   CeedElemRestrictionGetNumElements(elem_restr_q_sur, &loc_num_elem_sur);
922b916ea7SJeremy L Thompson   CeedVectorCreate(ceed, q_data_size_sur * loc_num_elem_sur * num_qpts_sur, &q_data_sur);
93368c645fSJames Wright 
94368c645fSJames Wright   // ---- CEED Operator
95368c645fSJames Wright   // ----- CEED Operator for Setup (geometric factors)
96368c645fSJames Wright   CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur);
972b916ea7SJeremy L Thompson   CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, ceed_data->basis_x_sur, CEED_VECTOR_ACTIVE);
982b916ea7SJeremy L Thompson   CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x_sur, CEED_VECTOR_NONE);
992b916ea7SJeremy L Thompson   CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
100368c645fSJames Wright 
101368c645fSJames Wright   // ----- CEED Operator for Physics
102368c645fSJames Wright   CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc);
1032b916ea7SJeremy L Thompson   CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
1042b916ea7SJeremy L Thompson   CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
1052b916ea7SJeremy L Thompson   CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur);
1062b916ea7SJeremy L Thompson   CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord);
1072b916ea7SJeremy L Thompson   CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
1082b916ea7SJeremy L Thompson   if (elem_restr_jd_i_sur) CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur);
109368c645fSJames Wright 
110368c645fSJames Wright   if (qf_apply_bc_jacobian) {
1112b916ea7SJeremy L Thompson     CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian);
1122b916ea7SJeremy L Thompson     CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
1132b916ea7SJeremy L Thompson     CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
1142b916ea7SJeremy L Thompson     CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur);
1152b916ea7SJeremy L Thompson     CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord);
1162b916ea7SJeremy L Thompson     CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur);
1172b916ea7SJeremy L Thompson     CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
118368c645fSJames Wright   }
119368c645fSJames Wright 
120368c645fSJames Wright   // ----- Apply CEED operator for Setup
1212b916ea7SJeremy L Thompson   CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE);
122368c645fSJames Wright 
123368c645fSJames Wright   // ----- Apply Sub-Operator for Physics
124368c645fSJames Wright   CeedCompositeOperatorAddSub(*op_apply, op_apply_bc);
1252b916ea7SJeremy L Thompson   if (op_apply_bc_jacobian) CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_bc_jacobian);
126368c645fSJames Wright 
127368c645fSJames Wright   // ----- Cleanup
128368c645fSJames Wright   CeedVectorDestroy(&q_data_sur);
129368c645fSJames Wright   CeedVectorDestroy(&jac_data_sur);
130368c645fSJames Wright   CeedElemRestrictionDestroy(&elem_restr_q_sur);
131368c645fSJames Wright   CeedElemRestrictionDestroy(&elem_restr_x_sur);
132368c645fSJames Wright   CeedElemRestrictionDestroy(&elem_restr_qd_i_sur);
133368c645fSJames Wright   CeedElemRestrictionDestroy(&elem_restr_jd_i_sur);
134368c645fSJames Wright   CeedOperatorDestroy(&op_setup_sur);
135368c645fSJames Wright   CeedOperatorDestroy(&op_apply_bc);
136368c645fSJames Wright   CeedOperatorDestroy(&op_apply_bc_jacobian);
137368c645fSJames Wright 
138d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
139368c645fSJames Wright }
140368c645fSJames Wright 
141a515125bSLeila Ghaffari // Utility function to create CEED Composite Operator for the entire domain
1422b916ea7SJeremy L Thompson PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol,
1432b916ea7SJeremy L Thompson                                        CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur,
1442b916ea7SJeremy L Thompson                                        CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) {
145a515125bSLeila Ghaffari   DMLabel domain_label;
146a515125bSLeila Ghaffari   PetscFunctionBeginUser;
147a515125bSLeila Ghaffari 
148a515125bSLeila Ghaffari   // Create Composite Operaters
149a515125bSLeila Ghaffari   CeedCompositeOperatorCreate(ceed, op_apply);
1502b916ea7SJeremy L Thompson   if (op_apply_ijacobian) CeedCompositeOperatorCreate(ceed, op_apply_ijacobian);
151a515125bSLeila Ghaffari 
152a515125bSLeila Ghaffari   // --Apply Sub-Operator for the volume
153a515125bSLeila Ghaffari   CeedCompositeOperatorAddSub(*op_apply, op_apply_vol);
1542b916ea7SJeremy L Thompson   if (op_apply_ijacobian) CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_ijacobian_vol);
155a515125bSLeila Ghaffari 
156a515125bSLeila Ghaffari   // -- Create Sub-Operator for in/outflow BCs
157bb8a0c61SJames Wright   if (phys->has_neumann || 1) {
158a515125bSLeila Ghaffari     // --- Setup
1592b916ea7SJeremy L Thompson     PetscCall(DMGetLabel(dm, "Face Sets", &domain_label));
160a515125bSLeila Ghaffari 
161002797a3SLeila Ghaffari     // --- Create Sub-Operator for inflow boundaries
162002797a3SLeila Ghaffari     for (CeedInt i = 0; i < bc->num_inflow; i++) {
1632b916ea7SJeremy L Thompson       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->inflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
1642b916ea7SJeremy L Thompson                                  ceed_data->qf_apply_inflow, ceed_data->qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian));
165a515125bSLeila Ghaffari     }
166002797a3SLeila Ghaffari     // --- Create Sub-Operator for outflow boundaries
167002797a3SLeila Ghaffari     for (CeedInt i = 0; i < bc->num_outflow; i++) {
1682b916ea7SJeremy L Thompson       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
1692b916ea7SJeremy L Thompson                                  ceed_data->qf_apply_outflow, ceed_data->qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian));
1702556a851SJed Brown     }
171df55ba5fSJames Wright     // --- Create Sub-Operator for freestream boundaries
172df55ba5fSJames Wright     for (CeedInt i = 0; i < bc->num_freestream; i++) {
1732b916ea7SJeremy L Thompson       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
1742b916ea7SJeremy L Thompson                                  ceed_data->qf_apply_freestream, ceed_data->qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian));
175002797a3SLeila Ghaffari     }
176002797a3SLeila Ghaffari   }
17792ada588SJames Wright 
17892ada588SJames Wright   // ----- Get Context Labels for Operator
179a5f46a7bSJeremy L Thompson   CeedOperatorGetContextFieldLabel(*op_apply, "solution time", &phys->solution_time_label);
180a5f46a7bSJeremy L Thompson   CeedOperatorGetContextFieldLabel(*op_apply, "timestep size", &phys->timestep_size_label);
18192ada588SJames Wright 
182d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
183a515125bSLeila Ghaffari }
184a515125bSLeila Ghaffari 
1852b916ea7SJeremy L Thompson PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur,
1862b916ea7SJeremy L Thompson                                  PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian,
18725988f00SJames Wright                                  CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) {
18825988f00SJames Wright   PetscFunctionBeginUser;
18925988f00SJames Wright 
19025988f00SJames Wright   if (apply_bc.qfunction) {
19125988f00SJames Wright     CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc);
19225988f00SJames Wright     CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context);
19325988f00SJames Wright     CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP);
19425988f00SJames Wright     CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD);
19525988f00SJames Wright     CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
19625988f00SJames Wright     CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP);
19725988f00SJames Wright     CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP);
1982b916ea7SJeremy L Thompson     if (jac_data_size_sur) CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE);
19925988f00SJames Wright   }
20025988f00SJames Wright   if (apply_bc_jacobian.qfunction) {
2012b916ea7SJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian);
20225988f00SJames Wright     CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context);
20325988f00SJames Wright     CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP);
20425988f00SJames Wright     CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD);
20525988f00SJames Wright     CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
20625988f00SJames Wright     CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP);
2072b916ea7SJeremy L Thompson     CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE);
20825988f00SJames Wright     CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP);
20925988f00SJames Wright   }
210d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
21125988f00SJames Wright }
21225988f00SJames Wright 
213*67263decSKenneth E. Jansen // -----------------------------------------------------------------------------
214*67263decSKenneth E. Jansen // Convert DM field to DS field
215*67263decSKenneth E. Jansen // -----------------------------------------------------------------------------
216*67263decSKenneth E. Jansen PetscErrorCode DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, PetscInt *ds_field) {
217*67263decSKenneth E. Jansen   PetscDS         ds;
218*67263decSKenneth E. Jansen   IS              field_is;
219*67263decSKenneth E. Jansen   const PetscInt *fields;
220*67263decSKenneth E. Jansen   PetscInt        num_fields;
221*67263decSKenneth E. Jansen 
222*67263decSKenneth E. Jansen   PetscFunctionBeginUser;
223*67263decSKenneth E. Jansen 
224*67263decSKenneth E. Jansen   // Translate dm_field to ds_field
225*67263decSKenneth E. Jansen   PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds, NULL));
226*67263decSKenneth E. Jansen   PetscCall(ISGetIndices(field_is, &fields));
227*67263decSKenneth E. Jansen   PetscCall(ISGetSize(field_is, &num_fields));
228*67263decSKenneth E. Jansen   for (PetscInt i = 0; i < num_fields; i++) {
229*67263decSKenneth E. Jansen     if (dm_field == fields[i]) {
230*67263decSKenneth E. Jansen       *ds_field = i;
231*67263decSKenneth E. Jansen       break;
232*67263decSKenneth E. Jansen     }
233*67263decSKenneth E. Jansen   }
234*67263decSKenneth E. Jansen   PetscCall(ISRestoreIndices(field_is, &fields));
235*67263decSKenneth E. Jansen 
236*67263decSKenneth E. Jansen   if (*ds_field == -1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field);
237*67263decSKenneth E. Jansen 
238*67263decSKenneth E. Jansen   PetscFunctionReturn(PETSC_SUCCESS);
239*67263decSKenneth E. Jansen }
240*67263decSKenneth E. Jansen 
241*67263decSKenneth E. Jansen // -----------------------------------------------------------------------------
242*67263decSKenneth E. Jansen // Utility function - convert from DMPolytopeType to CeedElemTopology
243*67263decSKenneth E. Jansen // -----------------------------------------------------------------------------
244*67263decSKenneth E. Jansen CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) {
245*67263decSKenneth E. Jansen   switch (cell_type) {
246*67263decSKenneth E. Jansen     case DM_POLYTOPE_TRIANGLE:
247*67263decSKenneth E. Jansen       return CEED_TOPOLOGY_TRIANGLE;
248*67263decSKenneth E. Jansen     case DM_POLYTOPE_QUADRILATERAL:
249*67263decSKenneth E. Jansen       return CEED_TOPOLOGY_QUAD;
250*67263decSKenneth E. Jansen     case DM_POLYTOPE_TETRAHEDRON:
251*67263decSKenneth E. Jansen       return CEED_TOPOLOGY_TET;
252*67263decSKenneth E. Jansen     case DM_POLYTOPE_HEXAHEDRON:
253*67263decSKenneth E. Jansen       return CEED_TOPOLOGY_HEX;
254*67263decSKenneth E. Jansen     default:
255*67263decSKenneth E. Jansen       return 0;
256*67263decSKenneth E. Jansen   }
257*67263decSKenneth E. Jansen }
258*67263decSKenneth E. Jansen 
259*67263decSKenneth E. Jansen // -----------------------------------------------------------------------------
260*67263decSKenneth E. Jansen // Create libCEED Basis from PetscTabulation
261*67263decSKenneth E. Jansen // -----------------------------------------------------------------------------
262*67263decSKenneth E. Jansen PetscErrorCode BasisCreateFromTabulation(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt face, PetscFE fe,
263*67263decSKenneth E. Jansen                                          PetscTabulation basis_tabulation, PetscQuadrature quadrature, CeedBasis *basis) {
264*67263decSKenneth E. Jansen   PetscInt           first_point;
265*67263decSKenneth E. Jansen   PetscInt           ids[1] = {label_value};
266*67263decSKenneth E. Jansen   DMLabel            depth_label;
267*67263decSKenneth E. Jansen   DMPolytopeType     cell_type;
268*67263decSKenneth E. Jansen   CeedElemTopology   elem_topo;
269*67263decSKenneth E. Jansen   PetscScalar       *q_points, *interp, *grad;
270*67263decSKenneth E. Jansen   const PetscScalar *q_weights;
271*67263decSKenneth E. Jansen   PetscDualSpace     dual_space;
272*67263decSKenneth E. Jansen   PetscInt           num_dual_basis_vectors;
273*67263decSKenneth E. Jansen   PetscInt           dim, num_comp, P, Q;
274*67263decSKenneth E. Jansen 
275*67263decSKenneth E. Jansen   PetscFunctionBeginUser;
276*67263decSKenneth E. Jansen 
277*67263decSKenneth E. Jansen   // General basis information
278*67263decSKenneth E. Jansen   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
279*67263decSKenneth E. Jansen   PetscCall(PetscFEGetNumComponents(fe, &num_comp));
280*67263decSKenneth E. Jansen   PetscCall(PetscFEGetDualSpace(fe, &dual_space));
281*67263decSKenneth E. Jansen   PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
282*67263decSKenneth E. Jansen   P = num_dual_basis_vectors / num_comp;
283*67263decSKenneth E. Jansen 
284*67263decSKenneth E. Jansen   // Use depth label if no domain label present
285*67263decSKenneth E. Jansen   if (!domain_label) {
286*67263decSKenneth E. Jansen     PetscInt depth;
287*67263decSKenneth E. Jansen 
288*67263decSKenneth E. Jansen     PetscCall(DMPlexGetDepth(dm, &depth));
289*67263decSKenneth E. Jansen     PetscCall(DMPlexGetDepthLabel(dm, &depth_label));
290*67263decSKenneth E. Jansen     ids[0] = depth - height;
291*67263decSKenneth E. Jansen   }
292*67263decSKenneth E. Jansen 
293*67263decSKenneth E. Jansen   // Get cell interp, grad, and quadrature data
294*67263decSKenneth E. Jansen   PetscCall(DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label, 1, ids, height, &first_point, NULL));
295*67263decSKenneth E. Jansen   PetscCall(DMPlexGetCellType(dm, first_point, &cell_type));
296*67263decSKenneth E. Jansen   elem_topo = ElemTopologyP2C(cell_type);
297*67263decSKenneth E. Jansen   if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMPlex topology not supported");
298*67263decSKenneth E. Jansen   {
299*67263decSKenneth E. Jansen     size_t             q_points_size;
300*67263decSKenneth E. Jansen     const PetscScalar *q_points_petsc;
301*67263decSKenneth E. Jansen     PetscInt           q_dim;
302*67263decSKenneth E. Jansen 
303*67263decSKenneth E. Jansen     PetscCall(PetscQuadratureGetData(quadrature, &q_dim, NULL, &Q, &q_points_petsc, &q_weights));
304*67263decSKenneth E. Jansen     q_points_size = Q * dim * sizeof(CeedScalar);
305*67263decSKenneth E. Jansen     PetscCall(PetscCalloc(q_points_size, &q_points));
306*67263decSKenneth E. Jansen     for (PetscInt q = 0; q < Q; q++) {
307*67263decSKenneth E. Jansen       for (PetscInt d = 0; d < q_dim; d++) q_points[q * dim + d] = q_points_petsc[q * q_dim + d];
308*67263decSKenneth E. Jansen     }
309*67263decSKenneth E. Jansen   }
310*67263decSKenneth E. Jansen 
311*67263decSKenneth E. Jansen   {  // Convert to libCEED orientation
312*67263decSKenneth E. Jansen     PetscBool       is_simplex  = PETSC_FALSE;
313*67263decSKenneth E. Jansen     IS              permutation = NULL;
314*67263decSKenneth E. Jansen     const PetscInt *permutation_indices;
315*67263decSKenneth E. Jansen 
316*67263decSKenneth E. Jansen     PetscCall(DMPlexIsSimplex(dm, &is_simplex));
317*67263decSKenneth E. Jansen     if (!is_simplex) {
318*67263decSKenneth E. Jansen       PetscSection section;
319*67263decSKenneth E. Jansen 
320*67263decSKenneth E. Jansen       // -- Get permutation
321*67263decSKenneth E. Jansen       PetscCall(DMGetLocalSection(dm, &section));
322*67263decSKenneth E. Jansen       PetscCall(PetscSectionGetClosurePermutation(section, (PetscObject)dm, dim, num_comp * P, &permutation));
323*67263decSKenneth E. Jansen       PetscCall(ISGetIndices(permutation, &permutation_indices));
324*67263decSKenneth E. Jansen     }
325*67263decSKenneth E. Jansen 
326*67263decSKenneth E. Jansen     // -- Copy interp, grad matrices
327*67263decSKenneth E. Jansen     PetscCall(PetscCalloc(P * Q * sizeof(CeedScalar), &interp));
328*67263decSKenneth E. Jansen     PetscCall(PetscCalloc(P * Q * dim * sizeof(CeedScalar), &grad));
329*67263decSKenneth E. Jansen     const CeedInt c = 0;
330*67263decSKenneth E. Jansen     for (CeedInt q = 0; q < Q; q++) {
331*67263decSKenneth E. Jansen       for (CeedInt p_ceed = 0; p_ceed < P; p_ceed++) {
332*67263decSKenneth E. Jansen         CeedInt p_petsc = is_simplex ? (p_ceed * num_comp) : permutation_indices[p_ceed * num_comp];
333*67263decSKenneth E. Jansen 
334*67263decSKenneth E. Jansen         interp[q * P + p_ceed] = basis_tabulation->T[0][((face * Q + q) * P * num_comp + p_petsc) * num_comp + c];
335*67263decSKenneth E. Jansen         for (CeedInt d = 0; d < dim; d++) {
336*67263decSKenneth 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];
337*67263decSKenneth E. Jansen         }
338*67263decSKenneth E. Jansen       }
339*67263decSKenneth E. Jansen     }
340*67263decSKenneth E. Jansen 
341*67263decSKenneth E. Jansen     // -- Cleanup
342*67263decSKenneth E. Jansen     if (permutation) PetscCall(ISRestoreIndices(permutation, &permutation_indices));
343*67263decSKenneth E. Jansen     PetscCall(ISDestroy(&permutation));
344*67263decSKenneth E. Jansen   }
345*67263decSKenneth E. Jansen 
346*67263decSKenneth E. Jansen   // Finally, create libCEED basis
347*67263decSKenneth E. Jansen   CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points, q_weights, basis);
348*67263decSKenneth E. Jansen   PetscCall(PetscFree(q_points));
349*67263decSKenneth E. Jansen   PetscCall(PetscFree(interp));
350*67263decSKenneth E. Jansen   PetscCall(PetscFree(grad));
351*67263decSKenneth E. Jansen 
352*67263decSKenneth E. Jansen   PetscFunctionReturn(PETSC_SUCCESS);
353*67263decSKenneth E. Jansen }
354*67263decSKenneth E. Jansen 
355*67263decSKenneth E. Jansen // -----------------------------------------------------------------------------
356*67263decSKenneth E. Jansen // Get CEED Basis from DMPlex
357*67263decSKenneth E. Jansen // -----------------------------------------------------------------------------
358*67263decSKenneth E. Jansen PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis) {
359*67263decSKenneth E. Jansen   PetscDS         ds;
360*67263decSKenneth E. Jansen   PetscFE         fe;
361*67263decSKenneth E. Jansen   PetscQuadrature quadrature;
362*67263decSKenneth E. Jansen   PetscBool       is_simplex = PETSC_TRUE;
363*67263decSKenneth E. Jansen   PetscInt        ds_field   = -1;
364*67263decSKenneth E. Jansen 
365*67263decSKenneth E. Jansen   PetscFunctionBeginUser;
366*67263decSKenneth E. Jansen 
367*67263decSKenneth E. Jansen   // Get element information
368*67263decSKenneth E. Jansen   PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL));
369*67263decSKenneth E. Jansen   PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field));
370*67263decSKenneth E. Jansen   PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe));
371*67263decSKenneth E. Jansen   PetscCall(PetscFEGetHeightSubspace(fe, height, &fe));
372*67263decSKenneth E. Jansen   PetscCall(PetscFEGetQuadrature(fe, &quadrature));
373*67263decSKenneth E. Jansen 
374*67263decSKenneth E. Jansen   // Check if simplex or tensor-product mesh
375*67263decSKenneth E. Jansen   PetscCall(DMPlexIsSimplex(dm, &is_simplex));
376*67263decSKenneth E. Jansen 
377*67263decSKenneth E. Jansen   // Build libCEED basis
378*67263decSKenneth E. Jansen   if (is_simplex) {
379*67263decSKenneth E. Jansen     PetscTabulation basis_tabulation;
380*67263decSKenneth E. Jansen     PetscInt        num_derivatives = 1, face = 0;
381*67263decSKenneth E. Jansen 
382*67263decSKenneth E. Jansen     PetscCall(PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation));
383*67263decSKenneth E. Jansen     PetscCall(BasisCreateFromTabulation(ceed, dm, domain_label, label_value, height, face, fe, basis_tabulation, quadrature, basis));
384*67263decSKenneth E. Jansen   } else {
385*67263decSKenneth E. Jansen     PetscDualSpace dual_space;
386*67263decSKenneth E. Jansen     PetscInt       num_dual_basis_vectors;
387*67263decSKenneth E. Jansen     PetscInt       dim, num_comp, P, Q;
388*67263decSKenneth E. Jansen 
389*67263decSKenneth E. Jansen     PetscCall(PetscFEGetSpatialDimension(fe, &dim));
390*67263decSKenneth E. Jansen     PetscCall(PetscFEGetNumComponents(fe, &num_comp));
391*67263decSKenneth E. Jansen     PetscCall(PetscFEGetDualSpace(fe, &dual_space));
392*67263decSKenneth E. Jansen     PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
393*67263decSKenneth E. Jansen     P = num_dual_basis_vectors / num_comp;
394*67263decSKenneth E. Jansen     PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL));
395*67263decSKenneth E. Jansen 
396*67263decSKenneth E. Jansen     CeedInt P_1d = (CeedInt)round(pow(P, 1.0 / dim));
397*67263decSKenneth E. Jansen     CeedInt Q_1d = (CeedInt)round(pow(Q, 1.0 / dim));
398*67263decSKenneth E. Jansen 
399*67263decSKenneth E. Jansen     CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d, CEED_GAUSS, basis);
400*67263decSKenneth E. Jansen   }
401*67263decSKenneth E. Jansen 
402*67263decSKenneth E. Jansen   PetscFunctionReturn(PETSC_SUCCESS);
403*67263decSKenneth E. Jansen }
404*67263decSKenneth E. Jansen 
4052b916ea7SJeremy L Thompson PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc) {
406a515125bSLeila Ghaffari   PetscFunctionBeginUser;
407a515125bSLeila Ghaffari 
408a515125bSLeila Ghaffari   // *****************************************************************************
409a515125bSLeila Ghaffari   // Set up CEED objects for the interior domain (volume)
410a515125bSLeila Ghaffari   // *****************************************************************************
411a515125bSLeila Ghaffari   const PetscInt num_comp_q = 5;
412*67263decSKenneth 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;
413752f40e3SJed Brown   CeedElemRestriction elem_restr_jd_i;
414752f40e3SJed Brown   CeedVector          jac_data;
415*67263decSKenneth E. Jansen   CeedInt             num_qpts;
416a515125bSLeila Ghaffari 
417a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
418a515125bSLeila Ghaffari   // CEED Bases
419a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
420*67263decSKenneth E. Jansen   DM dm_coord;
421*67263decSKenneth E. Jansen   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
422*67263decSKenneth E. Jansen 
423*67263decSKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, &ceed_data->basis_q));
424*67263decSKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, dm_coord, 0, 0, 0, 0, &ceed_data->basis_x));
425*67263decSKenneth E. Jansen   PetscCall(CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &ceed_data->basis_xc));
426*67263decSKenneth E. Jansen   CeedBasisGetNumQuadraturePoints(ceed_data->basis_q, &num_qpts);
427a515125bSLeila Ghaffari 
428a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
429a515125bSLeila Ghaffari   // CEED Restrictions
430a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
431a515125bSLeila Ghaffari   // -- Create restriction
432*67263decSKenneth 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,
4332b916ea7SJeremy L Thompson                                     &ceed_data->elem_restr_qd_i));
434752f40e3SJed Brown 
435*67263decSKenneth E. Jansen   PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, 0, num_qpts, jac_data_size_vol, NULL, NULL, &elem_restr_jd_i));
436a515125bSLeila Ghaffari   // -- Create E vectors
437a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL);
4382b916ea7SJeremy L Thompson   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL);
439a515125bSLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL);
440a515125bSLeila Ghaffari 
441a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
442a515125bSLeila Ghaffari   // CEED QFunctions
443a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
444a515125bSLeila Ghaffari   // -- Create QFunction for quadrature data
4452b916ea7SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, problem->setup_vol.qfunction_loc, &ceed_data->qf_setup_vol);
44615a3537eSJed Brown   if (problem->setup_vol.qfunction_context) {
4472b916ea7SJeremy L Thompson     CeedQFunctionSetContext(ceed_data->qf_setup_vol, problem->setup_vol.qfunction_context);
44815a3537eSJed Brown   }
4492b916ea7SJeremy L Thompson   CeedQFunctionAddInput(ceed_data->qf_setup_vol, "dx", num_comp_x * dim, CEED_EVAL_GRAD);
450a515125bSLeila Ghaffari   CeedQFunctionAddInput(ceed_data->qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT);
4512b916ea7SJeremy L Thompson   CeedQFunctionAddOutput(ceed_data->qf_setup_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE);
452a515125bSLeila Ghaffari 
453a515125bSLeila Ghaffari   // -- Create QFunction for ICs
4542b916ea7SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &ceed_data->qf_ics);
45515a3537eSJed Brown   CeedQFunctionSetContext(ceed_data->qf_ics, problem->ics.qfunction_context);
456a515125bSLeila Ghaffari   CeedQFunctionAddInput(ceed_data->qf_ics, "x", num_comp_x, CEED_EVAL_INTERP);
4574f0244d1SJeremy L Thompson   CeedQFunctionAddInput(ceed_data->qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD);
458a515125bSLeila Ghaffari   CeedQFunctionAddOutput(ceed_data->qf_ics, "q0", num_comp_q, CEED_EVAL_NONE);
459a515125bSLeila Ghaffari 
460a515125bSLeila Ghaffari   // -- Create QFunction for RHS
4619785fe93SJed Brown   if (problem->apply_vol_rhs.qfunction) {
4622b916ea7SJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &ceed_data->qf_rhs_vol);
4632b916ea7SJeremy L Thompson     CeedQFunctionSetContext(ceed_data->qf_rhs_vol, problem->apply_vol_rhs.qfunction_context);
464a515125bSLeila Ghaffari     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP);
4652b916ea7SJeremy L Thompson     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD);
4662b916ea7SJeremy L Thompson     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE);
467a515125bSLeila Ghaffari     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP);
4682b916ea7SJeremy L Thompson     CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP);
4692b916ea7SJeremy L Thompson     CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD);
470a515125bSLeila Ghaffari   }
471a515125bSLeila Ghaffari 
472a515125bSLeila Ghaffari   // -- Create QFunction for IFunction
4739785fe93SJed Brown   if (problem->apply_vol_ifunction.qfunction) {
4742b916ea7SJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc,
4752b916ea7SJeremy L Thompson                                 &ceed_data->qf_ifunction_vol);
4762b916ea7SJeremy L Thompson     CeedQFunctionSetContext(ceed_data->qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context);
4772b916ea7SJeremy L Thompson     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP);
4782b916ea7SJeremy L Thompson     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD);
4792b916ea7SJeremy L Thompson     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP);
4802b916ea7SJeremy L Thompson     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE);
4812b916ea7SJeremy L Thompson     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP);
4822b916ea7SJeremy L Thompson     CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP);
4832b916ea7SJeremy L Thompson     CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD);
4842b916ea7SJeremy L Thompson     CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE);
485a515125bSLeila Ghaffari   }
486a515125bSLeila Ghaffari 
487f0b65372SJed Brown   CeedQFunction qf_ijacobian_vol = NULL;
488f0b65372SJed Brown   if (problem->apply_vol_ijacobian.qfunction) {
4892b916ea7SJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc, &qf_ijacobian_vol);
4902b916ea7SJeremy L Thompson     CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context);
4912b916ea7SJeremy L Thompson     CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP);
4922b916ea7SJeremy L Thompson     CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD);
4932b916ea7SJeremy L Thompson     CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE);
4942b916ea7SJeremy L Thompson     CeedQFunctionAddInput(qf_ijacobian_vol, "x", num_comp_x, CEED_EVAL_INTERP);
4952b916ea7SJeremy L Thompson     CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE);
4962b916ea7SJeremy L Thompson     CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP);
4972b916ea7SJeremy L Thompson     CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD);
498f0b65372SJed Brown   }
499f0b65372SJed Brown 
500a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
501a515125bSLeila Ghaffari   // Element coordinates
502a515125bSLeila Ghaffari   // ---------------------------------------------------------------------------
503a515125bSLeila Ghaffari   // -- Create CEED vector
5042b916ea7SJeremy L Thompson   CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL);
505a515125bSLeila Ghaffari 
506a515125bSLeila Ghaffari   // -- Copy PETSc vector in CEED vector
507a515125bSLeila Ghaffari   Vec X_loc;
508cac8aa24SJed Brown   {
509cac8aa24SJed Brown     DM cdm;
5102b916ea7SJeremy L Thompson     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
5112b916ea7SJeremy L Thompson     if (cdm) {
5122b916ea7SJeremy L Thompson       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
5132b916ea7SJeremy L Thompson     } else {
5142b916ea7SJeremy L Thompson       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
515cac8aa24SJed Brown     }
5162b916ea7SJeremy L Thompson   }
5172b916ea7SJeremy L Thompson   PetscCall(VecScale(X_loc, problem->dm_scale));
518502d3feeSJames Wright   PetscCall(VecCopyP2C(X_loc, ceed_data->x_coord));
519a515125bSLeila Ghaffari 
520a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
521a515125bSLeila Ghaffari   // CEED vectors
522a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
523a515125bSLeila Ghaffari   // -- Create CEED vector for geometric data
5240b0430b7SJames Wright   CeedElemRestrictionCreateVector(ceed_data->elem_restr_qd_i, &ceed_data->q_data, NULL);
525752f40e3SJed Brown   CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL);
5260b0430b7SJames Wright 
527a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
528a515125bSLeila Ghaffari   // CEED Operators
529a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
530a515125bSLeila Ghaffari   // -- Create CEED operator for quadrature data
5312b916ea7SJeremy L Thompson   CeedOperatorCreate(ceed, ceed_data->qf_setup_vol, NULL, NULL, &ceed_data->op_setup_vol);
5322b916ea7SJeremy L Thompson   CeedOperatorSetField(ceed_data->op_setup_vol, "dx", ceed_data->elem_restr_x, ceed_data->basis_x, CEED_VECTOR_ACTIVE);
5332b916ea7SJeremy L Thompson   CeedOperatorSetField(ceed_data->op_setup_vol, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE);
5342b916ea7SJeremy L Thompson   CeedOperatorSetField(ceed_data->op_setup_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
535a515125bSLeila Ghaffari 
536a515125bSLeila Ghaffari   // -- Create CEED operator for ICs
5378f18bb8bSJames Wright   CeedOperator op_ics;
5388f18bb8bSJames Wright   CeedOperatorCreate(ceed, ceed_data->qf_ics, NULL, NULL, &op_ics);
5398f18bb8bSJames Wright   CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE);
5404f0244d1SJeremy L Thompson   CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE);
5418f18bb8bSJames Wright   CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
5428f18bb8bSJames Wright   CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label);
5438f18bb8bSJames Wright   PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx));
5448f18bb8bSJames Wright   CeedOperatorDestroy(&op_ics);
545a515125bSLeila Ghaffari 
546a515125bSLeila Ghaffari   // Create CEED operator for RHS
547a515125bSLeila Ghaffari   if (ceed_data->qf_rhs_vol) {
548a515125bSLeila Ghaffari     CeedOperator op;
549a515125bSLeila Ghaffari     CeedOperatorCreate(ceed, ceed_data->qf_rhs_vol, NULL, NULL, &op);
5502b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5512b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5522b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
5532b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
5542b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5552b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
556a515125bSLeila Ghaffari     user->op_rhs_vol = op;
557a515125bSLeila Ghaffari   }
558a515125bSLeila Ghaffari 
559a515125bSLeila Ghaffari   // -- CEED operator for IFunction
560a515125bSLeila Ghaffari   if (ceed_data->qf_ifunction_vol) {
561a515125bSLeila Ghaffari     CeedOperator op;
562a515125bSLeila Ghaffari     CeedOperatorCreate(ceed, ceed_data->qf_ifunction_vol, NULL, NULL, &op);
5632b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5642b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5652b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed);
5662b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
5672b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
5682b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5692b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5702b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data);
571752f40e3SJed Brown 
572a515125bSLeila Ghaffari     user->op_ifunction_vol = op;
573a515125bSLeila Ghaffari   }
574a515125bSLeila Ghaffari 
575f0b65372SJed Brown   CeedOperator op_ijacobian_vol = NULL;
576f0b65372SJed Brown   if (qf_ijacobian_vol) {
577f0b65372SJed Brown     CeedOperator op;
578f0b65372SJed Brown     CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op);
5792b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5802b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5812b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
5822b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
5832b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data);
5842b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5852b916ea7SJeremy L Thompson     CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
586f0b65372SJed Brown     op_ijacobian_vol = op;
587f0b65372SJed Brown     CeedQFunctionDestroy(&qf_ijacobian_vol);
588f0b65372SJed Brown   }
589f0b65372SJed Brown 
590a515125bSLeila Ghaffari   // *****************************************************************************
591a515125bSLeila Ghaffari   // Set up CEED objects for the exterior domain (surface)
592a515125bSLeila Ghaffari   // *****************************************************************************
5932b916ea7SJeremy L Thompson   CeedInt       height = 1, dim_sur = dim - height, P_sur = app_ctx->degree + 1, Q_sur = P_sur + app_ctx->q_extra;
5942b916ea7SJeremy L Thompson   const CeedInt q_data_size_sur = problem->q_data_size_sur, jac_data_size_sur = problem->jac_data_size_sur;
595a515125bSLeila Ghaffari 
596a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
597a515125bSLeila Ghaffari   // CEED Bases
598a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
599*67263decSKenneth E. Jansen 
600*67263decSKenneth E. Jansen   DMLabel  label   = 0;
601*67263decSKenneth E. Jansen   PetscInt face_id = 0;
602*67263decSKenneth E. Jansen   PetscInt field   = 0;  // Still want the normal, default field
603*67263decSKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, dm, label, face_id, height, field, &ceed_data->basis_q_sur));
604*67263decSKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, face_id, height, field, &ceed_data->basis_x_sur));
605*67263decSKenneth E. Jansen   PetscCall(CeedBasisCreateProjection(ceed_data->basis_x_sur, ceed_data->basis_q_sur, &ceed_data->basis_xc_sur));
606a515125bSLeila Ghaffari 
607a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
608a515125bSLeila Ghaffari   // CEED QFunctions
609a515125bSLeila Ghaffari   // -----------------------------------------------------------------------------
610a515125bSLeila Ghaffari   // -- Create QFunction for quadrature data
6112b916ea7SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, problem->setup_sur.qfunction_loc, &ceed_data->qf_setup_sur);
61215a3537eSJed Brown   if (problem->setup_sur.qfunction_context) {
6132b916ea7SJeremy L Thompson     CeedQFunctionSetContext(ceed_data->qf_setup_sur, problem->setup_sur.qfunction_context);
61415a3537eSJed Brown   }
6152b916ea7SJeremy L Thompson   CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x * dim_sur, CEED_EVAL_GRAD);
616a515125bSLeila Ghaffari   CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT);
6172b916ea7SJeremy L Thompson   CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
618a515125bSLeila Ghaffari 
6192b916ea7SJeremy L Thompson   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow,
6202b916ea7SJeremy L Thompson                               problem->apply_inflow_jacobian, &ceed_data->qf_apply_inflow, &ceed_data->qf_apply_inflow_jacobian));
6212b916ea7SJeremy L Thompson   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow,
6222b916ea7SJeremy L Thompson                               problem->apply_outflow_jacobian, &ceed_data->qf_apply_outflow, &ceed_data->qf_apply_outflow_jacobian));
6232b916ea7SJeremy L Thompson   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream,
6242b916ea7SJeremy L Thompson                               problem->apply_freestream_jacobian, &ceed_data->qf_apply_freestream, &ceed_data->qf_apply_freestream_jacobian));
625a515125bSLeila Ghaffari 
626a515125bSLeila Ghaffari   // *****************************************************************************
627a515125bSLeila Ghaffari   // CEED Operator Apply
628a515125bSLeila Ghaffari   // *****************************************************************************
629a515125bSLeila Ghaffari   // -- Apply CEED Operator for the geometric data
6302b916ea7SJeremy L Thompson   CeedOperatorApply(ceed_data->op_setup_vol, ceed_data->x_coord, ceed_data->q_data, CEED_REQUEST_IMMEDIATE);
631a515125bSLeila Ghaffari 
632a515125bSLeila Ghaffari   // -- Create and apply CEED Composite Operator for the entire domain
633a515125bSLeila Ghaffari   if (!user->phys->implicit) {  // RHS
634da5fe0e4SJames Wright     CeedOperator op_rhs;
635da5fe0e4SJames 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,
636da5fe0e4SJames Wright                                       NULL));
637da5fe0e4SJames Wright     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx));
638da5fe0e4SJames Wright     CeedOperatorDestroy(&op_rhs);
639a515125bSLeila Ghaffari   } else {  // IFunction
6402b916ea7SJeremy L Thompson     PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_ifunction_vol, op_ijacobian_vol, height, P_sur, Q_sur,
6412b916ea7SJeremy L Thompson                                       q_data_size_sur, jac_data_size_sur, &user->op_ifunction, op_ijacobian_vol ? &user->op_ijacobian : NULL));
642f0b65372SJed Brown     if (user->op_ijacobian) {
643a5f46a7bSJeremy L Thompson       CeedOperatorGetContextFieldLabel(user->op_ijacobian, "ijacobian time shift", &user->phys->ijacobian_time_shift_label);
644f0b65372SJed Brown     }
645*67263decSKenneth E. Jansen     if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc, q_data_size_sur));
64691933550SJames Wright     if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SGS_DD_ModelSetup(ceed, user, ceed_data, problem));
6476d0190e2SJames Wright   }
64891933550SJames Wright 
64991933550SJames Wright   if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem));
65088b07121SJames Wright   if (app_ctx->diff_filter_monitor) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem));
65191933550SJames Wright 
652752f40e3SJed Brown   CeedElemRestrictionDestroy(&elem_restr_jd_i);
653e22e0f1cSLeila Ghaffari   CeedOperatorDestroy(&op_ijacobian_vol);
654752f40e3SJed Brown   CeedVectorDestroy(&jac_data);
655d949ddfcSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
656a515125bSLeila Ghaffari }
657