xref: /libCEED/examples/petsc/src/swarmutils.c (revision 4d00b080eb3f95d2e04e55c0ff369c5c847bb288)
15aed82e4SJeremy L Thompson // Copyright (c) 2017-2024, Lawrence Livermore National Security, LLC and other CEED contributors.
298285ab4SZach Atkins // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
398285ab4SZach Atkins //
498285ab4SZach Atkins // SPDX-License-Identifier: BSD-2-Clause
598285ab4SZach Atkins //
698285ab4SZach Atkins // This file is part of CEED:  http://github.com/ceed
798285ab4SZach Atkins 
8b6972d74SZach Atkins #include "../include/swarmutils.h"
978f7fce3SZach Atkins #include "../include/matops.h"
10b6972d74SZach Atkins #include "../qfunctions/swarm/swarmmass.h"
11b6972d74SZach Atkins 
12b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
13b6972d74SZach Atkins // Context utilities
14b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
15b6972d74SZach Atkins PetscErrorCode DMSwarmCeedContextCreate(DM dm_swarm, const char *ceed_resource, DMSwarmCeedContext *ctx) {
16b6972d74SZach Atkins   DM                  dm_mesh, dm_coord;
17b6972d74SZach Atkins   CeedElemRestriction elem_restr_u_mesh, elem_restr_x_mesh, elem_restr_x_points, elem_restr_u_points, elem_restr_q_data_points;
18b6972d74SZach Atkins   CeedBasis           basis_u, basis_x;
19b6972d74SZach Atkins   CeedVector          x_ref_points, q_data_points;
20b6972d74SZach Atkins   CeedInt             num_comp;
21b6972d74SZach Atkins 
22b6972d74SZach Atkins   PetscFunctionBeginUser;
23b6972d74SZach Atkins   PetscCall(PetscNew(ctx));
24b6972d74SZach Atkins   PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
25b6972d74SZach Atkins   PetscCall(DMGetCoordinateDM(dm_mesh, &dm_coord));
26b6972d74SZach Atkins 
27b6972d74SZach Atkins   CeedInit(ceed_resource, &(*ctx)->ceed);
28b6972d74SZach Atkins   // Background mesh objects
29b6972d74SZach Atkins   {
30b6972d74SZach Atkins     BPData bp_data = {.q_mode = CEED_GAUSS};
31b6972d74SZach Atkins 
32b6972d74SZach Atkins     PetscCall(CreateBasisFromPlex((*ctx)->ceed, dm_mesh, NULL, 0, 0, 0, bp_data, &basis_u));
33b6972d74SZach Atkins     PetscCall(CreateBasisFromPlex((*ctx)->ceed, dm_coord, NULL, 0, 0, 0, bp_data, &basis_x));
34b6972d74SZach Atkins     PetscCall(CreateRestrictionFromPlex((*ctx)->ceed, dm_mesh, 0, NULL, 0, &elem_restr_u_mesh));
35b6972d74SZach Atkins     PetscCall(CreateRestrictionFromPlex((*ctx)->ceed, dm_coord, 0, NULL, 0, &elem_restr_x_mesh));
36b6972d74SZach Atkins 
37b6972d74SZach Atkins     // -- Mesh vectors
38b6972d74SZach Atkins     CeedElemRestrictionCreateVector(elem_restr_u_mesh, &(*ctx)->u_mesh, NULL);
39b6972d74SZach Atkins     CeedElemRestrictionCreateVector(elem_restr_u_mesh, &(*ctx)->v_mesh, NULL);
40b6972d74SZach Atkins   }
41b6972d74SZach Atkins   // Swarm objects
42b6972d74SZach Atkins   {
43b6972d74SZach Atkins     PetscInt        dim;
44b6972d74SZach Atkins     const PetscInt *cell_points;
45b6972d74SZach Atkins     IS              is_points;
46b6972d74SZach Atkins     Vec             X_ref;
47b6972d74SZach Atkins     CeedInt         num_elem;
48b6972d74SZach Atkins 
49b6972d74SZach Atkins     PetscCall(DMSwarmCreateReferenceCoordinates(dm_swarm, &is_points, &X_ref));
50b6972d74SZach Atkins     PetscCall(DMGetDimension(dm_mesh, &dim));
51b6972d74SZach Atkins     CeedElemRestrictionGetNumElements(elem_restr_u_mesh, &num_elem);
52b6972d74SZach Atkins     CeedElemRestrictionGetNumComponents(elem_restr_u_mesh, &num_comp);
53b6972d74SZach Atkins 
54b6972d74SZach Atkins     PetscCall(ISGetIndices(is_points, &cell_points));
55b6972d74SZach Atkins     PetscInt num_points = cell_points[num_elem + 1] - num_elem - 2;
56b6972d74SZach Atkins     CeedInt  offsets[num_elem + 1 + num_points];
57b6972d74SZach Atkins 
58b6972d74SZach Atkins     for (PetscInt i = 0; i < num_elem + 1; i++) offsets[i] = cell_points[i + 1] - 1;
59b6972d74SZach Atkins     for (PetscInt i = num_elem + 1; i < num_points + num_elem + 1; i++) offsets[i] = cell_points[i + 1];
60b6972d74SZach Atkins     PetscCall(ISRestoreIndices(is_points, &cell_points));
61b6972d74SZach Atkins 
62b6972d74SZach Atkins     // -- Points restrictions
63b6972d74SZach Atkins     CeedElemRestrictionCreateAtPoints((*ctx)->ceed, num_elem, num_points, num_comp, num_points * num_comp, CEED_MEM_HOST, CEED_COPY_VALUES, offsets,
64b6972d74SZach Atkins                                       &elem_restr_u_points);
65b6972d74SZach Atkins     CeedElemRestrictionCreateAtPoints((*ctx)->ceed, num_elem, num_points, dim, num_points * dim, CEED_MEM_HOST, CEED_COPY_VALUES, offsets,
66b6972d74SZach Atkins                                       &elem_restr_x_points);
67b6972d74SZach Atkins     CeedElemRestrictionCreateAtPoints((*ctx)->ceed, num_elem, num_points, 1, num_points, CEED_MEM_HOST, CEED_COPY_VALUES, offsets,
68b6972d74SZach Atkins                                       &elem_restr_q_data_points);
69b6972d74SZach Atkins 
70b6972d74SZach Atkins     // -- Points vectors
71b6972d74SZach Atkins     CeedElemRestrictionCreateVector(elem_restr_u_points, &(*ctx)->u_points, NULL);
72b6972d74SZach Atkins     CeedElemRestrictionCreateVector(elem_restr_q_data_points, &q_data_points, NULL);
73b6972d74SZach Atkins 
74b6972d74SZach Atkins     // -- Ref coordinates
75b6972d74SZach Atkins     {
76b6972d74SZach Atkins       PetscMemType       X_mem_type;
77b6972d74SZach Atkins       const PetscScalar *x;
78b6972d74SZach Atkins 
79b6972d74SZach Atkins       CeedVectorCreate((*ctx)->ceed, num_points * dim, &x_ref_points);
80b6972d74SZach Atkins 
81b6972d74SZach Atkins       PetscCall(VecGetArrayReadAndMemType(X_ref, (const PetscScalar **)&x, &X_mem_type));
82b6972d74SZach Atkins       CeedVectorSetArray(x_ref_points, MemTypeP2C(X_mem_type), CEED_COPY_VALUES, (CeedScalar *)x);
83b6972d74SZach Atkins       PetscCall(VecRestoreArrayReadAndMemType(X_ref, (const PetscScalar **)&x));
84b6972d74SZach Atkins     }
85b6972d74SZach Atkins 
86b6972d74SZach Atkins     // Create Q data
87b6972d74SZach Atkins     {
88b6972d74SZach Atkins       CeedQFunction qf_setup;
89b6972d74SZach Atkins       CeedOperator  op_setup;
90b6972d74SZach Atkins       CeedVector    x_coord;
91b6972d74SZach Atkins 
92b6972d74SZach Atkins       {
93b6972d74SZach Atkins         Vec                X_loc;
94*4d00b080SJeremy L Thompson         PetscInt           len;
95b6972d74SZach Atkins         const PetscScalar *x;
96b6972d74SZach Atkins 
97b6972d74SZach Atkins         PetscCall(DMGetCoordinatesLocal(dm_mesh, &X_loc));
98b6972d74SZach Atkins         PetscCall(VecGetLocalSize(X_loc, &len));
99b6972d74SZach Atkins         CeedVectorCreate((*ctx)->ceed, len, &x_coord);
100b6972d74SZach Atkins 
101b6972d74SZach Atkins         PetscCall(VecGetArrayRead(X_loc, &x));
102b6972d74SZach Atkins         CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, (CeedScalar *)x);
103b6972d74SZach Atkins         PetscCall(VecRestoreArrayRead(X_loc, &x));
104b6972d74SZach Atkins       }
105b6972d74SZach Atkins 
106b6972d74SZach Atkins       // Setup geometric scaling
107b6972d74SZach Atkins       CeedQFunctionCreateInterior((*ctx)->ceed, 1, SetupMass, SetupMass_loc, &qf_setup);
108b6972d74SZach Atkins       CeedQFunctionAddInput(qf_setup, "x", dim * dim, CEED_EVAL_GRAD);
109b6972d74SZach Atkins       CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT);
110b6972d74SZach Atkins       CeedQFunctionAddOutput(qf_setup, "rho", 1, CEED_EVAL_NONE);
111b6972d74SZach Atkins 
112b6972d74SZach Atkins       CeedOperatorCreateAtPoints((*ctx)->ceed, qf_setup, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup);
113b6972d74SZach Atkins       CeedOperatorSetField(op_setup, "x", elem_restr_x_mesh, basis_x, CEED_VECTOR_ACTIVE);
114b6972d74SZach Atkins       CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE);
115b6972d74SZach Atkins       CeedOperatorSetField(op_setup, "rho", elem_restr_q_data_points, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
116b6972d74SZach Atkins       CeedOperatorAtPointsSetPoints(op_setup, elem_restr_x_points, x_ref_points);
117b6972d74SZach Atkins 
118b6972d74SZach Atkins       CeedOperatorApply(op_setup, x_coord, q_data_points, CEED_REQUEST_IMMEDIATE);
119b6972d74SZach Atkins 
120b6972d74SZach Atkins       // Cleanup
121b6972d74SZach Atkins       CeedVectorDestroy(&x_coord);
122b6972d74SZach Atkins       CeedQFunctionDestroy(&qf_setup);
123b6972d74SZach Atkins       CeedOperatorDestroy(&op_setup);
124b6972d74SZach Atkins     }
125b6972d74SZach Atkins 
126b6972d74SZach Atkins     // -- Cleanup
127b6972d74SZach Atkins     PetscCall(ISDestroy(&is_points));
128b6972d74SZach Atkins     PetscCall(VecDestroy(&X_ref));
129b6972d74SZach Atkins   }
130b6972d74SZach Atkins 
131b6972d74SZach Atkins   PetscCall(DMSetApplicationContext(dm_mesh, (void *)(*ctx)));
132b6972d74SZach Atkins 
133b6972d74SZach Atkins   // Create operators
134b6972d74SZach Atkins   // Mesh to points interpolation operator
135b6972d74SZach Atkins   {
136b6972d74SZach Atkins     CeedQFunction qf_mesh_to_points;
137b6972d74SZach Atkins 
138b6972d74SZach Atkins     // -- Create operator
139b6972d74SZach Atkins     CeedQFunctionCreateIdentity((*ctx)->ceed, num_comp, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_mesh_to_points);
140b6972d74SZach Atkins 
141b6972d74SZach Atkins     CeedOperatorCreateAtPoints((*ctx)->ceed, qf_mesh_to_points, NULL, NULL, &(*ctx)->op_mesh_to_points);
142b6972d74SZach Atkins     CeedOperatorSetField((*ctx)->op_mesh_to_points, "input", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE);
143b6972d74SZach Atkins     CeedOperatorSetField((*ctx)->op_mesh_to_points, "output", elem_restr_u_points, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
144b6972d74SZach Atkins     CeedOperatorAtPointsSetPoints((*ctx)->op_mesh_to_points, elem_restr_x_points, x_ref_points);
145b6972d74SZach Atkins 
146b6972d74SZach Atkins     // -- Cleanup
147b6972d74SZach Atkins     CeedQFunctionDestroy(&qf_mesh_to_points);
148b6972d74SZach Atkins   }
149b6972d74SZach Atkins 
150b6972d74SZach Atkins   // RHS operator
151b6972d74SZach Atkins   {
152b6972d74SZach Atkins     CeedQFunction        qf_pts_to_mesh;
153b6972d74SZach Atkins     CeedQFunctionContext qf_ctx;
154b6972d74SZach Atkins 
155b6972d74SZach Atkins     // -- Mass QFunction
156b6972d74SZach Atkins     CeedQFunctionCreateInterior((*ctx)->ceed, 1, Mass, Mass_loc, &qf_pts_to_mesh);
157b6972d74SZach Atkins     CeedQFunctionAddInput(qf_pts_to_mesh, "q data", 1, CEED_EVAL_NONE);
158b6972d74SZach Atkins     CeedQFunctionAddInput(qf_pts_to_mesh, "u", num_comp, CEED_EVAL_NONE);
159b6972d74SZach Atkins     CeedQFunctionAddOutput(qf_pts_to_mesh, "v", num_comp, CEED_EVAL_INTERP);
160b6972d74SZach Atkins 
161b6972d74SZach Atkins     // -- QFunction context
162b6972d74SZach Atkins     CeedQFunctionContextCreate((*ctx)->ceed, &qf_ctx);
163b6972d74SZach Atkins     CeedQFunctionContextSetData(qf_ctx, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof(num_comp), &num_comp);
164b6972d74SZach Atkins     CeedQFunctionSetContext(qf_pts_to_mesh, qf_ctx);
165b6972d74SZach Atkins 
166b6972d74SZach Atkins     // -- Mass Operator
167b6972d74SZach Atkins     CeedOperatorCreateAtPoints((*ctx)->ceed, qf_pts_to_mesh, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &(*ctx)->op_points_to_mesh);
168b6972d74SZach Atkins     CeedOperatorSetField((*ctx)->op_points_to_mesh, "q data", elem_restr_q_data_points, CEED_BASIS_NONE, q_data_points);
169b6972d74SZach Atkins     CeedOperatorSetField((*ctx)->op_points_to_mesh, "u", elem_restr_u_points, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
170b6972d74SZach Atkins     CeedOperatorSetField((*ctx)->op_points_to_mesh, "v", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE);
171b6972d74SZach Atkins     CeedOperatorAtPointsSetPoints((*ctx)->op_points_to_mesh, elem_restr_x_points, x_ref_points);
172b6972d74SZach Atkins 
173b6972d74SZach Atkins     // -- Cleanup
174b6972d74SZach Atkins     CeedQFunctionContextDestroy(&qf_ctx);
175b6972d74SZach Atkins     CeedQFunctionDestroy(&qf_pts_to_mesh);
176b6972d74SZach Atkins   }
177b6972d74SZach Atkins 
178b6972d74SZach Atkins   // Mass operator
179b6972d74SZach Atkins   {
180b6972d74SZach Atkins     CeedQFunction        qf_mass;
181b6972d74SZach Atkins     CeedQFunctionContext ctx_mass;
182b6972d74SZach Atkins 
183b6972d74SZach Atkins     // -- Mass QFunction
184b6972d74SZach Atkins     CeedQFunctionCreateInterior((*ctx)->ceed, 1, Mass, Mass_loc, &qf_mass);
185b6972d74SZach Atkins     CeedQFunctionAddInput(qf_mass, "q data", 1, CEED_EVAL_NONE);
186b6972d74SZach Atkins     CeedQFunctionAddInput(qf_mass, "u", num_comp, CEED_EVAL_INTERP);
187b6972d74SZach Atkins     CeedQFunctionAddOutput(qf_mass, "v", num_comp, CEED_EVAL_INTERP);
188b6972d74SZach Atkins 
189b6972d74SZach Atkins     // -- QFunction context
190b6972d74SZach Atkins     CeedQFunctionContextCreate((*ctx)->ceed, &ctx_mass);
191b6972d74SZach Atkins     CeedQFunctionContextSetData(ctx_mass, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof(num_comp), &num_comp);
192b6972d74SZach Atkins     CeedQFunctionSetContext(qf_mass, ctx_mass);
193b6972d74SZach Atkins 
194b6972d74SZach Atkins     // -- Mass Operator
195b6972d74SZach Atkins     CeedOperatorCreateAtPoints((*ctx)->ceed, qf_mass, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &(*ctx)->op_mass);
196b6972d74SZach Atkins     CeedOperatorSetField((*ctx)->op_mass, "q data", elem_restr_q_data_points, CEED_BASIS_NONE, q_data_points);
197b6972d74SZach Atkins     CeedOperatorSetField((*ctx)->op_mass, "u", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE);
198b6972d74SZach Atkins     CeedOperatorSetField((*ctx)->op_mass, "v", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE);
199b6972d74SZach Atkins     CeedOperatorAtPointsSetPoints((*ctx)->op_mass, elem_restr_x_points, x_ref_points);
200b6972d74SZach Atkins 
201b6972d74SZach Atkins     // -- Cleanup
202b6972d74SZach Atkins     CeedQFunctionContextDestroy(&ctx_mass);
203b6972d74SZach Atkins     CeedQFunctionDestroy(&qf_mass);
204b6972d74SZach Atkins   }
205b6972d74SZach Atkins 
206b6972d74SZach Atkins   // Cleanup
207b6972d74SZach Atkins   CeedElemRestrictionDestroy(&elem_restr_u_mesh);
208b6972d74SZach Atkins   CeedElemRestrictionDestroy(&elem_restr_x_mesh);
209b6972d74SZach Atkins   CeedElemRestrictionDestroy(&elem_restr_u_points);
210b6972d74SZach Atkins   CeedElemRestrictionDestroy(&elem_restr_x_points);
211b6972d74SZach Atkins   CeedElemRestrictionDestroy(&elem_restr_q_data_points);
212b6972d74SZach Atkins   CeedBasisDestroy(&basis_u);
213b6972d74SZach Atkins   CeedBasisDestroy(&basis_x);
214b6972d74SZach Atkins   CeedVectorDestroy(&x_ref_points);
215b6972d74SZach Atkins   CeedVectorDestroy(&q_data_points);
216b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
217b6972d74SZach Atkins }
218b6972d74SZach Atkins 
219b6972d74SZach Atkins PetscErrorCode DMSwarmCeedContextDestroy(DMSwarmCeedContext *ctx) {
220b6972d74SZach Atkins   PetscFunctionBeginUser;
221b6972d74SZach Atkins   CeedDestroy(&(*ctx)->ceed);
222b6972d74SZach Atkins   CeedVectorDestroy(&(*ctx)->u_mesh);
223b6972d74SZach Atkins   CeedVectorDestroy(&(*ctx)->v_mesh);
224b6972d74SZach Atkins   CeedVectorDestroy(&(*ctx)->u_points);
225b6972d74SZach Atkins   CeedOperatorDestroy(&(*ctx)->op_mesh_to_points);
226b6972d74SZach Atkins   CeedOperatorDestroy(&(*ctx)->op_points_to_mesh);
227b6972d74SZach Atkins   CeedOperatorDestroy(&(*ctx)->op_mass);
228b6972d74SZach Atkins   PetscCall(PetscFree(*ctx));
229b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
230b6972d74SZach Atkins }
231b6972d74SZach Atkins 
232b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
233b6972d74SZach Atkins // PETSc-libCEED memory space utilities
234b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
235b6972d74SZach Atkins PetscErrorCode DMSwarmPICFieldP2C(DM dm_swarm, const char *field, CeedVector x_ceed) {
236b6972d74SZach Atkins   PetscScalar *x;
237b6972d74SZach Atkins 
238b6972d74SZach Atkins   PetscFunctionBeginUser;
239b6972d74SZach Atkins   PetscCall(DMSwarmGetField(dm_swarm, field, NULL, NULL, (void **)&x));
240b6972d74SZach Atkins   CeedVectorSetArray(x_ceed, CEED_MEM_HOST, CEED_USE_POINTER, (CeedScalar *)x);
241b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
242b6972d74SZach Atkins }
243b6972d74SZach Atkins 
244b6972d74SZach Atkins PetscErrorCode DMSwarmPICFieldC2P(DM dm_swarm, const char *field, CeedVector x_ceed) {
245b6972d74SZach Atkins   PetscScalar *x;
246b6972d74SZach Atkins 
247b6972d74SZach Atkins   PetscFunctionBeginUser;
248b6972d74SZach Atkins   CeedVectorTakeArray(x_ceed, CEED_MEM_HOST, (CeedScalar **)&x);
249b6972d74SZach Atkins   PetscCall(DMSwarmRestoreField(dm_swarm, field, NULL, NULL, (void **)&x));
250b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
251b6972d74SZach Atkins }
252b6972d74SZach Atkins 
253b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
254b6972d74SZach Atkins // Swarm point location utility
255b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
256b6972d74SZach Atkins PetscErrorCode DMSwarmInitalizePointLocations(DM dm_swarm, PointSwarmType point_swarm_type, PetscInt num_points, PetscInt num_points_per_cell) {
257b6972d74SZach Atkins   PetscFunctionBeginUser;
258b6972d74SZach Atkins   switch (point_swarm_type) {
259b6972d74SZach Atkins     case SWARM_GAUSS:
260b6972d74SZach Atkins     case SWARM_UNIFORM: {
261b6972d74SZach Atkins       // -- Set gauss or uniform point locations in each cell
262b6972d74SZach Atkins       PetscInt    num_points_per_cell_1d = round(cbrt(num_points_per_cell * 1.0)), dim = 3;
263b6972d74SZach Atkins       PetscScalar point_coords[num_points_per_cell * 3];
264b6972d74SZach Atkins       CeedScalar  points_1d[num_points_per_cell_1d], weights_1d[num_points_per_cell_1d];
265b6972d74SZach Atkins 
266b6972d74SZach Atkins       if (point_swarm_type == SWARM_GAUSS) {
267b6972d74SZach Atkins         PetscCall(CeedGaussQuadrature(num_points_per_cell_1d, points_1d, weights_1d));
268b6972d74SZach Atkins       } else {
26990dde50eSZach Atkins         for (PetscInt i = 0; i < num_points_per_cell_1d; i++) points_1d[i] = 2.0 * (PetscReal)(i + 0.5) / (PetscReal)num_points_per_cell_1d - 1;
270b6972d74SZach Atkins       }
271b6972d74SZach Atkins       for (PetscInt i = 0; i < num_points_per_cell_1d; i++) {
272b6972d74SZach Atkins         for (PetscInt j = 0; j < num_points_per_cell_1d; j++) {
273b6972d74SZach Atkins           for (PetscInt k = 0; k < num_points_per_cell_1d; k++) {
274b6972d74SZach Atkins             PetscInt p = (i * num_points_per_cell_1d + j) * num_points_per_cell_1d + k;
275b6972d74SZach Atkins 
276b6972d74SZach Atkins             point_coords[p * dim + 0] = points_1d[i];
277b6972d74SZach Atkins             point_coords[p * dim + 1] = points_1d[j];
278b6972d74SZach Atkins             point_coords[p * dim + 2] = points_1d[k];
279b6972d74SZach Atkins           }
280b6972d74SZach Atkins         }
281b6972d74SZach Atkins       }
282b6972d74SZach Atkins       PetscCall(DMSwarmSetPointCoordinatesCellwise(dm_swarm, num_points_per_cell_1d * num_points_per_cell_1d * num_points_per_cell_1d, point_coords));
283b6972d74SZach Atkins     } break;
284b6972d74SZach Atkins     case SWARM_CELL_RANDOM: {
28590dde50eSZach Atkins       DM dm_mesh;
286b6972d74SZach Atkins 
28790dde50eSZach Atkins       PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
28890dde50eSZach Atkins       // DMSwarmSetPointCoordinatesRandom expects local coordinates to be set up to ensure call is non-collective
28990dde50eSZach Atkins       PetscCall(DMGetCoordinatesLocalSetUp(dm_mesh));
29090dde50eSZach Atkins       PetscCall(DMSwarmSetPointCoordinatesRandom(dm_swarm, num_points_per_cell));
291b6972d74SZach Atkins     } break;
292b6972d74SZach Atkins     case SWARM_SINUSOIDAL: {
293b6972d74SZach Atkins       // -- Set points distributed per sinusoidal functions
294b6972d74SZach Atkins       PetscInt     dim = 3;
295b6972d74SZach Atkins       PetscScalar *point_coords;
296b6972d74SZach Atkins       DM           dm_mesh;
297b6972d74SZach Atkins 
298b6972d74SZach Atkins       PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
299b6972d74SZach Atkins       PetscCall(DMGetDimension(dm_mesh, &dim));
300b6972d74SZach Atkins 
301b6972d74SZach Atkins       PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords));
302b6972d74SZach Atkins       for (PetscInt p = 0; p < num_points; p++) {
303b6972d74SZach Atkins         point_coords[p * dim + 0] = -PetscCosReal((PetscReal)(p + 1) / (PetscReal)(num_points + 1) * PETSC_PI);
304b6972d74SZach Atkins         if (dim > 1) point_coords[p * dim + 1] = -PetscSinReal((PetscReal)(p + 1) / (PetscReal)(num_points + 1) * PETSC_PI);
305b6972d74SZach Atkins         if (dim > 2) point_coords[p * dim + 2] = PetscSinReal((PetscReal)(p + 1) / (PetscReal)(num_points + 1) * PETSC_PI);
306b6972d74SZach Atkins       }
307b6972d74SZach Atkins       PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords));
308b6972d74SZach Atkins     } break;
309b6972d74SZach Atkins   }
310b6972d74SZach Atkins   PetscCall(DMSwarmMigrate(dm_swarm, PETSC_TRUE));
311b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
312b6972d74SZach Atkins }
313b6972d74SZach Atkins 
314b6972d74SZach Atkins /*@C
315b6972d74SZach Atkins   DMSwarmCreateReferenceCoordinates - Compute the cell reference coordinates for local DMSwarm points.
316b6972d74SZach Atkins 
317b6972d74SZach Atkins   Collective
318b6972d74SZach Atkins 
319b6972d74SZach Atkins   Input Parameter:
320b6972d74SZach Atkins . dm_swarm  - the `DMSwarm`
321b6972d74SZach Atkins 
322b6972d74SZach Atkins   Output Parameters:
323b6972d74SZach Atkins + is_points    - The IS object for indexing into points per cell
324b6972d74SZach Atkins - X_points_ref - Vec holding the cell reference coordinates for local DMSwarm points
325b6972d74SZach Atkins 
326b6972d74SZach Atkins The index set contains ranges of indices for each local cell. This range contains the indices of every point in the cell.
327b6972d74SZach Atkins 
328b6972d74SZach Atkins ```
329b6972d74SZach Atkins total_num_cells
330b6972d74SZach Atkins cell_0_start_index
331b6972d74SZach Atkins cell_1_start_index
332b6972d74SZach Atkins cell_2_start_index
333b6972d74SZach Atkins ...
334b6972d74SZach Atkins cell_n_start_index
335b6972d74SZach Atkins cell_n_stop_index
336b6972d74SZach Atkins cell_0_point_0
337b6972d74SZach Atkins cell_0_point_0
338b6972d74SZach Atkins ...
339b6972d74SZach Atkins cell_n_point_m
340b6972d74SZach Atkins ```
341b6972d74SZach Atkins 
342b6972d74SZach Atkins   Level: beginner
343b6972d74SZach Atkins 
344b6972d74SZach Atkins .seealso: `DMSwarm`
345b6972d74SZach Atkins @*/
346b6972d74SZach Atkins PetscErrorCode DMSwarmCreateReferenceCoordinates(DM dm_swarm, IS *is_points, Vec *X_points_ref) {
347b6972d74SZach Atkins   PetscInt           cell_start, cell_end, num_cells_local, dim, num_points_local, *cell_points, points_offset;
348b6972d74SZach Atkins   PetscScalar       *coords_points_ref;
349b6972d74SZach Atkins   const PetscScalar *coords_points_true;
350b6972d74SZach Atkins   DM                 dm_mesh;
351b6972d74SZach Atkins 
352b6972d74SZach Atkins   PetscFunctionBeginUser;
353b6972d74SZach Atkins   PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
354b6972d74SZach Atkins 
355b6972d74SZach Atkins   // Create vector to hold reference coordinates
356b6972d74SZach Atkins   {
357b6972d74SZach Atkins     Vec X_points_true;
358b6972d74SZach Atkins 
359b6972d74SZach Atkins     PetscCall(DMSwarmCreateLocalVectorFromField(dm_swarm, DMSwarmPICField_coor, &X_points_true));
360b6972d74SZach Atkins     PetscCall(VecDuplicate(X_points_true, X_points_ref));
361b6972d74SZach Atkins     PetscCall(DMSwarmDestroyLocalVectorFromField(dm_swarm, DMSwarmPICField_coor, &X_points_true));
362b6972d74SZach Atkins   }
363b6972d74SZach Atkins 
364b6972d74SZach Atkins   // Allocate index set array
365b6972d74SZach Atkins   PetscCall(DMPlexGetHeightStratum(dm_mesh, 0, &cell_start, &cell_end));
366b6972d74SZach Atkins   num_cells_local = cell_end - cell_start;
367b6972d74SZach Atkins   points_offset   = num_cells_local + 2;
368b6972d74SZach Atkins   PetscCall(VecGetLocalSize(*X_points_ref, &num_points_local));
369b6972d74SZach Atkins   PetscCall(DMGetDimension(dm_mesh, &dim));
370b6972d74SZach Atkins   num_points_local /= dim;
371b6972d74SZach Atkins   PetscCall(PetscMalloc1(num_points_local + num_cells_local + 2, &cell_points));
372b6972d74SZach Atkins   cell_points[0] = num_cells_local;
373b6972d74SZach Atkins 
374b6972d74SZach Atkins   // Get reference coordinates for each swarm point wrt the elements in the background mesh
375b6972d74SZach Atkins   PetscCall(DMSwarmSortGetAccess(dm_swarm));
376b6972d74SZach Atkins   PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points_true));
377b6972d74SZach Atkins   PetscCall(VecGetArray(*X_points_ref, &coords_points_ref));
378b6972d74SZach Atkins   for (PetscInt cell = cell_start, num_points_processed = 0; cell < cell_end; cell++) {
379b6972d74SZach Atkins     PetscInt *points_in_cell, num_points_in_cell, local_cell = cell - cell_start;
380b6972d74SZach Atkins     PetscReal v[3], J[9], invJ[9], detJ, v0_ref[3] = {-1.0, -1.0, -1.0};
381b6972d74SZach Atkins 
382b6972d74SZach Atkins     PetscCall(DMSwarmSortGetPointsPerCell(dm_swarm, cell, &num_points_in_cell, &points_in_cell));
383b6972d74SZach Atkins     // -- Reference coordinates for swarm points in background mesh element
384b6972d74SZach Atkins     PetscCall(DMPlexComputeCellGeometryFEM(dm_mesh, cell, NULL, v, J, invJ, &detJ));
385b6972d74SZach Atkins     cell_points[local_cell + 1] = num_points_processed + points_offset;
386b6972d74SZach Atkins     for (PetscInt p = 0; p < num_points_in_cell; p++) {
387b6972d74SZach Atkins       const CeedInt point = points_in_cell[p];
388b6972d74SZach Atkins 
389b6972d74SZach Atkins       cell_points[points_offset + (num_points_processed++)] = point;
390b6972d74SZach Atkins       CoordinatesRealToRef(dim, dim, v0_ref, v, invJ, &coords_points_true[point * dim], &coords_points_ref[point * dim]);
391b6972d74SZach Atkins     }
392b6972d74SZach Atkins 
393b6972d74SZach Atkins     // -- Cleanup
394b6972d74SZach Atkins     PetscCall(PetscFree(points_in_cell));
395b6972d74SZach Atkins   }
396b6972d74SZach Atkins   cell_points[points_offset - 1] = num_points_local + points_offset;
397b6972d74SZach Atkins 
398b6972d74SZach Atkins   // Cleanup
399b6972d74SZach Atkins   PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points_true));
400b6972d74SZach Atkins   PetscCall(VecRestoreArray(*X_points_ref, &coords_points_ref));
401b6972d74SZach Atkins   PetscCall(DMSwarmSortRestoreAccess(dm_swarm));
402b6972d74SZach Atkins 
403b6972d74SZach Atkins   // Create index set
404b6972d74SZach Atkins   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, num_points_local + points_offset, cell_points, PETSC_OWN_POINTER, is_points));
405b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
406b6972d74SZach Atkins }
407b6972d74SZach Atkins 
408b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
409b6972d74SZach Atkins // RHS for Swarm to Mesh projection
410b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
41178f7fce3SZach Atkins PetscErrorCode DMSwarmCreateProjectionRHS(DM dm_swarm, const char *field, Vec U_points, Vec B_mesh) {
41278f7fce3SZach Atkins   PetscMemType       B_mem_type, U_mem_type;
413b6972d74SZach Atkins   DM                 dm_mesh;
414b6972d74SZach Atkins   Vec                B_mesh_loc;
41578f7fce3SZach Atkins   PetscBool          has_u_points;
416b6972d74SZach Atkins   DMSwarmCeedContext swarm_ceed_context;
417b6972d74SZach Atkins 
418b6972d74SZach Atkins   PetscFunctionBeginUser;
419b6972d74SZach Atkins   // Get mesh DM
420b6972d74SZach Atkins   PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
421b6972d74SZach Atkins   PetscCall(DMGetApplicationContext(dm_mesh, (void *)&swarm_ceed_context));
422b6972d74SZach Atkins 
42378f7fce3SZach Atkins   // Get swarm access
42478f7fce3SZach Atkins   has_u_points = !!U_points;
42578f7fce3SZach Atkins   if (!has_u_points) {
42678f7fce3SZach Atkins     PetscCall(DMSwarmSortGetAccess(dm_swarm));
42778f7fce3SZach Atkins     PetscCall(DMSwarmCreateLocalVectorFromField(dm_swarm, field, &U_points));
42878f7fce3SZach Atkins   }
42978f7fce3SZach Atkins 
430b6972d74SZach Atkins   // Get mesh values
43178f7fce3SZach Atkins   PetscCall(VecReadP2C(U_points, &U_mem_type, swarm_ceed_context->u_points));
432b6972d74SZach Atkins   PetscCall(DMGetLocalVector(dm_mesh, &B_mesh_loc));
433b6972d74SZach Atkins   PetscCall(VecZeroEntries(B_mesh_loc));
434b6972d74SZach Atkins   PetscCall(VecP2C(B_mesh_loc, &B_mem_type, swarm_ceed_context->v_mesh));
435b6972d74SZach Atkins 
436b6972d74SZach Atkins   // Interpolate field from swarm points to mesh
437b6972d74SZach Atkins   CeedOperatorApply(swarm_ceed_context->op_points_to_mesh, swarm_ceed_context->u_points, swarm_ceed_context->v_mesh, CEED_REQUEST_IMMEDIATE);
438b6972d74SZach Atkins 
439b6972d74SZach Atkins   // Restore PETSc Vecs and Local to Global
44078f7fce3SZach Atkins   PetscCall(VecReadC2P(swarm_ceed_context->u_points, U_mem_type, U_points));
441b6972d74SZach Atkins   PetscCall(VecC2P(swarm_ceed_context->v_mesh, B_mem_type, B_mesh_loc));
442b6972d74SZach Atkins   PetscCall(VecZeroEntries(B_mesh));
443b6972d74SZach Atkins   PetscCall(DMLocalToGlobal(dm_mesh, B_mesh_loc, ADD_VALUES, B_mesh));
444b6972d74SZach Atkins 
44578f7fce3SZach Atkins   // Restore swarm access
44678f7fce3SZach Atkins   if (!has_u_points) {
44778f7fce3SZach Atkins     PetscCall(DMSwarmDestroyLocalVectorFromField(dm_swarm, field, &U_points));
448b6972d74SZach Atkins     PetscCall(DMSwarmSortRestoreAccess(dm_swarm));
44978f7fce3SZach Atkins   }
45078f7fce3SZach Atkins 
45178f7fce3SZach Atkins   // Cleanup
452b6972d74SZach Atkins   PetscCall(DMRestoreLocalVector(dm_mesh, &B_mesh_loc));
453b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
454b6972d74SZach Atkins }
455b6972d74SZach Atkins 
456b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
457b6972d74SZach Atkins // Swarm "mass matrix"
458b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
459b6972d74SZach Atkins PetscErrorCode MatMult_SwarmMass(Mat A, Vec U_mesh, Vec V_mesh) {
460b6972d74SZach Atkins   PetscMemType       U_mem_type, V_mem_type;
461b6972d74SZach Atkins   DM                 dm_mesh;
462b6972d74SZach Atkins   Vec                U_mesh_loc, V_mesh_loc;
463b6972d74SZach Atkins   DMSwarmCeedContext swarm_ceed_context;
464b6972d74SZach Atkins 
465b6972d74SZach Atkins   PetscFunctionBeginUser;
466b6972d74SZach Atkins   // Get mesh DM
467b6972d74SZach Atkins   PetscCall(MatGetDM(A, &dm_mesh));
468b6972d74SZach Atkins   PetscCall(DMGetApplicationContext(dm_mesh, (void *)&swarm_ceed_context));
469b6972d74SZach Atkins 
470b6972d74SZach Atkins   // Global to Local and get PETSc Vec access
471b6972d74SZach Atkins   PetscCall(DMGetLocalVector(dm_mesh, &U_mesh_loc));
472b6972d74SZach Atkins   PetscCall(VecZeroEntries(U_mesh_loc));
473b6972d74SZach Atkins   PetscCall(DMGlobalToLocal(dm_mesh, U_mesh, INSERT_VALUES, U_mesh_loc));
474b6972d74SZach Atkins   PetscCall(VecReadP2C(U_mesh_loc, &U_mem_type, swarm_ceed_context->u_mesh));
475b6972d74SZach Atkins   PetscCall(DMGetLocalVector(dm_mesh, &V_mesh_loc));
476b6972d74SZach Atkins   PetscCall(VecZeroEntries(V_mesh_loc));
477b6972d74SZach Atkins   PetscCall(VecP2C(V_mesh_loc, &V_mem_type, swarm_ceed_context->v_mesh));
478b6972d74SZach Atkins 
479b6972d74SZach Atkins   // Apply swarm mass operator
480b6972d74SZach Atkins   CeedOperatorApply(swarm_ceed_context->op_mass, swarm_ceed_context->u_mesh, swarm_ceed_context->v_mesh, CEED_REQUEST_IMMEDIATE);
481b6972d74SZach Atkins 
482b6972d74SZach Atkins   // Restore PETSc Vecs and Local to Global
483b6972d74SZach Atkins   PetscCall(VecReadC2P(swarm_ceed_context->u_mesh, U_mem_type, U_mesh_loc));
484b6972d74SZach Atkins   PetscCall(VecC2P(swarm_ceed_context->v_mesh, V_mem_type, V_mesh_loc));
485b6972d74SZach Atkins   PetscCall(VecZeroEntries(V_mesh));
486b6972d74SZach Atkins   PetscCall(DMLocalToGlobal(dm_mesh, V_mesh_loc, ADD_VALUES, V_mesh));
487b6972d74SZach Atkins 
488b6972d74SZach Atkins   // Cleanup
489b6972d74SZach Atkins   PetscCall(DMRestoreLocalVector(dm_mesh, &U_mesh_loc));
490b6972d74SZach Atkins   PetscCall(DMRestoreLocalVector(dm_mesh, &V_mesh_loc));
491b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
492b6972d74SZach Atkins }
493b6972d74SZach Atkins 
494b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
495b6972d74SZach Atkins // Swarm to mesh projection
496b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
49778f7fce3SZach Atkins PetscErrorCode DMSwarmProjectFromSwarmToCells(DM dm_swarm, const char *field, Vec U_points, Vec U_mesh) {
498b6972d74SZach Atkins   PetscBool          test_mode;
499b6972d74SZach Atkins   Vec                B_mesh;
500b6972d74SZach Atkins   Mat                M;
501b6972d74SZach Atkins   KSP                ksp;
502b6972d74SZach Atkins   DM                 dm_mesh;
503b6972d74SZach Atkins   DMSwarmCeedContext swarm_ceed_context;
504b6972d74SZach Atkins   MPI_Comm           comm;
505b6972d74SZach Atkins 
506b6972d74SZach Atkins   PetscFunctionBeginUser;
507b6972d74SZach Atkins   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swarm-to-Mesh Projection Options", NULL);
508b6972d74SZach Atkins   PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL));
509b6972d74SZach Atkins   PetscOptionsEnd();
510b6972d74SZach Atkins 
511b6972d74SZach Atkins   comm = PetscObjectComm((PetscObject)dm_swarm);
512b6972d74SZach Atkins   PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
513b6972d74SZach Atkins   PetscCall(DMGetApplicationContext(dm_mesh, (void *)&swarm_ceed_context));
514b6972d74SZach Atkins   PetscCall(VecDuplicate(U_mesh, &B_mesh));
515b6972d74SZach Atkins 
516b6972d74SZach Atkins   // Setup "mass matrix"
517b6972d74SZach Atkins   {
518b6972d74SZach Atkins     PetscInt l_size, g_size;
519b6972d74SZach Atkins 
520b6972d74SZach Atkins     PetscCall(VecGetLocalSize(U_mesh, &l_size));
521b6972d74SZach Atkins     PetscCall(VecGetSize(U_mesh, &g_size));
522b6972d74SZach Atkins     PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, swarm_ceed_context, &M));
523b6972d74SZach Atkins     PetscCall(MatSetDM(M, dm_mesh));
524b6972d74SZach Atkins     PetscCall(MatShellSetOperation(M, MATOP_MULT, (void (*)(void))MatMult_SwarmMass));
525b6972d74SZach Atkins   }
526b6972d74SZach Atkins 
527b6972d74SZach Atkins   // Setup KSP
528b6972d74SZach Atkins   {
529b6972d74SZach Atkins     PC pc;
530b6972d74SZach Atkins 
531b6972d74SZach Atkins     PetscCall(KSPCreate(comm, &ksp));
532b6972d74SZach Atkins     PetscCall(KSPGetPC(ksp, &pc));
533b6972d74SZach Atkins     PetscCall(PCSetType(pc, PCJACOBI));
534b6972d74SZach Atkins     PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
535b6972d74SZach Atkins     PetscCall(KSPSetType(ksp, KSPCG));
536b6972d74SZach Atkins     PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
537b6972d74SZach Atkins     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
538b6972d74SZach Atkins     PetscCall(KSPSetOperators(ksp, M, M));
539b6972d74SZach Atkins     PetscCall(KSPSetFromOptions(ksp));
540b6972d74SZach Atkins     PetscCall(PetscObjectSetName((PetscObject)ksp, "Swarm-to-Mesh Projection"));
541b6972d74SZach Atkins     PetscCall(KSPViewFromOptions(ksp, NULL, "-ksp_projection_view"));
542b6972d74SZach Atkins   }
543b6972d74SZach Atkins 
544b6972d74SZach Atkins   // Setup RHS
54578f7fce3SZach Atkins   PetscCall(DMSwarmCreateProjectionRHS(dm_swarm, field, U_points, B_mesh));
546b6972d74SZach Atkins 
547b6972d74SZach Atkins   // Solve
548b6972d74SZach Atkins   PetscCall(VecZeroEntries(U_mesh));
549b6972d74SZach Atkins   PetscCall(KSPSolve(ksp, B_mesh, U_mesh));
550b6972d74SZach Atkins 
551b6972d74SZach Atkins   // KSP summary
552b6972d74SZach Atkins   {
553b6972d74SZach Atkins     KSPType            ksp_type;
554b6972d74SZach Atkins     KSPConvergedReason reason;
555b6972d74SZach Atkins     PetscReal          rnorm;
556b6972d74SZach Atkins     PetscInt           its;
557b6972d74SZach Atkins     PetscCall(KSPGetType(ksp, &ksp_type));
558b6972d74SZach Atkins     PetscCall(KSPGetConvergedReason(ksp, &reason));
559b6972d74SZach Atkins     PetscCall(KSPGetIterationNumber(ksp, &its));
560b6972d74SZach Atkins     PetscCall(KSPGetResidualNorm(ksp, &rnorm));
561b6972d74SZach Atkins 
562b6972d74SZach Atkins     if (!test_mode || reason < 0 || rnorm > 1e-8) {
563b6972d74SZach Atkins       PetscCall(PetscPrintf(comm,
564b6972d74SZach Atkins                             "Swarm-to-Mesh Projection KSP Solve:\n"
565b6972d74SZach Atkins                             "  KSP type: %s\n"
566b6972d74SZach Atkins                             "  KSP convergence: %s\n"
567b6972d74SZach Atkins                             "  Total KSP iterations: %" PetscInt_FMT "\n"
568b6972d74SZach Atkins                             "  Final rnorm: %e\n",
569b6972d74SZach Atkins                             ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
570b6972d74SZach Atkins     }
571b6972d74SZach Atkins   }
572b6972d74SZach Atkins 
573b6972d74SZach Atkins   // Optional viewing
574b6972d74SZach Atkins   PetscCall(KSPViewFromOptions(ksp, NULL, "-ksp_view"));
575b6972d74SZach Atkins 
576b6972d74SZach Atkins   // Cleanup
577b6972d74SZach Atkins   PetscCall(VecDestroy(&B_mesh));
578b6972d74SZach Atkins   PetscCall(MatDestroy(&M));
579b6972d74SZach Atkins   PetscCall(KSPDestroy(&ksp));
580b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
581b6972d74SZach Atkins }
58278f7fce3SZach Atkins 
58378f7fce3SZach Atkins // ------------------------------------------------------------------------------------------------
58478f7fce3SZach Atkins // BP setup
58578f7fce3SZach Atkins // ------------------------------------------------------------------------------------------------
58678f7fce3SZach Atkins PetscErrorCode SetupProblemSwarm(DM dm_swarm, Ceed ceed, BPData bp_data, CeedData data, PetscBool setup_rhs, Vec rhs, Vec target) {
58778f7fce3SZach Atkins   DM                  dm_mesh, dm_coord;
58878f7fce3SZach Atkins   CeedElemRestriction elem_restr_u_mesh, elem_restr_x_mesh, elem_restr_x_points, elem_restr_u_points, elem_restr_q_data_points;
58978f7fce3SZach Atkins   CeedBasis           basis_u, basis_x;
59078f7fce3SZach Atkins   CeedVector          x_coord, x_ref_points, q_data_points;
591*4d00b080SJeremy L Thompson   PetscInt            X_loc_len, dim;
592*4d00b080SJeremy L Thompson   CeedInt             num_comp, q_data_size = bp_data.q_data_size;
59378f7fce3SZach Atkins   CeedScalar          R = 1;                         // radius of the sphere
59478f7fce3SZach Atkins   CeedScalar          l = 1.0 / PetscSqrtReal(3.0);  // half edge of the inscribed cube
59578f7fce3SZach Atkins   Vec                 X_loc;
59678f7fce3SZach Atkins   PetscMemType        X_mem_type;
59778f7fce3SZach Atkins 
59878f7fce3SZach Atkins   PetscFunctionBeginUser;
59978f7fce3SZach Atkins   PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
60078f7fce3SZach Atkins   PetscCall(DMGetCoordinateDM(dm_mesh, &dm_coord));
60178f7fce3SZach Atkins 
60278f7fce3SZach Atkins   // Get coordinates
60378f7fce3SZach Atkins   PetscCall(DMGetCoordinatesLocal(dm_mesh, &X_loc));
60478f7fce3SZach Atkins   PetscCall(VecGetLocalSize(X_loc, &X_loc_len));
60578f7fce3SZach Atkins   CeedVectorCreate(ceed, X_loc_len, &x_coord);
60678f7fce3SZach Atkins   PetscCall(VecReadP2C(X_loc, &X_mem_type, x_coord));
60778f7fce3SZach Atkins 
60878f7fce3SZach Atkins   // Background mesh objects
60978f7fce3SZach Atkins   PetscCall(CreateBasisFromPlex(ceed, dm_mesh, NULL, 0, 0, 0, bp_data, &basis_u));
61078f7fce3SZach Atkins   PetscCall(CreateBasisFromPlex(ceed, dm_coord, NULL, 0, 0, 0, bp_data, &basis_x));
61178f7fce3SZach Atkins   PetscCall(CreateRestrictionFromPlex(ceed, dm_mesh, 0, NULL, 0, &elem_restr_u_mesh));
61278f7fce3SZach Atkins   PetscCall(CreateRestrictionFromPlex(ceed, dm_coord, 0, NULL, 0, &elem_restr_x_mesh));
61378f7fce3SZach Atkins 
61478f7fce3SZach Atkins   CeedElemRestrictionCreateVector(elem_restr_u_mesh, &data->x_ceed, NULL);
61578f7fce3SZach Atkins   CeedElemRestrictionCreateVector(elem_restr_u_mesh, &data->y_ceed, NULL);
61678f7fce3SZach Atkins 
61778f7fce3SZach Atkins   // Swarm objects
61878f7fce3SZach Atkins   {
61978f7fce3SZach Atkins     const PetscInt *cell_points;
62078f7fce3SZach Atkins     IS              is_points;
62178f7fce3SZach Atkins     Vec             X_ref;
62278f7fce3SZach Atkins     CeedInt         num_elem;
62378f7fce3SZach Atkins 
62478f7fce3SZach Atkins     PetscCall(DMSwarmCreateReferenceCoordinates(dm_swarm, &is_points, &X_ref));
62578f7fce3SZach Atkins     PetscCall(DMGetDimension(dm_mesh, &dim));
62678f7fce3SZach Atkins     CeedElemRestrictionGetNumElements(elem_restr_u_mesh, &num_elem);
62778f7fce3SZach Atkins     CeedElemRestrictionGetNumComponents(elem_restr_u_mesh, &num_comp);
62878f7fce3SZach Atkins 
62978f7fce3SZach Atkins     PetscCall(ISGetIndices(is_points, &cell_points));
63078f7fce3SZach Atkins     PetscInt num_points = cell_points[num_elem + 1] - num_elem - 2;
63178f7fce3SZach Atkins     CeedInt  offsets[num_elem + 1 + num_points];
63278f7fce3SZach Atkins 
63378f7fce3SZach Atkins     for (PetscInt i = 0; i < num_elem + 1; i++) offsets[i] = cell_points[i + 1] - 1;
63478f7fce3SZach Atkins     for (PetscInt i = num_elem + 1; i < num_points + num_elem + 1; i++) offsets[i] = cell_points[i + 1];
63578f7fce3SZach Atkins     PetscCall(ISRestoreIndices(is_points, &cell_points));
63678f7fce3SZach Atkins 
63778f7fce3SZach Atkins     // -- Points restrictions
63878f7fce3SZach Atkins     CeedElemRestrictionCreateAtPoints(ceed, num_elem, num_points, num_comp, num_points * num_comp, CEED_MEM_HOST, CEED_COPY_VALUES, offsets,
63978f7fce3SZach Atkins                                       &elem_restr_u_points);
64078f7fce3SZach Atkins     CeedElemRestrictionCreateAtPoints(ceed, num_elem, num_points, dim, num_points * dim, CEED_MEM_HOST, CEED_COPY_VALUES, offsets,
64178f7fce3SZach Atkins                                       &elem_restr_x_points);
64278f7fce3SZach Atkins     CeedElemRestrictionCreateAtPoints(ceed, num_elem, num_points, q_data_size, num_points * q_data_size, CEED_MEM_HOST, CEED_COPY_VALUES, offsets,
64378f7fce3SZach Atkins                                       &elem_restr_q_data_points);
64478f7fce3SZach Atkins 
64578f7fce3SZach Atkins     // -- Points vectors
64678f7fce3SZach Atkins     CeedElemRestrictionCreateVector(elem_restr_q_data_points, &q_data_points, NULL);
64778f7fce3SZach Atkins 
64878f7fce3SZach Atkins     // -- Ref coordinates
64978f7fce3SZach Atkins     {
65078f7fce3SZach Atkins       PetscMemType       X_mem_type;
65178f7fce3SZach Atkins       const PetscScalar *x;
65278f7fce3SZach Atkins 
65378f7fce3SZach Atkins       CeedVectorCreate(ceed, num_points * dim, &x_ref_points);
65478f7fce3SZach Atkins 
65578f7fce3SZach Atkins       PetscCall(VecGetArrayReadAndMemType(X_ref, (const PetscScalar **)&x, &X_mem_type));
65678f7fce3SZach Atkins       CeedVectorSetArray(x_ref_points, MemTypeP2C(X_mem_type), CEED_COPY_VALUES, (CeedScalar *)x);
65778f7fce3SZach Atkins       PetscCall(VecRestoreArrayReadAndMemType(X_ref, (const PetscScalar **)&x));
65878f7fce3SZach Atkins     }
65978f7fce3SZach Atkins 
66078f7fce3SZach Atkins     // Create Q data
66178f7fce3SZach Atkins     {
66278f7fce3SZach Atkins       CeedQFunction qf_setup;
66378f7fce3SZach Atkins       CeedOperator  op_setup;
66478f7fce3SZach Atkins 
66578f7fce3SZach Atkins       // Setup geometric scaling
66678f7fce3SZach Atkins       CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_geo, bp_data.setup_geo_loc, &qf_setup);
66778f7fce3SZach Atkins       CeedQFunctionAddInput(qf_setup, "x", dim, CEED_EVAL_INTERP);
66878f7fce3SZach Atkins       CeedQFunctionAddInput(qf_setup, "dx", dim * dim, CEED_EVAL_GRAD);
66978f7fce3SZach Atkins       CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT);
67078f7fce3SZach Atkins       CeedQFunctionAddOutput(qf_setup, "qdata", q_data_size, CEED_EVAL_NONE);
67178f7fce3SZach Atkins 
67278f7fce3SZach Atkins       CeedOperatorCreateAtPoints(ceed, qf_setup, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup);
67378f7fce3SZach Atkins       CeedOperatorSetField(op_setup, "x", elem_restr_x_mesh, basis_x, CEED_VECTOR_ACTIVE);
67478f7fce3SZach Atkins       CeedOperatorSetField(op_setup, "dx", elem_restr_x_mesh, basis_x, CEED_VECTOR_ACTIVE);
67578f7fce3SZach Atkins       CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE);
67678f7fce3SZach Atkins       CeedOperatorSetField(op_setup, "qdata", elem_restr_q_data_points, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
67778f7fce3SZach Atkins       CeedOperatorAtPointsSetPoints(op_setup, elem_restr_x_points, x_ref_points);
67878f7fce3SZach Atkins 
67978f7fce3SZach Atkins       CeedOperatorApply(op_setup, x_coord, q_data_points, CEED_REQUEST_IMMEDIATE);
68078f7fce3SZach Atkins 
68178f7fce3SZach Atkins       // Cleanup
68278f7fce3SZach Atkins       CeedQFunctionDestroy(&qf_setup);
68378f7fce3SZach Atkins       CeedOperatorDestroy(&op_setup);
68478f7fce3SZach Atkins     }
68578f7fce3SZach Atkins 
68678f7fce3SZach Atkins     // Cleanup
68778f7fce3SZach Atkins     PetscCall(ISDestroy(&is_points));
68878f7fce3SZach Atkins     PetscCall(VecDestroy(&X_ref));
68978f7fce3SZach Atkins   }
69078f7fce3SZach Atkins 
69178f7fce3SZach Atkins   // Set up PDE operator
69278f7fce3SZach Atkins 
69378f7fce3SZach Atkins   CeedQFunction qf_apply;
69478f7fce3SZach Atkins   CeedOperator  op_apply;
69578f7fce3SZach Atkins   CeedInt       in_scale  = bp_data.in_mode == CEED_EVAL_GRAD ? dim : 1;
69678f7fce3SZach Atkins   CeedInt       out_scale = bp_data.out_mode == CEED_EVAL_GRAD ? dim : 1;
69778f7fce3SZach Atkins 
69878f7fce3SZach Atkins   CeedQFunctionCreateInterior(ceed, 1, bp_data.apply, bp_data.apply_loc, &qf_apply);
69978f7fce3SZach Atkins   CeedQFunctionAddInput(qf_apply, "u", num_comp * in_scale, bp_data.in_mode);
70078f7fce3SZach Atkins   CeedQFunctionAddInput(qf_apply, "qdata", q_data_size, CEED_EVAL_NONE);
70178f7fce3SZach Atkins   CeedQFunctionAddOutput(qf_apply, "v", num_comp * out_scale, bp_data.out_mode);
70278f7fce3SZach Atkins 
70378f7fce3SZach Atkins   // Create the mass or diff operator
70478f7fce3SZach Atkins   CeedOperatorCreateAtPoints(ceed, qf_apply, NULL, NULL, &op_apply);
70578f7fce3SZach Atkins   CeedOperatorSetField(op_apply, "u", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE);
70678f7fce3SZach Atkins   CeedOperatorSetField(op_apply, "qdata", elem_restr_q_data_points, CEED_BASIS_NONE, q_data_points);
70778f7fce3SZach Atkins   CeedOperatorSetField(op_apply, "v", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE);
70878f7fce3SZach Atkins   CeedOperatorAtPointsSetPoints(op_apply, elem_restr_x_points, x_ref_points);
70978f7fce3SZach Atkins 
71078f7fce3SZach Atkins   // Set up RHS if needed
71178f7fce3SZach Atkins   if (setup_rhs) {
71278f7fce3SZach Atkins     CeedQFunction qf_setup_rhs;
71378f7fce3SZach Atkins     CeedOperator  op_setup_rhs;
71478f7fce3SZach Atkins     Vec           rhs_loc;
71578f7fce3SZach Atkins     CeedVector    rhs_ceed, target_ceed;
71678f7fce3SZach Atkins     PetscMemType  rhs_mem_type, target_mem_type;
71778f7fce3SZach Atkins 
71878f7fce3SZach Atkins     // Create RHS vector
71978f7fce3SZach Atkins     PetscCall(DMCreateLocalVector(dm_mesh, &rhs_loc));
72078f7fce3SZach Atkins 
72178f7fce3SZach Atkins     CeedElemRestrictionCreateVector(elem_restr_u_points, &target_ceed, NULL);
72278f7fce3SZach Atkins     CeedElemRestrictionCreateVector(elem_restr_u_mesh, &rhs_ceed, NULL);
72378f7fce3SZach Atkins 
72478f7fce3SZach Atkins     // Create the q-function that sets up the RHS and true solution
72578f7fce3SZach Atkins     CeedQFunctionCreateInterior(ceed, 1, bp_data.setup_rhs, bp_data.setup_rhs_loc, &qf_setup_rhs);
72678f7fce3SZach Atkins     CeedQFunctionAddInput(qf_setup_rhs, "x", dim, CEED_EVAL_INTERP);
72778f7fce3SZach Atkins     CeedQFunctionAddInput(qf_setup_rhs, "qdata", q_data_size, CEED_EVAL_NONE);
72878f7fce3SZach Atkins     CeedQFunctionAddOutput(qf_setup_rhs, "true solution", num_comp, CEED_EVAL_NONE);
72978f7fce3SZach Atkins     CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp, CEED_EVAL_INTERP);
73078f7fce3SZach Atkins 
73178f7fce3SZach Atkins     // Create the operator that builds the RHS and true solution
73278f7fce3SZach Atkins     CeedOperatorCreateAtPoints(ceed, qf_setup_rhs, NULL, NULL, &op_setup_rhs);
73378f7fce3SZach Atkins     CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x_mesh, basis_x, CEED_VECTOR_ACTIVE);
73478f7fce3SZach Atkins     CeedOperatorSetField(op_setup_rhs, "qdata", elem_restr_q_data_points, CEED_BASIS_NONE, q_data_points);
73578f7fce3SZach Atkins     CeedOperatorSetField(op_setup_rhs, "true solution", elem_restr_u_points, CEED_BASIS_NONE, target_ceed);
73678f7fce3SZach Atkins     CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE);
73778f7fce3SZach Atkins     CeedOperatorAtPointsSetPoints(op_setup_rhs, elem_restr_x_points, x_ref_points);
73878f7fce3SZach Atkins 
73978f7fce3SZach Atkins     // Set up the libCEED context
74078f7fce3SZach Atkins     CeedQFunctionContext ctx_rhs_setup;
74178f7fce3SZach Atkins     CeedQFunctionContextCreate(ceed, &ctx_rhs_setup);
74278f7fce3SZach Atkins     CeedScalar rhs_setup_data[2] = {R, l};
74378f7fce3SZach Atkins     CeedQFunctionContextSetData(ctx_rhs_setup, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof rhs_setup_data, &rhs_setup_data);
74478f7fce3SZach Atkins     CeedQFunctionSetContext(qf_setup_rhs, ctx_rhs_setup);
74578f7fce3SZach Atkins     CeedQFunctionContextDestroy(&ctx_rhs_setup);
74678f7fce3SZach Atkins 
74778f7fce3SZach Atkins     // Setup RHS and target
74878f7fce3SZach Atkins     PetscCall(VecP2C(rhs_loc, &rhs_mem_type, rhs_ceed));
74978f7fce3SZach Atkins     PetscCall(VecP2C(target, &target_mem_type, target_ceed));
75078f7fce3SZach Atkins     CeedOperatorApply(op_setup_rhs, x_coord, rhs_ceed, CEED_REQUEST_IMMEDIATE);
75178f7fce3SZach Atkins     PetscCall(VecC2P(rhs_ceed, rhs_mem_type, rhs_loc));
75278f7fce3SZach Atkins     PetscCall(VecC2P(target_ceed, target_mem_type, target));
75378f7fce3SZach Atkins 
75478f7fce3SZach Atkins     // Local-to-global
75578f7fce3SZach Atkins     PetscCall(VecZeroEntries(rhs));
75678f7fce3SZach Atkins     PetscCall(DMLocalToGlobal(dm_mesh, rhs_loc, ADD_VALUES, rhs));
75778f7fce3SZach Atkins 
75878f7fce3SZach Atkins     PetscCall(VecViewFromOptions(rhs, NULL, "-rhs_view"));
75978f7fce3SZach Atkins 
76078f7fce3SZach Atkins     // Cleanup
76178f7fce3SZach Atkins     PetscCall(DMRestoreLocalVector(dm_mesh, &rhs_loc));
76278f7fce3SZach Atkins     CeedVectorDestroy(&rhs_ceed);
76378f7fce3SZach Atkins     CeedVectorDestroy(&target_ceed);
76478f7fce3SZach Atkins     CeedQFunctionDestroy(&qf_setup_rhs);
76578f7fce3SZach Atkins     CeedOperatorDestroy(&op_setup_rhs);
76678f7fce3SZach Atkins   }
76778f7fce3SZach Atkins 
76878f7fce3SZach Atkins   // Save libCEED data
76978f7fce3SZach Atkins   data->basis_x         = basis_x;
77078f7fce3SZach Atkins   data->basis_u         = basis_u;
77178f7fce3SZach Atkins   data->elem_restr_x    = elem_restr_x_mesh;
77278f7fce3SZach Atkins   data->elem_restr_u    = elem_restr_u_mesh;
77378f7fce3SZach Atkins   data->elem_restr_u_i  = elem_restr_u_points;
77478f7fce3SZach Atkins   data->elem_restr_qd_i = elem_restr_q_data_points;
77578f7fce3SZach Atkins   data->qf_apply        = qf_apply;
77678f7fce3SZach Atkins   data->op_apply        = op_apply;
77778f7fce3SZach Atkins   data->q_data          = q_data_points;
77878f7fce3SZach Atkins   data->q_data_size     = q_data_size;
77978f7fce3SZach Atkins 
78078f7fce3SZach Atkins   // Cleanup
78178f7fce3SZach Atkins   PetscCall(VecReadC2P(x_coord, X_mem_type, X_loc));
78278f7fce3SZach Atkins   CeedVectorDestroy(&x_coord);
78378f7fce3SZach Atkins   CeedVectorDestroy(&x_ref_points);
78478f7fce3SZach Atkins   CeedElemRestrictionDestroy(&elem_restr_x_points);
78578f7fce3SZach Atkins 
78678f7fce3SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
78778f7fce3SZach Atkins }
788