xref: /libCEED/examples/petsc/src/libceedsetup.c (revision 509d4af65d23546c690c9766d8b29e47dc3b3afb)
1 // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 #include "../include/libceedsetup.h"
9 
10 #include <stdio.h>
11 
12 #include "../include/libceedsetup.h"
13 #include "../include/petscutils.h"
14 
15 // -----------------------------------------------------------------------------
16 // Destroy libCEED operator objects
17 // -----------------------------------------------------------------------------
18 PetscErrorCode CeedDataDestroy(CeedInt i, CeedData data) {
19   PetscFunctionBeginUser;
20   CeedVectorDestroy(&data->q_data);
21   CeedVectorDestroy(&data->x_ceed);
22   CeedVectorDestroy(&data->y_ceed);
23   CeedBasisDestroy(&data->basis_x);
24   CeedBasisDestroy(&data->basis_u);
25   CeedElemRestrictionDestroy(&data->elem_restr_u);
26   CeedElemRestrictionDestroy(&data->elem_restr_x);
27   CeedElemRestrictionDestroy(&data->elem_restr_u_i);
28   CeedElemRestrictionDestroy(&data->elem_restr_qd_i);
29   CeedQFunctionDestroy(&data->qf_apply);
30   CeedOperatorDestroy(&data->op_apply);
31   if (i > 0) {
32     CeedOperatorDestroy(&data->op_prolong);
33     CeedOperatorDestroy(&data->op_restrict);
34   }
35   PetscCall(PetscFree(data));
36   PetscFunctionReturn(PETSC_SUCCESS);
37 };
38 
39 // -----------------------------------------------------------------------------
40 // Set up libCEED for a given degree
41 // -----------------------------------------------------------------------------
42 PetscErrorCode SetupLibceedByDegree(DM dm, Ceed ceed, CeedInt degree, CeedInt topo_dim, CeedInt q_extra, PetscInt num_comp_x, PetscInt num_comp_u,
43                                     PetscInt g_size, PetscInt xl_size, BPData bp_data, CeedData data, PetscBool setup_rhs, CeedVector rhs_ceed,
44                                     CeedVector *target) {
45   DM                  dm_coord;
46   Vec                 coords;
47   const PetscScalar  *coord_array;
48   CeedBasis           basis_x, basis_u;
49   CeedElemRestriction elem_restr_x, elem_restr_u, elem_restr_u_i, elem_restr_qd_i;
50   CeedQFunction       qf_setup_geo, qf_apply;
51   CeedOperator        op_setup_geo, op_apply;
52   CeedVector          x_coord, q_data, x_ceed, y_ceed;
53   PetscInt            c_start, c_end, num_elem;
54   CeedInt             num_qpts, q_data_size = bp_data.q_data_size;
55   CeedScalar          R = 1;                         // radius of the sphere
56   CeedScalar          l = 1.0 / PetscSqrtReal(3.0);  // half edge of the inscribed cube
57 
58   PetscFunctionBeginUser;
59   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
60 
61   // CEED bases
62   PetscCall(CreateBasisFromPlex(ceed, dm_coord, 0, 0, 0, 0, bp_data, &basis_x));
63   PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, bp_data, &basis_u));
64 
65   // CEED restrictions
66   PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, 0, 0, 0, &elem_restr_x));
67   PetscCall(CreateRestrictionFromPlex(ceed, dm, 0, 0, 0, &elem_restr_u));
68 
69   PetscCall(DMPlexGetHeightStratum(dm, 0, &c_start, &c_end));
70   num_elem = c_end - c_start;
71   CeedBasisGetNumQuadraturePoints(basis_u, &num_qpts);
72 
73   CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, num_comp_u, num_comp_u * num_elem * num_qpts, CEED_STRIDES_BACKEND, &elem_restr_u_i);
74   CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, q_data_size, q_data_size * num_elem * num_qpts, CEED_STRIDES_BACKEND, &elem_restr_qd_i);
75 
76   // Element coordinates
77   PetscCall(DMGetCoordinatesLocal(dm, &coords));
78   PetscCall(VecGetArrayRead(coords, &coord_array));
79 
80   CeedElemRestrictionCreateVector(elem_restr_x, &x_coord, NULL);
81   CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, (PetscScalar *)coord_array);
82   PetscCall(VecRestoreArrayRead(coords, &coord_array));
83 
84   // Create the persistent vectors that will be needed in setup and apply
85   CeedVectorCreate(ceed, q_data_size * num_elem * num_qpts, &q_data);
86   CeedVectorCreate(ceed, xl_size, &x_ceed);
87   CeedVectorCreate(ceed, xl_size, &y_ceed);
88 
89   // Create the QFunction that builds the context data
90   CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_geo, bp_data.setup_geo_loc, &qf_setup_geo);
91   CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP);
92   CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x * topo_dim, CEED_EVAL_GRAD);
93   CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT);
94   CeedQFunctionAddOutput(qf_setup_geo, "qdata", q_data_size, CEED_EVAL_NONE);
95 
96   // Create the operator that builds the quadrature data
97   CeedOperatorCreate(ceed, qf_setup_geo, NULL, NULL, &op_setup_geo);
98   CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
99   CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
100   CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE);
101   CeedOperatorSetField(op_setup_geo, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
102 
103   // Setup q_data
104   CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE);
105 
106   // Set up PDE operator
107   CeedInt in_scale  = bp_data.in_mode == CEED_EVAL_GRAD ? topo_dim : 1;
108   CeedInt out_scale = bp_data.out_mode == CEED_EVAL_GRAD ? topo_dim : 1;
109   CeedQFunctionCreateInterior(ceed, 1, bp_data.apply, bp_data.apply_loc, &qf_apply);
110   CeedQFunctionAddInput(qf_apply, "u", num_comp_u * in_scale, bp_data.in_mode);
111   CeedQFunctionAddInput(qf_apply, "qdata", q_data_size, CEED_EVAL_NONE);
112   CeedQFunctionAddOutput(qf_apply, "v", num_comp_u * out_scale, bp_data.out_mode);
113 
114   // Create the mass or diff operator
115   CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply);
116   CeedOperatorSetField(op_apply, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
117   CeedOperatorSetField(op_apply, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
118   CeedOperatorSetField(op_apply, "v", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
119 
120   // Set up RHS if needed
121   if (setup_rhs) {
122     CeedQFunction qf_setup_rhs;
123     CeedOperator  op_setup_rhs;
124     CeedVectorCreate(ceed, num_elem * num_qpts * num_comp_u, target);
125     // Create the q-function that sets up the RHS and true solution
126     CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_rhs, bp_data.setup_rhs_loc, &qf_setup_rhs);
127     CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP);
128     CeedQFunctionAddInput(qf_setup_rhs, "qdata", q_data_size, CEED_EVAL_NONE);
129     CeedQFunctionAddOutput(qf_setup_rhs, "true solution", num_comp_u, CEED_EVAL_NONE);
130     CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP);
131 
132     // Create the operator that builds the RHS and true solution
133     CeedOperatorCreate(ceed, qf_setup_rhs, NULL, NULL, &op_setup_rhs);
134     CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
135     CeedOperatorSetField(op_setup_rhs, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
136     CeedOperatorSetField(op_setup_rhs, "true solution", elem_restr_u_i, CEED_BASIS_NONE, *target);
137     CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
138 
139     // Set up the libCEED context
140     CeedQFunctionContext ctx_rhs_setup;
141     CeedQFunctionContextCreate(ceed, &ctx_rhs_setup);
142     CeedScalar rhs_setup_data[2] = {R, l};
143     CeedQFunctionContextSetData(ctx_rhs_setup, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof rhs_setup_data, &rhs_setup_data);
144     CeedQFunctionSetContext(qf_setup_rhs, ctx_rhs_setup);
145     CeedQFunctionContextDestroy(&ctx_rhs_setup);
146 
147     // Setup RHS and target
148     CeedOperatorApply(op_setup_rhs, x_coord, rhs_ceed, CEED_REQUEST_IMMEDIATE);
149 
150     // Cleanup
151     CeedQFunctionDestroy(&qf_setup_rhs);
152     CeedOperatorDestroy(&op_setup_rhs);
153   }
154 
155   // Cleanup
156   CeedQFunctionDestroy(&qf_setup_geo);
157   CeedOperatorDestroy(&op_setup_geo);
158   CeedVectorDestroy(&x_coord);
159 
160   // Save libCEED data required for level
161   data->basis_x         = basis_x;
162   data->basis_u         = basis_u;
163   data->elem_restr_x    = elem_restr_x;
164   data->elem_restr_u    = elem_restr_u;
165   data->elem_restr_u_i  = elem_restr_u_i;
166   data->elem_restr_qd_i = elem_restr_qd_i;
167   data->qf_apply        = qf_apply;
168   data->op_apply        = op_apply;
169   data->q_data          = q_data;
170   data->x_ceed          = x_ceed;
171   data->y_ceed          = y_ceed;
172   data->q_data_size     = q_data_size;
173   PetscFunctionReturn(PETSC_SUCCESS);
174 };
175 
176 // -----------------------------------------------------------------------------
177 // Setup libCEED level transfer operator objects
178 // -----------------------------------------------------------------------------
179 PetscErrorCode CeedLevelTransferSetup(DM dm, Ceed ceed, CeedInt level, CeedInt num_comp_u, CeedData *data, BPData bp_data, Vec fine_mult) {
180   PetscFunctionBeginUser;
181   // Restriction - Fine to corse
182   CeedOperator op_restrict;
183   // Interpolation - Corse to fine
184   CeedOperator op_prolong;
185   // Coarse grid operator
186   CeedOperator op_apply;
187   // Basis
188   CeedBasis basis_u;
189   PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, bp_data, &basis_u));
190 
191   // ---------------------------------------------------------------------------
192   // Coarse Grid, Prolongation, and Restriction Operators
193   // ---------------------------------------------------------------------------
194   // Create the Operators that compute the prolongation and
195   //   restriction between the p-multigrid levels and the coarse grid eval.
196   // ---------------------------------------------------------------------------
197   // Place in libCEED array
198   PetscMemType m_mem_type;
199   PetscCall(VecReadP2C(fine_mult, &m_mem_type, data[level]->x_ceed));
200 
201   CeedOperatorMultigridLevelCreate(data[level]->op_apply, data[level]->x_ceed, data[level - 1]->elem_restr_u, basis_u, &op_apply, &op_prolong,
202                                    &op_restrict);
203 
204   // Restore PETSc vector
205   PetscCall(VecReadC2P(data[level]->x_ceed, m_mem_type, fine_mult));
206   PetscCall(VecZeroEntries(fine_mult));
207   // -- Save libCEED data
208   data[level - 1]->op_apply = op_apply;
209   data[level]->op_prolong   = op_prolong;
210   data[level]->op_restrict  = op_restrict;
211 
212   CeedBasisDestroy(&basis_u);
213   PetscFunctionReturn(PETSC_SUCCESS);
214 };
215 
216 PetscErrorCode SetupErrorOperator(DM dm, Ceed ceed, BPData bp_data, CeedInt topo_dim, PetscInt num_comp_x, PetscInt num_comp_u,
217                                   CeedOperator *op_error) {
218   DM                  dm_coord;
219   Vec                 coords;
220   const PetscScalar  *coord_array;
221   CeedBasis           basis_x, basis_u;
222   CeedElemRestriction elem_restr_x, elem_restr_u, elem_restr_u_i, elem_restr_qd_i;
223   CeedQFunction       qf_setup_geo, qf_setup_rhs, qf_error;
224   CeedOperator        op_setup_geo, op_setup_rhs;
225   CeedVector          x_coord, q_data, target, rhs;
226   PetscInt            c_start, c_end, num_elem;
227   CeedInt             num_qpts, q_data_size = bp_data.q_data_size;
228   CeedScalar          R = 1;                         // radius of the sphere
229   CeedScalar          l = 1.0 / PetscSqrtReal(3.0);  // half edge of the inscribed cube
230 
231   PetscFunctionBeginUser;
232   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
233 
234   // CEED bases
235   PetscCall(CreateBasisFromPlex(ceed, dm_coord, 0, 0, 0, 0, bp_data, &basis_x));
236   PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, bp_data, &basis_u));
237 
238   // CEED restrictions
239   PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, 0, 0, 0, &elem_restr_x));
240   PetscCall(CreateRestrictionFromPlex(ceed, dm, 0, 0, 0, &elem_restr_u));
241 
242   PetscCall(DMPlexGetHeightStratum(dm, 0, &c_start, &c_end));
243   num_elem = c_end - c_start;
244   CeedBasisGetNumQuadraturePoints(basis_u, &num_qpts);
245 
246   CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, num_comp_u, num_comp_u * num_elem * num_qpts, CEED_STRIDES_BACKEND, &elem_restr_u_i);
247   CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, q_data_size, q_data_size * num_elem * num_qpts, CEED_STRIDES_BACKEND, &elem_restr_qd_i);
248 
249   // Element coordinates
250   PetscCall(DMGetCoordinatesLocal(dm, &coords));
251   PetscCall(VecGetArrayRead(coords, &coord_array));
252 
253   CeedElemRestrictionCreateVector(elem_restr_x, &x_coord, NULL);
254   CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, (PetscScalar *)coord_array);
255   PetscCall(VecRestoreArrayRead(coords, &coord_array));
256 
257   // Create the persistent vectors that will be needed in setup and apply
258   CeedVectorCreate(ceed, q_data_size * num_elem * num_qpts, &q_data);
259 
260   // Create the QFunction that builds the context data
261   CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_geo, bp_data.setup_geo_loc, &qf_setup_geo);
262   CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP);
263   CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x * topo_dim, CEED_EVAL_GRAD);
264   CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT);
265   CeedQFunctionAddOutput(qf_setup_geo, "qdata", q_data_size, CEED_EVAL_NONE);
266 
267   // Create the operator that builds the quadrature data
268   CeedOperatorCreate(ceed, qf_setup_geo, NULL, NULL, &op_setup_geo);
269   CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
270   CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
271   CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE);
272   CeedOperatorSetField(op_setup_geo, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
273 
274   // Setup q_data
275   CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE);
276 
277   // Set up target vector
278   CeedElemRestrictionCreateVector(elem_restr_u, &rhs, NULL);
279   CeedVectorCreate(ceed, num_elem * num_qpts * num_comp_u, &target);
280   // Create the q-function that sets up the RHS and true solution
281   CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_rhs, bp_data.setup_rhs_loc, &qf_setup_rhs);
282   CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP);
283   CeedQFunctionAddInput(qf_setup_rhs, "qdata", q_data_size, CEED_EVAL_NONE);
284   CeedQFunctionAddOutput(qf_setup_rhs, "true solution", num_comp_u, CEED_EVAL_NONE);
285   CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP);
286 
287   // Create the operator that builds the RHS and true solution
288   CeedOperatorCreate(ceed, qf_setup_rhs, NULL, NULL, &op_setup_rhs);
289   CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
290   CeedOperatorSetField(op_setup_rhs, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
291   CeedOperatorSetField(op_setup_rhs, "true solution", elem_restr_u_i, CEED_BASIS_NONE, target);
292   CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
293 
294   // Set up the libCEED context
295   CeedQFunctionContext ctx_rhs_setup;
296   CeedQFunctionContextCreate(ceed, &ctx_rhs_setup);
297   CeedScalar rhs_setup_data[2] = {R, l};
298   CeedQFunctionContextSetData(ctx_rhs_setup, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof rhs_setup_data, &rhs_setup_data);
299   CeedQFunctionSetContext(qf_setup_rhs, ctx_rhs_setup);
300   CeedQFunctionContextDestroy(&ctx_rhs_setup);
301 
302   // Setup RHS and target
303   CeedOperatorApply(op_setup_rhs, x_coord, rhs, CEED_REQUEST_IMMEDIATE);
304 
305   // Set up error operator
306   // Create the error QFunction
307   CeedQFunctionCreateInterior(ceed, 1, bp_data.error, bp_data.error_loc, &qf_error);
308   CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
309   CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
310   CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE);
311   CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_INTERP);
312 
313   // Create the error operator
314   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, op_error);
315   CeedOperatorSetField(*op_error, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
316   CeedOperatorSetField(*op_error, "true_soln", elem_restr_u_i, CEED_BASIS_NONE, target);
317   CeedOperatorSetField(*op_error, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
318   CeedOperatorSetField(*op_error, "error", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
319 
320   // Cleanup
321   CeedQFunctionDestroy(&qf_setup_rhs);
322   CeedOperatorDestroy(&op_setup_rhs);
323   CeedQFunctionDestroy(&qf_setup_geo);
324   CeedOperatorDestroy(&op_setup_geo);
325   CeedQFunctionDestroy(&qf_error);
326   CeedVectorDestroy(&x_coord);
327   CeedVectorDestroy(&rhs);
328   CeedVectorDestroy(&target);
329   CeedVectorDestroy(&q_data);
330   CeedElemRestrictionDestroy(&elem_restr_u_i);
331   CeedElemRestrictionDestroy(&elem_restr_qd_i);
332   CeedElemRestrictionDestroy(&elem_restr_x);
333   CeedElemRestrictionDestroy(&elem_restr_u);
334   CeedBasisDestroy(&basis_x);
335   CeedBasisDestroy(&basis_u);
336 
337   PetscFunctionReturn(PETSC_SUCCESS);
338 }
339