xref: /libCEED/examples/petsc/src/libceedsetup.c (revision 4dbe2ad5260ad67d4b6134091b73ad6d31c08219)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
298285ab4SZach Atkins // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
398285ab4SZach Atkins //
498285ab4SZach Atkins // SPDX-License-Identifier: BSD-2-Clause
598285ab4SZach Atkins //
698285ab4SZach Atkins // This file is part of CEED:  http://github.com/ceed
798285ab4SZach Atkins 
8e83e87a5Sjeremylt #include "../include/libceedsetup.h"
92b730f8bSJeremy L Thompson 
102b730f8bSJeremy L Thompson #include <stdio.h>
112b730f8bSJeremy L Thompson 
1278f7fce3SZach Atkins #include "../include/libceedsetup.h"
13e83e87a5Sjeremylt #include "../include/petscutils.h"
14e83e87a5Sjeremylt 
15e83e87a5Sjeremylt // -----------------------------------------------------------------------------
16e83e87a5Sjeremylt // Destroy libCEED operator objects
17e83e87a5Sjeremylt // -----------------------------------------------------------------------------
18e83e87a5Sjeremylt PetscErrorCode CeedDataDestroy(CeedInt i, CeedData data) {
195dfaedb8SJed Brown   PetscFunctionBeginUser;
209b072555Sjeremylt   CeedVectorDestroy(&data->q_data);
219b072555Sjeremylt   CeedVectorDestroy(&data->x_ceed);
229b072555Sjeremylt   CeedVectorDestroy(&data->y_ceed);
239b072555Sjeremylt   CeedBasisDestroy(&data->basis_x);
249b072555Sjeremylt   CeedBasisDestroy(&data->basis_u);
259b072555Sjeremylt   CeedElemRestrictionDestroy(&data->elem_restr_u);
269b072555Sjeremylt   CeedElemRestrictionDestroy(&data->elem_restr_x);
279b072555Sjeremylt   CeedElemRestrictionDestroy(&data->elem_restr_u_i);
289b072555Sjeremylt   CeedElemRestrictionDestroy(&data->elem_restr_qd_i);
299b072555Sjeremylt   CeedQFunctionDestroy(&data->qf_apply);
309b072555Sjeremylt   CeedOperatorDestroy(&data->op_apply);
31e83e87a5Sjeremylt   if (i > 0) {
329b072555Sjeremylt     CeedOperatorDestroy(&data->op_prolong);
339b072555Sjeremylt     CeedOperatorDestroy(&data->op_restrict);
34e83e87a5Sjeremylt   }
352b730f8bSJeremy L Thompson   PetscCall(PetscFree(data));
36ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
37e83e87a5Sjeremylt };
38e83e87a5Sjeremylt 
39e83e87a5Sjeremylt // -----------------------------------------------------------------------------
40e83e87a5Sjeremylt // Set up libCEED for a given degree
41e83e87a5Sjeremylt // -----------------------------------------------------------------------------
422b730f8bSJeremy L Thompson PetscErrorCode SetupLibceedByDegree(DM dm, Ceed ceed, CeedInt degree, CeedInt topo_dim, CeedInt q_extra, PetscInt num_comp_x, PetscInt num_comp_u,
43*4dbe2ad5SJeremy L Thompson                                     PetscInt g_size, PetscInt xl_size, BPData bp_data, CeedData data, PetscBool setup_rhs, PetscBool is_fine_level,
44*4dbe2ad5SJeremy L Thompson                                     CeedVector rhs_ceed, CeedVector *target) {
459b072555Sjeremylt   DM                  dm_coord;
46e83e87a5Sjeremylt   Vec                 coords;
479b072555Sjeremylt   const PetscScalar  *coord_array;
489b072555Sjeremylt   CeedBasis           basis_x, basis_u;
499b072555Sjeremylt   CeedElemRestriction elem_restr_x, elem_restr_u, elem_restr_u_i, elem_restr_qd_i;
50*4dbe2ad5SJeremy L Thompson   CeedQFunction       qf_setup_geo = NULL, qf_apply = NULL;
519b072555Sjeremylt   CeedOperator        op_setup_geo, op_apply;
529b072555Sjeremylt   CeedVector          x_coord, q_data, x_ceed, y_ceed;
534d00b080SJeremy L Thompson   PetscInt            c_start, c_end, num_elem;
544d00b080SJeremy L Thompson   CeedInt             num_qpts, q_data_size = bp_data.q_data_size;
552b730f8bSJeremy L Thompson   CeedScalar          R = 1;                         // radius of the sphere
562b730f8bSJeremy L Thompson   CeedScalar          l = 1.0 / PetscSqrtReal(3.0);  // half edge of the inscribed cube
57e83e87a5Sjeremylt 
585dfaedb8SJed Brown   PetscFunctionBeginUser;
592b730f8bSJeremy L Thompson   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
60129d8736Srezgarshakeri 
61129d8736Srezgarshakeri   // CEED bases
622b730f8bSJeremy L Thompson   PetscCall(CreateBasisFromPlex(ceed, dm_coord, 0, 0, 0, 0, bp_data, &basis_x));
632b730f8bSJeremy L Thompson   PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, bp_data, &basis_u));
64129d8736Srezgarshakeri 
65129d8736Srezgarshakeri   // CEED restrictions
662b730f8bSJeremy L Thompson   PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, 0, 0, 0, &elem_restr_x));
672b730f8bSJeremy L Thompson   PetscCall(CreateRestrictionFromPlex(ceed, dm, 0, 0, 0, &elem_restr_u));
68e83e87a5Sjeremylt 
692b730f8bSJeremy L Thompson   PetscCall(DMPlexGetHeightStratum(dm, 0, &c_start, &c_end));
709b072555Sjeremylt   num_elem = c_end - c_start;
71129d8736Srezgarshakeri   CeedBasisGetNumQuadraturePoints(basis_u, &num_qpts);
72e83e87a5Sjeremylt 
732b730f8bSJeremy L Thompson   CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, num_comp_u, num_comp_u * num_elem * num_qpts, CEED_STRIDES_BACKEND, &elem_restr_u_i);
742b730f8bSJeremy L Thompson   CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, q_data_size, q_data_size * num_elem * num_qpts, CEED_STRIDES_BACKEND, &elem_restr_qd_i);
75e83e87a5Sjeremylt 
76e83e87a5Sjeremylt   // Element coordinates
772b730f8bSJeremy L Thompson   PetscCall(DMGetCoordinatesLocal(dm, &coords));
782b730f8bSJeremy L Thompson   PetscCall(VecGetArrayRead(coords, &coord_array));
79e83e87a5Sjeremylt 
809b072555Sjeremylt   CeedElemRestrictionCreateVector(elem_restr_x, &x_coord, NULL);
812b730f8bSJeremy L Thompson   CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, (PetscScalar *)coord_array);
822b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayRead(coords, &coord_array));
83e83e87a5Sjeremylt 
84e83e87a5Sjeremylt   // Create the persistent vectors that will be needed in setup and apply
859b072555Sjeremylt   CeedVectorCreate(ceed, q_data_size * num_elem * num_qpts, &q_data);
869b072555Sjeremylt   CeedVectorCreate(ceed, xl_size, &x_ceed);
879b072555Sjeremylt   CeedVectorCreate(ceed, xl_size, &y_ceed);
88e83e87a5Sjeremylt 
89*4dbe2ad5SJeremy L Thompson   if (is_fine_level) {
90e83e87a5Sjeremylt     // Create the QFunction that builds the context data
912b730f8bSJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_geo, bp_data.setup_geo_loc, &qf_setup_geo);
929b072555Sjeremylt     CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP);
939b072555Sjeremylt     CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x * topo_dim, CEED_EVAL_GRAD);
949b072555Sjeremylt     CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT);
95a61c78d6SJeremy L Thompson     CeedQFunctionAddOutput(qf_setup_geo, "qdata", q_data_size, CEED_EVAL_NONE);
96e83e87a5Sjeremylt 
97e83e87a5Sjeremylt     // Create the operator that builds the quadrature data
989b072555Sjeremylt     CeedOperatorCreate(ceed, qf_setup_geo, NULL, NULL, &op_setup_geo);
992b730f8bSJeremy L Thompson     CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
1002b730f8bSJeremy L Thompson     CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
1012b730f8bSJeremy L Thompson     CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE);
102356036faSJeremy L Thompson     CeedOperatorSetField(op_setup_geo, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
103e83e87a5Sjeremylt 
1049b072555Sjeremylt     // Setup q_data
1059b072555Sjeremylt     CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE);
106e83e87a5Sjeremylt 
107e83e87a5Sjeremylt     // Set up PDE operator
1089b072555Sjeremylt     CeedInt in_scale  = bp_data.in_mode == CEED_EVAL_GRAD ? topo_dim : 1;
1099b072555Sjeremylt     CeedInt out_scale = bp_data.out_mode == CEED_EVAL_GRAD ? topo_dim : 1;
1102b730f8bSJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, bp_data.apply, bp_data.apply_loc, &qf_apply);
1119b072555Sjeremylt     CeedQFunctionAddInput(qf_apply, "u", num_comp_u * in_scale, bp_data.in_mode);
112a61c78d6SJeremy L Thompson     CeedQFunctionAddInput(qf_apply, "qdata", q_data_size, CEED_EVAL_NONE);
1139b072555Sjeremylt     CeedQFunctionAddOutput(qf_apply, "v", num_comp_u * out_scale, bp_data.out_mode);
114e83e87a5Sjeremylt 
115e83e87a5Sjeremylt     // Create the mass or diff operator
1169b072555Sjeremylt     CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply);
1179b072555Sjeremylt     CeedOperatorSetField(op_apply, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
118356036faSJeremy L Thompson     CeedOperatorSetField(op_apply, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
1199b072555Sjeremylt     CeedOperatorSetField(op_apply, "v", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
120e83e87a5Sjeremylt 
121*4dbe2ad5SJeremy L Thompson     // Cleanup
122*4dbe2ad5SJeremy L Thompson     CeedQFunctionDestroy(&qf_setup_geo);
123*4dbe2ad5SJeremy L Thompson     CeedOperatorDestroy(&op_setup_geo);
124*4dbe2ad5SJeremy L Thompson   }
125*4dbe2ad5SJeremy L Thompson 
126e83e87a5Sjeremylt   // Set up RHS if needed
127e83e87a5Sjeremylt   if (setup_rhs) {
1289b072555Sjeremylt     CeedQFunction qf_setup_rhs;
1299b072555Sjeremylt     CeedOperator  op_setup_rhs;
1309b072555Sjeremylt     CeedVectorCreate(ceed, num_elem * num_qpts * num_comp_u, target);
131e83e87a5Sjeremylt     // Create the q-function that sets up the RHS and true solution
1322b730f8bSJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_rhs, bp_data.setup_rhs_loc, &qf_setup_rhs);
1339b072555Sjeremylt     CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP);
134a61c78d6SJeremy L Thompson     CeedQFunctionAddInput(qf_setup_rhs, "qdata", q_data_size, CEED_EVAL_NONE);
1352b730f8bSJeremy L Thompson     CeedQFunctionAddOutput(qf_setup_rhs, "true solution", num_comp_u, CEED_EVAL_NONE);
1369b072555Sjeremylt     CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP);
137e83e87a5Sjeremylt 
138e83e87a5Sjeremylt     // Create the operator that builds the RHS and true solution
1399b072555Sjeremylt     CeedOperatorCreate(ceed, qf_setup_rhs, NULL, NULL, &op_setup_rhs);
1402b730f8bSJeremy L Thompson     CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
141356036faSJeremy L Thompson     CeedOperatorSetField(op_setup_rhs, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
142356036faSJeremy L Thompson     CeedOperatorSetField(op_setup_rhs, "true solution", elem_restr_u_i, CEED_BASIS_NONE, *target);
1432b730f8bSJeremy L Thompson     CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
144e83e87a5Sjeremylt 
145e83e87a5Sjeremylt     // Set up the libCEED context
1469b072555Sjeremylt     CeedQFunctionContext ctx_rhs_setup;
1479b072555Sjeremylt     CeedQFunctionContextCreate(ceed, &ctx_rhs_setup);
1489b072555Sjeremylt     CeedScalar rhs_setup_data[2] = {R, l};
1492b730f8bSJeremy L Thompson     CeedQFunctionContextSetData(ctx_rhs_setup, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof rhs_setup_data, &rhs_setup_data);
1509b072555Sjeremylt     CeedQFunctionSetContext(qf_setup_rhs, ctx_rhs_setup);
1519b072555Sjeremylt     CeedQFunctionContextDestroy(&ctx_rhs_setup);
152e83e87a5Sjeremylt 
153e83e87a5Sjeremylt     // Setup RHS and target
1549b072555Sjeremylt     CeedOperatorApply(op_setup_rhs, x_coord, rhs_ceed, CEED_REQUEST_IMMEDIATE);
155e83e87a5Sjeremylt 
156e83e87a5Sjeremylt     // Cleanup
1579b072555Sjeremylt     CeedQFunctionDestroy(&qf_setup_rhs);
1589b072555Sjeremylt     CeedOperatorDestroy(&op_setup_rhs);
159e83e87a5Sjeremylt   }
160e83e87a5Sjeremylt   // Cleanup
1619b072555Sjeremylt   CeedVectorDestroy(&x_coord);
162e83e87a5Sjeremylt 
163e83e87a5Sjeremylt   // Save libCEED data required for level
1642b730f8bSJeremy L Thompson   data->basis_x         = basis_x;
1652b730f8bSJeremy L Thompson   data->basis_u         = basis_u;
1669b072555Sjeremylt   data->elem_restr_x    = elem_restr_x;
1679b072555Sjeremylt   data->elem_restr_u    = elem_restr_u;
1689b072555Sjeremylt   data->elem_restr_u_i  = elem_restr_u_i;
1699b072555Sjeremylt   data->elem_restr_qd_i = elem_restr_qd_i;
1709b072555Sjeremylt   data->qf_apply        = qf_apply;
1719b072555Sjeremylt   data->op_apply        = op_apply;
1729b072555Sjeremylt   data->q_data          = q_data;
1739b072555Sjeremylt   data->x_ceed          = x_ceed;
1749b072555Sjeremylt   data->y_ceed          = y_ceed;
175d4d45553Srezgarshakeri   data->q_data_size     = q_data_size;
176ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
177e83e87a5Sjeremylt };
178e83e87a5Sjeremylt 
179e83e87a5Sjeremylt // -----------------------------------------------------------------------------
180e83e87a5Sjeremylt // Setup libCEED level transfer operator objects
181e83e87a5Sjeremylt // -----------------------------------------------------------------------------
1822b730f8bSJeremy L Thompson PetscErrorCode CeedLevelTransferSetup(DM dm, Ceed ceed, CeedInt level, CeedInt num_comp_u, CeedData *data, BPData bp_data, Vec fine_mult) {
183751eb813Srezgarshakeri   PetscFunctionBeginUser;
184e83e87a5Sjeremylt   // Restriction - Fine to corse
1859b072555Sjeremylt   CeedOperator op_restrict;
186e83e87a5Sjeremylt   // Interpolation - Corse to fine
1879b072555Sjeremylt   CeedOperator op_prolong;
1881c9a79dbSrezgarshakeri   // Coarse grid operator
1891c9a79dbSrezgarshakeri   CeedOperator op_apply;
1901c9a79dbSrezgarshakeri   // Basis
1911c9a79dbSrezgarshakeri   CeedBasis basis_u;
1922b730f8bSJeremy L Thompson   PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, bp_data, &basis_u));
193e83e87a5Sjeremylt 
1941c9a79dbSrezgarshakeri   // ---------------------------------------------------------------------------
1951c9a79dbSrezgarshakeri   // Coarse Grid, Prolongation, and Restriction Operators
1961c9a79dbSrezgarshakeri   // ---------------------------------------------------------------------------
1971c9a79dbSrezgarshakeri   // Create the Operators that compute the prolongation and
1981c9a79dbSrezgarshakeri   //   restriction between the p-multigrid levels and the coarse grid eval.
1991c9a79dbSrezgarshakeri   // ---------------------------------------------------------------------------
2001c9a79dbSrezgarshakeri   // Place in libCEED array
2011c9a79dbSrezgarshakeri   PetscMemType m_mem_type;
202179e5961SZach Atkins   PetscCall(VecReadP2C(fine_mult, &m_mem_type, data[level]->x_ceed));
203e83e87a5Sjeremylt 
2042b730f8bSJeremy L Thompson   CeedOperatorMultigridLevelCreate(data[level]->op_apply, data[level]->x_ceed, data[level - 1]->elem_restr_u, basis_u, &op_apply, &op_prolong,
2052b730f8bSJeremy L Thompson                                    &op_restrict);
206e83e87a5Sjeremylt 
2071c9a79dbSrezgarshakeri   // Restore PETSc vector
208179e5961SZach Atkins   PetscCall(VecReadC2P(data[level]->x_ceed, m_mem_type, fine_mult));
2092b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(fine_mult));
2101c9a79dbSrezgarshakeri   // -- Save libCEED data
2111c9a79dbSrezgarshakeri   data[level - 1]->op_apply = op_apply;
2121c9a79dbSrezgarshakeri   data[level]->op_prolong   = op_prolong;
2131c9a79dbSrezgarshakeri   data[level]->op_restrict  = op_restrict;
2141c9a79dbSrezgarshakeri 
2151c9a79dbSrezgarshakeri   CeedBasisDestroy(&basis_u);
216ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
217e83e87a5Sjeremylt };
218e83e87a5Sjeremylt 
21978f7fce3SZach Atkins PetscErrorCode SetupErrorOperator(DM dm, Ceed ceed, BPData bp_data, CeedInt topo_dim, PetscInt num_comp_x, PetscInt num_comp_u,
22078f7fce3SZach Atkins                                   CeedOperator *op_error) {
22178f7fce3SZach Atkins   DM                  dm_coord;
22278f7fce3SZach Atkins   Vec                 coords;
22378f7fce3SZach Atkins   const PetscScalar  *coord_array;
22478f7fce3SZach Atkins   CeedBasis           basis_x, basis_u;
22578f7fce3SZach Atkins   CeedElemRestriction elem_restr_x, elem_restr_u, elem_restr_u_i, elem_restr_qd_i;
22678f7fce3SZach Atkins   CeedQFunction       qf_setup_geo, qf_setup_rhs, qf_error;
22778f7fce3SZach Atkins   CeedOperator        op_setup_geo, op_setup_rhs;
22878f7fce3SZach Atkins   CeedVector          x_coord, q_data, target, rhs;
2294d00b080SJeremy L Thompson   PetscInt            c_start, c_end, num_elem;
2304d00b080SJeremy L Thompson   CeedInt             num_qpts, q_data_size = bp_data.q_data_size;
23178f7fce3SZach Atkins   CeedScalar          R = 1;                         // radius of the sphere
23278f7fce3SZach Atkins   CeedScalar          l = 1.0 / PetscSqrtReal(3.0);  // half edge of the inscribed cube
23378f7fce3SZach Atkins 
23478f7fce3SZach Atkins   PetscFunctionBeginUser;
23578f7fce3SZach Atkins   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
23678f7fce3SZach Atkins 
23778f7fce3SZach Atkins   // CEED bases
23878f7fce3SZach Atkins   PetscCall(CreateBasisFromPlex(ceed, dm_coord, 0, 0, 0, 0, bp_data, &basis_x));
23978f7fce3SZach Atkins   PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, bp_data, &basis_u));
24078f7fce3SZach Atkins 
24178f7fce3SZach Atkins   // CEED restrictions
24278f7fce3SZach Atkins   PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, 0, 0, 0, &elem_restr_x));
24378f7fce3SZach Atkins   PetscCall(CreateRestrictionFromPlex(ceed, dm, 0, 0, 0, &elem_restr_u));
24478f7fce3SZach Atkins 
24578f7fce3SZach Atkins   PetscCall(DMPlexGetHeightStratum(dm, 0, &c_start, &c_end));
24678f7fce3SZach Atkins   num_elem = c_end - c_start;
24778f7fce3SZach Atkins   CeedBasisGetNumQuadraturePoints(basis_u, &num_qpts);
24878f7fce3SZach Atkins 
24978f7fce3SZach Atkins   CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, num_comp_u, num_comp_u * num_elem * num_qpts, CEED_STRIDES_BACKEND, &elem_restr_u_i);
25078f7fce3SZach Atkins   CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, q_data_size, q_data_size * num_elem * num_qpts, CEED_STRIDES_BACKEND, &elem_restr_qd_i);
25178f7fce3SZach Atkins 
25278f7fce3SZach Atkins   // Element coordinates
25378f7fce3SZach Atkins   PetscCall(DMGetCoordinatesLocal(dm, &coords));
25478f7fce3SZach Atkins   PetscCall(VecGetArrayRead(coords, &coord_array));
25578f7fce3SZach Atkins 
25678f7fce3SZach Atkins   CeedElemRestrictionCreateVector(elem_restr_x, &x_coord, NULL);
25778f7fce3SZach Atkins   CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, (PetscScalar *)coord_array);
25878f7fce3SZach Atkins   PetscCall(VecRestoreArrayRead(coords, &coord_array));
25978f7fce3SZach Atkins 
26078f7fce3SZach Atkins   // Create the persistent vectors that will be needed in setup and apply
26178f7fce3SZach Atkins   CeedVectorCreate(ceed, q_data_size * num_elem * num_qpts, &q_data);
26278f7fce3SZach Atkins 
26378f7fce3SZach Atkins   // Create the QFunction that builds the context data
26478f7fce3SZach Atkins   CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_geo, bp_data.setup_geo_loc, &qf_setup_geo);
26578f7fce3SZach Atkins   CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP);
26678f7fce3SZach Atkins   CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x * topo_dim, CEED_EVAL_GRAD);
26778f7fce3SZach Atkins   CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT);
26878f7fce3SZach Atkins   CeedQFunctionAddOutput(qf_setup_geo, "qdata", q_data_size, CEED_EVAL_NONE);
26978f7fce3SZach Atkins 
27078f7fce3SZach Atkins   // Create the operator that builds the quadrature data
27178f7fce3SZach Atkins   CeedOperatorCreate(ceed, qf_setup_geo, NULL, NULL, &op_setup_geo);
27278f7fce3SZach Atkins   CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
27378f7fce3SZach Atkins   CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
27478f7fce3SZach Atkins   CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE);
27578f7fce3SZach Atkins   CeedOperatorSetField(op_setup_geo, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
27678f7fce3SZach Atkins 
27778f7fce3SZach Atkins   // Setup q_data
27878f7fce3SZach Atkins   CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE);
27978f7fce3SZach Atkins 
28078f7fce3SZach Atkins   // Set up target vector
28178f7fce3SZach Atkins   CeedElemRestrictionCreateVector(elem_restr_u, &rhs, NULL);
28278f7fce3SZach Atkins   CeedVectorCreate(ceed, num_elem * num_qpts * num_comp_u, &target);
28378f7fce3SZach Atkins   // Create the q-function that sets up the RHS and true solution
28478f7fce3SZach Atkins   CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_rhs, bp_data.setup_rhs_loc, &qf_setup_rhs);
28578f7fce3SZach Atkins   CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP);
28678f7fce3SZach Atkins   CeedQFunctionAddInput(qf_setup_rhs, "qdata", q_data_size, CEED_EVAL_NONE);
28778f7fce3SZach Atkins   CeedQFunctionAddOutput(qf_setup_rhs, "true solution", num_comp_u, CEED_EVAL_NONE);
28878f7fce3SZach Atkins   CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP);
28978f7fce3SZach Atkins 
29078f7fce3SZach Atkins   // Create the operator that builds the RHS and true solution
29178f7fce3SZach Atkins   CeedOperatorCreate(ceed, qf_setup_rhs, NULL, NULL, &op_setup_rhs);
29278f7fce3SZach Atkins   CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
29378f7fce3SZach Atkins   CeedOperatorSetField(op_setup_rhs, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
29478f7fce3SZach Atkins   CeedOperatorSetField(op_setup_rhs, "true solution", elem_restr_u_i, CEED_BASIS_NONE, target);
29578f7fce3SZach Atkins   CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
29678f7fce3SZach Atkins 
29778f7fce3SZach Atkins   // Set up the libCEED context
29878f7fce3SZach Atkins   CeedQFunctionContext ctx_rhs_setup;
29978f7fce3SZach Atkins   CeedQFunctionContextCreate(ceed, &ctx_rhs_setup);
30078f7fce3SZach Atkins   CeedScalar rhs_setup_data[2] = {R, l};
30178f7fce3SZach Atkins   CeedQFunctionContextSetData(ctx_rhs_setup, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof rhs_setup_data, &rhs_setup_data);
30278f7fce3SZach Atkins   CeedQFunctionSetContext(qf_setup_rhs, ctx_rhs_setup);
30378f7fce3SZach Atkins   CeedQFunctionContextDestroy(&ctx_rhs_setup);
30478f7fce3SZach Atkins 
30578f7fce3SZach Atkins   // Setup RHS and target
30678f7fce3SZach Atkins   CeedOperatorApply(op_setup_rhs, x_coord, rhs, CEED_REQUEST_IMMEDIATE);
30778f7fce3SZach Atkins 
30878f7fce3SZach Atkins   // Set up error operator
30978f7fce3SZach Atkins   // Create the error QFunction
31078f7fce3SZach Atkins   CeedQFunctionCreateInterior(ceed, 1, bp_data.error, bp_data.error_loc, &qf_error);
31178f7fce3SZach Atkins   CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
31278f7fce3SZach Atkins   CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
31378f7fce3SZach Atkins   CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE);
31478f7fce3SZach Atkins   CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_INTERP);
31578f7fce3SZach Atkins 
31678f7fce3SZach Atkins   // Create the error operator
31778f7fce3SZach Atkins   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, op_error);
31878f7fce3SZach Atkins   CeedOperatorSetField(*op_error, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
31978f7fce3SZach Atkins   CeedOperatorSetField(*op_error, "true_soln", elem_restr_u_i, CEED_BASIS_NONE, target);
32078f7fce3SZach Atkins   CeedOperatorSetField(*op_error, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
32178f7fce3SZach Atkins   CeedOperatorSetField(*op_error, "error", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
32278f7fce3SZach Atkins 
32378f7fce3SZach Atkins   // Cleanup
32478f7fce3SZach Atkins   CeedQFunctionDestroy(&qf_setup_rhs);
32578f7fce3SZach Atkins   CeedOperatorDestroy(&op_setup_rhs);
32678f7fce3SZach Atkins   CeedQFunctionDestroy(&qf_setup_geo);
32778f7fce3SZach Atkins   CeedOperatorDestroy(&op_setup_geo);
32878f7fce3SZach Atkins   CeedQFunctionDestroy(&qf_error);
32978f7fce3SZach Atkins   CeedVectorDestroy(&x_coord);
33078f7fce3SZach Atkins   CeedVectorDestroy(&rhs);
33178f7fce3SZach Atkins   CeedVectorDestroy(&target);
33278f7fce3SZach Atkins   CeedVectorDestroy(&q_data);
33378f7fce3SZach Atkins   CeedElemRestrictionDestroy(&elem_restr_u_i);
33478f7fce3SZach Atkins   CeedElemRestrictionDestroy(&elem_restr_qd_i);
33578f7fce3SZach Atkins   CeedElemRestrictionDestroy(&elem_restr_x);
33678f7fce3SZach Atkins   CeedElemRestrictionDestroy(&elem_restr_u);
33778f7fce3SZach Atkins   CeedBasisDestroy(&basis_x);
33878f7fce3SZach Atkins   CeedBasisDestroy(&basis_u);
33978f7fce3SZach Atkins 
34078f7fce3SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
34178f7fce3SZach Atkins }
342