xref: /libCEED/examples/fluids/src/setuplibceed.c (revision a424bcd0dd58ea8f8b80ddde5af211268f0b358a)
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 
140814d5a7SKenneth 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));
28*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreate(ceed, num_elem, elem_size, num_comp, 1, num_dof, CEED_MEM_HOST, CEED_COPY_VALUES,
29*a424bcd0SJames Wright                                                 elem_restr_offsets_ceed, elem_restr));
30457e73b2SJames Wright   PetscCall(PetscFree(elem_restr_offsets_ceed));
3177841947SLeila Ghaffari 
32ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
3377841947SLeila Ghaffari }
3477841947SLeila Ghaffari 
3577841947SLeila Ghaffari // Utility function to get Ceed Restriction for each domain
36431cd09aSJames Wright PetscErrorCode GetRestrictionForDomain(Ceed ceed, DM dm, CeedInt height, DMLabel domain_label, PetscInt label_value, PetscInt dm_field, CeedInt Q,
37431cd09aSJames Wright                                        CeedInt q_data_size, CeedElemRestriction *elem_restr_q, CeedElemRestriction *elem_restr_x,
38431cd09aSJames Wright                                        CeedElemRestriction *elem_restr_qd_i) {
3977841947SLeila Ghaffari   DM                  dm_coord;
400814d5a7SKenneth E. Jansen   CeedInt             loc_num_elem;
41457e73b2SJames Wright   PetscInt            dim;
422534dcc8SJed Brown   CeedElemRestriction elem_restr_tmp;
4377841947SLeila Ghaffari   PetscFunctionBeginUser;
4477841947SLeila Ghaffari 
452b730f8bSJeremy L Thompson   PetscCall(DMGetDimension(dm, &dim));
4677841947SLeila Ghaffari   dim -= height;
47431cd09aSJames Wright   PetscCall(CreateRestrictionFromPlex(ceed, dm, height, domain_label, label_value, dm_field, &elem_restr_tmp));
482534dcc8SJed Brown   if (elem_restr_q) *elem_restr_q = elem_restr_tmp;
492534dcc8SJed Brown   if (elem_restr_x) {
502b730f8bSJeremy L Thompson     PetscCall(DMGetCellCoordinateDM(dm, &dm_coord));
513796c488SJed Brown     if (!dm_coord) {
522b730f8bSJeremy L Thompson       PetscCall(DMGetCoordinateDM(dm, &dm_coord));
533796c488SJed Brown     }
542b730f8bSJeremy L Thompson     PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL));
55431cd09aSJames Wright     PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, height, domain_label, label_value, dm_field, elem_restr_x));
562534dcc8SJed Brown   }
572534dcc8SJed Brown   if (elem_restr_qd_i) {
58*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_tmp, &loc_num_elem));
59*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateStrided(ceed, loc_num_elem, Q, q_data_size, q_data_size * loc_num_elem * Q, CEED_STRIDES_BACKEND,
60*a424bcd0SJames Wright                                                          elem_restr_qd_i));
612534dcc8SJed Brown   }
62*a424bcd0SJames Wright   if (!elem_restr_q) PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_tmp));
63ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
6477841947SLeila Ghaffari }
6577841947SLeila Ghaffari 
66431cd09aSJames Wright PetscErrorCode AddBCSubOperator(Ceed ceed, DM dm, CeedData ceed_data, DMLabel domain_label, PetscInt label_value, CeedInt height, CeedInt Q_sur,
672b730f8bSJeremy L Thompson                                 CeedInt q_data_size_sur, CeedInt jac_data_size_sur, CeedQFunction qf_apply_bc, CeedQFunction qf_apply_bc_jacobian,
689eb9c072SJames Wright                                 CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) {
699eb9c072SJames Wright   CeedVector          q_data_sur, jac_data_sur;
702b730f8bSJeremy L Thompson   CeedOperator        op_setup_sur, op_apply_bc, op_apply_bc_jacobian = NULL;
712b730f8bSJeremy L Thompson   CeedElemRestriction elem_restr_x_sur, elem_restr_q_sur, elem_restr_qd_i_sur, elem_restr_jd_i_sur;
729eb9c072SJames Wright   CeedInt             num_qpts_sur;
739eb9c072SJames Wright   PetscFunctionBeginUser;
749eb9c072SJames Wright 
759eb9c072SJames Wright   // --- Get number of quadrature points for the boundaries
76*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(ceed_data->basis_q_sur, &num_qpts_sur));
779eb9c072SJames Wright 
789eb9c072SJames Wright   // ---- CEED Restriction
790814d5a7SKenneth E. Jansen   PetscCall(GetRestrictionForDomain(ceed, dm, height, domain_label, label_value, 0, num_qpts_sur, q_data_size_sur, &elem_restr_q_sur,
800814d5a7SKenneth E. Jansen                                     &elem_restr_x_sur, &elem_restr_qd_i_sur));
819eb9c072SJames Wright   if (jac_data_size_sur > 0) {
829eb9c072SJames Wright     // State-dependent data will be passed from residual to Jacobian. This will be collocated.
830814d5a7SKenneth E. Jansen     PetscCall(
840814d5a7SKenneth E. Jansen         GetRestrictionForDomain(ceed, dm, height, domain_label, label_value, 0, num_qpts_sur, jac_data_size_sur, NULL, NULL, &elem_restr_jd_i_sur));
85*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i_sur, &jac_data_sur, NULL));
869eb9c072SJames Wright   } else {
879eb9c072SJames Wright     elem_restr_jd_i_sur = NULL;
889eb9c072SJames Wright     jac_data_sur        = NULL;
899eb9c072SJames Wright   }
909eb9c072SJames Wright 
919eb9c072SJames Wright   // ---- CEED Vector
92457e73b2SJames Wright   CeedInt loc_num_elem_sur;
93*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionGetNumElements(elem_restr_q_sur, &loc_num_elem_sur));
94*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorCreate(ceed, q_data_size_sur * loc_num_elem_sur * num_qpts_sur, &q_data_sur));
959eb9c072SJames Wright 
969eb9c072SJames Wright   // ---- CEED Operator
979eb9c072SJames Wright   // ----- CEED Operator for Setup (geometric factors)
98*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_sur, NULL, NULL, &op_setup_sur));
99*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "dx", elem_restr_x_sur, ceed_data->basis_x_sur, CEED_VECTOR_ACTIVE));
100*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x_sur, CEED_VECTOR_NONE));
101*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_setup_sur, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
1029eb9c072SJames Wright 
1039eb9c072SJames Wright   // ----- CEED Operator for Physics
104*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc, NULL, NULL, &op_apply_bc));
105*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
106*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "Grad_q", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
107*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur));
108*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord));
109*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
110*a424bcd0SJames Wright   if (elem_restr_jd_i_sur)
111*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur));
1129eb9c072SJames Wright 
1139eb9c072SJames Wright   if (qf_apply_bc_jacobian) {
114*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_apply_bc_jacobian, NULL, NULL, &op_apply_bc_jacobian));
115*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
116*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "Grad_dq", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
117*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "surface qdata", elem_restr_qd_i_sur, CEED_BASIS_COLLOCATED, q_data_sur));
118*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "x", elem_restr_x_sur, ceed_data->basis_x_sur, ceed_data->x_coord));
119*a424bcd0SJames Wright     PetscCallCeed(ceed,
120*a424bcd0SJames Wright                   CeedOperatorSetField(op_apply_bc_jacobian, "surface jacobian data", elem_restr_jd_i_sur, CEED_BASIS_COLLOCATED, jac_data_sur));
121*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op_apply_bc_jacobian, "v", elem_restr_q_sur, ceed_data->basis_q_sur, CEED_VECTOR_ACTIVE));
1229eb9c072SJames Wright   }
1239eb9c072SJames Wright 
1249eb9c072SJames Wright   // ----- Apply CEED operator for Setup
125*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorApply(op_setup_sur, ceed_data->x_coord, q_data_sur, CEED_REQUEST_IMMEDIATE));
1269eb9c072SJames Wright 
1279eb9c072SJames Wright   // ----- Apply Sub-Operator for Physics
128*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply, op_apply_bc));
129*a424bcd0SJames Wright   if (op_apply_bc_jacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_bc_jacobian));
1309eb9c072SJames Wright 
1319eb9c072SJames Wright   // ----- Cleanup
132*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&q_data_sur));
133*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data_sur));
134*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_q_sur));
135*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_x_sur));
136*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_qd_i_sur));
137*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i_sur));
138*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_setup_sur));
139*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc));
140*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_apply_bc_jacobian));
1419eb9c072SJames Wright 
142ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
1439eb9c072SJames Wright }
1449eb9c072SJames Wright 
14577841947SLeila Ghaffari // Utility function to create CEED Composite Operator for the entire domain
1462b730f8bSJeremy L Thompson PetscErrorCode CreateOperatorForDomain(Ceed ceed, DM dm, SimpleBC bc, CeedData ceed_data, Physics phys, CeedOperator op_apply_vol,
1472b730f8bSJeremy L Thompson                                        CeedOperator op_apply_ijacobian_vol, CeedInt height, CeedInt P_sur, CeedInt Q_sur, CeedInt q_data_size_sur,
1482b730f8bSJeremy L Thompson                                        CeedInt jac_data_size_sur, CeedOperator *op_apply, CeedOperator *op_apply_ijacobian) {
14977841947SLeila Ghaffari   DMLabel domain_label;
15077841947SLeila Ghaffari   PetscFunctionBeginUser;
15177841947SLeila Ghaffari 
15277841947SLeila Ghaffari   // Create Composite Operaters
153*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_apply));
154*a424bcd0SJames Wright   if (op_apply_ijacobian) PetscCallCeed(ceed, CeedCompositeOperatorCreate(ceed, op_apply_ijacobian));
15577841947SLeila Ghaffari 
15677841947SLeila Ghaffari   // --Apply Sub-Operator for the volume
157*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply, op_apply_vol));
158*a424bcd0SJames Wright   if (op_apply_ijacobian) PetscCallCeed(ceed, CeedCompositeOperatorAddSub(*op_apply_ijacobian, op_apply_ijacobian_vol));
15977841947SLeila Ghaffari 
16077841947SLeila Ghaffari   // -- Create Sub-Operator for in/outflow BCs
16188626eedSJames Wright   if (phys->has_neumann || 1) {
16277841947SLeila Ghaffari     // --- Setup
1632b730f8bSJeremy L Thompson     PetscCall(DMGetLabel(dm, "Face Sets", &domain_label));
16477841947SLeila Ghaffari 
1652fe7aee7SLeila Ghaffari     // --- Create Sub-Operator for inflow boundaries
1662fe7aee7SLeila Ghaffari     for (CeedInt i = 0; i < bc->num_inflow; i++) {
1672b730f8bSJeremy L Thompson       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->inflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
1682b730f8bSJeremy L Thompson                                  ceed_data->qf_apply_inflow, ceed_data->qf_apply_inflow_jacobian, op_apply, op_apply_ijacobian));
16977841947SLeila Ghaffari     }
1702fe7aee7SLeila Ghaffari     // --- Create Sub-Operator for outflow boundaries
1712fe7aee7SLeila Ghaffari     for (CeedInt i = 0; i < bc->num_outflow; i++) {
1722b730f8bSJeremy L Thompson       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->outflows[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
1732b730f8bSJeremy L Thompson                                  ceed_data->qf_apply_outflow, ceed_data->qf_apply_outflow_jacobian, op_apply, op_apply_ijacobian));
17439c69132SJed Brown     }
175f4ca79c2SJames Wright     // --- Create Sub-Operator for freestream boundaries
176f4ca79c2SJames Wright     for (CeedInt i = 0; i < bc->num_freestream; i++) {
1772b730f8bSJeremy L Thompson       PetscCall(AddBCSubOperator(ceed, dm, ceed_data, domain_label, bc->freestreams[i], height, Q_sur, q_data_size_sur, jac_data_size_sur,
1782b730f8bSJeremy L Thompson                                  ceed_data->qf_apply_freestream, ceed_data->qf_apply_freestream_jacobian, op_apply, op_apply_ijacobian));
1792fe7aee7SLeila Ghaffari     }
1802fe7aee7SLeila Ghaffari   }
18111436a05SJames Wright 
18211436a05SJames Wright   // ----- Get Context Labels for Operator
183*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(*op_apply, "solution time", &phys->solution_time_label));
184*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(*op_apply, "timestep size", &phys->timestep_size_label));
18511436a05SJames Wright 
186ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
18777841947SLeila Ghaffari }
18877841947SLeila Ghaffari 
1892b730f8bSJeremy L Thompson PetscErrorCode SetupBCQFunctions(Ceed ceed, PetscInt dim_sur, PetscInt num_comp_x, PetscInt num_comp_q, PetscInt q_data_size_sur,
1902b730f8bSJeremy L Thompson                                  PetscInt jac_data_size_sur, ProblemQFunctionSpec apply_bc, ProblemQFunctionSpec apply_bc_jacobian,
1915b881d11SJames Wright                                  CeedQFunction *qf_apply_bc, CeedQFunction *qf_apply_bc_jacobian) {
1925b881d11SJames Wright   PetscFunctionBeginUser;
1935b881d11SJames Wright 
1945b881d11SJames Wright   if (apply_bc.qfunction) {
195*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc.qfunction, apply_bc.qfunction_loc, qf_apply_bc));
196*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc, apply_bc.qfunction_context));
197*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "q", num_comp_q, CEED_EVAL_INTERP));
198*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "Grad_q", num_comp_q * dim_sur, CEED_EVAL_GRAD));
199*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
200*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc, "x", num_comp_x, CEED_EVAL_INTERP));
201*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "v", num_comp_q, CEED_EVAL_INTERP));
202*a424bcd0SJames Wright     if (jac_data_size_sur) PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
2035b881d11SJames Wright   }
2045b881d11SJames Wright   if (apply_bc_jacobian.qfunction) {
205*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, apply_bc_jacobian.qfunction, apply_bc_jacobian.qfunction_loc, qf_apply_bc_jacobian));
206*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(*qf_apply_bc_jacobian, apply_bc_jacobian.qfunction_context));
207*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "dq", num_comp_q, CEED_EVAL_INTERP));
208*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "Grad_dq", num_comp_q * dim_sur, CEED_EVAL_GRAD));
209*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
210*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "x", num_comp_x, CEED_EVAL_INTERP));
211*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(*qf_apply_bc_jacobian, "surface jacobian data", jac_data_size_sur, CEED_EVAL_NONE));
212*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(*qf_apply_bc_jacobian, "v", num_comp_q, CEED_EVAL_INTERP));
2135b881d11SJames Wright   }
214ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
2155b881d11SJames Wright }
2165b881d11SJames Wright 
2170814d5a7SKenneth E. Jansen // -----------------------------------------------------------------------------
2180814d5a7SKenneth E. Jansen // Convert DM field to DS field
2190814d5a7SKenneth E. Jansen // -----------------------------------------------------------------------------
2200814d5a7SKenneth E. Jansen PetscErrorCode DMFieldToDSField(DM dm, DMLabel domain_label, PetscInt dm_field, PetscInt *ds_field) {
2210814d5a7SKenneth E. Jansen   PetscDS         ds;
2220814d5a7SKenneth E. Jansen   IS              field_is;
2230814d5a7SKenneth E. Jansen   const PetscInt *fields;
2240814d5a7SKenneth E. Jansen   PetscInt        num_fields;
2250814d5a7SKenneth E. Jansen 
2260814d5a7SKenneth E. Jansen   PetscFunctionBeginUser;
2270814d5a7SKenneth E. Jansen 
2280814d5a7SKenneth E. Jansen   // Translate dm_field to ds_field
2290814d5a7SKenneth E. Jansen   PetscCall(DMGetRegionDS(dm, domain_label, &field_is, &ds, NULL));
2300814d5a7SKenneth E. Jansen   PetscCall(ISGetIndices(field_is, &fields));
2310814d5a7SKenneth E. Jansen   PetscCall(ISGetSize(field_is, &num_fields));
2320814d5a7SKenneth E. Jansen   for (PetscInt i = 0; i < num_fields; i++) {
2330814d5a7SKenneth E. Jansen     if (dm_field == fields[i]) {
2340814d5a7SKenneth E. Jansen       *ds_field = i;
2350814d5a7SKenneth E. Jansen       break;
2360814d5a7SKenneth E. Jansen     }
2370814d5a7SKenneth E. Jansen   }
2380814d5a7SKenneth E. Jansen   PetscCall(ISRestoreIndices(field_is, &fields));
2390814d5a7SKenneth E. Jansen 
2400814d5a7SKenneth E. Jansen   if (*ds_field == -1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Could not find dm_field %" PetscInt_FMT " in DS", dm_field);
2410814d5a7SKenneth E. Jansen 
2420814d5a7SKenneth E. Jansen   PetscFunctionReturn(PETSC_SUCCESS);
2430814d5a7SKenneth E. Jansen }
2440814d5a7SKenneth E. Jansen 
2450814d5a7SKenneth E. Jansen // -----------------------------------------------------------------------------
2460814d5a7SKenneth E. Jansen // Utility function - convert from DMPolytopeType to CeedElemTopology
2470814d5a7SKenneth E. Jansen // -----------------------------------------------------------------------------
2480814d5a7SKenneth E. Jansen CeedElemTopology ElemTopologyP2C(DMPolytopeType cell_type) {
2490814d5a7SKenneth E. Jansen   switch (cell_type) {
2500814d5a7SKenneth E. Jansen     case DM_POLYTOPE_TRIANGLE:
2510814d5a7SKenneth E. Jansen       return CEED_TOPOLOGY_TRIANGLE;
2520814d5a7SKenneth E. Jansen     case DM_POLYTOPE_QUADRILATERAL:
2530814d5a7SKenneth E. Jansen       return CEED_TOPOLOGY_QUAD;
2540814d5a7SKenneth E. Jansen     case DM_POLYTOPE_TETRAHEDRON:
2550814d5a7SKenneth E. Jansen       return CEED_TOPOLOGY_TET;
2560814d5a7SKenneth E. Jansen     case DM_POLYTOPE_HEXAHEDRON:
2570814d5a7SKenneth E. Jansen       return CEED_TOPOLOGY_HEX;
2580814d5a7SKenneth E. Jansen     default:
2590814d5a7SKenneth E. Jansen       return 0;
2600814d5a7SKenneth E. Jansen   }
2610814d5a7SKenneth E. Jansen }
2620814d5a7SKenneth E. Jansen 
2630814d5a7SKenneth E. Jansen // -----------------------------------------------------------------------------
2640814d5a7SKenneth E. Jansen // Create libCEED Basis from PetscTabulation
2650814d5a7SKenneth E. Jansen // -----------------------------------------------------------------------------
2660814d5a7SKenneth E. Jansen PetscErrorCode BasisCreateFromTabulation(Ceed ceed, DM dm, DMLabel domain_label, PetscInt label_value, PetscInt height, PetscInt face, PetscFE fe,
2670814d5a7SKenneth E. Jansen                                          PetscTabulation basis_tabulation, PetscQuadrature quadrature, CeedBasis *basis) {
2680814d5a7SKenneth E. Jansen   PetscInt           first_point;
2690814d5a7SKenneth E. Jansen   PetscInt           ids[1] = {label_value};
2700814d5a7SKenneth E. Jansen   DMLabel            depth_label;
2710814d5a7SKenneth E. Jansen   DMPolytopeType     cell_type;
2720814d5a7SKenneth E. Jansen   CeedElemTopology   elem_topo;
2730814d5a7SKenneth E. Jansen   PetscScalar       *q_points, *interp, *grad;
2740814d5a7SKenneth E. Jansen   const PetscScalar *q_weights;
2750814d5a7SKenneth E. Jansen   PetscDualSpace     dual_space;
2760814d5a7SKenneth E. Jansen   PetscInt           num_dual_basis_vectors;
2770814d5a7SKenneth E. Jansen   PetscInt           dim, num_comp, P, Q;
2780814d5a7SKenneth E. Jansen 
2790814d5a7SKenneth E. Jansen   PetscFunctionBeginUser;
2800814d5a7SKenneth E. Jansen 
2810814d5a7SKenneth E. Jansen   // General basis information
2820814d5a7SKenneth E. Jansen   PetscCall(PetscFEGetSpatialDimension(fe, &dim));
2830814d5a7SKenneth E. Jansen   PetscCall(PetscFEGetNumComponents(fe, &num_comp));
2840814d5a7SKenneth E. Jansen   PetscCall(PetscFEGetDualSpace(fe, &dual_space));
2850814d5a7SKenneth E. Jansen   PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
2860814d5a7SKenneth E. Jansen   P = num_dual_basis_vectors / num_comp;
2870814d5a7SKenneth E. Jansen 
2880814d5a7SKenneth E. Jansen   // Use depth label if no domain label present
2890814d5a7SKenneth E. Jansen   if (!domain_label) {
2900814d5a7SKenneth E. Jansen     PetscInt depth;
2910814d5a7SKenneth E. Jansen 
2920814d5a7SKenneth E. Jansen     PetscCall(DMPlexGetDepth(dm, &depth));
2930814d5a7SKenneth E. Jansen     PetscCall(DMPlexGetDepthLabel(dm, &depth_label));
2940814d5a7SKenneth E. Jansen     ids[0] = depth - height;
2950814d5a7SKenneth E. Jansen   }
2960814d5a7SKenneth E. Jansen 
2970814d5a7SKenneth E. Jansen   // Get cell interp, grad, and quadrature data
2980814d5a7SKenneth E. Jansen   PetscCall(DMGetFirstLabeledPoint(dm, dm, domain_label ? domain_label : depth_label, 1, ids, height, &first_point, NULL));
2990814d5a7SKenneth E. Jansen   PetscCall(DMPlexGetCellType(dm, first_point, &cell_type));
3000814d5a7SKenneth E. Jansen   elem_topo = ElemTopologyP2C(cell_type);
3010814d5a7SKenneth E. Jansen   if (!elem_topo) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMPlex topology not supported");
3020814d5a7SKenneth E. Jansen   {
3030814d5a7SKenneth E. Jansen     size_t             q_points_size;
3040814d5a7SKenneth E. Jansen     const PetscScalar *q_points_petsc;
3050814d5a7SKenneth E. Jansen     PetscInt           q_dim;
3060814d5a7SKenneth E. Jansen 
3070814d5a7SKenneth E. Jansen     PetscCall(PetscQuadratureGetData(quadrature, &q_dim, NULL, &Q, &q_points_petsc, &q_weights));
3080814d5a7SKenneth E. Jansen     q_points_size = Q * dim * sizeof(CeedScalar);
3090814d5a7SKenneth E. Jansen     PetscCall(PetscCalloc(q_points_size, &q_points));
3100814d5a7SKenneth E. Jansen     for (PetscInt q = 0; q < Q; q++) {
3110814d5a7SKenneth E. Jansen       for (PetscInt d = 0; d < q_dim; d++) q_points[q * dim + d] = q_points_petsc[q * q_dim + d];
3120814d5a7SKenneth E. Jansen     }
3130814d5a7SKenneth E. Jansen   }
3140814d5a7SKenneth E. Jansen 
3150814d5a7SKenneth E. Jansen   {  // Convert to libCEED orientation
3160814d5a7SKenneth E. Jansen     PetscBool       is_simplex  = PETSC_FALSE;
3170814d5a7SKenneth E. Jansen     IS              permutation = NULL;
3180814d5a7SKenneth E. Jansen     const PetscInt *permutation_indices;
3190814d5a7SKenneth E. Jansen 
3200814d5a7SKenneth E. Jansen     PetscCall(DMPlexIsSimplex(dm, &is_simplex));
3210814d5a7SKenneth E. Jansen     if (!is_simplex) {
3220814d5a7SKenneth E. Jansen       PetscSection section;
3230814d5a7SKenneth E. Jansen 
3240814d5a7SKenneth E. Jansen       // -- Get permutation
3250814d5a7SKenneth E. Jansen       PetscCall(DMGetLocalSection(dm, &section));
3260814d5a7SKenneth E. Jansen       PetscCall(PetscSectionGetClosurePermutation(section, (PetscObject)dm, dim, num_comp * P, &permutation));
3270814d5a7SKenneth E. Jansen       PetscCall(ISGetIndices(permutation, &permutation_indices));
3280814d5a7SKenneth E. Jansen     }
3290814d5a7SKenneth E. Jansen 
3300814d5a7SKenneth E. Jansen     // -- Copy interp, grad matrices
3310814d5a7SKenneth E. Jansen     PetscCall(PetscCalloc(P * Q * sizeof(CeedScalar), &interp));
3320814d5a7SKenneth E. Jansen     PetscCall(PetscCalloc(P * Q * dim * sizeof(CeedScalar), &grad));
3330814d5a7SKenneth E. Jansen     const CeedInt c = 0;
3340814d5a7SKenneth E. Jansen     for (CeedInt q = 0; q < Q; q++) {
3350814d5a7SKenneth E. Jansen       for (CeedInt p_ceed = 0; p_ceed < P; p_ceed++) {
3360814d5a7SKenneth E. Jansen         CeedInt p_petsc = is_simplex ? (p_ceed * num_comp) : permutation_indices[p_ceed * num_comp];
3370814d5a7SKenneth E. Jansen 
3380814d5a7SKenneth E. Jansen         interp[q * P + p_ceed] = basis_tabulation->T[0][((face * Q + q) * P * num_comp + p_petsc) * num_comp + c];
3390814d5a7SKenneth E. Jansen         for (CeedInt d = 0; d < dim; d++) {
3400814d5a7SKenneth 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];
3410814d5a7SKenneth E. Jansen         }
3420814d5a7SKenneth E. Jansen       }
3430814d5a7SKenneth E. Jansen     }
3440814d5a7SKenneth E. Jansen 
3450814d5a7SKenneth E. Jansen     // -- Cleanup
3460814d5a7SKenneth E. Jansen     if (permutation) PetscCall(ISRestoreIndices(permutation, &permutation_indices));
3470814d5a7SKenneth E. Jansen     PetscCall(ISDestroy(&permutation));
3480814d5a7SKenneth E. Jansen   }
3490814d5a7SKenneth E. Jansen 
3500814d5a7SKenneth E. Jansen   // Finally, create libCEED basis
351*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisCreateH1(ceed, elem_topo, num_comp, P, Q, interp, grad, q_points, q_weights, basis));
3520814d5a7SKenneth E. Jansen   PetscCall(PetscFree(q_points));
3530814d5a7SKenneth E. Jansen   PetscCall(PetscFree(interp));
3540814d5a7SKenneth E. Jansen   PetscCall(PetscFree(grad));
3550814d5a7SKenneth E. Jansen 
3560814d5a7SKenneth E. Jansen   PetscFunctionReturn(PETSC_SUCCESS);
3570814d5a7SKenneth E. Jansen }
3580814d5a7SKenneth E. Jansen 
3590814d5a7SKenneth E. Jansen // -----------------------------------------------------------------------------
3600814d5a7SKenneth E. Jansen // Get CEED Basis from DMPlex
3610814d5a7SKenneth E. Jansen // -----------------------------------------------------------------------------
3620814d5a7SKenneth E. Jansen PetscErrorCode CreateBasisFromPlex(Ceed ceed, DM dm, DMLabel domain_label, CeedInt label_value, CeedInt height, CeedInt dm_field, CeedBasis *basis) {
3630814d5a7SKenneth E. Jansen   PetscDS         ds;
3640814d5a7SKenneth E. Jansen   PetscFE         fe;
3650814d5a7SKenneth E. Jansen   PetscQuadrature quadrature;
3660814d5a7SKenneth E. Jansen   PetscBool       is_simplex = PETSC_TRUE;
3670814d5a7SKenneth E. Jansen   PetscInt        ds_field   = -1;
3680814d5a7SKenneth E. Jansen 
3690814d5a7SKenneth E. Jansen   PetscFunctionBeginUser;
3700814d5a7SKenneth E. Jansen 
3710814d5a7SKenneth E. Jansen   // Get element information
3720814d5a7SKenneth E. Jansen   PetscCall(DMGetRegionDS(dm, domain_label, NULL, &ds, NULL));
3730814d5a7SKenneth E. Jansen   PetscCall(DMFieldToDSField(dm, domain_label, dm_field, &ds_field));
3740814d5a7SKenneth E. Jansen   PetscCall(PetscDSGetDiscretization(ds, ds_field, (PetscObject *)&fe));
3750814d5a7SKenneth E. Jansen   PetscCall(PetscFEGetHeightSubspace(fe, height, &fe));
3760814d5a7SKenneth E. Jansen   PetscCall(PetscFEGetQuadrature(fe, &quadrature));
3770814d5a7SKenneth E. Jansen 
3780814d5a7SKenneth E. Jansen   // Check if simplex or tensor-product mesh
3790814d5a7SKenneth E. Jansen   PetscCall(DMPlexIsSimplex(dm, &is_simplex));
3800814d5a7SKenneth E. Jansen 
3810814d5a7SKenneth E. Jansen   // Build libCEED basis
3820814d5a7SKenneth E. Jansen   if (is_simplex) {
3830814d5a7SKenneth E. Jansen     PetscTabulation basis_tabulation;
3840814d5a7SKenneth E. Jansen     PetscInt        num_derivatives = 1, face = 0;
3850814d5a7SKenneth E. Jansen 
3860814d5a7SKenneth E. Jansen     PetscCall(PetscFEGetCellTabulation(fe, num_derivatives, &basis_tabulation));
3870814d5a7SKenneth E. Jansen     PetscCall(BasisCreateFromTabulation(ceed, dm, domain_label, label_value, height, face, fe, basis_tabulation, quadrature, basis));
3880814d5a7SKenneth E. Jansen   } else {
3890814d5a7SKenneth E. Jansen     PetscDualSpace dual_space;
3900814d5a7SKenneth E. Jansen     PetscInt       num_dual_basis_vectors;
3910814d5a7SKenneth E. Jansen     PetscInt       dim, num_comp, P, Q;
3920814d5a7SKenneth E. Jansen 
3930814d5a7SKenneth E. Jansen     PetscCall(PetscFEGetSpatialDimension(fe, &dim));
3940814d5a7SKenneth E. Jansen     PetscCall(PetscFEGetNumComponents(fe, &num_comp));
3950814d5a7SKenneth E. Jansen     PetscCall(PetscFEGetDualSpace(fe, &dual_space));
3960814d5a7SKenneth E. Jansen     PetscCall(PetscDualSpaceGetDimension(dual_space, &num_dual_basis_vectors));
3970814d5a7SKenneth E. Jansen     P = num_dual_basis_vectors / num_comp;
3980814d5a7SKenneth E. Jansen     PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &Q, NULL, NULL));
3990814d5a7SKenneth E. Jansen 
4000814d5a7SKenneth E. Jansen     CeedInt P_1d = (CeedInt)round(pow(P, 1.0 / dim));
4010814d5a7SKenneth E. Jansen     CeedInt Q_1d = (CeedInt)round(pow(Q, 1.0 / dim));
4020814d5a7SKenneth E. Jansen 
403*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp, P_1d, Q_1d, CEED_GAUSS, basis));
4040814d5a7SKenneth E. Jansen   }
4050814d5a7SKenneth E. Jansen 
4060814d5a7SKenneth E. Jansen   PetscFunctionReturn(PETSC_SUCCESS);
4070814d5a7SKenneth E. Jansen }
4080814d5a7SKenneth E. Jansen 
4092b730f8bSJeremy L Thompson PetscErrorCode SetupLibceed(Ceed ceed, CeedData ceed_data, DM dm, User user, AppCtx app_ctx, ProblemData *problem, SimpleBC bc) {
41077841947SLeila Ghaffari   PetscFunctionBeginUser;
41177841947SLeila Ghaffari 
41277841947SLeila Ghaffari   // *****************************************************************************
41377841947SLeila Ghaffari   // Set up CEED objects for the interior domain (volume)
41477841947SLeila Ghaffari   // *****************************************************************************
41577841947SLeila Ghaffari   const PetscInt num_comp_q = 5;
4160814d5a7SKenneth 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;
417a3ae0734SJed Brown   CeedElemRestriction elem_restr_jd_i;
418a3ae0734SJed Brown   CeedVector          jac_data;
4190814d5a7SKenneth E. Jansen   CeedInt             num_qpts;
42077841947SLeila Ghaffari 
42177841947SLeila Ghaffari   // -----------------------------------------------------------------------------
42277841947SLeila Ghaffari   // CEED Bases
42377841947SLeila Ghaffari   // -----------------------------------------------------------------------------
4240814d5a7SKenneth E. Jansen   DM dm_coord;
4250814d5a7SKenneth E. Jansen   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
4260814d5a7SKenneth E. Jansen 
4270814d5a7SKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, &ceed_data->basis_q));
4280814d5a7SKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, dm_coord, 0, 0, 0, 0, &ceed_data->basis_x));
429*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x, ceed_data->basis_q, &ceed_data->basis_xc));
430*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisGetNumQuadraturePoints(ceed_data->basis_q, &num_qpts));
43177841947SLeila Ghaffari 
43277841947SLeila Ghaffari   // -----------------------------------------------------------------------------
43377841947SLeila Ghaffari   // CEED Restrictions
43477841947SLeila Ghaffari   // -----------------------------------------------------------------------------
43577841947SLeila Ghaffari   // -- Create restriction
4360814d5a7SKenneth 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,
4372b730f8bSJeremy L Thompson                                     &ceed_data->elem_restr_qd_i));
438a3ae0734SJed Brown 
4390814d5a7SKenneth E. Jansen   PetscCall(GetRestrictionForDomain(ceed, dm, 0, 0, 0, 0, num_qpts, jac_data_size_vol, NULL, NULL, &elem_restr_jd_i));
44077841947SLeila Ghaffari   // -- Create E vectors
441*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_ceed, NULL));
442*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->q_dot_ceed, NULL));
443*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_q, &user->g_ceed, NULL));
44477841947SLeila Ghaffari 
44577841947SLeila Ghaffari   // -----------------------------------------------------------------------------
44677841947SLeila Ghaffari   // CEED QFunctions
44777841947SLeila Ghaffari   // -----------------------------------------------------------------------------
44877841947SLeila Ghaffari   // -- Create QFunction for quadrature data
449*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->setup_vol.qfunction, problem->setup_vol.qfunction_loc, &ceed_data->qf_setup_vol));
450841e4c73SJed Brown   if (problem->setup_vol.qfunction_context) {
451*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_setup_vol, problem->setup_vol.qfunction_context));
452841e4c73SJed Brown   }
453*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_vol, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
454*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_vol, "weight", 1, CEED_EVAL_WEIGHT));
455*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_setup_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
45677841947SLeila Ghaffari 
45777841947SLeila Ghaffari   // -- Create QFunction for ICs
458*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->ics.qfunction, problem->ics.qfunction_loc, &ceed_data->qf_ics));
459*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_ics, problem->ics.qfunction_context));
460*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ics, "x", num_comp_x, CEED_EVAL_INTERP));
461*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ics, "dx", num_comp_x * dim, CEED_EVAL_GRAD));
462*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_ics, "q0", num_comp_q, CEED_EVAL_NONE));
46377841947SLeila Ghaffari 
46477841947SLeila Ghaffari   // -- Create QFunction for RHS
46591e5af17SJed Brown   if (problem->apply_vol_rhs.qfunction) {
466*a424bcd0SJames Wright     PetscCallCeed(
467*a424bcd0SJames Wright         ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_rhs.qfunction, problem->apply_vol_rhs.qfunction_loc, &ceed_data->qf_rhs_vol));
468*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_rhs_vol, problem->apply_vol_rhs.qfunction_context));
469*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "q", num_comp_q, CEED_EVAL_INTERP));
470*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
471*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
472*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_rhs_vol, "x", num_comp_x, CEED_EVAL_INTERP));
473*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "v", num_comp_q, CEED_EVAL_INTERP));
474*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_rhs_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
47577841947SLeila Ghaffari   }
47677841947SLeila Ghaffari 
47777841947SLeila Ghaffari   // -- Create QFunction for IFunction
47891e5af17SJed Brown   if (problem->apply_vol_ifunction.qfunction) {
479*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ifunction.qfunction, problem->apply_vol_ifunction.qfunction_loc,
480*a424bcd0SJames Wright                                                     &ceed_data->qf_ifunction_vol));
481*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_ifunction_vol, problem->apply_vol_ifunction.qfunction_context));
482*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q", num_comp_q, CEED_EVAL_INTERP));
483*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "Grad_q", num_comp_q * dim, CEED_EVAL_GRAD));
484*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "q dot", num_comp_q, CEED_EVAL_INTERP));
485*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
486*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_ifunction_vol, "x", num_comp_x, CEED_EVAL_INTERP));
487*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "v", num_comp_q, CEED_EVAL_INTERP));
488*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
489*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_ifunction_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
49077841947SLeila Ghaffari   }
49177841947SLeila Ghaffari 
492e334ad8fSJed Brown   CeedQFunction qf_ijacobian_vol = NULL;
493e334ad8fSJed Brown   if (problem->apply_vol_ijacobian.qfunction) {
494*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->apply_vol_ijacobian.qfunction, problem->apply_vol_ijacobian.qfunction_loc,
495*a424bcd0SJames Wright                                                     &qf_ijacobian_vol));
496*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(qf_ijacobian_vol, problem->apply_vol_ijacobian.qfunction_context));
497*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "dq", num_comp_q, CEED_EVAL_INTERP));
498*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "Grad_dq", num_comp_q * dim, CEED_EVAL_GRAD));
499*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "qdata", q_data_size_vol, CEED_EVAL_NONE));
500*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "x", num_comp_x, CEED_EVAL_INTERP));
501*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddInput(qf_ijacobian_vol, "jac_data", jac_data_size_vol, CEED_EVAL_NONE));
502*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "v", num_comp_q, CEED_EVAL_INTERP));
503*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionAddOutput(qf_ijacobian_vol, "Grad_v", num_comp_q * dim, CEED_EVAL_GRAD));
504e334ad8fSJed Brown   }
505e334ad8fSJed Brown 
50677841947SLeila Ghaffari   // ---------------------------------------------------------------------------
50777841947SLeila Ghaffari   // Element coordinates
50877841947SLeila Ghaffari   // ---------------------------------------------------------------------------
50977841947SLeila Ghaffari   // -- Create CEED vector
510*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_x, &ceed_data->x_coord, NULL));
51177841947SLeila Ghaffari 
51277841947SLeila Ghaffari   // -- Copy PETSc vector in CEED vector
51377841947SLeila Ghaffari   Vec X_loc;
5143796c488SJed Brown   {
5153796c488SJed Brown     DM cdm;
5162b730f8bSJeremy L Thompson     PetscCall(DMGetCellCoordinateDM(dm, &cdm));
5172b730f8bSJeremy L Thompson     if (cdm) {
5182b730f8bSJeremy L Thompson       PetscCall(DMGetCellCoordinatesLocal(dm, &X_loc));
5192b730f8bSJeremy L Thompson     } else {
5202b730f8bSJeremy L Thompson       PetscCall(DMGetCoordinatesLocal(dm, &X_loc));
5213796c488SJed Brown     }
5222b730f8bSJeremy L Thompson   }
5232b730f8bSJeremy L Thompson   PetscCall(VecScale(X_loc, problem->dm_scale));
524fe69b334SJames Wright   PetscCall(VecCopyP2C(X_loc, ceed_data->x_coord));
52577841947SLeila Ghaffari 
52677841947SLeila Ghaffari   // -----------------------------------------------------------------------------
52777841947SLeila Ghaffari   // CEED vectors
52877841947SLeila Ghaffari   // -----------------------------------------------------------------------------
52977841947SLeila Ghaffari   // -- Create CEED vector for geometric data
530*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(ceed_data->elem_restr_qd_i, &ceed_data->q_data, NULL));
531*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionCreateVector(elem_restr_jd_i, &jac_data, NULL));
532d52d2babSJames Wright 
53377841947SLeila Ghaffari   // -----------------------------------------------------------------------------
53477841947SLeila Ghaffari   // CEED Operators
53577841947SLeila Ghaffari   // -----------------------------------------------------------------------------
53677841947SLeila Ghaffari   // -- Create CEED operator for quadrature data
537*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_setup_vol, NULL, NULL, &ceed_data->op_setup_vol));
538*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(ceed_data->op_setup_vol, "dx", ceed_data->elem_restr_x, ceed_data->basis_x, CEED_VECTOR_ACTIVE));
539*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(ceed_data->op_setup_vol, "weight", CEED_ELEMRESTRICTION_NONE, ceed_data->basis_x, CEED_VECTOR_NONE));
540*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(ceed_data->op_setup_vol, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
54177841947SLeila Ghaffari 
54277841947SLeila Ghaffari   // -- Create CEED operator for ICs
5435263e9c6SJames Wright   CeedOperator op_ics;
544*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_ics, NULL, NULL, &op_ics));
545*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "x", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE));
546*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "dx", ceed_data->elem_restr_x, ceed_data->basis_xc, CEED_VECTOR_ACTIVE));
547*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorSetField(op_ics, "q0", ceed_data->elem_restr_q, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE));
548*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(op_ics, "evaluation time", &user->phys->ics_time_label));
5495263e9c6SJames Wright   PetscCall(OperatorApplyContextCreate(NULL, dm, user->ceed, op_ics, ceed_data->x_coord, NULL, NULL, user->Q_loc, &ceed_data->op_ics_ctx));
550*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ics));
55177841947SLeila Ghaffari 
55277841947SLeila Ghaffari   // Create CEED operator for RHS
55377841947SLeila Ghaffari   if (ceed_data->qf_rhs_vol) {
55477841947SLeila Ghaffari     CeedOperator op;
555*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_rhs_vol, NULL, NULL, &op));
556*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
557*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
558*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data));
559*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
560*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
561*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
56277841947SLeila Ghaffari     user->op_rhs_vol = op;
56377841947SLeila Ghaffari   }
56477841947SLeila Ghaffari 
56577841947SLeila Ghaffari   // -- CEED operator for IFunction
56677841947SLeila Ghaffari   if (ceed_data->qf_ifunction_vol) {
56777841947SLeila Ghaffari     CeedOperator op;
568*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, ceed_data->qf_ifunction_vol, NULL, NULL, &op));
569*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
570*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_q", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
571*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "q dot", ceed_data->elem_restr_q, ceed_data->basis_q, user->q_dot_ceed));
572*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data));
573*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
574*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
575*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
576*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data));
577a3ae0734SJed Brown 
57877841947SLeila Ghaffari     user->op_ifunction_vol = op;
57977841947SLeila Ghaffari   }
58077841947SLeila Ghaffari 
581e334ad8fSJed Brown   CeedOperator op_ijacobian_vol = NULL;
582e334ad8fSJed Brown   if (qf_ijacobian_vol) {
583e334ad8fSJed Brown     CeedOperator op;
584*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorCreate(ceed, qf_ijacobian_vol, NULL, NULL, &op));
585*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
586*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_dq", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
587*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data));
588*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "x", ceed_data->elem_restr_x, ceed_data->basis_x, ceed_data->x_coord));
589*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "jac_data", elem_restr_jd_i, CEED_BASIS_COLLOCATED, jac_data));
590*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
591*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorSetField(op, "Grad_v", ceed_data->elem_restr_q, ceed_data->basis_q, CEED_VECTOR_ACTIVE));
592e334ad8fSJed Brown     op_ijacobian_vol = op;
593*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionDestroy(&qf_ijacobian_vol));
594e334ad8fSJed Brown   }
595e334ad8fSJed Brown 
59677841947SLeila Ghaffari   // *****************************************************************************
59777841947SLeila Ghaffari   // Set up CEED objects for the exterior domain (surface)
59877841947SLeila Ghaffari   // *****************************************************************************
5992b730f8bSJeremy L Thompson   CeedInt       height = 1, dim_sur = dim - height, P_sur = app_ctx->degree + 1, Q_sur = P_sur + app_ctx->q_extra;
6002b730f8bSJeremy L Thompson   const CeedInt q_data_size_sur = problem->q_data_size_sur, jac_data_size_sur = problem->jac_data_size_sur;
60177841947SLeila Ghaffari 
60277841947SLeila Ghaffari   // -----------------------------------------------------------------------------
60377841947SLeila Ghaffari   // CEED Bases
60477841947SLeila Ghaffari   // -----------------------------------------------------------------------------
6050814d5a7SKenneth E. Jansen 
6060814d5a7SKenneth E. Jansen   DMLabel  label   = 0;
6070814d5a7SKenneth E. Jansen   PetscInt face_id = 0;
6080814d5a7SKenneth E. Jansen   PetscInt field   = 0;  // Still want the normal, default field
6090814d5a7SKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, dm, label, face_id, height, field, &ceed_data->basis_q_sur));
6100814d5a7SKenneth E. Jansen   PetscCall(CreateBasisFromPlex(ceed, dm_coord, label, face_id, height, field, &ceed_data->basis_x_sur));
611*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedBasisCreateProjection(ceed_data->basis_x_sur, ceed_data->basis_q_sur, &ceed_data->basis_xc_sur));
61277841947SLeila Ghaffari 
61377841947SLeila Ghaffari   // -----------------------------------------------------------------------------
61477841947SLeila Ghaffari   // CEED QFunctions
61577841947SLeila Ghaffari   // -----------------------------------------------------------------------------
61677841947SLeila Ghaffari   // -- Create QFunction for quadrature data
617*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionCreateInterior(ceed, 1, problem->setup_sur.qfunction, problem->setup_sur.qfunction_loc, &ceed_data->qf_setup_sur));
618841e4c73SJed Brown   if (problem->setup_sur.qfunction_context) {
619*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedQFunctionSetContext(ceed_data->qf_setup_sur, problem->setup_sur.qfunction_context));
620841e4c73SJed Brown   }
621*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "dx", num_comp_x * dim_sur, CEED_EVAL_GRAD));
622*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddInput(ceed_data->qf_setup_sur, "weight", 1, CEED_EVAL_WEIGHT));
623*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedQFunctionAddOutput(ceed_data->qf_setup_sur, "surface qdata", q_data_size_sur, CEED_EVAL_NONE));
62477841947SLeila Ghaffari 
6252b730f8bSJeremy L Thompson   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_inflow,
6262b730f8bSJeremy L Thompson                               problem->apply_inflow_jacobian, &ceed_data->qf_apply_inflow, &ceed_data->qf_apply_inflow_jacobian));
6272b730f8bSJeremy L Thompson   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_outflow,
6282b730f8bSJeremy L Thompson                               problem->apply_outflow_jacobian, &ceed_data->qf_apply_outflow, &ceed_data->qf_apply_outflow_jacobian));
6292b730f8bSJeremy L Thompson   PetscCall(SetupBCQFunctions(ceed, dim_sur, num_comp_x, num_comp_q, q_data_size_sur, jac_data_size_sur, problem->apply_freestream,
6302b730f8bSJeremy L Thompson                               problem->apply_freestream_jacobian, &ceed_data->qf_apply_freestream, &ceed_data->qf_apply_freestream_jacobian));
63177841947SLeila Ghaffari 
63277841947SLeila Ghaffari   // *****************************************************************************
63377841947SLeila Ghaffari   // CEED Operator Apply
63477841947SLeila Ghaffari   // *****************************************************************************
63577841947SLeila Ghaffari   // -- Apply CEED Operator for the geometric data
636*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorApply(ceed_data->op_setup_vol, ceed_data->x_coord, ceed_data->q_data, CEED_REQUEST_IMMEDIATE));
63777841947SLeila Ghaffari 
63877841947SLeila Ghaffari   // -- Create and apply CEED Composite Operator for the entire domain
63977841947SLeila Ghaffari   if (!user->phys->implicit) {  // RHS
6409ad5e8e4SJames Wright     CeedOperator op_rhs;
6419ad5e8e4SJames 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,
6429ad5e8e4SJames Wright                                       NULL));
6439ad5e8e4SJames Wright     PetscCall(OperatorApplyContextCreate(dm, dm, ceed, op_rhs, user->q_ceed, user->g_ceed, user->Q_loc, NULL, &user->op_rhs_ctx));
644*a424bcd0SJames Wright     PetscCallCeed(ceed, CeedOperatorDestroy(&op_rhs));
64577841947SLeila Ghaffari   } else {  // IFunction
6462b730f8bSJeremy L Thompson     PetscCall(CreateOperatorForDomain(ceed, dm, bc, ceed_data, user->phys, user->op_ifunction_vol, op_ijacobian_vol, height, P_sur, Q_sur,
6472b730f8bSJeremy L Thompson                                       q_data_size_sur, jac_data_size_sur, &user->op_ifunction, op_ijacobian_vol ? &user->op_ijacobian : NULL));
648e334ad8fSJed Brown     if (user->op_ijacobian) {
649*a424bcd0SJames Wright       PetscCallCeed(ceed, CeedOperatorGetContextFieldLabel(user->op_ijacobian, "ijacobian time shift", &user->phys->ijacobian_time_shift_label));
650e334ad8fSJed Brown     }
6510814d5a7SKenneth E. Jansen     if (problem->use_strong_bc_ceed) PetscCall(SetupStrongBC_Ceed(ceed, ceed_data, dm, user, problem, bc, q_data_size_sur));
652f5452247SJames Wright     if (app_ctx->sgs_model_type == SGS_MODEL_DATA_DRIVEN) PetscCall(SGS_DD_ModelSetup(ceed, user, ceed_data, problem));
653dada6cc0SJames Wright   }
654f5452247SJames Wright 
655f5452247SJames Wright   if (app_ctx->turb_spanstats_enable) PetscCall(TurbulenceStatisticsSetup(ceed, user, ceed_data, problem));
6564e9802d1SJames Wright   if (app_ctx->diff_filter_monitor) PetscCall(DifferentialFilterSetup(ceed, user, ceed_data, problem));
657f5452247SJames Wright 
658*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedElemRestrictionDestroy(&elem_restr_jd_i));
659*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedOperatorDestroy(&op_ijacobian_vol));
660*a424bcd0SJames Wright   PetscCallCeed(ceed, CeedVectorDestroy(&jac_data));
661ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
66277841947SLeila Ghaffari }
663