xref: /libCEED/examples/petsc/src/libceedsetup.c (revision a61c78d6a6d5ea69db49949746e6dc59b544c365)
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 
119b072555Sjeremylt   CeedVectorDestroy(&data->q_data);
129b072555Sjeremylt   CeedVectorDestroy(&data->x_ceed);
139b072555Sjeremylt   CeedVectorDestroy(&data->y_ceed);
149b072555Sjeremylt   CeedBasisDestroy(&data->basis_x);
159b072555Sjeremylt   CeedBasisDestroy(&data->basis_u);
169b072555Sjeremylt   CeedElemRestrictionDestroy(&data->elem_restr_u);
179b072555Sjeremylt   CeedElemRestrictionDestroy(&data->elem_restr_x);
189b072555Sjeremylt   CeedElemRestrictionDestroy(&data->elem_restr_u_i);
199b072555Sjeremylt   CeedElemRestrictionDestroy(&data->elem_restr_qd_i);
209b072555Sjeremylt   CeedQFunctionDestroy(&data->qf_apply);
219b072555Sjeremylt   CeedOperatorDestroy(&data->op_apply);
22e83e87a5Sjeremylt   if (i > 0) {
239b072555Sjeremylt     CeedOperatorDestroy(&data->op_prolong);
249b072555Sjeremylt     CeedBasisDestroy(&data->basis_c_to_f);
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;
519b072555Sjeremylt   CeedInt P, Q, 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 
56e83e87a5Sjeremylt   // CEED bases
57e83e87a5Sjeremylt   P = degree + 1;
589b072555Sjeremylt   Q = P + q_extra;
599b072555Sjeremylt   CeedBasisCreateTensorH1Lagrange(ceed, topo_dim, num_comp_u, P, Q,
609b072555Sjeremylt                                   bp_data.q_mode,
619b072555Sjeremylt                                   &basis_u);
629b072555Sjeremylt   CeedBasisCreateTensorH1Lagrange(ceed, topo_dim, num_comp_x, 2, Q,
639b072555Sjeremylt                                   bp_data.q_mode,
649b072555Sjeremylt                                   &basis_x);
659b072555Sjeremylt   CeedBasisGetNumQuadraturePoints(basis_u, &num_qpts);
66e83e87a5Sjeremylt 
67e83e87a5Sjeremylt   // CEED restrictions
687ed3e4cdSJeremy L Thompson   ierr = DMSetCoordinateDim(dm, topo_dim); CHKERRQ(ierr);
699b072555Sjeremylt   ierr = DMGetCoordinateDM(dm, &dm_coord); CHKERRQ(ierr);
709b072555Sjeremylt   ierr = DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL);
71e83e87a5Sjeremylt   CHKERRQ(ierr);
727ed3e4cdSJeremy L Thompson   ierr = CreateRestrictionFromPlex(ceed, dm_coord, 0, 0, 0, &elem_restr_x);
73e83e87a5Sjeremylt   CHKERRQ(ierr);
747ed3e4cdSJeremy L Thompson   ierr = CreateRestrictionFromPlex(ceed, dm, 0, 0, 0, &elem_restr_u);
75e83e87a5Sjeremylt   CHKERRQ(ierr);
76e83e87a5Sjeremylt 
779b072555Sjeremylt   ierr = DMPlexGetHeightStratum(dm, 0, &c_start, &c_end); CHKERRQ(ierr);
789b072555Sjeremylt   num_elem = c_end - c_start;
79e83e87a5Sjeremylt 
809b072555Sjeremylt   CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, num_comp_u,
819b072555Sjeremylt                                    num_comp_u*num_elem*num_qpts,
829b072555Sjeremylt                                    CEED_STRIDES_BACKEND, &elem_restr_u_i);
839b072555Sjeremylt   CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, q_data_size,
849b072555Sjeremylt                                    q_data_size*num_elem*num_qpts,
859b072555Sjeremylt                                    CEED_STRIDES_BACKEND, &elem_restr_qd_i);
86e83e87a5Sjeremylt 
87e83e87a5Sjeremylt   // Element coordinates
88e83e87a5Sjeremylt   ierr = DMGetCoordinatesLocal(dm, &coords); CHKERRQ(ierr);
899b072555Sjeremylt   ierr = VecGetArrayRead(coords, &coord_array); CHKERRQ(ierr);
90e83e87a5Sjeremylt 
919b072555Sjeremylt   CeedElemRestrictionCreateVector(elem_restr_x, &x_coord, NULL);
929b072555Sjeremylt   CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_COPY_VALUES,
939b072555Sjeremylt                      (PetscScalar *)coord_array);
949b072555Sjeremylt   ierr = VecRestoreArrayRead(coords, &coord_array);
95e83e87a5Sjeremylt 
96e83e87a5Sjeremylt   // Create the persistent vectors that will be needed in setup and apply
979b072555Sjeremylt   CeedVectorCreate(ceed, q_data_size*num_elem*num_qpts, &q_data);
989b072555Sjeremylt   CeedVectorCreate(ceed, xl_size, &x_ceed);
999b072555Sjeremylt   CeedVectorCreate(ceed, xl_size, &y_ceed);
100e83e87a5Sjeremylt 
101e83e87a5Sjeremylt   // Create the QFunction that builds the context data
1029b072555Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_geo, bp_data.setup_geo_loc,
1039b072555Sjeremylt                               &qf_setup_geo);
1049b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP);
1059b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x*topo_dim, CEED_EVAL_GRAD);
1069b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT);
107*a61c78d6SJeremy L Thompson   CeedQFunctionAddOutput(qf_setup_geo, "qdata", q_data_size, CEED_EVAL_NONE);
108e83e87a5Sjeremylt 
109e83e87a5Sjeremylt   // Create the operator that builds the quadrature data
1109b072555Sjeremylt   CeedOperatorCreate(ceed, qf_setup_geo, NULL, NULL, &op_setup_geo);
1119b072555Sjeremylt   CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x,
1129b072555Sjeremylt                        CEED_VECTOR_ACTIVE);
1139b072555Sjeremylt   CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x,
1149b072555Sjeremylt                        CEED_VECTOR_ACTIVE);
1159b072555Sjeremylt   CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x,
116e83e87a5Sjeremylt                        CEED_VECTOR_NONE);
117*a61c78d6SJeremy L Thompson   CeedOperatorSetField(op_setup_geo, "qdata", elem_restr_qd_i,
118e83e87a5Sjeremylt                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
119e83e87a5Sjeremylt 
1209b072555Sjeremylt   // Setup q_data
1219b072555Sjeremylt   CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE);
122e83e87a5Sjeremylt 
123e83e87a5Sjeremylt   // Set up PDE operator
1249b072555Sjeremylt   CeedInt in_scale = bp_data.in_mode == CEED_EVAL_GRAD ? topo_dim : 1;
1259b072555Sjeremylt   CeedInt out_scale = bp_data.out_mode == CEED_EVAL_GRAD ? topo_dim : 1;
1269b072555Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bp_data.apply, bp_data.apply_loc,
1279b072555Sjeremylt                               &qf_apply);
1289b072555Sjeremylt   CeedQFunctionAddInput(qf_apply, "u", num_comp_u*in_scale, bp_data.in_mode);
129*a61c78d6SJeremy L Thompson   CeedQFunctionAddInput(qf_apply, "qdata", q_data_size, CEED_EVAL_NONE);
1309b072555Sjeremylt   CeedQFunctionAddOutput(qf_apply, "v", num_comp_u*out_scale, bp_data.out_mode);
131e83e87a5Sjeremylt 
132e83e87a5Sjeremylt   // Create the mass or diff operator
1339b072555Sjeremylt   CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply);
1349b072555Sjeremylt   CeedOperatorSetField(op_apply, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
135*a61c78d6SJeremy L Thompson   CeedOperatorSetField(op_apply, "qdata", elem_restr_qd_i, CEED_BASIS_COLLOCATED,
1369b072555Sjeremylt                        q_data);
1379b072555Sjeremylt   CeedOperatorSetField(op_apply, "v", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
138e83e87a5Sjeremylt 
139e83e87a5Sjeremylt   // Set up RHS if needed
140e83e87a5Sjeremylt   if (setup_rhs) {
1419b072555Sjeremylt     CeedQFunction qf_setup_rhs;
1429b072555Sjeremylt     CeedOperator op_setup_rhs;
1439b072555Sjeremylt     CeedVectorCreate(ceed, num_elem*num_qpts*num_comp_u, target);
144e83e87a5Sjeremylt 
145e83e87a5Sjeremylt     // Create the q-function that sets up the RHS and true solution
1469b072555Sjeremylt     CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_rhs, bp_data.setup_rhs_loc,
1479b072555Sjeremylt                                 &qf_setup_rhs);
1489b072555Sjeremylt     CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP);
149*a61c78d6SJeremy L Thompson     CeedQFunctionAddInput(qf_setup_rhs, "qdata", q_data_size, CEED_EVAL_NONE);
150*a61c78d6SJeremy L Thompson     CeedQFunctionAddOutput(qf_setup_rhs, "true solution", num_comp_u,
151*a61c78d6SJeremy L Thompson                            CEED_EVAL_NONE);
1529b072555Sjeremylt     CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP);
153e83e87a5Sjeremylt 
154e83e87a5Sjeremylt     // Create the operator that builds the RHS and true solution
1559b072555Sjeremylt     CeedOperatorCreate(ceed, qf_setup_rhs, NULL, NULL, &op_setup_rhs);
1569b072555Sjeremylt     CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x,
1579b072555Sjeremylt                          CEED_VECTOR_ACTIVE);
158*a61c78d6SJeremy L Thompson     CeedOperatorSetField(op_setup_rhs, "qdata", elem_restr_qd_i,
159*a61c78d6SJeremy L Thompson                          CEED_BASIS_COLLOCATED, q_data);
160*a61c78d6SJeremy L Thompson     CeedOperatorSetField(op_setup_rhs, "true solution", elem_restr_u_i,
161e83e87a5Sjeremylt                          CEED_BASIS_COLLOCATED, *target);
1629b072555Sjeremylt     CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u,
163e83e87a5Sjeremylt                          CEED_VECTOR_ACTIVE);
164e83e87a5Sjeremylt 
165e83e87a5Sjeremylt     // Set up the libCEED context
1669b072555Sjeremylt     CeedQFunctionContext ctx_rhs_setup;
1679b072555Sjeremylt     CeedQFunctionContextCreate(ceed, &ctx_rhs_setup);
1689b072555Sjeremylt     CeedScalar rhs_setup_data[2] = {R, l};
1699b072555Sjeremylt     CeedQFunctionContextSetData(ctx_rhs_setup, CEED_MEM_HOST, CEED_COPY_VALUES,
1709b072555Sjeremylt                                 sizeof rhs_setup_data, &rhs_setup_data);
1719b072555Sjeremylt     CeedQFunctionSetContext(qf_setup_rhs, ctx_rhs_setup);
1729b072555Sjeremylt     CeedQFunctionContextDestroy(&ctx_rhs_setup);
173e83e87a5Sjeremylt 
174e83e87a5Sjeremylt     // Setup RHS and target
1759b072555Sjeremylt     CeedOperatorApply(op_setup_rhs, x_coord, rhs_ceed, CEED_REQUEST_IMMEDIATE);
176e83e87a5Sjeremylt 
177e83e87a5Sjeremylt     // Cleanup
1789b072555Sjeremylt     CeedQFunctionDestroy(&qf_setup_rhs);
1799b072555Sjeremylt     CeedOperatorDestroy(&op_setup_rhs);
180e83e87a5Sjeremylt   }
181e83e87a5Sjeremylt 
182e83e87a5Sjeremylt   // Cleanup
1839b072555Sjeremylt   CeedQFunctionDestroy(&qf_setup_geo);
1849b072555Sjeremylt   CeedOperatorDestroy(&op_setup_geo);
1859b072555Sjeremylt   CeedVectorDestroy(&x_coord);
186e83e87a5Sjeremylt 
187e83e87a5Sjeremylt   // Save libCEED data required for level
1889b072555Sjeremylt   data->basis_x = basis_x; data->basis_u = basis_u;
1899b072555Sjeremylt   data->elem_restr_x = elem_restr_x;
1909b072555Sjeremylt   data->elem_restr_u = elem_restr_u;
1919b072555Sjeremylt   data->elem_restr_u_i = elem_restr_u_i;
1929b072555Sjeremylt   data->elem_restr_qd_i = elem_restr_qd_i;
1939b072555Sjeremylt   data->qf_apply = qf_apply;
1949b072555Sjeremylt   data->op_apply = op_apply;
1959b072555Sjeremylt   data->q_data = q_data;
1969b072555Sjeremylt   data->x_ceed = x_ceed;
1979b072555Sjeremylt   data->y_ceed = y_ceed;
198e83e87a5Sjeremylt 
199e83e87a5Sjeremylt   PetscFunctionReturn(0);
200e83e87a5Sjeremylt };
201e83e87a5Sjeremylt 
202e83e87a5Sjeremylt // -----------------------------------------------------------------------------
203e83e87a5Sjeremylt // Setup libCEED level transfer operator objects
204e83e87a5Sjeremylt // -----------------------------------------------------------------------------
2059b072555Sjeremylt PetscErrorCode CeedLevelTransferSetup(Ceed ceed, CeedInt num_levels,
2069b072555Sjeremylt                                       CeedInt num_comp_u, CeedData *data,
2079b072555Sjeremylt                                       CeedInt *level_degrees,
2089b072555Sjeremylt                                       CeedQFunction qf_restrict, CeedQFunction qf_prolong) {
2099b072555Sjeremylt   // Return early if num_levels=1
2109b072555Sjeremylt   if (num_levels == 1)
211e83e87a5Sjeremylt     PetscFunctionReturn(0);
212e83e87a5Sjeremylt 
213e83e87a5Sjeremylt   // Set up each level
2149b072555Sjeremylt   for (CeedInt i=1; i<num_levels; i++) {
215e83e87a5Sjeremylt     // P coarse and P fine
2169b072555Sjeremylt     CeedInt Pc = level_degrees[i-1] + 1;
2179b072555Sjeremylt     CeedInt Pf = level_degrees[i] + 1;
218e83e87a5Sjeremylt 
219e83e87a5Sjeremylt     // Restriction - Fine to corse
2209b072555Sjeremylt     CeedBasis basis_c_to_f;
2219b072555Sjeremylt     CeedOperator op_restrict;
222e83e87a5Sjeremylt 
223e83e87a5Sjeremylt     // Basis
2249b072555Sjeremylt     CeedBasisCreateTensorH1Lagrange(ceed, 3, num_comp_u, Pc, Pf,
2259b072555Sjeremylt                                     CEED_GAUSS_LOBATTO, &basis_c_to_f);
226e83e87a5Sjeremylt 
227e83e87a5Sjeremylt     // Create the restriction operator
2289b072555Sjeremylt     CeedOperatorCreate(ceed, qf_restrict, CEED_QFUNCTION_NONE,
2299b072555Sjeremylt                        CEED_QFUNCTION_NONE, &op_restrict);
2309b072555Sjeremylt     CeedOperatorSetField(op_restrict, "input", data[i]->elem_restr_u,
231e83e87a5Sjeremylt                          CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
2329b072555Sjeremylt     CeedOperatorSetField(op_restrict, "output", data[i-1]->elem_restr_u,
2339b072555Sjeremylt                          basis_c_to_f, CEED_VECTOR_ACTIVE);
234e83e87a5Sjeremylt 
235e83e87a5Sjeremylt     // Save libCEED data required for level
2369b072555Sjeremylt     data[i]->basis_c_to_f = basis_c_to_f;
2379b072555Sjeremylt     data[i]->op_restrict = op_restrict;
238e83e87a5Sjeremylt 
239e83e87a5Sjeremylt     // Interpolation - Corse to fine
2409b072555Sjeremylt     CeedOperator op_prolong;
241e83e87a5Sjeremylt 
242e83e87a5Sjeremylt     // Create the prolongation operator
2439b072555Sjeremylt     CeedOperatorCreate(ceed, qf_prolong, CEED_QFUNCTION_NONE,
2449b072555Sjeremylt                        CEED_QFUNCTION_NONE, &op_prolong);
2459b072555Sjeremylt     CeedOperatorSetField(op_prolong, "input", data[i-1]->elem_restr_u,
2469b072555Sjeremylt                          basis_c_to_f, CEED_VECTOR_ACTIVE);
2479b072555Sjeremylt     CeedOperatorSetField(op_prolong, "output", data[i]->elem_restr_u,
248e83e87a5Sjeremylt                          CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
249e83e87a5Sjeremylt 
250e83e87a5Sjeremylt     // Save libCEED data required for level
2519b072555Sjeremylt     data[i]->op_prolong = op_prolong;
252e83e87a5Sjeremylt   }
253e83e87a5Sjeremylt 
254e83e87a5Sjeremylt   PetscFunctionReturn(0);
255e83e87a5Sjeremylt };
256e83e87a5Sjeremylt 
257e83e87a5Sjeremylt // -----------------------------------------------------------------------------
258