xref: /libCEED/examples/petsc/src/libceedsetup.c (revision 9b072555b57804a6f4e0fc2b1ad83be89838f0e5)
1e83e87a5Sjeremylt #include "../include/libceedsetup.h"
2e83e87a5Sjeremylt #include "../include/petscutils.h"
3e83e87a5Sjeremylt #include <stdio.h>
4e83e87a5Sjeremylt 
5e83e87a5Sjeremylt // -----------------------------------------------------------------------------
6e83e87a5Sjeremylt // Destroy libCEED operator objects
7e83e87a5Sjeremylt // -----------------------------------------------------------------------------
8e83e87a5Sjeremylt PetscErrorCode CeedDataDestroy(CeedInt i, CeedData data) {
9e83e87a5Sjeremylt   int ierr;
10e83e87a5Sjeremylt 
11*9b072555Sjeremylt   CeedVectorDestroy(&data->q_data);
12*9b072555Sjeremylt   CeedVectorDestroy(&data->x_ceed);
13*9b072555Sjeremylt   CeedVectorDestroy(&data->y_ceed);
14*9b072555Sjeremylt   CeedBasisDestroy(&data->basis_x);
15*9b072555Sjeremylt   CeedBasisDestroy(&data->basis_u);
16*9b072555Sjeremylt   CeedElemRestrictionDestroy(&data->elem_restr_u);
17*9b072555Sjeremylt   CeedElemRestrictionDestroy(&data->elem_restr_x);
18*9b072555Sjeremylt   CeedElemRestrictionDestroy(&data->elem_restr_u_i);
19*9b072555Sjeremylt   CeedElemRestrictionDestroy(&data->elem_restr_qd_i);
20*9b072555Sjeremylt   CeedQFunctionDestroy(&data->qf_apply);
21*9b072555Sjeremylt   CeedOperatorDestroy(&data->op_apply);
22e83e87a5Sjeremylt   if (i > 0) {
23*9b072555Sjeremylt     CeedOperatorDestroy(&data->op_prolong);
24*9b072555Sjeremylt     CeedBasisDestroy(&data->basis_c_to_f);
25*9b072555Sjeremylt     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,
36*9b072555Sjeremylt                                     CeedInt topo_dim, CeedInt q_extra,
37*9b072555Sjeremylt                                     PetscInt num_comp_x, PetscInt num_comp_u,
38*9b072555Sjeremylt                                     PetscInt g_size, PetscInt xl_size,
39*9b072555Sjeremylt                                     BPData bp_data, CeedData data,
40*9b072555Sjeremylt                                     PetscBool setup_rhs, CeedVector rhs_ceed,
41e83e87a5Sjeremylt                                     CeedVector *target) {
42e83e87a5Sjeremylt   int ierr;
43*9b072555Sjeremylt   DM dm_coord;
44e83e87a5Sjeremylt   Vec coords;
45*9b072555Sjeremylt   const PetscScalar *coord_array;
46*9b072555Sjeremylt   CeedBasis basis_x, basis_u;
47*9b072555Sjeremylt   CeedElemRestriction elem_restr_x, elem_restr_u, elem_restr_u_i, elem_restr_qd_i;
48*9b072555Sjeremylt   CeedQFunction qf_setup_geo, qf_apply;
49*9b072555Sjeremylt   CeedOperator op_setup_geo, op_apply;
50*9b072555Sjeremylt   CeedVector x_coord, q_data, x_ceed, y_ceed;
51*9b072555Sjeremylt   CeedInt P, Q, num_qpts, c_start, c_end, num_elem,
52*9b072555Sjeremylt           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;
58*9b072555Sjeremylt   Q = P + q_extra;
59*9b072555Sjeremylt   CeedBasisCreateTensorH1Lagrange(ceed, topo_dim, num_comp_u, P, Q,
60*9b072555Sjeremylt                                   bp_data.q_mode,
61*9b072555Sjeremylt                                   &basis_u);
62*9b072555Sjeremylt   CeedBasisCreateTensorH1Lagrange(ceed, topo_dim, num_comp_x, 2, Q,
63*9b072555Sjeremylt                                   bp_data.q_mode,
64*9b072555Sjeremylt                                   &basis_x);
65*9b072555Sjeremylt   CeedBasisGetNumQuadraturePoints(basis_u, &num_qpts);
66e83e87a5Sjeremylt 
67e83e87a5Sjeremylt   // CEED restrictions
68*9b072555Sjeremylt   ierr = DMGetCoordinateDM(dm, &dm_coord); CHKERRQ(ierr);
69*9b072555Sjeremylt   ierr = DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL);
70e83e87a5Sjeremylt   CHKERRQ(ierr);
71*9b072555Sjeremylt   ierr = CreateRestrictionFromPlex(ceed, dm_coord, 2, topo_dim, 0, 0, 0,
72*9b072555Sjeremylt                                    &elem_restr_x);
73e83e87a5Sjeremylt   CHKERRQ(ierr);
74*9b072555Sjeremylt   ierr = CreateRestrictionFromPlex(ceed, dm, P, topo_dim, 0, 0, 0, &elem_restr_u);
75e83e87a5Sjeremylt   CHKERRQ(ierr);
76e83e87a5Sjeremylt 
77*9b072555Sjeremylt   ierr = DMPlexGetHeightStratum(dm, 0, &c_start, &c_end); CHKERRQ(ierr);
78*9b072555Sjeremylt   num_elem = c_end - c_start;
79e83e87a5Sjeremylt 
80*9b072555Sjeremylt   CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, num_comp_u,
81*9b072555Sjeremylt                                    num_comp_u*num_elem*num_qpts,
82*9b072555Sjeremylt                                    CEED_STRIDES_BACKEND, &elem_restr_u_i);
83*9b072555Sjeremylt   CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, q_data_size,
84*9b072555Sjeremylt                                    q_data_size*num_elem*num_qpts,
85*9b072555Sjeremylt                                    CEED_STRIDES_BACKEND, &elem_restr_qd_i);
86e83e87a5Sjeremylt 
87e83e87a5Sjeremylt   // Element coordinates
88e83e87a5Sjeremylt   ierr = DMGetCoordinatesLocal(dm, &coords); CHKERRQ(ierr);
89*9b072555Sjeremylt   ierr = VecGetArrayRead(coords, &coord_array); CHKERRQ(ierr);
90e83e87a5Sjeremylt 
91*9b072555Sjeremylt   CeedElemRestrictionCreateVector(elem_restr_x, &x_coord, NULL);
92*9b072555Sjeremylt   CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_COPY_VALUES,
93*9b072555Sjeremylt                      (PetscScalar *)coord_array);
94*9b072555Sjeremylt   ierr = VecRestoreArrayRead(coords, &coord_array);
95e83e87a5Sjeremylt 
96e83e87a5Sjeremylt   // Create the persistent vectors that will be needed in setup and apply
97*9b072555Sjeremylt   CeedVectorCreate(ceed, q_data_size*num_elem*num_qpts, &q_data);
98*9b072555Sjeremylt   CeedVectorCreate(ceed, xl_size, &x_ceed);
99*9b072555Sjeremylt   CeedVectorCreate(ceed, xl_size, &y_ceed);
100e83e87a5Sjeremylt 
101e83e87a5Sjeremylt   // Create the QFunction that builds the context data
102*9b072555Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_geo, bp_data.setup_geo_loc,
103*9b072555Sjeremylt                               &qf_setup_geo);
104*9b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP);
105*9b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x*topo_dim, CEED_EVAL_GRAD);
106*9b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT);
107*9b072555Sjeremylt   CeedQFunctionAddOutput(qf_setup_geo, "q_data", q_data_size, CEED_EVAL_NONE);
108e83e87a5Sjeremylt 
109e83e87a5Sjeremylt   // Create the operator that builds the quadrature data
110*9b072555Sjeremylt   CeedOperatorCreate(ceed, qf_setup_geo, NULL, NULL, &op_setup_geo);
111*9b072555Sjeremylt   CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x,
112*9b072555Sjeremylt                        CEED_VECTOR_ACTIVE);
113*9b072555Sjeremylt   CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x,
114*9b072555Sjeremylt                        CEED_VECTOR_ACTIVE);
115*9b072555Sjeremylt   CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x,
116e83e87a5Sjeremylt                        CEED_VECTOR_NONE);
117*9b072555Sjeremylt   CeedOperatorSetField(op_setup_geo, "q_data", elem_restr_qd_i,
118e83e87a5Sjeremylt                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
119e83e87a5Sjeremylt 
120*9b072555Sjeremylt   // Setup q_data
121*9b072555Sjeremylt   CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE);
122e83e87a5Sjeremylt 
123e83e87a5Sjeremylt   // Set up PDE operator
124*9b072555Sjeremylt   CeedInt in_scale = bp_data.in_mode == CEED_EVAL_GRAD ? topo_dim : 1;
125*9b072555Sjeremylt   CeedInt out_scale = bp_data.out_mode == CEED_EVAL_GRAD ? topo_dim : 1;
126*9b072555Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bp_data.apply, bp_data.apply_loc,
127*9b072555Sjeremylt                               &qf_apply);
128*9b072555Sjeremylt   CeedQFunctionAddInput(qf_apply, "u", num_comp_u*in_scale, bp_data.in_mode);
129*9b072555Sjeremylt   CeedQFunctionAddInput(qf_apply, "q_data", q_data_size, CEED_EVAL_NONE);
130*9b072555Sjeremylt   CeedQFunctionAddOutput(qf_apply, "v", num_comp_u*out_scale, bp_data.out_mode);
131e83e87a5Sjeremylt 
132e83e87a5Sjeremylt   // Create the mass or diff operator
133*9b072555Sjeremylt   CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply);
134*9b072555Sjeremylt   CeedOperatorSetField(op_apply, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
135*9b072555Sjeremylt   CeedOperatorSetField(op_apply, "q_data", elem_restr_qd_i, CEED_BASIS_COLLOCATED,
136*9b072555Sjeremylt                        q_data);
137*9b072555Sjeremylt   CeedOperatorSetField(op_apply, "v", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
138e83e87a5Sjeremylt 
139e83e87a5Sjeremylt   // Set up RHS if needed
140e83e87a5Sjeremylt   if (setup_rhs) {
141*9b072555Sjeremylt     CeedQFunction qf_setup_rhs;
142*9b072555Sjeremylt     CeedOperator op_setup_rhs;
143*9b072555Sjeremylt     CeedVectorCreate(ceed, num_elem*num_qpts*num_comp_u, target);
144e83e87a5Sjeremylt 
145e83e87a5Sjeremylt     // Create the q-function that sets up the RHS and true solution
146*9b072555Sjeremylt     CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_rhs, bp_data.setup_rhs_loc,
147*9b072555Sjeremylt                                 &qf_setup_rhs);
148*9b072555Sjeremylt     CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP);
149*9b072555Sjeremylt     CeedQFunctionAddInput(qf_setup_rhs, "q_data", q_data_size, CEED_EVAL_NONE);
150*9b072555Sjeremylt     CeedQFunctionAddOutput(qf_setup_rhs, "true_soln", num_comp_u, CEED_EVAL_NONE);
151*9b072555Sjeremylt     CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP);
152e83e87a5Sjeremylt 
153e83e87a5Sjeremylt     // Create the operator that builds the RHS and true solution
154*9b072555Sjeremylt     CeedOperatorCreate(ceed, qf_setup_rhs, NULL, NULL, &op_setup_rhs);
155*9b072555Sjeremylt     CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x,
156*9b072555Sjeremylt                          CEED_VECTOR_ACTIVE);
157*9b072555Sjeremylt     CeedOperatorSetField(op_setup_rhs, "q_data", elem_restr_qd_i,
158*9b072555Sjeremylt                          CEED_BASIS_COLLOCATED,
159*9b072555Sjeremylt                          q_data);
160*9b072555Sjeremylt     CeedOperatorSetField(op_setup_rhs, "true_soln", elem_restr_u_i,
161e83e87a5Sjeremylt                          CEED_BASIS_COLLOCATED, *target);
162*9b072555Sjeremylt     CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u,
163e83e87a5Sjeremylt                          CEED_VECTOR_ACTIVE);
164e83e87a5Sjeremylt 
165e83e87a5Sjeremylt     // Set up the libCEED context
166*9b072555Sjeremylt     CeedQFunctionContext ctx_rhs_setup;
167*9b072555Sjeremylt     CeedQFunctionContextCreate(ceed, &ctx_rhs_setup);
168*9b072555Sjeremylt     CeedScalar rhs_setup_data[2] = {R, l};
169*9b072555Sjeremylt     CeedQFunctionContextSetData(ctx_rhs_setup, CEED_MEM_HOST, CEED_COPY_VALUES,
170*9b072555Sjeremylt                                 sizeof rhs_setup_data, &rhs_setup_data);
171*9b072555Sjeremylt     CeedQFunctionSetContext(qf_setup_rhs, ctx_rhs_setup);
172*9b072555Sjeremylt     CeedQFunctionContextDestroy(&ctx_rhs_setup);
173e83e87a5Sjeremylt 
174e83e87a5Sjeremylt     // Setup RHS and target
175*9b072555Sjeremylt     CeedOperatorApply(op_setup_rhs, x_coord, rhs_ceed, CEED_REQUEST_IMMEDIATE);
176e83e87a5Sjeremylt 
177e83e87a5Sjeremylt     // Cleanup
178*9b072555Sjeremylt     CeedQFunctionDestroy(&qf_setup_rhs);
179*9b072555Sjeremylt     CeedOperatorDestroy(&op_setup_rhs);
180e83e87a5Sjeremylt   }
181e83e87a5Sjeremylt 
182e83e87a5Sjeremylt   // Cleanup
183*9b072555Sjeremylt   CeedQFunctionDestroy(&qf_setup_geo);
184*9b072555Sjeremylt   CeedOperatorDestroy(&op_setup_geo);
185*9b072555Sjeremylt   CeedVectorDestroy(&x_coord);
186e83e87a5Sjeremylt 
187e83e87a5Sjeremylt   // Save libCEED data required for level
188*9b072555Sjeremylt   data->basis_x = basis_x; data->basis_u = basis_u;
189*9b072555Sjeremylt   data->elem_restr_x = elem_restr_x;
190*9b072555Sjeremylt   data->elem_restr_u = elem_restr_u;
191*9b072555Sjeremylt   data->elem_restr_u_i = elem_restr_u_i;
192*9b072555Sjeremylt   data->elem_restr_qd_i = elem_restr_qd_i;
193*9b072555Sjeremylt   data->qf_apply = qf_apply;
194*9b072555Sjeremylt   data->op_apply = op_apply;
195*9b072555Sjeremylt   data->q_data = q_data;
196*9b072555Sjeremylt   data->x_ceed = x_ceed;
197*9b072555Sjeremylt   data->y_ceed = y_ceed;
198e83e87a5Sjeremylt 
199e83e87a5Sjeremylt   PetscFunctionReturn(0);
200e83e87a5Sjeremylt };
201e83e87a5Sjeremylt 
202e83e87a5Sjeremylt // -----------------------------------------------------------------------------
203e83e87a5Sjeremylt // Setup libCEED level transfer operator objects
204e83e87a5Sjeremylt // -----------------------------------------------------------------------------
205*9b072555Sjeremylt PetscErrorCode CeedLevelTransferSetup(Ceed ceed, CeedInt num_levels,
206*9b072555Sjeremylt                                       CeedInt num_comp_u, CeedData *data,
207*9b072555Sjeremylt                                       CeedInt *level_degrees,
208*9b072555Sjeremylt                                       CeedQFunction qf_restrict, CeedQFunction qf_prolong) {
209*9b072555Sjeremylt   // Return early if num_levels=1
210*9b072555Sjeremylt   if (num_levels == 1)
211e83e87a5Sjeremylt     PetscFunctionReturn(0);
212e83e87a5Sjeremylt 
213e83e87a5Sjeremylt   // Set up each level
214*9b072555Sjeremylt   for (CeedInt i=1; i<num_levels; i++) {
215e83e87a5Sjeremylt     // P coarse and P fine
216*9b072555Sjeremylt     CeedInt Pc = level_degrees[i-1] + 1;
217*9b072555Sjeremylt     CeedInt Pf = level_degrees[i] + 1;
218e83e87a5Sjeremylt 
219e83e87a5Sjeremylt     // Restriction - Fine to corse
220*9b072555Sjeremylt     CeedBasis basis_c_to_f;
221*9b072555Sjeremylt     CeedOperator op_restrict;
222e83e87a5Sjeremylt 
223e83e87a5Sjeremylt     // Basis
224*9b072555Sjeremylt     CeedBasisCreateTensorH1Lagrange(ceed, 3, num_comp_u, Pc, Pf,
225*9b072555Sjeremylt                                     CEED_GAUSS_LOBATTO, &basis_c_to_f);
226e83e87a5Sjeremylt 
227e83e87a5Sjeremylt     // Create the restriction operator
228*9b072555Sjeremylt     CeedOperatorCreate(ceed, qf_restrict, CEED_QFUNCTION_NONE,
229*9b072555Sjeremylt                        CEED_QFUNCTION_NONE, &op_restrict);
230*9b072555Sjeremylt     CeedOperatorSetField(op_restrict, "input", data[i]->elem_restr_u,
231e83e87a5Sjeremylt                          CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
232*9b072555Sjeremylt     CeedOperatorSetField(op_restrict, "output", data[i-1]->elem_restr_u,
233*9b072555Sjeremylt                          basis_c_to_f, CEED_VECTOR_ACTIVE);
234e83e87a5Sjeremylt 
235e83e87a5Sjeremylt     // Save libCEED data required for level
236*9b072555Sjeremylt     data[i]->basis_c_to_f = basis_c_to_f;
237*9b072555Sjeremylt     data[i]->op_restrict = op_restrict;
238e83e87a5Sjeremylt 
239e83e87a5Sjeremylt     // Interpolation - Corse to fine
240*9b072555Sjeremylt     CeedOperator op_prolong;
241e83e87a5Sjeremylt 
242e83e87a5Sjeremylt     // Create the prolongation operator
243*9b072555Sjeremylt     CeedOperatorCreate(ceed, qf_prolong, CEED_QFUNCTION_NONE,
244*9b072555Sjeremylt                        CEED_QFUNCTION_NONE, &op_prolong);
245*9b072555Sjeremylt     CeedOperatorSetField(op_prolong, "input", data[i-1]->elem_restr_u,
246*9b072555Sjeremylt                          basis_c_to_f, CEED_VECTOR_ACTIVE);
247*9b072555Sjeremylt     CeedOperatorSetField(op_prolong, "output", data[i]->elem_restr_u,
248e83e87a5Sjeremylt                          CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
249e83e87a5Sjeremylt 
250e83e87a5Sjeremylt     // Save libCEED data required for level
251*9b072555Sjeremylt     data[i]->op_prolong = op_prolong;
252e83e87a5Sjeremylt   }
253e83e87a5Sjeremylt 
254e83e87a5Sjeremylt   PetscFunctionReturn(0);
255e83e87a5Sjeremylt };
256e83e87a5Sjeremylt 
257e83e87a5Sjeremylt // -----------------------------------------------------------------------------
258