1636cccdbSjeremylt #include <stdio.h> 2e83e87a5Sjeremylt #include "../include/libceedsetup.h" 3e83e87a5Sjeremylt #include "../include/petscutils.h" 4e83e87a5Sjeremylt 5e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 6e83e87a5Sjeremylt // Destroy libCEED operator objects 7e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 8e83e87a5Sjeremylt PetscErrorCode CeedDataDestroy(CeedInt i, CeedData data) { 9e83e87a5Sjeremylt int ierr; 10e83e87a5Sjeremylt 115dfaedb8SJed Brown PetscFunctionBeginUser; 129b072555Sjeremylt CeedVectorDestroy(&data->q_data); 139b072555Sjeremylt CeedVectorDestroy(&data->x_ceed); 149b072555Sjeremylt CeedVectorDestroy(&data->y_ceed); 159b072555Sjeremylt CeedBasisDestroy(&data->basis_x); 169b072555Sjeremylt CeedBasisDestroy(&data->basis_u); 179b072555Sjeremylt CeedElemRestrictionDestroy(&data->elem_restr_u); 189b072555Sjeremylt CeedElemRestrictionDestroy(&data->elem_restr_x); 199b072555Sjeremylt CeedElemRestrictionDestroy(&data->elem_restr_u_i); 209b072555Sjeremylt CeedElemRestrictionDestroy(&data->elem_restr_qd_i); 219b072555Sjeremylt CeedQFunctionDestroy(&data->qf_apply); 229b072555Sjeremylt CeedOperatorDestroy(&data->op_apply); 23e83e87a5Sjeremylt if (i > 0) { 249b072555Sjeremylt CeedOperatorDestroy(&data->op_prolong); 259b072555Sjeremylt CeedOperatorDestroy(&data->op_restrict); 26e83e87a5Sjeremylt } 27e83e87a5Sjeremylt ierr = PetscFree(data); CHKERRQ(ierr); 28e83e87a5Sjeremylt 29e83e87a5Sjeremylt PetscFunctionReturn(0); 30e83e87a5Sjeremylt }; 31e83e87a5Sjeremylt 32e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 33e83e87a5Sjeremylt // Set up libCEED for a given degree 34e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 35e83e87a5Sjeremylt PetscErrorCode SetupLibceedByDegree(DM dm, Ceed ceed, CeedInt degree, 369b072555Sjeremylt CeedInt topo_dim, CeedInt q_extra, 379b072555Sjeremylt PetscInt num_comp_x, PetscInt num_comp_u, 389b072555Sjeremylt PetscInt g_size, PetscInt xl_size, 399b072555Sjeremylt BPData bp_data, CeedData data, 409b072555Sjeremylt PetscBool setup_rhs, CeedVector rhs_ceed, 41e83e87a5Sjeremylt CeedVector *target) { 42e83e87a5Sjeremylt int ierr; 439b072555Sjeremylt DM dm_coord; 44e83e87a5Sjeremylt Vec coords; 459b072555Sjeremylt const PetscScalar *coord_array; 469b072555Sjeremylt CeedBasis basis_x, basis_u; 479b072555Sjeremylt CeedElemRestriction elem_restr_x, elem_restr_u, elem_restr_u_i, elem_restr_qd_i; 489b072555Sjeremylt CeedQFunction qf_setup_geo, qf_apply; 499b072555Sjeremylt CeedOperator op_setup_geo, op_apply; 509b072555Sjeremylt CeedVector x_coord, q_data, x_ceed, y_ceed; 51129d8736Srezgarshakeri CeedInt num_qpts, c_start, c_end, num_elem, 529b072555Sjeremylt q_data_size = bp_data.q_data_size; 53e83e87a5Sjeremylt CeedScalar R = 1, // radius of the sphere 54e83e87a5Sjeremylt l = 1.0/PetscSqrtReal(3.0); // half edge of the inscribed cube 55e83e87a5Sjeremylt 565dfaedb8SJed Brown PetscFunctionBeginUser; 579b072555Sjeremylt ierr = DMGetCoordinateDM(dm, &dm_coord); CHKERRQ(ierr); 58129d8736Srezgarshakeri 59129d8736Srezgarshakeri // CEED bases 60*d4d45553Srezgarshakeri ierr = CreateBasisFromPlex(ceed, dm_coord, 0, 0, 0, 0, bp_data, &basis_x); 61*d4d45553Srezgarshakeri CHKERRQ(ierr); 62*d4d45553Srezgarshakeri ierr = CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, bp_data, &basis_u); 63*d4d45553Srezgarshakeri CHKERRQ(ierr); 64129d8736Srezgarshakeri 65129d8736Srezgarshakeri // CEED restrictions 667ed3e4cdSJeremy L Thompson ierr = CreateRestrictionFromPlex(ceed, dm_coord, 0, 0, 0, &elem_restr_x); 67e83e87a5Sjeremylt CHKERRQ(ierr); 687ed3e4cdSJeremy L Thompson ierr = CreateRestrictionFromPlex(ceed, dm, 0, 0, 0, &elem_restr_u); 69e83e87a5Sjeremylt CHKERRQ(ierr); 70e83e87a5Sjeremylt 719b072555Sjeremylt ierr = DMPlexGetHeightStratum(dm, 0, &c_start, &c_end); CHKERRQ(ierr); 729b072555Sjeremylt num_elem = c_end - c_start; 73129d8736Srezgarshakeri CeedBasisGetNumQuadraturePoints(basis_u, &num_qpts); 74e83e87a5Sjeremylt 759b072555Sjeremylt CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, num_comp_u, 769b072555Sjeremylt num_comp_u*num_elem*num_qpts, 779b072555Sjeremylt CEED_STRIDES_BACKEND, &elem_restr_u_i); 789b072555Sjeremylt CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, q_data_size, 799b072555Sjeremylt q_data_size*num_elem*num_qpts, 809b072555Sjeremylt CEED_STRIDES_BACKEND, &elem_restr_qd_i); 81e83e87a5Sjeremylt 82e83e87a5Sjeremylt // Element coordinates 83e83e87a5Sjeremylt ierr = DMGetCoordinatesLocal(dm, &coords); CHKERRQ(ierr); 849b072555Sjeremylt ierr = VecGetArrayRead(coords, &coord_array); CHKERRQ(ierr); 85e83e87a5Sjeremylt 869b072555Sjeremylt CeedElemRestrictionCreateVector(elem_restr_x, &x_coord, NULL); 879b072555Sjeremylt CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, 889b072555Sjeremylt (PetscScalar *)coord_array); 899b072555Sjeremylt ierr = VecRestoreArrayRead(coords, &coord_array); 90e83e87a5Sjeremylt 91e83e87a5Sjeremylt // Create the persistent vectors that will be needed in setup and apply 929b072555Sjeremylt CeedVectorCreate(ceed, q_data_size*num_elem*num_qpts, &q_data); 939b072555Sjeremylt CeedVectorCreate(ceed, xl_size, &x_ceed); 949b072555Sjeremylt CeedVectorCreate(ceed, xl_size, &y_ceed); 95e83e87a5Sjeremylt 96e83e87a5Sjeremylt // Create the QFunction that builds the context data 979b072555Sjeremylt CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_geo, bp_data.setup_geo_loc, 989b072555Sjeremylt &qf_setup_geo); 999b072555Sjeremylt CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP); 1009b072555Sjeremylt CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x*topo_dim, CEED_EVAL_GRAD); 1019b072555Sjeremylt CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT); 102a61c78d6SJeremy L Thompson CeedQFunctionAddOutput(qf_setup_geo, "qdata", q_data_size, CEED_EVAL_NONE); 103e83e87a5Sjeremylt 104e83e87a5Sjeremylt // Create the operator that builds the quadrature data 1059b072555Sjeremylt CeedOperatorCreate(ceed, qf_setup_geo, NULL, NULL, &op_setup_geo); 1069b072555Sjeremylt CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x, 1079b072555Sjeremylt CEED_VECTOR_ACTIVE); 1089b072555Sjeremylt CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x, 1099b072555Sjeremylt CEED_VECTOR_ACTIVE); 1109b072555Sjeremylt CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, 111e83e87a5Sjeremylt CEED_VECTOR_NONE); 112a61c78d6SJeremy L Thompson CeedOperatorSetField(op_setup_geo, "qdata", elem_restr_qd_i, 113e83e87a5Sjeremylt CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 114e83e87a5Sjeremylt 1159b072555Sjeremylt // Setup q_data 1169b072555Sjeremylt CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE); 117e83e87a5Sjeremylt 118e83e87a5Sjeremylt // Set up PDE operator 1199b072555Sjeremylt CeedInt in_scale = bp_data.in_mode == CEED_EVAL_GRAD ? topo_dim : 1; 1209b072555Sjeremylt CeedInt out_scale = bp_data.out_mode == CEED_EVAL_GRAD ? topo_dim : 1; 1219b072555Sjeremylt CeedQFunctionCreateInterior(ceed, 1, bp_data.apply, bp_data.apply_loc, 1229b072555Sjeremylt &qf_apply); 1239b072555Sjeremylt CeedQFunctionAddInput(qf_apply, "u", num_comp_u*in_scale, bp_data.in_mode); 124a61c78d6SJeremy L Thompson CeedQFunctionAddInput(qf_apply, "qdata", q_data_size, CEED_EVAL_NONE); 1259b072555Sjeremylt CeedQFunctionAddOutput(qf_apply, "v", num_comp_u*out_scale, bp_data.out_mode); 126e83e87a5Sjeremylt 127e83e87a5Sjeremylt // Create the mass or diff operator 1289b072555Sjeremylt CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply); 1299b072555Sjeremylt CeedOperatorSetField(op_apply, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE); 130a61c78d6SJeremy L Thompson CeedOperatorSetField(op_apply, "qdata", elem_restr_qd_i, CEED_BASIS_COLLOCATED, 1319b072555Sjeremylt q_data); 1329b072555Sjeremylt CeedOperatorSetField(op_apply, "v", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE); 133e83e87a5Sjeremylt 134e83e87a5Sjeremylt // Set up RHS if needed 135e83e87a5Sjeremylt if (setup_rhs) { 1369b072555Sjeremylt CeedQFunction qf_setup_rhs; 1379b072555Sjeremylt CeedOperator op_setup_rhs; 1389b072555Sjeremylt CeedVectorCreate(ceed, num_elem*num_qpts*num_comp_u, target); 139e83e87a5Sjeremylt // Create the q-function that sets up the RHS and true solution 1409b072555Sjeremylt CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_rhs, bp_data.setup_rhs_loc, 1419b072555Sjeremylt &qf_setup_rhs); 1429b072555Sjeremylt CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP); 143a61c78d6SJeremy L Thompson CeedQFunctionAddInput(qf_setup_rhs, "qdata", q_data_size, CEED_EVAL_NONE); 144a61c78d6SJeremy L Thompson CeedQFunctionAddOutput(qf_setup_rhs, "true solution", num_comp_u, 145a61c78d6SJeremy L Thompson CEED_EVAL_NONE); 1469b072555Sjeremylt CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP); 147e83e87a5Sjeremylt 148e83e87a5Sjeremylt // Create the operator that builds the RHS and true solution 1499b072555Sjeremylt CeedOperatorCreate(ceed, qf_setup_rhs, NULL, NULL, &op_setup_rhs); 1509b072555Sjeremylt CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x, 1519b072555Sjeremylt CEED_VECTOR_ACTIVE); 152a61c78d6SJeremy L Thompson CeedOperatorSetField(op_setup_rhs, "qdata", elem_restr_qd_i, 153a61c78d6SJeremy L Thompson CEED_BASIS_COLLOCATED, q_data); 154a61c78d6SJeremy L Thompson CeedOperatorSetField(op_setup_rhs, "true solution", elem_restr_u_i, 155e83e87a5Sjeremylt CEED_BASIS_COLLOCATED, *target); 1569b072555Sjeremylt CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u, 157e83e87a5Sjeremylt CEED_VECTOR_ACTIVE); 158e83e87a5Sjeremylt 159e83e87a5Sjeremylt // Set up the libCEED context 1609b072555Sjeremylt CeedQFunctionContext ctx_rhs_setup; 1619b072555Sjeremylt CeedQFunctionContextCreate(ceed, &ctx_rhs_setup); 1629b072555Sjeremylt CeedScalar rhs_setup_data[2] = {R, l}; 1639b072555Sjeremylt CeedQFunctionContextSetData(ctx_rhs_setup, CEED_MEM_HOST, CEED_COPY_VALUES, 1649b072555Sjeremylt sizeof rhs_setup_data, &rhs_setup_data); 1659b072555Sjeremylt CeedQFunctionSetContext(qf_setup_rhs, ctx_rhs_setup); 1669b072555Sjeremylt CeedQFunctionContextDestroy(&ctx_rhs_setup); 167e83e87a5Sjeremylt 168e83e87a5Sjeremylt // Setup RHS and target 1699b072555Sjeremylt CeedOperatorApply(op_setup_rhs, x_coord, rhs_ceed, CEED_REQUEST_IMMEDIATE); 170e83e87a5Sjeremylt 171e83e87a5Sjeremylt // Cleanup 1729b072555Sjeremylt CeedQFunctionDestroy(&qf_setup_rhs); 1739b072555Sjeremylt CeedOperatorDestroy(&op_setup_rhs); 174e83e87a5Sjeremylt } 175e83e87a5Sjeremylt 176e83e87a5Sjeremylt // Cleanup 1779b072555Sjeremylt CeedQFunctionDestroy(&qf_setup_geo); 1789b072555Sjeremylt CeedOperatorDestroy(&op_setup_geo); 1799b072555Sjeremylt CeedVectorDestroy(&x_coord); 180e83e87a5Sjeremylt 181e83e87a5Sjeremylt // Save libCEED data required for level 1829b072555Sjeremylt data->basis_x = basis_x; data->basis_u = basis_u; 1839b072555Sjeremylt data->elem_restr_x = elem_restr_x; 1849b072555Sjeremylt data->elem_restr_u = elem_restr_u; 1859b072555Sjeremylt data->elem_restr_u_i = elem_restr_u_i; 1869b072555Sjeremylt data->elem_restr_qd_i = elem_restr_qd_i; 1879b072555Sjeremylt data->qf_apply = qf_apply; 1889b072555Sjeremylt data->op_apply = op_apply; 1899b072555Sjeremylt data->q_data = q_data; 1909b072555Sjeremylt data->x_ceed = x_ceed; 1919b072555Sjeremylt data->y_ceed = y_ceed; 192*d4d45553Srezgarshakeri data->q_data_size = q_data_size; 193e83e87a5Sjeremylt 194e83e87a5Sjeremylt PetscFunctionReturn(0); 195e83e87a5Sjeremylt }; 196e83e87a5Sjeremylt 197e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 198e83e87a5Sjeremylt // Setup libCEED level transfer operator objects 199e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 2001c9a79dbSrezgarshakeri PetscErrorCode CeedLevelTransferSetup(DM dm, Ceed ceed, CeedInt level, 2019b072555Sjeremylt CeedInt num_comp_u, CeedData *data, 202*d4d45553Srezgarshakeri BPData bp_data, Vec fine_mult) { 2031c9a79dbSrezgarshakeri int ierr; 204e83e87a5Sjeremylt 205751eb813Srezgarshakeri PetscFunctionBeginUser; 206e83e87a5Sjeremylt // Restriction - Fine to corse 2079b072555Sjeremylt CeedOperator op_restrict; 208e83e87a5Sjeremylt // Interpolation - Corse to fine 2099b072555Sjeremylt CeedOperator op_prolong; 2101c9a79dbSrezgarshakeri // Coarse grid operator 2111c9a79dbSrezgarshakeri CeedOperator op_apply; 2121c9a79dbSrezgarshakeri // Basis 2131c9a79dbSrezgarshakeri CeedBasis basis_u; 214*d4d45553Srezgarshakeri ierr = CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, bp_data, &basis_u); 2151c9a79dbSrezgarshakeri CHKERRQ(ierr); 216e83e87a5Sjeremylt 2171c9a79dbSrezgarshakeri // --------------------------------------------------------------------------- 2181c9a79dbSrezgarshakeri // Coarse Grid, Prolongation, and Restriction Operators 2191c9a79dbSrezgarshakeri // --------------------------------------------------------------------------- 2201c9a79dbSrezgarshakeri // Create the Operators that compute the prolongation and 2211c9a79dbSrezgarshakeri // restriction between the p-multigrid levels and the coarse grid eval. 2221c9a79dbSrezgarshakeri // --------------------------------------------------------------------------- 2231c9a79dbSrezgarshakeri // Place in libCEED array 2241c9a79dbSrezgarshakeri const PetscScalar *m; 2251c9a79dbSrezgarshakeri PetscMemType m_mem_type; 2261c9a79dbSrezgarshakeri ierr = VecGetArrayReadAndMemType(fine_mult, &m, &m_mem_type); 2271c9a79dbSrezgarshakeri CHKERRQ(ierr); 2281c9a79dbSrezgarshakeri CeedVectorSetArray(data[level]->x_ceed, MemTypeP2C(m_mem_type), 2291c9a79dbSrezgarshakeri CEED_USE_POINTER, (CeedScalar *)m); 230e83e87a5Sjeremylt 2311c9a79dbSrezgarshakeri CeedOperatorMultigridLevelCreate(data[level]->op_apply, data[level]->x_ceed, 2321c9a79dbSrezgarshakeri data[level-1]->elem_restr_u, basis_u, 2331c9a79dbSrezgarshakeri &op_apply, &op_prolong, &op_restrict); 234e83e87a5Sjeremylt 2351c9a79dbSrezgarshakeri // Restore PETSc vector 2361c9a79dbSrezgarshakeri CeedVectorTakeArray(data[level]->x_ceed, MemTypeP2C(m_mem_type), 2371c9a79dbSrezgarshakeri (CeedScalar **)&m); 2381c9a79dbSrezgarshakeri ierr = VecRestoreArrayReadAndMemType(fine_mult, &m); CHKERRQ(ierr); 2391c9a79dbSrezgarshakeri ierr = VecZeroEntries(fine_mult); CHKERRQ(ierr); 2401c9a79dbSrezgarshakeri // -- Save libCEED data 2411c9a79dbSrezgarshakeri data[level-1]->op_apply = op_apply; 2421c9a79dbSrezgarshakeri data[level]->op_prolong = op_prolong; 2431c9a79dbSrezgarshakeri data[level]->op_restrict = op_restrict; 2441c9a79dbSrezgarshakeri 2451c9a79dbSrezgarshakeri CeedBasisDestroy(&basis_u); 246e83e87a5Sjeremylt PetscFunctionReturn(0); 247e83e87a5Sjeremylt }; 248e83e87a5Sjeremylt 249e83e87a5Sjeremylt // ----------------------------------------------------------------------------- 250