xref: /libCEED/examples/petsc/src/libceedsetup.c (revision 98285ab464d104dd6040959f61a83e9969073ceb)
1*98285ab4SZach Atkins // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2*98285ab4SZach Atkins // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3*98285ab4SZach Atkins //
4*98285ab4SZach Atkins // SPDX-License-Identifier: BSD-2-Clause
5*98285ab4SZach Atkins //
6*98285ab4SZach Atkins // This file is part of CEED:  http://github.com/ceed
7*98285ab4SZach Atkins 
8e83e87a5Sjeremylt #include "../include/libceedsetup.h"
92b730f8bSJeremy L Thompson 
102b730f8bSJeremy L Thompson #include <stdio.h>
112b730f8bSJeremy L Thompson 
12e83e87a5Sjeremylt #include "../include/petscutils.h"
13e83e87a5Sjeremylt 
14e83e87a5Sjeremylt // -----------------------------------------------------------------------------
15e83e87a5Sjeremylt // Destroy libCEED operator objects
16e83e87a5Sjeremylt // -----------------------------------------------------------------------------
17e83e87a5Sjeremylt PetscErrorCode CeedDataDestroy(CeedInt i, CeedData data) {
185dfaedb8SJed Brown   PetscFunctionBeginUser;
199b072555Sjeremylt   CeedVectorDestroy(&data->q_data);
209b072555Sjeremylt   CeedVectorDestroy(&data->x_ceed);
219b072555Sjeremylt   CeedVectorDestroy(&data->y_ceed);
229b072555Sjeremylt   CeedBasisDestroy(&data->basis_x);
239b072555Sjeremylt   CeedBasisDestroy(&data->basis_u);
249b072555Sjeremylt   CeedElemRestrictionDestroy(&data->elem_restr_u);
259b072555Sjeremylt   CeedElemRestrictionDestroy(&data->elem_restr_x);
269b072555Sjeremylt   CeedElemRestrictionDestroy(&data->elem_restr_u_i);
279b072555Sjeremylt   CeedElemRestrictionDestroy(&data->elem_restr_qd_i);
289b072555Sjeremylt   CeedQFunctionDestroy(&data->qf_apply);
299b072555Sjeremylt   CeedOperatorDestroy(&data->op_apply);
30e83e87a5Sjeremylt   if (i > 0) {
319b072555Sjeremylt     CeedOperatorDestroy(&data->op_prolong);
329b072555Sjeremylt     CeedOperatorDestroy(&data->op_restrict);
33e83e87a5Sjeremylt   }
342b730f8bSJeremy L Thompson   PetscCall(PetscFree(data));
35ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
36e83e87a5Sjeremylt };
37e83e87a5Sjeremylt 
38e83e87a5Sjeremylt // -----------------------------------------------------------------------------
39e83e87a5Sjeremylt // Set up libCEED for a given degree
40e83e87a5Sjeremylt // -----------------------------------------------------------------------------
412b730f8bSJeremy L Thompson PetscErrorCode SetupLibceedByDegree(DM dm, Ceed ceed, CeedInt degree, CeedInt topo_dim, CeedInt q_extra, PetscInt num_comp_x, PetscInt num_comp_u,
422b730f8bSJeremy L Thompson                                     PetscInt g_size, PetscInt xl_size, BPData bp_data, CeedData data, PetscBool setup_rhs, CeedVector rhs_ceed,
43e83e87a5Sjeremylt                                     CeedVector *target) {
449b072555Sjeremylt   DM                  dm_coord;
45e83e87a5Sjeremylt   Vec                 coords;
469b072555Sjeremylt   const PetscScalar  *coord_array;
479b072555Sjeremylt   CeedBasis           basis_x, basis_u;
489b072555Sjeremylt   CeedElemRestriction elem_restr_x, elem_restr_u, elem_restr_u_i, elem_restr_qd_i;
499b072555Sjeremylt   CeedQFunction       qf_setup_geo, qf_apply;
509b072555Sjeremylt   CeedOperator        op_setup_geo, op_apply;
519b072555Sjeremylt   CeedVector          x_coord, q_data, x_ceed, y_ceed;
522b730f8bSJeremy L Thompson   CeedInt             num_qpts, c_start, c_end, num_elem, q_data_size = bp_data.q_data_size;
532b730f8bSJeremy L Thompson   CeedScalar          R = 1;                         // radius of the sphere
542b730f8bSJeremy L Thompson   CeedScalar          l = 1.0 / PetscSqrtReal(3.0);  // half edge of the inscribed cube
55e83e87a5Sjeremylt 
565dfaedb8SJed Brown   PetscFunctionBeginUser;
572b730f8bSJeremy L Thompson   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
58129d8736Srezgarshakeri 
59129d8736Srezgarshakeri   // CEED bases
602b730f8bSJeremy L Thompson   PetscCall(CreateBasisFromPlex(ceed, dm_coord, 0, 0, 0, 0, bp_data, &basis_x));
612b730f8bSJeremy L Thompson   PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, bp_data, &basis_u));
62129d8736Srezgarshakeri 
63129d8736Srezgarshakeri   // CEED restrictions
642b730f8bSJeremy L Thompson   PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, 0, 0, 0, &elem_restr_x));
652b730f8bSJeremy L Thompson   PetscCall(CreateRestrictionFromPlex(ceed, dm, 0, 0, 0, &elem_restr_u));
66e83e87a5Sjeremylt 
672b730f8bSJeremy L Thompson   PetscCall(DMPlexGetHeightStratum(dm, 0, &c_start, &c_end));
689b072555Sjeremylt   num_elem = c_end - c_start;
69129d8736Srezgarshakeri   CeedBasisGetNumQuadraturePoints(basis_u, &num_qpts);
70e83e87a5Sjeremylt 
712b730f8bSJeremy 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);
722b730f8bSJeremy 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);
73e83e87a5Sjeremylt 
74e83e87a5Sjeremylt   // Element coordinates
752b730f8bSJeremy L Thompson   PetscCall(DMGetCoordinatesLocal(dm, &coords));
762b730f8bSJeremy L Thompson   PetscCall(VecGetArrayRead(coords, &coord_array));
77e83e87a5Sjeremylt 
789b072555Sjeremylt   CeedElemRestrictionCreateVector(elem_restr_x, &x_coord, NULL);
792b730f8bSJeremy L Thompson   CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, (PetscScalar *)coord_array);
802b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayRead(coords, &coord_array));
81e83e87a5Sjeremylt 
82e83e87a5Sjeremylt   // Create the persistent vectors that will be needed in setup and apply
839b072555Sjeremylt   CeedVectorCreate(ceed, q_data_size * num_elem * num_qpts, &q_data);
849b072555Sjeremylt   CeedVectorCreate(ceed, xl_size, &x_ceed);
859b072555Sjeremylt   CeedVectorCreate(ceed, xl_size, &y_ceed);
86e83e87a5Sjeremylt 
87e83e87a5Sjeremylt   // Create the QFunction that builds the context data
882b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_geo, bp_data.setup_geo_loc, &qf_setup_geo);
899b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP);
909b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x * topo_dim, CEED_EVAL_GRAD);
919b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT);
92a61c78d6SJeremy L Thompson   CeedQFunctionAddOutput(qf_setup_geo, "qdata", q_data_size, CEED_EVAL_NONE);
93e83e87a5Sjeremylt 
94e83e87a5Sjeremylt   // Create the operator that builds the quadrature data
959b072555Sjeremylt   CeedOperatorCreate(ceed, qf_setup_geo, NULL, NULL, &op_setup_geo);
962b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
972b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
982b730f8bSJeremy L Thompson   CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE);
99356036faSJeremy L Thompson   CeedOperatorSetField(op_setup_geo, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
100e83e87a5Sjeremylt 
1019b072555Sjeremylt   // Setup q_data
1029b072555Sjeremylt   CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE);
103e83e87a5Sjeremylt 
104e83e87a5Sjeremylt   // Set up PDE operator
1059b072555Sjeremylt   CeedInt in_scale  = bp_data.in_mode == CEED_EVAL_GRAD ? topo_dim : 1;
1069b072555Sjeremylt   CeedInt out_scale = bp_data.out_mode == CEED_EVAL_GRAD ? topo_dim : 1;
1072b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bp_data.apply, bp_data.apply_loc, &qf_apply);
1089b072555Sjeremylt   CeedQFunctionAddInput(qf_apply, "u", num_comp_u * in_scale, bp_data.in_mode);
109a61c78d6SJeremy L Thompson   CeedQFunctionAddInput(qf_apply, "qdata", q_data_size, CEED_EVAL_NONE);
1109b072555Sjeremylt   CeedQFunctionAddOutput(qf_apply, "v", num_comp_u * out_scale, bp_data.out_mode);
111e83e87a5Sjeremylt 
112e83e87a5Sjeremylt   // Create the mass or diff operator
1139b072555Sjeremylt   CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply);
1149b072555Sjeremylt   CeedOperatorSetField(op_apply, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
115356036faSJeremy L Thompson   CeedOperatorSetField(op_apply, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
1169b072555Sjeremylt   CeedOperatorSetField(op_apply, "v", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
117e83e87a5Sjeremylt 
118e83e87a5Sjeremylt   // Set up RHS if needed
119e83e87a5Sjeremylt   if (setup_rhs) {
1209b072555Sjeremylt     CeedQFunction qf_setup_rhs;
1219b072555Sjeremylt     CeedOperator  op_setup_rhs;
1229b072555Sjeremylt     CeedVectorCreate(ceed, num_elem * num_qpts * num_comp_u, target);
123e83e87a5Sjeremylt     // Create the q-function that sets up the RHS and true solution
1242b730f8bSJeremy L Thompson     CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_rhs, bp_data.setup_rhs_loc, &qf_setup_rhs);
1259b072555Sjeremylt     CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP);
126a61c78d6SJeremy L Thompson     CeedQFunctionAddInput(qf_setup_rhs, "qdata", q_data_size, CEED_EVAL_NONE);
1272b730f8bSJeremy L Thompson     CeedQFunctionAddOutput(qf_setup_rhs, "true solution", num_comp_u, CEED_EVAL_NONE);
1289b072555Sjeremylt     CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP);
129e83e87a5Sjeremylt 
130e83e87a5Sjeremylt     // Create the operator that builds the RHS and true solution
1319b072555Sjeremylt     CeedOperatorCreate(ceed, qf_setup_rhs, NULL, NULL, &op_setup_rhs);
1322b730f8bSJeremy L Thompson     CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
133356036faSJeremy L Thompson     CeedOperatorSetField(op_setup_rhs, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
134356036faSJeremy L Thompson     CeedOperatorSetField(op_setup_rhs, "true solution", elem_restr_u_i, CEED_BASIS_NONE, *target);
1352b730f8bSJeremy L Thompson     CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
136e83e87a5Sjeremylt 
137e83e87a5Sjeremylt     // Set up the libCEED context
1389b072555Sjeremylt     CeedQFunctionContext ctx_rhs_setup;
1399b072555Sjeremylt     CeedQFunctionContextCreate(ceed, &ctx_rhs_setup);
1409b072555Sjeremylt     CeedScalar rhs_setup_data[2] = {R, l};
1412b730f8bSJeremy L Thompson     CeedQFunctionContextSetData(ctx_rhs_setup, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof rhs_setup_data, &rhs_setup_data);
1429b072555Sjeremylt     CeedQFunctionSetContext(qf_setup_rhs, ctx_rhs_setup);
1439b072555Sjeremylt     CeedQFunctionContextDestroy(&ctx_rhs_setup);
144e83e87a5Sjeremylt 
145e83e87a5Sjeremylt     // Setup RHS and target
1469b072555Sjeremylt     CeedOperatorApply(op_setup_rhs, x_coord, rhs_ceed, CEED_REQUEST_IMMEDIATE);
147e83e87a5Sjeremylt 
148e83e87a5Sjeremylt     // Cleanup
1499b072555Sjeremylt     CeedQFunctionDestroy(&qf_setup_rhs);
1509b072555Sjeremylt     CeedOperatorDestroy(&op_setup_rhs);
151e83e87a5Sjeremylt   }
152e83e87a5Sjeremylt 
153e83e87a5Sjeremylt   // Cleanup
1549b072555Sjeremylt   CeedQFunctionDestroy(&qf_setup_geo);
1559b072555Sjeremylt   CeedOperatorDestroy(&op_setup_geo);
1569b072555Sjeremylt   CeedVectorDestroy(&x_coord);
157e83e87a5Sjeremylt 
158e83e87a5Sjeremylt   // Save libCEED data required for level
1592b730f8bSJeremy L Thompson   data->basis_x         = basis_x;
1602b730f8bSJeremy L Thompson   data->basis_u         = basis_u;
1619b072555Sjeremylt   data->elem_restr_x    = elem_restr_x;
1629b072555Sjeremylt   data->elem_restr_u    = elem_restr_u;
1639b072555Sjeremylt   data->elem_restr_u_i  = elem_restr_u_i;
1649b072555Sjeremylt   data->elem_restr_qd_i = elem_restr_qd_i;
1659b072555Sjeremylt   data->qf_apply        = qf_apply;
1669b072555Sjeremylt   data->op_apply        = op_apply;
1679b072555Sjeremylt   data->q_data          = q_data;
1689b072555Sjeremylt   data->x_ceed          = x_ceed;
1699b072555Sjeremylt   data->y_ceed          = y_ceed;
170d4d45553Srezgarshakeri   data->q_data_size     = q_data_size;
171ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
172e83e87a5Sjeremylt };
173e83e87a5Sjeremylt 
174e83e87a5Sjeremylt // -----------------------------------------------------------------------------
175e83e87a5Sjeremylt // Setup libCEED level transfer operator objects
176e83e87a5Sjeremylt // -----------------------------------------------------------------------------
1772b730f8bSJeremy L Thompson PetscErrorCode CeedLevelTransferSetup(DM dm, Ceed ceed, CeedInt level, CeedInt num_comp_u, CeedData *data, BPData bp_data, Vec fine_mult) {
178751eb813Srezgarshakeri   PetscFunctionBeginUser;
179e83e87a5Sjeremylt   // Restriction - Fine to corse
1809b072555Sjeremylt   CeedOperator op_restrict;
181e83e87a5Sjeremylt   // Interpolation - Corse to fine
1829b072555Sjeremylt   CeedOperator op_prolong;
1831c9a79dbSrezgarshakeri   // Coarse grid operator
1841c9a79dbSrezgarshakeri   CeedOperator op_apply;
1851c9a79dbSrezgarshakeri   // Basis
1861c9a79dbSrezgarshakeri   CeedBasis basis_u;
1872b730f8bSJeremy L Thompson   PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, bp_data, &basis_u));
188e83e87a5Sjeremylt 
1891c9a79dbSrezgarshakeri   // ---------------------------------------------------------------------------
1901c9a79dbSrezgarshakeri   // Coarse Grid, Prolongation, and Restriction Operators
1911c9a79dbSrezgarshakeri   // ---------------------------------------------------------------------------
1921c9a79dbSrezgarshakeri   // Create the Operators that compute the prolongation and
1931c9a79dbSrezgarshakeri   //   restriction between the p-multigrid levels and the coarse grid eval.
1941c9a79dbSrezgarshakeri   // ---------------------------------------------------------------------------
1951c9a79dbSrezgarshakeri   // Place in libCEED array
1961c9a79dbSrezgarshakeri   PetscMemType m_mem_type;
197179e5961SZach Atkins   PetscCall(VecReadP2C(fine_mult, &m_mem_type, data[level]->x_ceed));
198e83e87a5Sjeremylt 
1992b730f8bSJeremy L Thompson   CeedOperatorMultigridLevelCreate(data[level]->op_apply, data[level]->x_ceed, data[level - 1]->elem_restr_u, basis_u, &op_apply, &op_prolong,
2002b730f8bSJeremy L Thompson                                    &op_restrict);
201e83e87a5Sjeremylt 
2021c9a79dbSrezgarshakeri   // Restore PETSc vector
203179e5961SZach Atkins   PetscCall(VecReadC2P(data[level]->x_ceed, m_mem_type, fine_mult));
2042b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(fine_mult));
2051c9a79dbSrezgarshakeri   // -- Save libCEED data
2061c9a79dbSrezgarshakeri   data[level - 1]->op_apply = op_apply;
2071c9a79dbSrezgarshakeri   data[level]->op_prolong   = op_prolong;
2081c9a79dbSrezgarshakeri   data[level]->op_restrict  = op_restrict;
2091c9a79dbSrezgarshakeri 
2101c9a79dbSrezgarshakeri   CeedBasisDestroy(&basis_u);
211ee4ca9cbSJames Wright   PetscFunctionReturn(PETSC_SUCCESS);
212e83e87a5Sjeremylt };
213e83e87a5Sjeremylt 
214e83e87a5Sjeremylt // -----------------------------------------------------------------------------
215