xref: /libCEED/examples/fluids/src/setuplibceed.c (revision 0814d5a7c7108099a9a6106ea0d2f7f6b264fad8)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
377841947SLeila Ghaffari //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
577841947SLeila Ghaffari //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
777841947SLeila Ghaffari 
877841947SLeila Ghaffari /// @file
977841947SLeila Ghaffari /// Setup libCEED for Navier-Stokes example using PETSc
1077841947SLeila Ghaffari 
1149aac155SJeremy L Thompson #include <ceed.h>
1249aac155SJeremy L Thompson #include <petscdmplex.h>
1349aac155SJeremy L Thompson 
14*0814d5a7SKenneth E. Jansen #include <petscds.h>
1577841947SLeila Ghaffari #include "../navierstokes.h"
1677841947SLeila Ghaffari 
1777841947SLeila Ghaffari // Utility function to create local CEED restriction
18431cd09aSJames Wright PetscErrorCode CreateRestrictionFromPlex(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, CeedInt label_value, PetscInt dm_field,
19431cd09aSJames Wright                                          CeedElemRestriction *elem_restr) {
20457e73b2SJames Wright   PetscInt num_elem, elem_size, num_dof, num_comp, *elem_restr_offsets_petsc;
21457e73b2SJames Wright   CeedInt *elem_restr_offsets_ceed;
227ed3e4cdSJeremy L Thompson 
2377841947SLeila Ghaffari   PetscFunctionBeginUser;
24457e73b2SJames Wright   PetscCall(
25457e73b2SJames Wright       DMPlexGetLocalOffsets(dm, domain_label, label_value, height, dm_field, &num_elem, &elem_size, &num_comp, &num_dof, &elem_restr_offsets_petsc));
2677841947SLeila Ghaffari 
27457e73b2SJames Wright   PetscCall(IntArrayP2C(num_elem * elem_size, &elem_restr_offsets_petsc, &elem_restr_offsets_ceed));
28457e73b2SJames Wright   CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES, elem_restr_offsets_ceed, elem_restr);
29457e73b2SJames Wright   PetscCall(PetscFree(elem_restr_offsets_ceed));
3077841947SLeila Ghaffari 
31ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3277841947SLeila Ghaffari }
3377841947SLeila Ghaffari 
3477841947SLeila Ghaffari // Utility function to get Ceed Restriction for each domain
35431cd09aSJames Wright PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, PetscInt label_value, PetscInt dm_field, CeedInt Q,
36431cd09aSJames Wright                                        CeedInt q_data_size, CeedElemRestriction *elem_restr_q, CeedElemRestriction *elem_restr_x,
37431cd09aSJames Wright                                        CeedElemRestriction *elem_restr_qd_i) {
3877841947SLeila Ghaffari   DM                  dm_coord;
39*0814d5a7SKenneth E. Jansen   CeedInt             loc_num_elem;
40457e73b2SJames Wright   PetscInt            dim;
412534dcc8SJed Brown   CeedElemRestriction elem_restr_tmp;
4277841947SLeila Ghaffari   PetscFunctionBeginUser;
4377841947SLeila Ghaffari 
442b730f8bSJeremy L Thompson   PetscCall(DMGetDimension(dm, &dim));
4577841947SLeila Ghaffari   dim -= height;
46431cd09aSJames Wright   PetscCall(CreateRestrictionFromPlex(ceed, dm, height, domain_label, label_value, dm_field, &elem_restr_tmp));
472534dcc8SJed Brown   if (elem_restr_q) *elem_restr_q = elem_restr_tmp;
482534dcc8SJed Brown   if (elem_restr_x) {
492b730f8bSJeremy L Thompson     PetscCall(DMGetCellCoordinateDM(dm, &dm_coord));
503796c488SJed Brown     if (!dm_coord) {
512b730f8bSJeremy L Thompson       PetscCall(DMGetCoordinateDM(dm, &dm_coord));
523796c488SJed Brown     }
532b730f8bSJeremy L Thompson     PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL));
54431cd09aSJames Wright     PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, height, domain_label, label_value, dm_field, elem_restr_x));
552534dcc8SJed Brown   }
562534dcc8SJed Brown   if (elem_restr_qd_i) {
572534dcc8SJed Brown     CeedElemRestrictionGetNumElements(elem_restr_tmp, &loc_num_elem);
58*0814d5a7SKenneth 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);
592534dcc8SJed Brown   }
602534dcc8SJed Brown   if (!elem_restr_q) CeedElemRestrictionDestroy(&elem_restr_tmp);
61ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
6277841947SLeila Ghaffari }
6377841947SLeila Ghaffari 
64431cd09aSJames Wright PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height, CeedInt Q_sur,
652b730f8bSJeremy L Thompson                                 CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian,
669eb9c072SJames Wright                                 CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) {
679eb9c072SJames Wright   CeedVector          q_data_sur, jac_data_sur;
682b730f8bSJeremy L Thompson   CeedOperator        op_setup_sur, op_apply_bc, op_apply_bc_jacobian = NULL;
692b730f8bSJeremy L Thompson   CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur;
709eb9c072SJames Wright   CeedInt             num_qpts_sur;
719eb9c072SJames Wright   PetscFunctionBeginUser;
729eb9c072SJames Wright 
739eb9c072SJames Wright   // --- Get number of quadrature points for the boundaries
749eb9c072SJames Wright   CeedBasisGetNumQuadraturePoints(ceed_data->basis_q_sur, &num_qpts_sur);
759eb9c072SJames Wright 
769eb9c072SJames Wright   // ---- CEED Restriction
77*0814d5a7SKenneth E. Jansen   PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, label_value, 0, num_qpts_sur, q_data_size_sur, &elem_restr_q_sur,
78*0814d5a7SKenneth E. Jansen                                     &elem_restr_x_sur, &elem_restr_qd_i_sur));
799eb9c072SJames Wright   if (jac_data_size_sur > 0) {
809eb9c072SJames Wright     // State-dependent data will be passed from residual to Jacobian. This will be collocated.
81*0814d5a7SKenneth E. Jansen     PetscCall(
82*0814d5a7SKenneth 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));
839eb9c072SJames Wright     CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL);
849eb9c072SJames Wright   } else {
859eb9c072SJames Wright     elem_restr_jd_i_sur = NULL;
869eb9c072SJames Wright     jac_data_sur        = NULL;
879eb9c072SJames Wright   }
889eb9c072SJames Wright 
899eb9c072SJames Wright   // ---- CEED Vector
90457e73b2SJames Wright   CeedInt loc_num_elem_sur;
919eb9c072SJames Wright   CeedElemRestrictionGetNumElements(elem_restr_q_sur, &loc_num_elem_sur);
922b730f8bSJeremy L Thompson   CeedVectorCreate(ceed, q_data_size_sur * loc_num_elem_sur * num_qpts_sur, &q_data_sur);
939eb9c072SJames Wright 
949eb9c072SJames Wright   // ---- CEED Operator
959eb9c072SJames Wright   // ----- CEED Operator for Setup (geometric factors)
969eb9c072SJames Wright   CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur);
972b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, ceed_data->basis_x_sur, CEED_VECTOR_ACTIVE);
982b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x_sur, CEED_VECTOR_NONE);
992b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
1009eb9c072SJames Wright 
1019eb9c072SJames Wright   // ----- CEED Operator for Physics
1029eb9c072SJames Wright   CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc);
1032b730f8bSJeremy L Thompson   CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
1042b730f8bSJeremy L Thompson   CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
1052b730f8bSJeremy L Thompson   CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur);
1062b730f8bSJeremy L Thompson   CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord);
1072b730f8bSJeremy L Thompson   CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
1082b730f8bSJeremy 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);
1099eb9c072SJames Wright 
1109eb9c072SJames Wright   if (qf_apply_bc_jacobian) {
1112b730f8bSJeremy L Thompson     CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian);
1122b730f8bSJeremy L Thompson     CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
1132b730f8bSJeremy L Thompson     CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
1142b730f8bSJeremy L Thompson     CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur);
1152b730f8bSJeremy L Thompson     CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord);
1162b730f8bSJeremy L Thompson     CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur);
1172b730f8bSJeremy L Thompson     CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE);
1189eb9c072SJames Wright   }
1199eb9c072SJames Wright 
1209eb9c072SJames Wright   // ----- Apply CEED operator for Setup
1212b730f8bSJeremy L Thompson   CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE);
1229eb9c072SJames Wright 
1239eb9c072SJames Wright   // ----- Apply Sub-Operator for Physics
1249eb9c072SJames Wright   CeedCompositeOperatorAddSub(*op_apply, op_apply_bc);
1252b730f8bSJeremy L Thompson   if (op_apply_bc_jacobian) CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_bc_jacobian);
1269eb9c072SJames Wright 
1279eb9c072SJames Wright   // ----- Cleanup
1289eb9c072SJames Wright   CeedVectorDestroy(&q_data_sur);
1299eb9c072SJames Wright   CeedVectorDestroy(&jac_data_sur);
1309eb9c072SJames Wright   CeedElemRestrictionDestroy(&elem_restr_q_sur);
1319eb9c072SJames Wright   CeedElemRestrictionDestroy(&elem_restr_x_sur);
1329eb9c072SJames Wright   CeedElemRestrictionDestroy(&elem_restr_qd_i_sur);
1339eb9c072SJames Wright   CeedElemRestrictionDestroy(&elem_restr_jd_i_sur);
1349eb9c072SJames Wright   CeedOperatorDestroy(&op_setup_sur);
1359eb9c072SJames Wright   CeedOperatorDestroy(&op_apply_bc);
1369eb9c072SJames Wright   CeedOperatorDestroy(&op_apply_bc_jacobian);
1379eb9c072SJames Wright 
138ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1399eb9c072SJames Wright }
1409eb9c072SJames Wright 
14177841947SLeila Ghaffari // Utility function to create CEED Composite Operator for the entire domain
1422b730f8bSJeremy L Thompson PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol,
1432b730f8bSJeremy L Thompson                                        CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur,
1442b730f8bSJeremy L Thompson                                        CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) {
14577841947SLeila Ghaffari   DMLabel domain_label;
14677841947SLeila Ghaffari   PetscFunctionBeginUser;
14777841947SLeila Ghaffari 
14877841947SLeila Ghaffari   // Create Composite Operaters
14977841947SLeila Ghaffari   CeedCompositeOperatorCreate(ceed, op_apply);
1502b730f8bSJeremy L Thompson   if (op_apply_ijacobian) CeedCompositeOperatorCreate(ceed, op_apply_ijacobian);
15177841947SLeila Ghaffari 
15277841947SLeila Ghaffari   // --Apply Sub-Operator for the volume
15377841947SLeila Ghaffari   CeedCompositeOperatorAddSub(*op_apply, op_apply_vol);
1542b730f8bSJeremy L Thompson   if (op_apply_ijacobian) CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_ijacobian_vol);
15577841947SLeila Ghaffari 
15677841947SLeila Ghaffari   // -- Create Sub-Operator for in/outflow BCs
15788626eedSJames Wright   if (phys->has_neumann || 1) {
15877841947SLeila Ghaffari     // --- Setup
1592b730f8bSJeremy L Thompson     PetscCall(DMGetLabel(dm, "Face Sets", &domain_label));
16077841947SLeila Ghaffari 
1612fe7aee7SLeila Ghaffari     // --- Create Sub-Operator for inflow boundaries
1622fe7aee7SLeila Ghaffari     for (CeedInt i = 0; i < bc->num_inflow; i++) {
1632b730f8bSJeremy L Thompson       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->inflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
1642b730f8bSJeremy L Thompson                                  ceed_data->qf_apply_inflow, ceed_data->qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian));
16577841947SLeila Ghaffari     }
1662fe7aee7SLeila Ghaffari     // --- Create Sub-Operator for outflow boundaries
1672fe7aee7SLeila Ghaffari     for (CeedInt i = 0; i < bc->num_outflow; i++) {
1682b730f8bSJeremy L Thompson       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
1692b730f8bSJeremy L Thompson                                  ceed_data->qf_apply_outflow, ceed_data->qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian));
17039c69132SJed Brown     }
171f4ca79c2SJames Wright     // --- Create Sub-Operator for freestream boundaries
172f4ca79c2SJames Wright     for (CeedInt i = 0; i < bc->num_freestream; i++) {
1732b730f8bSJeremy L Thompson       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
1742b730f8bSJeremy L Thompson                                  ceed_data->qf_apply_freestream, ceed_data->qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian));
1752fe7aee7SLeila Ghaffari     }
1762fe7aee7SLeila Ghaffari   }
17711436a05SJames Wright 
17811436a05SJames Wright   // ----- Get Context Labels for Operator
17917b0d5c6SJeremy L Thompson   CeedOperatorGetContextFieldLabel(*op_apply, "solution time", &phys->solution_time_label);
18017b0d5c6SJeremy L Thompson   CeedOperatorGetContextFieldLabel(*op_apply, "timestep size", &phys->timestep_size_label);
18111436a05SJames Wright 
182ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
18377841947SLeila Ghaffari }
18477841947SLeila Ghaffari 
1852b730f8bSJeremy L Thompson PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur,
1862b730f8bSJeremy L Thompson                                  PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian,
1875b881d11SJames Wright                                  CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) {
1885b881d11SJames Wright   PetscFunctionBeginUser;
1895b881d11SJames Wright 
1905b881d11SJames Wright   if (apply_bc.qfunction) {
1915b881d11SJames Wright     CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc);
1925b881d11SJames Wright     CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context);
1935b881d11SJames Wright     CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP);
1945b881d11SJames Wright     CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD);
1955b881d11SJames Wright     CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
1965b881d11SJames Wright     CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP);
1975b881d11SJames Wright     CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP);
1982b730f8bSJeremy L Thompson     if (jac_data_size_sur) CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE);
1995b881d11SJames Wright   }
2005b881d11SJames Wright   if (apply_bc_jacobian.qfunction) {
2012b730f8bSJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian);
2025b881d11SJames Wright     CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context);
2035b881d11SJames Wright     CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP);
2045b881d11SJames Wright     CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD);
2055b881d11SJames Wright     CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
2065b881d11SJames Wright     CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP);
2072b730f8bSJeremy L Thompson     CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE);
2085b881d11SJames Wright     CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP);
2095b881d11SJames Wright   }
210ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2115b881d11SJames Wright }
2125b881d11SJames Wright 
213*0814d5a7SKenneth E. Jansen // -----------------------------------------------------------------------------
214*0814d5a7SKenneth E. Jansen // Convert DM field to DS field
215*0814d5a7SKenneth E. Jansen // -----------------------------------------------------------------------------
216*0814d5a7SKenneth E. Jansen PetscErrorCode DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, PetscInt *ds_field) {
217*0814d5a7SKenneth E. Jansen   PetscDS         ds;
218*0814d5a7SKenneth E. Jansen   IS              field_is;
219*0814d5a7SKenneth E. Jansen   const PetscInt *fields;
220*0814d5a7SKenneth E. Jansen   PetscInt        num_fields;
221*0814d5a7SKenneth E. Jansen 
222*0814d5a7SKenneth E. Jansen   PetscFunctionBeginUser;
223*0814d5a7SKenneth E. Jansen 
224*0814d5a7SKenneth E. Jansen   // Translate dm_field to ds_field
225*0814d5a7SKenneth E. Jansen   PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds, NULL));
226*0814d5a7SKenneth E. Jansen   PetscCall(ISGetIndices(field_is, &fields));
227*0814d5a7SKenneth E. Jansen   PetscCall(ISGetSize(field_is, &num_fields));
228*0814d5a7SKenneth E. Jansen   for (PetscInt i = 0; i < num_fields; i++) {
229*0814d5a7SKenneth E. Jansen     if (dm_field == fields[i]) {
230*0814d5a7SKenneth E. Jansen       *ds_field = i;
231*0814d5a7SKenneth E. Jansen       break;
232*0814d5a7SKenneth E. Jansen     }
233*0814d5a7SKenneth E. Jansen   }
234*0814d5a7SKenneth E. Jansen   PetscCall(ISRestoreIndices(field_is, &fields));
235*0814d5a7SKenneth E. Jansen 
236*0814d5a7SKenneth 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*0814d5a7SKenneth E. Jansen 
238*0814d5a7SKenneth E. Jansen   PetscFunctionReturn(PETSC_SUCCESS);
239*0814d5a7SKenneth E. Jansen }
240*0814d5a7SKenneth E. Jansen 
241*0814d5a7SKenneth E. Jansen // -----------------------------------------------------------------------------
242*0814d5a7SKenneth E. Jansen // Utility function - convert from DMPolytopeType to CeedElemTopology
243*0814d5a7SKenneth E. Jansen // -----------------------------------------------------------------------------
244*0814d5a7SKenneth E. Jansen CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) {
245*0814d5a7SKenneth E. Jansen   switch (cell_type) {
246*0814d5a7SKenneth E. Jansen     case DM_POLYTOPE_TRIANGLE:
247*0814d5a7SKenneth E. Jansen       return CEED_TOPOLOGY_TRIANGLE;
248*0814d5a7SKenneth E. Jansen     case DM_POLYTOPE_QUADRILATERAL:
249*0814d5a7SKenneth E. Jansen       return CEED_TOPOLOGY_QUAD;
250*0814d5a7SKenneth E. Jansen     case DM_POLYTOPE_TETRAHEDRON:
251*0814d5a7SKenneth E. Jansen       return CEED_TOPOLOGY_TET;
252*0814d5a7SKenneth E. Jansen     case DM_POLYTOPE_HEXAHEDRON:
253*0814d5a7SKenneth E. Jansen       return CEED_TOPOLOGY_HEX;
254*0814d5a7SKenneth E. Jansen     default:
255*0814d5a7SKenneth E. Jansen       return 0;
256*0814d5a7SKenneth E. Jansen   }
257*0814d5a7SKenneth E. Jansen }
258*0814d5a7SKenneth E. Jansen 
259*0814d5a7SKenneth E. Jansen // -----------------------------------------------------------------------------
260*0814d5a7SKenneth E. Jansen // Create libCEED Basis from PetscTabulation
261*0814d5a7SKenneth E. Jansen // -----------------------------------------------------------------------------
262*0814d5a7SKenneth E. Jansen PetscErrorCode BasisCreateFromTabulation(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt face, PetscFE fe,
263*0814d5a7SKenneth E. Jansen                                          PetscTabulation basis_tabulation, PetscQuadrature quadrature, CeedBasis *basis) {
264*0814d5a7SKenneth E. Jansen   PetscInt           first_point;
265*0814d5a7SKenneth E. Jansen   PetscInt           ids[1] = {label_value};
266*0814d5a7SKenneth E. Jansen   DMLabel            depth_label;
267*0814d5a7SKenneth E. Jansen   DMPolytopeType     cell_type;
268*0814d5a7SKenneth E. Jansen   CeedElemTopology   elem_topo;
269*0814d5a7SKenneth E. Jansen   PetscScalar       *q_points, *interp, *grad;
270*0814d5a7SKenneth E. Jansen   const PetscScalar *q_weights;
271*0814d5a7SKenneth E. Jansen   PetscDualSpace     dual_space;
272*0814d5a7SKenneth E. Jansen   PetscInt           num_dual_basis_vectors;
273*0814d5a7SKenneth E. Jansen   PetscInt           dim, num_comp, P, Q;
274*0814d5a7SKenneth E. Jansen 
275*0814d5a7SKenneth E. Jansen   PetscFunctionBeginUser;
276*0814d5a7SKenneth E. Jansen 
277*0814d5a7SKenneth E. Jansen   // General basis information
278*0814d5a7SKenneth E. Jansen   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
279*0814d5a7SKenneth E. Jansen   PetscCall(PetscFEGetNumComponents(fe, &num_comp));
280*0814d5a7SKenneth E. Jansen   PetscCall(PetscFEGetDualSpace(fe, &dual_space));
281*0814d5a7SKenneth E. Jansen   PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
282*0814d5a7SKenneth E. Jansen   P = num_dual_basis_vectors / num_comp;
283*0814d5a7SKenneth E. Jansen 
284*0814d5a7SKenneth E. Jansen   // Use depth label if no domain label present
285*0814d5a7SKenneth E. Jansen   if (!domain_label) {
286*0814d5a7SKenneth E. Jansen     PetscInt depth;
287*0814d5a7SKenneth E. Jansen 
288*0814d5a7SKenneth E. Jansen     PetscCall(DMPlexGetDepth(dm, &depth));
289*0814d5a7SKenneth E. Jansen     PetscCall(DMPlexGetDepthLabel(dm, &depth_label));
290*0814d5a7SKenneth E. Jansen     ids[0] = depth - height;
291*0814d5a7SKenneth E. Jansen   }
292*0814d5a7SKenneth E. Jansen 
293*0814d5a7SKenneth E. Jansen   // Get cell interp, grad, and quadrature data
294*0814d5a7SKenneth E. Jansen   PetscCall(DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label, 1, ids, height, &first_point, NULL));
295*0814d5a7SKenneth E. Jansen   PetscCall(DMPlexGetCellType(dm, first_point, &cell_type));
296*0814d5a7SKenneth E. Jansen   elem_topo = ElemTopologyP2C(cell_type);
297*0814d5a7SKenneth E. Jansen   if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMPlex topology not supported");
298*0814d5a7SKenneth E. Jansen   {
299*0814d5a7SKenneth E. Jansen     size_t             q_points_size;
300*0814d5a7SKenneth E. Jansen     const PetscScalar *q_points_petsc;
301*0814d5a7SKenneth E. Jansen     PetscInt           q_dim;
302*0814d5a7SKenneth E. Jansen 
303*0814d5a7SKenneth E. Jansen     PetscCall(PetscQuadratureGetData(quadrature, &q_dim, NULL, &Q, &q_points_petsc, &q_weights));
304*0814d5a7SKenneth E. Jansen     q_points_size = Q * dim * sizeof(CeedScalar);
305*0814d5a7SKenneth E. Jansen     PetscCall(PetscCalloc(q_points_size, &q_points));
306*0814d5a7SKenneth E. Jansen     for (PetscInt q = 0; q < Q; q++) {
307*0814d5a7SKenneth E. Jansen       for (PetscInt d = 0; d < q_dim; d++) q_points[q * dim + d] = q_points_petsc[q * q_dim + d];
308*0814d5a7SKenneth E. Jansen     }
309*0814d5a7SKenneth E. Jansen   }
310*0814d5a7SKenneth E. Jansen 
311*0814d5a7SKenneth E. Jansen   {  // Convert to libCEED orientation
312*0814d5a7SKenneth E. Jansen     PetscBool       is_simplex  = PETSC_FALSE;
313*0814d5a7SKenneth E. Jansen     IS              permutation = NULL;
314*0814d5a7SKenneth E. Jansen     const PetscInt *permutation_indices;
315*0814d5a7SKenneth E. Jansen 
316*0814d5a7SKenneth E. Jansen     PetscCall(DMPlexIsSimplex(dm, &is_simplex));
317*0814d5a7SKenneth E. Jansen     if (!is_simplex) {
318*0814d5a7SKenneth E. Jansen       PetscSection section;
319*0814d5a7SKenneth E. Jansen 
320*0814d5a7SKenneth E. Jansen       // -- Get permutation
321*0814d5a7SKenneth E. Jansen       PetscCall(DMGetLocalSection(dm, &section));
322*0814d5a7SKenneth E. Jansen       PetscCall(PetscSectionGetClosurePermutation(section, (PetscObject)dm, dim, num_comp * P, &permutation));
323*0814d5a7SKenneth E. Jansen       PetscCall(ISGetIndices(permutation, &permutation_indices));
324*0814d5a7SKenneth E. Jansen     }
325*0814d5a7SKenneth E. Jansen 
326*0814d5a7SKenneth E. Jansen     // -- Copy interp, grad matrices
327*0814d5a7SKenneth E. Jansen     PetscCall(PetscCalloc(P * Q * sizeof(CeedScalar), &interp));
328*0814d5a7SKenneth E. Jansen     PetscCall(PetscCalloc(P * Q * dim * sizeof(CeedScalar), &grad));
329*0814d5a7SKenneth E. Jansen     const CeedInt c = 0;
330*0814d5a7SKenneth E. Jansen     for (CeedInt q = 0; q < Q; q++) {
331*0814d5a7SKenneth E. Jansen       for (CeedInt p_ceed = 0; p_ceed < P; p_ceed++) {
332*0814d5a7SKenneth E. Jansen         CeedInt p_petsc = is_simplex ? (p_ceed * num_comp) : permutation_indices[p_ceed * num_comp];
333*0814d5a7SKenneth E. Jansen 
334*0814d5a7SKenneth E. Jansen         interp[q * P + p_ceed] = basis_tabulation->T[0][((face * Q + q) * P * num_comp + p_petsc) * num_comp + c];
335*0814d5a7SKenneth E. Jansen         for (CeedInt d = 0; d < dim; d++) {
336*0814d5a7SKenneth 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*0814d5a7SKenneth E. Jansen         }
338*0814d5a7SKenneth E. Jansen       }
339*0814d5a7SKenneth E. Jansen     }
340*0814d5a7SKenneth E. Jansen 
341*0814d5a7SKenneth E. Jansen     // -- Cleanup
342*0814d5a7SKenneth E. Jansen     if (permutation) PetscCall(ISRestoreIndices(permutation, &permutation_indices));
343*0814d5a7SKenneth E. Jansen     PetscCall(ISDestroy(&permutation));
344*0814d5a7SKenneth E. Jansen   }
345*0814d5a7SKenneth E. Jansen 
346*0814d5a7SKenneth E. Jansen   // Finally, create libCEED basis
347*0814d5a7SKenneth E. Jansen   CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points, q_weights, basis);
348*0814d5a7SKenneth E. Jansen   PetscCall(PetscFree(q_points));
349*0814d5a7SKenneth E. Jansen   PetscCall(PetscFree(interp));
350*0814d5a7SKenneth E. Jansen   PetscCall(PetscFree(grad));
351*0814d5a7SKenneth E. Jansen 
352*0814d5a7SKenneth E. Jansen   PetscFunctionReturn(PETSC_SUCCESS);
353*0814d5a7SKenneth E. Jansen }
354*0814d5a7SKenneth E. Jansen 
355*0814d5a7SKenneth E. Jansen // -----------------------------------------------------------------------------
356*0814d5a7SKenneth E. Jansen // Get CEED Basis from DMPlex
357*0814d5a7SKenneth E. Jansen // -----------------------------------------------------------------------------
358*0814d5a7SKenneth E. Jansen PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis) {
359*0814d5a7SKenneth E. Jansen   PetscDS         ds;
360*0814d5a7SKenneth E. Jansen   PetscFE         fe;
361*0814d5a7SKenneth E. Jansen   PetscQuadrature quadrature;
362*0814d5a7SKenneth E. Jansen   PetscBool       is_simplex = PETSC_TRUE;
363*0814d5a7SKenneth E. Jansen   PetscInt        ds_field   = -1;
364*0814d5a7SKenneth E. Jansen 
365*0814d5a7SKenneth E. Jansen   PetscFunctionBeginUser;
366*0814d5a7SKenneth E. Jansen 
367*0814d5a7SKenneth E. Jansen   // Get element information
368*0814d5a7SKenneth E. Jansen   PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL));
369*0814d5a7SKenneth E. Jansen   PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field));
370*0814d5a7SKenneth E. Jansen   PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe));
371*0814d5a7SKenneth E. Jansen   PetscCall(PetscFEGetHeightSubspace(fe, height, &fe));
372*0814d5a7SKenneth E. Jansen   PetscCall(PetscFEGetQuadrature(fe, &quadrature));
373*0814d5a7SKenneth E. Jansen 
374*0814d5a7SKenneth E. Jansen   // Check if simplex or tensor-product mesh
375*0814d5a7SKenneth E. Jansen   PetscCall(DMPlexIsSimplex(dm, &is_simplex));
376*0814d5a7SKenneth E. Jansen 
377*0814d5a7SKenneth E. Jansen   // Build libCEED basis
378*0814d5a7SKenneth E. Jansen   if (is_simplex) {
379*0814d5a7SKenneth E. Jansen     PetscTabulation basis_tabulation;
380*0814d5a7SKenneth E. Jansen     PetscInt        num_derivatives = 1, face = 0;
381*0814d5a7SKenneth E. Jansen 
382*0814d5a7SKenneth E. Jansen     PetscCall(PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation));
383*0814d5a7SKenneth E. Jansen     PetscCall(BasisCreateFromTabulation(ceed, dm, domain_label, label_value, height, face, fe, basis_tabulation, quadrature, basis));
384*0814d5a7SKenneth E. Jansen   } else {
385*0814d5a7SKenneth E. Jansen     PetscDualSpace dual_space;
386*0814d5a7SKenneth E. Jansen     PetscInt       num_dual_basis_vectors;
387*0814d5a7SKenneth E. Jansen     PetscInt       dim, num_comp, P, Q;
388*0814d5a7SKenneth E. Jansen 
389*0814d5a7SKenneth E. Jansen     PetscCall(PetscFEGetSpatialDimension(fe, &dim));
390*0814d5a7SKenneth E. Jansen     PetscCall(PetscFEGetNumComponents(fe, &num_comp));
391*0814d5a7SKenneth E. Jansen     PetscCall(PetscFEGetDualSpace(fe, &dual_space));
392*0814d5a7SKenneth E. Jansen     PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
393*0814d5a7SKenneth E. Jansen     P = num_dual_basis_vectors / num_comp;
394*0814d5a7SKenneth E. Jansen     PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL));
395*0814d5a7SKenneth E. Jansen 
396*0814d5a7SKenneth E. Jansen     CeedInt P_1d = (CeedInt)round(pow(P, 1.0 / dim));
397*0814d5a7SKenneth E. Jansen     CeedInt Q_1d = (CeedInt)round(pow(Q, 1.0 / dim));
398*0814d5a7SKenneth E. Jansen 
399*0814d5a7SKenneth E. Jansen     CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d, CEED_GAUSS, basis);
400*0814d5a7SKenneth E. Jansen   }
401*0814d5a7SKenneth E. Jansen 
402*0814d5a7SKenneth E. Jansen   PetscFunctionReturn(PETSC_SUCCESS);
403*0814d5a7SKenneth E. Jansen }
404*0814d5a7SKenneth E. Jansen 
4052b730f8bSJeremy L Thompson PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc) {
40677841947SLeila Ghaffari   PetscFunctionBeginUser;
40777841947SLeila Ghaffari 
40877841947SLeila Ghaffari   // *****************************************************************************
40977841947SLeila Ghaffari   // Set up CEED objects for the interior domain (volume)
41077841947SLeila Ghaffari   // *****************************************************************************
41177841947SLeila Ghaffari   const PetscInt num_comp_q = 5;
412*0814d5a7SKenneth 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;
413a3ae0734SJed Brown   CeedElemRestriction elem_restr_jd_i;
414a3ae0734SJed Brown   CeedVector          jac_data;
415*0814d5a7SKenneth E. Jansen   CeedInt             num_qpts;
41677841947SLeila Ghaffari 
41777841947SLeila Ghaffari   // -----------------------------------------------------------------------------
41877841947SLeila Ghaffari   // CEED Bases
41977841947SLeila Ghaffari   // -----------------------------------------------------------------------------
420*0814d5a7SKenneth E. Jansen   DM dm_coord;
421*0814d5a7SKenneth E. Jansen   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
422*0814d5a7SKenneth E. Jansen 
423*0814d5a7SKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, &ceed_data->basis_q));
424*0814d5a7SKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, dm_coord, 0, 0, 0, 0, &ceed_data->basis_x));
425*0814d5a7SKenneth E. Jansen   PetscCall(CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &ceed_data->basis_xc));
426*0814d5a7SKenneth E. Jansen   CeedBasisGetNumQuadraturePoints(ceed_data->basis_q, &num_qpts);
42777841947SLeila Ghaffari 
42877841947SLeila Ghaffari   // -----------------------------------------------------------------------------
42977841947SLeila Ghaffari   // CEED Restrictions
43077841947SLeila Ghaffari   // -----------------------------------------------------------------------------
43177841947SLeila Ghaffari   // -- Create restriction
432*0814d5a7SKenneth 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,
4332b730f8bSJeremy L Thompson                                     &ceed_data->elem_restr_qd_i));
434a3ae0734SJed Brown 
435*0814d5a7SKenneth E. Jansen   PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, 0, num_qpts, jac_data_size_vol, NULL, NULL, &elem_restr_jd_i));
43677841947SLeila Ghaffari   // -- Create E vectors
43777841947SLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL);
4382b730f8bSJeremy L Thompson   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL);
43977841947SLeila Ghaffari   CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL);
44077841947SLeila Ghaffari 
44177841947SLeila Ghaffari   // -----------------------------------------------------------------------------
44277841947SLeila Ghaffari   // CEED QFunctions
44377841947SLeila Ghaffari   // -----------------------------------------------------------------------------
44477841947SLeila Ghaffari   // -- Create QFunction for quadrature data
4452b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, problem->setup_vol.qfunction_loc, &ceed_data->qf_setup_vol);
446841e4c73SJed Brown   if (problem->setup_vol.qfunction_context) {
4472b730f8bSJeremy L Thompson     CeedQFunctionSetContext(ceed_data->qf_setup_vol, problem->setup_vol.qfunction_context);
448841e4c73SJed Brown   }
4492b730f8bSJeremy L Thompson   CeedQFunctionAddInput(ceed_data->qf_setup_vol, "dx", num_comp_x * dim, CEED_EVAL_GRAD);
45077841947SLeila Ghaffari   CeedQFunctionAddInput(ceed_data->qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT);
4512b730f8bSJeremy L Thompson   CeedQFunctionAddOutput(ceed_data->qf_setup_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE);
45277841947SLeila Ghaffari 
45377841947SLeila Ghaffari   // -- Create QFunction for ICs
4542b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &ceed_data->qf_ics);
455841e4c73SJed Brown   CeedQFunctionSetContext(ceed_data->qf_ics, problem->ics.qfunction_context);
45677841947SLeila Ghaffari   CeedQFunctionAddInput(ceed_data->qf_ics, "x", num_comp_x, CEED_EVAL_INTERP);
4571c299e57SJeremy L Thompson   CeedQFunctionAddInput(ceed_data->qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD);
45877841947SLeila Ghaffari   CeedQFunctionAddOutput(ceed_data->qf_ics, "q0", num_comp_q, CEED_EVAL_NONE);
45977841947SLeila Ghaffari 
46077841947SLeila Ghaffari   // -- Create QFunction for RHS
46191e5af17SJed Brown   if (problem->apply_vol_rhs.qfunction) {
4622b730f8bSJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &ceed_data->qf_rhs_vol);
4632b730f8bSJeremy L Thompson     CeedQFunctionSetContext(ceed_data->qf_rhs_vol, problem->apply_vol_rhs.qfunction_context);
46477841947SLeila Ghaffari     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP);
4652b730f8bSJeremy L Thompson     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD);
4662b730f8bSJeremy L Thompson     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE);
46777841947SLeila Ghaffari     CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP);
4682b730f8bSJeremy L Thompson     CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP);
4692b730f8bSJeremy L Thompson     CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD);
47077841947SLeila Ghaffari   }
47177841947SLeila Ghaffari 
47277841947SLeila Ghaffari   // -- Create QFunction for IFunction
47391e5af17SJed Brown   if (problem->apply_vol_ifunction.qfunction) {
4742b730f8bSJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc,
4752b730f8bSJeremy L Thompson                                 &ceed_data->qf_ifunction_vol);
4762b730f8bSJeremy L Thompson     CeedQFunctionSetContext(ceed_data->qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context);
4772b730f8bSJeremy L Thompson     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP);
4782b730f8bSJeremy L Thompson     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD);
4792b730f8bSJeremy L Thompson     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP);
4802b730f8bSJeremy L Thompson     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE);
4812b730f8bSJeremy L Thompson     CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP);
4822b730f8bSJeremy L Thompson     CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP);
4832b730f8bSJeremy L Thompson     CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD);
4842b730f8bSJeremy L Thompson     CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE);
48577841947SLeila Ghaffari   }
48677841947SLeila Ghaffari 
487e334ad8fSJed Brown   CeedQFunction qf_ijacobian_vol = NULL;
488e334ad8fSJed Brown   if (problem->apply_vol_ijacobian.qfunction) {
4892b730f8bSJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc, &qf_ijacobian_vol);
4902b730f8bSJeremy L Thompson     CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context);
4912b730f8bSJeremy L Thompson     CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP);
4922b730f8bSJeremy L Thompson     CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD);
4932b730f8bSJeremy L Thompson     CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE);
4942b730f8bSJeremy L Thompson     CeedQFunctionAddInput(qf_ijacobian_vol, "x", num_comp_x, CEED_EVAL_INTERP);
4952b730f8bSJeremy L Thompson     CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE);
4962b730f8bSJeremy L Thompson     CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP);
4972b730f8bSJeremy L Thompson     CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD);
498e334ad8fSJed Brown   }
499e334ad8fSJed Brown 
50077841947SLeila Ghaffari   // ---------------------------------------------------------------------------
50177841947SLeila Ghaffari   // Element coordinates
50277841947SLeila Ghaffari   // ---------------------------------------------------------------------------
50377841947SLeila Ghaffari   // -- Create CEED vector
5042b730f8bSJeremy L Thompson   CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL);
50577841947SLeila Ghaffari 
50677841947SLeila Ghaffari   // -- Copy PETSc vector in CEED vector
50777841947SLeila Ghaffari   Vec X_loc;
5083796c488SJed Brown   {
5093796c488SJed Brown     DM cdm;
5102b730f8bSJeremy L Thompson     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
5112b730f8bSJeremy L Thompson     if (cdm) {
5122b730f8bSJeremy L Thompson       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
5132b730f8bSJeremy L Thompson     } else {
5142b730f8bSJeremy L Thompson       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
5153796c488SJed Brown     }
5162b730f8bSJeremy L Thompson   }
5172b730f8bSJeremy L Thompson   PetscCall(VecScale(X_loc, problem->dm_scale));
518fe69b334SJames Wright   PetscCall(VecCopyP2C(X_loc, ceed_data->x_coord));
51977841947SLeila Ghaffari 
52077841947SLeila Ghaffari   // -----------------------------------------------------------------------------
52177841947SLeila Ghaffari   // CEED vectors
52277841947SLeila Ghaffari   // -----------------------------------------------------------------------------
52377841947SLeila Ghaffari   // -- Create CEED vector for geometric data
524d52d2babSJames Wright   CeedElemRestrictionCreateVector(ceed_data->elem_restr_qd_i, &ceed_data->q_data, NULL);
525a3ae0734SJed Brown   CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL);
526d52d2babSJames Wright 
52777841947SLeila Ghaffari   // -----------------------------------------------------------------------------
52877841947SLeila Ghaffari   // CEED Operators
52977841947SLeila Ghaffari   // -----------------------------------------------------------------------------
53077841947SLeila Ghaffari   // -- Create CEED operator for quadrature data
5312b730f8bSJeremy L Thompson   CeedOperatorCreate(ceed, ceed_data->qf_setup_vol, NULL, NULL, &ceed_data->op_setup_vol);
5322b730f8bSJeremy L Thompson   CeedOperatorSetField(ceed_data->op_setup_vol, "dx", ceed_data->elem_restr_x, ceed_data->basis_x, CEED_VECTOR_ACTIVE);
5332b730f8bSJeremy L Thompson   CeedOperatorSetField(ceed_data->op_setup_vol, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE);
5342b730f8bSJeremy L Thompson   CeedOperatorSetField(ceed_data->op_setup_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
53577841947SLeila Ghaffari 
53677841947SLeila Ghaffari   // -- Create CEED operator for ICs
5375263e9c6SJames Wright   CeedOperator op_ics;
5385263e9c6SJames Wright   CeedOperatorCreate(ceed, ceed_data->qf_ics, NULL, NULL, &op_ics);
5395263e9c6SJames Wright   CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE);
5401c299e57SJeremy L Thompson   CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE);
5415263e9c6SJames Wright   CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
5425263e9c6SJames Wright   CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label);
5435263e9c6SJames Wright   PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx));
5445263e9c6SJames Wright   CeedOperatorDestroy(&op_ics);
54577841947SLeila Ghaffari 
54677841947SLeila Ghaffari   // Create CEED operator for RHS
54777841947SLeila Ghaffari   if (ceed_data->qf_rhs_vol) {
54877841947SLeila Ghaffari     CeedOperator op;
54977841947SLeila Ghaffari     CeedOperatorCreate(ceed, ceed_data->qf_rhs_vol, NULL, NULL, &op);
5502b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5512b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5522b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
5532b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
5542b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5552b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
55677841947SLeila Ghaffari     user->op_rhs_vol = op;
55777841947SLeila Ghaffari   }
55877841947SLeila Ghaffari 
55977841947SLeila Ghaffari   // -- CEED operator for IFunction
56077841947SLeila Ghaffari   if (ceed_data->qf_ifunction_vol) {
56177841947SLeila Ghaffari     CeedOperator op;
56277841947SLeila Ghaffari     CeedOperatorCreate(ceed, ceed_data->qf_ifunction_vol, NULL, NULL, &op);
5632b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5642b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5652b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed);
5662b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
5672b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
5682b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5692b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5702b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data);
571a3ae0734SJed Brown 
57277841947SLeila Ghaffari     user->op_ifunction_vol = op;
57377841947SLeila Ghaffari   }
57477841947SLeila Ghaffari 
575e334ad8fSJed Brown   CeedOperator op_ijacobian_vol = NULL;
576e334ad8fSJed Brown   if (qf_ijacobian_vol) {
577e334ad8fSJed Brown     CeedOperator op;
578e334ad8fSJed Brown     CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op);
5792b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5802b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5812b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data);
5822b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord);
5832b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data);
5842b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
5852b730f8bSJeremy L Thompson     CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE);
586e334ad8fSJed Brown     op_ijacobian_vol = op;
587e334ad8fSJed Brown     CeedQFunctionDestroy(&qf_ijacobian_vol);
588e334ad8fSJed Brown   }
589e334ad8fSJed Brown 
59077841947SLeila Ghaffari   // *****************************************************************************
59177841947SLeila Ghaffari   // Set up CEED objects for the exterior domain (surface)
59277841947SLeila Ghaffari   // *****************************************************************************
5932b730f8bSJeremy L Thompson   CeedInt       height = 1, dim_sur = dim - height, P_sur = app_ctx->degree + 1, Q_sur = P_sur + app_ctx->q_extra;
5942b730f8bSJeremy L Thompson   const CeedInt q_data_size_sur = problem->q_data_size_sur, jac_data_size_sur = problem->jac_data_size_sur;
59577841947SLeila Ghaffari 
59677841947SLeila Ghaffari   // -----------------------------------------------------------------------------
59777841947SLeila Ghaffari   // CEED Bases
59877841947SLeila Ghaffari   // -----------------------------------------------------------------------------
599*0814d5a7SKenneth E. Jansen 
600*0814d5a7SKenneth E. Jansen   DMLabel  label   = 0;
601*0814d5a7SKenneth E. Jansen   PetscInt face_id = 0;
602*0814d5a7SKenneth E. Jansen   PetscInt field   = 0;  // Still want the normal, default field
603*0814d5a7SKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, dm, label, face_id, height, field, &ceed_data->basis_q_sur));
604*0814d5a7SKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, face_id, height, field, &ceed_data->basis_x_sur));
605*0814d5a7SKenneth E. Jansen   PetscCall(CeedBasisCreateProjection(ceed_data->basis_x_sur, ceed_data->basis_q_sur, &ceed_data->basis_xc_sur));
60677841947SLeila Ghaffari 
60777841947SLeila Ghaffari   // -----------------------------------------------------------------------------
60877841947SLeila Ghaffari   // CEED QFunctions
60977841947SLeila Ghaffari   // -----------------------------------------------------------------------------
61077841947SLeila Ghaffari   // -- Create QFunction for quadrature data
6112b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, problem->setup_sur.qfunction_loc, &ceed_data->qf_setup_sur);
612841e4c73SJed Brown   if (problem->setup_sur.qfunction_context) {
6132b730f8bSJeremy L Thompson     CeedQFunctionSetContext(ceed_data->qf_setup_sur, problem->setup_sur.qfunction_context);
614841e4c73SJed Brown   }
6152b730f8bSJeremy L Thompson   CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x * dim_sur, CEED_EVAL_GRAD);
61677841947SLeila Ghaffari   CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT);
6172b730f8bSJeremy L Thompson   CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", q_data_size_sur, CEED_EVAL_NONE);
61877841947SLeila Ghaffari 
6192b730f8bSJeremy L Thompson   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow,
6202b730f8bSJeremy L Thompson                               problem->apply_inflow_jacobian, &ceed_data->qf_apply_inflow, &ceed_data->qf_apply_inflow_jacobian));
6212b730f8bSJeremy L Thompson   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow,
6222b730f8bSJeremy L Thompson                               problem->apply_outflow_jacobian, &ceed_data->qf_apply_outflow, &ceed_data->qf_apply_outflow_jacobian));
6232b730f8bSJeremy L Thompson   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream,
6242b730f8bSJeremy L Thompson                               problem->apply_freestream_jacobian, &ceed_data->qf_apply_freestream, &ceed_data->qf_apply_freestream_jacobian));
62577841947SLeila Ghaffari 
62677841947SLeila Ghaffari   // *****************************************************************************
62777841947SLeila Ghaffari   // CEED Operator Apply
62877841947SLeila Ghaffari   // *****************************************************************************
62977841947SLeila Ghaffari   // -- Apply CEED Operator for the geometric data
6302b730f8bSJeremy L Thompson   CeedOperatorApply(ceed_data->op_setup_vol, ceed_data->x_coord, ceed_data->q_data, CEED_REQUEST_IMMEDIATE);
63177841947SLeila Ghaffari 
63277841947SLeila Ghaffari   // -- Create and apply CEED Composite Operator for the entire domain
63377841947SLeila Ghaffari   if (!user->phys->implicit) {  // RHS
6349ad5e8e4SJames Wright     CeedOperator op_rhs;
6359ad5e8e4SJames 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,
6369ad5e8e4SJames Wright                                       NULL));
6379ad5e8e4SJames Wright     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx));
6389ad5e8e4SJames Wright     CeedOperatorDestroy(&op_rhs);
63977841947SLeila Ghaffari   } else {  // IFunction
6402b730f8bSJeremy L Thompson     PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_ifunction_vol, op_ijacobian_vol, height, P_sur, Q_sur,
6412b730f8bSJeremy L Thompson                                       q_data_size_sur, jac_data_size_sur, &user->op_ifunction, op_ijacobian_vol ? &user->op_ijacobian : NULL));
642e334ad8fSJed Brown     if (user->op_ijacobian) {
64317b0d5c6SJeremy L Thompson       CeedOperatorGetContextFieldLabel(user->op_ijacobian, "ijacobian time shift", &user->phys->ijacobian_time_shift_label);
644e334ad8fSJed Brown     }
645*0814d5a7SKenneth E. Jansen     if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc, q_data_size_sur));
646f5452247SJames Wright     if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SGS_DD_ModelSetup(ceed, user, ceed_data, problem));
647dada6cc0SJames Wright   }
648f5452247SJames Wright 
649f5452247SJames Wright   if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem));
6504e9802d1SJames Wright   if (app_ctx->diff_filter_monitor) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem));
651f5452247SJames Wright 
652a3ae0734SJed Brown   CeedElemRestrictionDestroy(&elem_restr_jd_i);
65322d320f5SLeila Ghaffari   CeedOperatorDestroy(&op_ijacobian_vol);
654a3ae0734SJed Brown   CeedVectorDestroy(&jac_data);
655ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
65677841947SLeila Ghaffari }
657