1e83e87a5Sjeremylt #include "../include/libceedsetup.h" 22b730f8bSJeremy L Thompson 32b730f8bSJeremy L Thompson #include <stdio.h> 42b730f8bSJeremy L Thompson 5e83e87a5Sjeremylt #include "../include/petscutils.h" 6e83e87a5Sjeremylt 7e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 8e83e87a5Sjeremylt // Destroy libCEED operator objects 9e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 10e83e87a5Sjeremylt PetscErrorCode CeedDataDestroy(CeedInt i, CeedData data) { 115dfaedb8SJed Brown PetscFunctionBeginUser; 122b730f8bSJeremy L Thompson 139b072555Sjeremylt CeedVectorDestroy(&data->q_data); 149b072555Sjeremylt CeedVectorDestroy(&data->x_ceed); 159b072555Sjeremylt CeedVectorDestroy(&data->y_ceed); 169b072555Sjeremylt CeedBasisDestroy(&data->basis_x); 179b072555Sjeremylt CeedBasisDestroy(&data->basis_u); 189b072555Sjeremylt CeedElemRestrictionDestroy(&data->elem_restr_u); 199b072555Sjeremylt CeedElemRestrictionDestroy(&data->elem_restr_x); 209b072555Sjeremylt CeedElemRestrictionDestroy(&data->elem_restr_u_i); 219b072555Sjeremylt CeedElemRestrictionDestroy(&data->elem_restr_qd_i); 229b072555Sjeremylt CeedQFunctionDestroy(&data->qf_apply); 239b072555Sjeremylt CeedOperatorDestroy(&data->op_apply); 24e83e87a5Sjeremylt if (i > 0) { 259b072555Sjeremylt CeedOperatorDestroy(&data->op_prolong); 269b072555Sjeremylt CeedOperatorDestroy(&data->op_restrict); 27e83e87a5Sjeremylt } 282b730f8bSJeremy L Thompson PetscCall(PetscFree(data)); 29e83e87a5Sjeremylt 30*ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 31e83e87a5Sjeremylt }; 32e83e87a5Sjeremylt 33e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 34e83e87a5Sjeremylt // Set up libCEED for a given degree 35e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 362b730f8bSJeremy L Thompson PetscErrorCode SetupLibceedByDegree(DM dm, Ceed ceed, CeedInt degree, CeedInt topo_dim, CeedInt q_extra, PetscInt num_comp_x, PetscInt num_comp_u, 372b730f8bSJeremy L Thompson PetscInt g_size, PetscInt xl_size, BPData bp_data, CeedData data, PetscBool setup_rhs, CeedVector rhs_ceed, 38e83e87a5Sjeremylt CeedVector *target) { 399b072555Sjeremylt DM dm_coord; 40e83e87a5Sjeremylt Vec coords; 419b072555Sjeremylt const PetscScalar *coord_array; 429b072555Sjeremylt CeedBasis basis_x, basis_u; 439b072555Sjeremylt CeedElemRestriction elem_restr_x, elem_restr_u, elem_restr_u_i, elem_restr_qd_i; 449b072555Sjeremylt CeedQFunction qf_setup_geo, qf_apply; 459b072555Sjeremylt CeedOperator op_setup_geo, op_apply; 469b072555Sjeremylt CeedVector x_coord, q_data, x_ceed, y_ceed; 472b730f8bSJeremy L Thompson CeedInt num_qpts, c_start, c_end, num_elem, q_data_size = bp_data.q_data_size; 482b730f8bSJeremy L Thompson CeedScalar R = 1; // radius of the sphere 492b730f8bSJeremy L Thompson CeedScalar l = 1.0 / PetscSqrtReal(3.0); // half edge of the inscribed cube 50e83e87a5Sjeremylt 515dfaedb8SJed Brown PetscFunctionBeginUser; 522b730f8bSJeremy L Thompson PetscCall(DMGetCoordinateDM(dm, &dm_coord)); 53129d8736Srezgarshakeri 54129d8736Srezgarshakeri // CEED bases 552b730f8bSJeremy L Thompson PetscCall(CreateBasisFromPlex(ceed, dm_coord, 0, 0, 0, 0, bp_data, &basis_x)); 562b730f8bSJeremy L Thompson PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, bp_data, &basis_u)); 57129d8736Srezgarshakeri 58129d8736Srezgarshakeri // CEED restrictions 592b730f8bSJeremy L Thompson PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, 0, 0, 0, &elem_restr_x)); 602b730f8bSJeremy L Thompson PetscCall(CreateRestrictionFromPlex(ceed, dm, 0, 0, 0, &elem_restr_u)); 61e83e87a5Sjeremylt 622b730f8bSJeremy L Thompson PetscCall(DMPlexGetHeightStratum(dm, 0, &c_start, &c_end)); 639b072555Sjeremylt num_elem = c_end - c_start; 64129d8736Srezgarshakeri CeedBasisGetNumQuadraturePoints(basis_u, &num_qpts); 65e83e87a5Sjeremylt 662b730f8bSJeremy 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); 672b730f8bSJeremy 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); 68e83e87a5Sjeremylt 69e83e87a5Sjeremylt // Element coordinates 702b730f8bSJeremy L Thompson PetscCall(DMGetCoordinatesLocal(dm, &coords)); 712b730f8bSJeremy L Thompson PetscCall(VecGetArrayRead(coords, &coord_array)); 72e83e87a5Sjeremylt 739b072555Sjeremylt CeedElemRestrictionCreateVector(elem_restr_x, &x_coord, NULL); 742b730f8bSJeremy L Thompson CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, (PetscScalar *)coord_array); 752b730f8bSJeremy L Thompson PetscCall(VecRestoreArrayRead(coords, &coord_array)); 76e83e87a5Sjeremylt 77e83e87a5Sjeremylt // Create the persistent vectors that will be needed in setup and apply 789b072555Sjeremylt CeedVectorCreate(ceed, q_data_size * num_elem * num_qpts, &q_data); 799b072555Sjeremylt CeedVectorCreate(ceed, xl_size, &x_ceed); 809b072555Sjeremylt CeedVectorCreate(ceed, xl_size, &y_ceed); 81e83e87a5Sjeremylt 82e83e87a5Sjeremylt // Create the QFunction that builds the context data 832b730f8bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_geo, bp_data.setup_geo_loc, &qf_setup_geo); 849b072555Sjeremylt CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP); 859b072555Sjeremylt CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x * topo_dim, CEED_EVAL_GRAD); 869b072555Sjeremylt CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT); 87a61c78d6SJeremy L Thompson CeedQFunctionAddOutput(qf_setup_geo, "qdata", q_data_size, CEED_EVAL_NONE); 88e83e87a5Sjeremylt 89e83e87a5Sjeremylt // Create the operator that builds the quadrature data 909b072555Sjeremylt CeedOperatorCreate(ceed, qf_setup_geo, NULL, NULL, &op_setup_geo); 912b730f8bSJeremy L Thompson CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE); 922b730f8bSJeremy L Thompson CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE); 932b730f8bSJeremy L Thompson CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE); 942b730f8bSJeremy L Thompson CeedOperatorSetField(op_setup_geo, "qdata", elem_restr_qd_i, CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 95e83e87a5Sjeremylt 969b072555Sjeremylt // Setup q_data 979b072555Sjeremylt CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE); 98e83e87a5Sjeremylt 99e83e87a5Sjeremylt // Set up PDE operator 1009b072555Sjeremylt CeedInt in_scale = bp_data.in_mode == CEED_EVAL_GRAD ? topo_dim : 1; 1019b072555Sjeremylt CeedInt out_scale = bp_data.out_mode == CEED_EVAL_GRAD ? topo_dim : 1; 1022b730f8bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, bp_data.apply, bp_data.apply_loc, &qf_apply); 1039b072555Sjeremylt CeedQFunctionAddInput(qf_apply, "u", num_comp_u * in_scale, bp_data.in_mode); 104a61c78d6SJeremy L Thompson CeedQFunctionAddInput(qf_apply, "qdata", q_data_size, CEED_EVAL_NONE); 1059b072555Sjeremylt CeedQFunctionAddOutput(qf_apply, "v", num_comp_u * out_scale, bp_data.out_mode); 106e83e87a5Sjeremylt 107e83e87a5Sjeremylt // Create the mass or diff operator 1089b072555Sjeremylt CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply); 1099b072555Sjeremylt CeedOperatorSetField(op_apply, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE); 1102b730f8bSJeremy L Thompson CeedOperatorSetField(op_apply, "qdata", elem_restr_qd_i, CEED_BASIS_COLLOCATED, q_data); 1119b072555Sjeremylt CeedOperatorSetField(op_apply, "v", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE); 112e83e87a5Sjeremylt 113e83e87a5Sjeremylt // Set up RHS if needed 114e83e87a5Sjeremylt if (setup_rhs) { 1159b072555Sjeremylt CeedQFunction qf_setup_rhs; 1169b072555Sjeremylt CeedOperator op_setup_rhs; 1179b072555Sjeremylt CeedVectorCreate(ceed, num_elem * num_qpts * num_comp_u, target); 118e83e87a5Sjeremylt // Create the q-function that sets up the RHS and true solution 1192b730f8bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_rhs, bp_data.setup_rhs_loc, &qf_setup_rhs); 1209b072555Sjeremylt CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP); 121a61c78d6SJeremy L Thompson CeedQFunctionAddInput(qf_setup_rhs, "qdata", q_data_size, CEED_EVAL_NONE); 1222b730f8bSJeremy L Thompson CeedQFunctionAddOutput(qf_setup_rhs, "true solution", num_comp_u, CEED_EVAL_NONE); 1239b072555Sjeremylt CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP); 124e83e87a5Sjeremylt 125e83e87a5Sjeremylt // Create the operator that builds the RHS and true solution 1269b072555Sjeremylt CeedOperatorCreate(ceed, qf_setup_rhs, NULL, NULL, &op_setup_rhs); 1272b730f8bSJeremy L Thompson CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE); 1282b730f8bSJeremy L Thompson CeedOperatorSetField(op_setup_rhs, "qdata", elem_restr_qd_i, CEED_BASIS_COLLOCATED, q_data); 1292b730f8bSJeremy L Thompson CeedOperatorSetField(op_setup_rhs, "true solution", elem_restr_u_i, CEED_BASIS_COLLOCATED, *target); 1302b730f8bSJeremy L Thompson CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE); 131e83e87a5Sjeremylt 132e83e87a5Sjeremylt // Set up the libCEED context 1339b072555Sjeremylt CeedQFunctionContext ctx_rhs_setup; 1349b072555Sjeremylt CeedQFunctionContextCreate(ceed, &ctx_rhs_setup); 1359b072555Sjeremylt CeedScalar rhs_setup_data[2] = {R, l}; 1362b730f8bSJeremy L Thompson CeedQFunctionContextSetData(ctx_rhs_setup, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof rhs_setup_data, &rhs_setup_data); 1379b072555Sjeremylt CeedQFunctionSetContext(qf_setup_rhs, ctx_rhs_setup); 1389b072555Sjeremylt CeedQFunctionContextDestroy(&ctx_rhs_setup); 139e83e87a5Sjeremylt 140e83e87a5Sjeremylt // Setup RHS and target 1419b072555Sjeremylt CeedOperatorApply(op_setup_rhs, x_coord, rhs_ceed, CEED_REQUEST_IMMEDIATE); 142e83e87a5Sjeremylt 143e83e87a5Sjeremylt // Cleanup 1449b072555Sjeremylt CeedQFunctionDestroy(&qf_setup_rhs); 1459b072555Sjeremylt CeedOperatorDestroy(&op_setup_rhs); 146e83e87a5Sjeremylt } 147e83e87a5Sjeremylt 148e83e87a5Sjeremylt // Cleanup 1499b072555Sjeremylt CeedQFunctionDestroy(&qf_setup_geo); 1509b072555Sjeremylt CeedOperatorDestroy(&op_setup_geo); 1519b072555Sjeremylt CeedVectorDestroy(&x_coord); 152e83e87a5Sjeremylt 153e83e87a5Sjeremylt // Save libCEED data required for level 1542b730f8bSJeremy L Thompson data->basis_x = basis_x; 1552b730f8bSJeremy L Thompson data->basis_u = basis_u; 1569b072555Sjeremylt data->elem_restr_x = elem_restr_x; 1579b072555Sjeremylt data->elem_restr_u = elem_restr_u; 1589b072555Sjeremylt data->elem_restr_u_i = elem_restr_u_i; 1599b072555Sjeremylt data->elem_restr_qd_i = elem_restr_qd_i; 1609b072555Sjeremylt data->qf_apply = qf_apply; 1619b072555Sjeremylt data->op_apply = op_apply; 1629b072555Sjeremylt data->q_data = q_data; 1639b072555Sjeremylt data->x_ceed = x_ceed; 1649b072555Sjeremylt data->y_ceed = y_ceed; 165d4d45553Srezgarshakeri data->q_data_size = q_data_size; 166e83e87a5Sjeremylt 167*ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 168e83e87a5Sjeremylt }; 169e83e87a5Sjeremylt 170e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 171e83e87a5Sjeremylt // Setup libCEED level transfer operator objects 172e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 1732b730f8bSJeremy L Thompson PetscErrorCode CeedLevelTransferSetup(DM dm, Ceed ceed, CeedInt level, CeedInt num_comp_u, CeedData *data, BPData bp_data, Vec fine_mult) { 174751eb813Srezgarshakeri PetscFunctionBeginUser; 175e83e87a5Sjeremylt // Restriction - Fine to corse 1769b072555Sjeremylt CeedOperator op_restrict; 177e83e87a5Sjeremylt // Interpolation - Corse to fine 1789b072555Sjeremylt CeedOperator op_prolong; 1791c9a79dbSrezgarshakeri // Coarse grid operator 1801c9a79dbSrezgarshakeri CeedOperator op_apply; 1811c9a79dbSrezgarshakeri // Basis 1821c9a79dbSrezgarshakeri CeedBasis basis_u; 1832b730f8bSJeremy L Thompson PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, bp_data, &basis_u)); 184e83e87a5Sjeremylt 1851c9a79dbSrezgarshakeri // --------------------------------------------------------------------------- 1861c9a79dbSrezgarshakeri // Coarse Grid, Prolongation, and Restriction Operators 1871c9a79dbSrezgarshakeri // --------------------------------------------------------------------------- 1881c9a79dbSrezgarshakeri // Create the Operators that compute the prolongation and 1891c9a79dbSrezgarshakeri // restriction between the p-multigrid levels and the coarse grid eval. 1901c9a79dbSrezgarshakeri // --------------------------------------------------------------------------- 1911c9a79dbSrezgarshakeri // Place in libCEED array 1921c9a79dbSrezgarshakeri const PetscScalar *m; 1931c9a79dbSrezgarshakeri PetscMemType m_mem_type; 1942b730f8bSJeremy L Thompson PetscCall(VecGetArrayReadAndMemType(fine_mult, &m, &m_mem_type)); 1952b730f8bSJeremy L Thompson CeedVectorSetArray(data[level]->x_ceed, MemTypeP2C(m_mem_type), CEED_USE_POINTER, (CeedScalar *)m); 196e83e87a5Sjeremylt 1972b730f8bSJeremy L Thompson CeedOperatorMultigridLevelCreate(data[level]->op_apply, data[level]->x_ceed, data[level - 1]->elem_restr_u, basis_u, &op_apply, &op_prolong, 1982b730f8bSJeremy L Thompson &op_restrict); 199e83e87a5Sjeremylt 2001c9a79dbSrezgarshakeri // Restore PETSc vector 2012b730f8bSJeremy L Thompson CeedVectorTakeArray(data[level]->x_ceed, MemTypeP2C(m_mem_type), (CeedScalar **)&m); 2022b730f8bSJeremy L Thompson PetscCall(VecRestoreArrayReadAndMemType(fine_mult, &m)); 2032b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(fine_mult)); 2041c9a79dbSrezgarshakeri // -- Save libCEED data 2051c9a79dbSrezgarshakeri data[level - 1]->op_apply = op_apply; 2061c9a79dbSrezgarshakeri data[level]->op_prolong = op_prolong; 2071c9a79dbSrezgarshakeri data[level]->op_restrict = op_restrict; 2081c9a79dbSrezgarshakeri 2091c9a79dbSrezgarshakeri CeedBasisDestroy(&basis_u); 210*ee4ca9cbSJames Wright PetscFunctionReturn(PETSC_SUCCESS); 211e83e87a5Sjeremylt }; 212e83e87a5Sjeremylt 213e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 214