xref: /libCEED/examples/petsc/src/libceedsetup.c (revision e03fef56705b317edc4a39dfee40c8366660a6d6)
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, PetscBool is_fine_level,
44                                     CeedVector rhs_ceed, 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 = NULL, qf_apply = NULL;
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   if (is_fine_level) {
90     // Create the QFunction that builds the context data
91     CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_geo, bp_data.setup_geo_loc, &qf_setup_geo);
92     CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP);
93     CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x * topo_dim, CEED_EVAL_GRAD);
94     CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT);
95     CeedQFunctionAddOutput(qf_setup_geo, "qdata", q_data_size, CEED_EVAL_NONE);
96 
97     // Create the operator that builds the quadrature data
98     CeedOperatorCreate(ceed, qf_setup_geo, NULL, NULL, &op_setup_geo);
99     CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
100     CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
101     CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE);
102     CeedOperatorSetField(op_setup_geo, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
103 
104     // Setup q_data
105     CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE);
106 
107     // Set up PDE operator
108     CeedInt in_scale  = bp_data.in_mode == CEED_EVAL_GRAD ? topo_dim : 1;
109     CeedInt out_scale = bp_data.out_mode == CEED_EVAL_GRAD ? topo_dim : 1;
110     CeedQFunctionCreateInterior(ceed, 1, bp_data.apply, bp_data.apply_loc, &qf_apply);
111     CeedQFunctionAddInput(qf_apply, "u", num_comp_u * in_scale, bp_data.in_mode);
112     CeedQFunctionAddInput(qf_apply, "qdata", q_data_size, CEED_EVAL_NONE);
113     CeedQFunctionAddOutput(qf_apply, "v", num_comp_u * out_scale, bp_data.out_mode);
114 
115     // Create the mass or diff operator
116     CeedOperatorCreate(ceed, qf_apply, NULL, NULL, &op_apply);
117     CeedOperatorSetField(op_apply, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
118     CeedOperatorSetField(op_apply, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
119     CeedOperatorSetField(op_apply, "v", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
120 
121     // Cleanup
122     CeedQFunctionDestroy(&qf_setup_geo);
123     CeedOperatorDestroy(&op_setup_geo);
124   }
125 
126   // Set up RHS if needed
127   if (setup_rhs) {
128     CeedQFunction qf_setup_rhs;
129     CeedOperator  op_setup_rhs;
130     CeedVectorCreate(ceed, num_elem * num_qpts * num_comp_u, target);
131     // Create the q-function that sets up the RHS and true solution
132     CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_rhs, bp_data.setup_rhs_loc, &qf_setup_rhs);
133     CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP);
134     CeedQFunctionAddInput(qf_setup_rhs, "qdata", q_data_size, CEED_EVAL_NONE);
135     CeedQFunctionAddOutput(qf_setup_rhs, "true solution", num_comp_u, CEED_EVAL_NONE);
136     CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP);
137 
138     // Create the operator that builds the RHS and true solution
139     CeedOperatorCreate(ceed, qf_setup_rhs, NULL, NULL, &op_setup_rhs);
140     CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
141     CeedOperatorSetField(op_setup_rhs, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
142     CeedOperatorSetField(op_setup_rhs, "true solution", elem_restr_u_i, CEED_BASIS_NONE, *target);
143     CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
144 
145     // Set up the libCEED context
146     CeedQFunctionContext ctx_rhs_setup;
147     CeedQFunctionContextCreate(ceed, &ctx_rhs_setup);
148     CeedScalar rhs_setup_data[2] = {R, l};
149     CeedQFunctionContextSetData(ctx_rhs_setup, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof rhs_setup_data, &rhs_setup_data);
150     CeedQFunctionSetContext(qf_setup_rhs, ctx_rhs_setup);
151     CeedQFunctionContextDestroy(&ctx_rhs_setup);
152 
153     // Setup RHS and target
154     CeedOperatorApply(op_setup_rhs, x_coord, rhs_ceed, CEED_REQUEST_IMMEDIATE);
155 
156     // Cleanup
157     CeedQFunctionDestroy(&qf_setup_rhs);
158     CeedOperatorDestroy(&op_setup_rhs);
159   }
160   // Cleanup
161   CeedVectorDestroy(&x_coord);
162 
163   // Save libCEED data required for level
164   data->basis_x         = basis_x;
165   data->basis_u         = basis_u;
166   data->elem_restr_x    = elem_restr_x;
167   data->elem_restr_u    = elem_restr_u;
168   data->elem_restr_u_i  = elem_restr_u_i;
169   data->elem_restr_qd_i = elem_restr_qd_i;
170   data->qf_apply        = qf_apply;
171   data->op_apply        = op_apply;
172   data->q_data          = q_data;
173   data->x_ceed          = x_ceed;
174   data->y_ceed          = y_ceed;
175   data->q_data_size     = q_data_size;
176   PetscFunctionReturn(PETSC_SUCCESS);
177 };
178 
179 // -----------------------------------------------------------------------------
180 // Setup libCEED level transfer operator objects
181 // -----------------------------------------------------------------------------
182 PetscErrorCode CeedLevelTransferSetup(DM dm, Ceed ceed, CeedInt level, CeedInt num_comp_u, CeedData *data, BPData bp_data, Vec fine_mult) {
183   PetscFunctionBeginUser;
184   // Restriction - Fine to corse
185   CeedOperator op_restrict;
186   // Interpolation - Corse to fine
187   CeedOperator op_prolong;
188   // Coarse grid operator
189   CeedOperator op_apply;
190   // Basis
191   CeedBasis basis_u;
192   PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, bp_data, &basis_u));
193 
194   // ---------------------------------------------------------------------------
195   // Coarse Grid, Prolongation, and Restriction Operators
196   // ---------------------------------------------------------------------------
197   // Create the Operators that compute the prolongation and
198   //   restriction between the p-multigrid levels and the coarse grid eval.
199   // ---------------------------------------------------------------------------
200   // Place in libCEED array
201   PetscMemType m_mem_type;
202   PetscCall(VecReadP2C(fine_mult, &m_mem_type, data[level]->x_ceed));
203 
204   CeedOperatorMultigridLevelCreate(data[level]->op_apply, data[level]->x_ceed, data[level - 1]->elem_restr_u, basis_u, &op_apply, &op_prolong,
205                                    &op_restrict);
206 
207   // Restore PETSc vector
208   PetscCall(VecReadC2P(data[level]->x_ceed, m_mem_type, fine_mult));
209   PetscCall(VecZeroEntries(fine_mult));
210   // -- Save libCEED data
211   data[level - 1]->op_apply = op_apply;
212   data[level]->op_prolong   = op_prolong;
213   data[level]->op_restrict  = op_restrict;
214 
215   CeedBasisDestroy(&basis_u);
216   PetscFunctionReturn(PETSC_SUCCESS);
217 };
218 
219 PetscErrorCode SetupErrorOperator(DM dm, Ceed ceed, BPData bp_data, CeedInt topo_dim, PetscInt num_comp_x, PetscInt num_comp_u,
220                                   CeedOperator *op_error) {
221   DM                  dm_coord;
222   Vec                 coords;
223   const PetscScalar  *coord_array;
224   CeedBasis           basis_x, basis_u;
225   CeedElemRestriction elem_restr_x, elem_restr_u, elem_restr_u_i, elem_restr_qd_i;
226   CeedQFunction       qf_setup_geo, qf_setup_rhs, qf_error;
227   CeedOperator        op_setup_geo, op_setup_rhs;
228   CeedVector          x_coord, q_data, target, rhs;
229   PetscInt            c_start, c_end, num_elem;
230   CeedInt             num_qpts, q_data_size = bp_data.q_data_size;
231   CeedScalar          R = 1;                         // radius of the sphere
232   CeedScalar          l = 1.0 / PetscSqrtReal(3.0);  // half edge of the inscribed cube
233 
234   PetscFunctionBeginUser;
235   PetscCall(DMGetCoordinateDM(dm, &dm_coord));
236 
237   // CEED bases
238   PetscCall(CreateBasisFromPlex(ceed, dm_coord, 0, 0, 0, 0, bp_data, &basis_x));
239   PetscCall(CreateBasisFromPlex(ceed, dm, 0, 0, 0, 0, bp_data, &basis_u));
240 
241   // CEED restrictions
242   PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, 0, 0, 0, &elem_restr_x));
243   PetscCall(CreateRestrictionFromPlex(ceed, dm, 0, 0, 0, &elem_restr_u));
244 
245   PetscCall(DMPlexGetHeightStratum(dm, 0, &c_start, &c_end));
246   num_elem = c_end - c_start;
247   CeedBasisGetNumQuadraturePoints(basis_u, &num_qpts);
248 
249   CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, num_comp_u, num_comp_u * num_elem * num_qpts, CEED_STRIDES_BACKEND, &elem_restr_u_i);
250   CeedElemRestrictionCreateStrided(ceed, num_elem, num_qpts, q_data_size, q_data_size * num_elem * num_qpts, CEED_STRIDES_BACKEND, &elem_restr_qd_i);
251 
252   // Element coordinates
253   PetscCall(DMGetCoordinatesLocal(dm, &coords));
254   PetscCall(VecGetArrayRead(coords, &coord_array));
255 
256   CeedElemRestrictionCreateVector(elem_restr_x, &x_coord, NULL);
257   CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, (PetscScalar *)coord_array);
258   PetscCall(VecRestoreArrayRead(coords, &coord_array));
259 
260   // Create the persistent vectors that will be needed in setup and apply
261   CeedVectorCreate(ceed, q_data_size * num_elem * num_qpts, &q_data);
262 
263   // Create the QFunction that builds the context data
264   CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_geo, bp_data.setup_geo_loc, &qf_setup_geo);
265   CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP);
266   CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x * topo_dim, CEED_EVAL_GRAD);
267   CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT);
268   CeedQFunctionAddOutput(qf_setup_geo, "qdata", q_data_size, CEED_EVAL_NONE);
269 
270   // Create the operator that builds the quadrature data
271   CeedOperatorCreate(ceed, qf_setup_geo, NULL, NULL, &op_setup_geo);
272   CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
273   CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
274   CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE);
275   CeedOperatorSetField(op_setup_geo, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
276 
277   // Setup q_data
278   CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE);
279 
280   // Set up target vector
281   CeedElemRestrictionCreateVector(elem_restr_u, &rhs, NULL);
282   CeedVectorCreate(ceed, num_elem * num_qpts * num_comp_u, &target);
283   // Create the q-function that sets up the RHS and true solution
284   CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_rhs, bp_data.setup_rhs_loc, &qf_setup_rhs);
285   CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP);
286   CeedQFunctionAddInput(qf_setup_rhs, "qdata", q_data_size, CEED_EVAL_NONE);
287   CeedQFunctionAddOutput(qf_setup_rhs, "true solution", num_comp_u, CEED_EVAL_NONE);
288   CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP);
289 
290   // Create the operator that builds the RHS and true solution
291   CeedOperatorCreate(ceed, qf_setup_rhs, NULL, NULL, &op_setup_rhs);
292   CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x, CEED_VECTOR_ACTIVE);
293   CeedOperatorSetField(op_setup_rhs, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
294   CeedOperatorSetField(op_setup_rhs, "true solution", elem_restr_u_i, CEED_BASIS_NONE, target);
295   CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
296 
297   // Set up the libCEED context
298   CeedQFunctionContext ctx_rhs_setup;
299   CeedQFunctionContextCreate(ceed, &ctx_rhs_setup);
300   CeedScalar rhs_setup_data[2] = {R, l};
301   CeedQFunctionContextSetData(ctx_rhs_setup, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof rhs_setup_data, &rhs_setup_data);
302   CeedQFunctionSetContext(qf_setup_rhs, ctx_rhs_setup);
303   CeedQFunctionContextDestroy(&ctx_rhs_setup);
304 
305   // Setup RHS and target
306   CeedOperatorApply(op_setup_rhs, x_coord, rhs, CEED_REQUEST_IMMEDIATE);
307 
308   // Set up error operator
309   // Create the error QFunction
310   CeedQFunctionCreateInterior(ceed, 1, bp_data.error, bp_data.error_loc, &qf_error);
311   CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
312   CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
313   CeedQFunctionAddInput(qf_error, "qdata", q_data_size, CEED_EVAL_NONE);
314   CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_INTERP);
315 
316   // Create the error operator
317   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, op_error);
318   CeedOperatorSetField(*op_error, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
319   CeedOperatorSetField(*op_error, "true_soln", elem_restr_u_i, CEED_BASIS_NONE, target);
320   CeedOperatorSetField(*op_error, "qdata", elem_restr_qd_i, CEED_BASIS_NONE, q_data);
321   CeedOperatorSetField(*op_error, "error", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
322 
323   // Cleanup
324   CeedQFunctionDestroy(&qf_setup_rhs);
325   CeedOperatorDestroy(&op_setup_rhs);
326   CeedQFunctionDestroy(&qf_setup_geo);
327   CeedOperatorDestroy(&op_setup_geo);
328   CeedQFunctionDestroy(&qf_error);
329   CeedVectorDestroy(&x_coord);
330   CeedVectorDestroy(&rhs);
331   CeedVectorDestroy(&target);
332   CeedVectorDestroy(&q_data);
333   CeedElemRestrictionDestroy(&elem_restr_u_i);
334   CeedElemRestrictionDestroy(&elem_restr_qd_i);
335   CeedElemRestrictionDestroy(&elem_restr_x);
336   CeedElemRestrictionDestroy(&elem_restr_u);
337   CeedBasisDestroy(&basis_x);
338   CeedBasisDestroy(&basis_u);
339 
340   PetscFunctionReturn(PETSC_SUCCESS);
341 }
342