xref: /libCEED/examples/petsc/src/swarmutils.c (revision b6972d7456611f84b0e462eb1490bcb662442e6a)
1*b6972d74SZach Atkins #include "../include/swarmutils.h"
2*b6972d74SZach Atkins #include "../qfunctions/swarm/swarmmass.h"
3*b6972d74SZach Atkins 
4*b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
5*b6972d74SZach Atkins // Context utilities
6*b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
7*b6972d74SZach Atkins PetscErrorCode DMSwarmCeedContextCreate(DM dm_swarm, const char *ceed_resource, DMSwarmCeedContext *ctx) {
8*b6972d74SZach Atkins   DM                  dm_mesh, dm_coord;
9*b6972d74SZach Atkins   CeedElemRestriction elem_restr_u_mesh, elem_restr_x_mesh, elem_restr_x_points, elem_restr_u_points, elem_restr_q_data_points;
10*b6972d74SZach Atkins   CeedBasis           basis_u, basis_x;
11*b6972d74SZach Atkins   CeedVector          x_ref_points, q_data_points;
12*b6972d74SZach Atkins   CeedInt             num_comp;
13*b6972d74SZach Atkins 
14*b6972d74SZach Atkins   PetscFunctionBeginUser;
15*b6972d74SZach Atkins   PetscCall(PetscNew(ctx));
16*b6972d74SZach Atkins   PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
17*b6972d74SZach Atkins   PetscCall(DMGetCoordinateDM(dm_mesh, &dm_coord));
18*b6972d74SZach Atkins 
19*b6972d74SZach Atkins   CeedInit(ceed_resource, &(*ctx)->ceed);
20*b6972d74SZach Atkins   // Background mesh objects
21*b6972d74SZach Atkins   {
22*b6972d74SZach Atkins     BPData bp_data = {.q_mode = CEED_GAUSS};
23*b6972d74SZach Atkins 
24*b6972d74SZach Atkins     PetscCall(CreateBasisFromPlex((*ctx)->ceed, dm_mesh, NULL, 0, 0, 0, bp_data, &basis_u));
25*b6972d74SZach Atkins     PetscCall(CreateBasisFromPlex((*ctx)->ceed, dm_coord, NULL, 0, 0, 0, bp_data, &basis_x));
26*b6972d74SZach Atkins     PetscCall(CreateRestrictionFromPlex((*ctx)->ceed, dm_mesh, 0, NULL, 0, &elem_restr_u_mesh));
27*b6972d74SZach Atkins     PetscCall(CreateRestrictionFromPlex((*ctx)->ceed, dm_coord, 0, NULL, 0, &elem_restr_x_mesh));
28*b6972d74SZach Atkins 
29*b6972d74SZach Atkins     // -- Mesh vectors
30*b6972d74SZach Atkins     CeedElemRestrictionCreateVector(elem_restr_u_mesh, &(*ctx)->u_mesh, NULL);
31*b6972d74SZach Atkins     CeedElemRestrictionCreateVector(elem_restr_u_mesh, &(*ctx)->v_mesh, NULL);
32*b6972d74SZach Atkins   }
33*b6972d74SZach Atkins   // Swarm objects
34*b6972d74SZach Atkins   {
35*b6972d74SZach Atkins     PetscInt        dim;
36*b6972d74SZach Atkins     const PetscInt *cell_points;
37*b6972d74SZach Atkins     IS              is_points;
38*b6972d74SZach Atkins     Vec             X_ref;
39*b6972d74SZach Atkins     CeedInt         num_elem;
40*b6972d74SZach Atkins 
41*b6972d74SZach Atkins     PetscCall(DMSwarmCreateReferenceCoordinates(dm_swarm, &is_points, &X_ref));
42*b6972d74SZach Atkins     PetscCall(DMGetDimension(dm_mesh, &dim));
43*b6972d74SZach Atkins     CeedElemRestrictionGetNumElements(elem_restr_u_mesh, &num_elem);
44*b6972d74SZach Atkins     CeedElemRestrictionGetNumComponents(elem_restr_u_mesh, &num_comp);
45*b6972d74SZach Atkins 
46*b6972d74SZach Atkins     PetscCall(ISGetIndices(is_points, &cell_points));
47*b6972d74SZach Atkins     PetscInt num_points = cell_points[num_elem + 1] - num_elem - 2;
48*b6972d74SZach Atkins     CeedInt  offsets[num_elem + 1 + num_points];
49*b6972d74SZach Atkins 
50*b6972d74SZach Atkins     for (PetscInt i = 0; i < num_elem + 1; i++) offsets[i] = cell_points[i + 1] - 1;
51*b6972d74SZach Atkins     for (PetscInt i = num_elem + 1; i < num_points + num_elem + 1; i++) offsets[i] = cell_points[i + 1];
52*b6972d74SZach Atkins     PetscCall(ISRestoreIndices(is_points, &cell_points));
53*b6972d74SZach Atkins 
54*b6972d74SZach Atkins     // -- Points restrictions
55*b6972d74SZach Atkins     CeedElemRestrictionCreateAtPoints((*ctx)->ceed, num_elem, num_points, num_comp, num_points * num_comp, CEED_MEM_HOST, CEED_COPY_VALUES, offsets,
56*b6972d74SZach Atkins                                       &elem_restr_u_points);
57*b6972d74SZach Atkins     CeedElemRestrictionCreateAtPoints((*ctx)->ceed, num_elem, num_points, dim, num_points * dim, CEED_MEM_HOST, CEED_COPY_VALUES, offsets,
58*b6972d74SZach Atkins                                       &elem_restr_x_points);
59*b6972d74SZach Atkins     CeedElemRestrictionCreateAtPoints((*ctx)->ceed, num_elem, num_points, 1, num_points, CEED_MEM_HOST, CEED_COPY_VALUES, offsets,
60*b6972d74SZach Atkins                                       &elem_restr_q_data_points);
61*b6972d74SZach Atkins 
62*b6972d74SZach Atkins     // -- Points vectors
63*b6972d74SZach Atkins     CeedElemRestrictionCreateVector(elem_restr_u_points, &(*ctx)->u_points, NULL);
64*b6972d74SZach Atkins     CeedElemRestrictionCreateVector(elem_restr_q_data_points, &q_data_points, NULL);
65*b6972d74SZach Atkins 
66*b6972d74SZach Atkins     // -- Ref coordinates
67*b6972d74SZach Atkins     {
68*b6972d74SZach Atkins       PetscMemType       X_mem_type;
69*b6972d74SZach Atkins       const PetscScalar *x;
70*b6972d74SZach Atkins 
71*b6972d74SZach Atkins       CeedVectorCreate((*ctx)->ceed, num_points * dim, &x_ref_points);
72*b6972d74SZach Atkins 
73*b6972d74SZach Atkins       PetscCall(VecGetArrayReadAndMemType(X_ref, (const PetscScalar **)&x, &X_mem_type));
74*b6972d74SZach Atkins       CeedVectorSetArray(x_ref_points, MemTypeP2C(X_mem_type), CEED_COPY_VALUES, (CeedScalar *)x);
75*b6972d74SZach Atkins       PetscCall(VecRestoreArrayReadAndMemType(X_ref, (const PetscScalar **)&x));
76*b6972d74SZach Atkins     }
77*b6972d74SZach Atkins 
78*b6972d74SZach Atkins     // Create Q data
79*b6972d74SZach Atkins     {
80*b6972d74SZach Atkins       CeedQFunction qf_setup;
81*b6972d74SZach Atkins       CeedOperator  op_setup;
82*b6972d74SZach Atkins       CeedVector    x_coord;
83*b6972d74SZach Atkins 
84*b6972d74SZach Atkins       {
85*b6972d74SZach Atkins         Vec                X_loc;
86*b6972d74SZach Atkins         CeedInt            len;
87*b6972d74SZach Atkins         const PetscScalar *x;
88*b6972d74SZach Atkins 
89*b6972d74SZach Atkins         PetscCall(DMGetCoordinatesLocal(dm_mesh, &X_loc));
90*b6972d74SZach Atkins         PetscCall(VecGetLocalSize(X_loc, &len));
91*b6972d74SZach Atkins         CeedVectorCreate((*ctx)->ceed, len, &x_coord);
92*b6972d74SZach Atkins 
93*b6972d74SZach Atkins         PetscCall(VecGetArrayRead(X_loc, &x));
94*b6972d74SZach Atkins         CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_COPY_VALUES, (CeedScalar *)x);
95*b6972d74SZach Atkins         PetscCall(VecRestoreArrayRead(X_loc, &x));
96*b6972d74SZach Atkins       }
97*b6972d74SZach Atkins 
98*b6972d74SZach Atkins       // Setup geometric scaling
99*b6972d74SZach Atkins       CeedQFunctionCreateInterior((*ctx)->ceed, 1, SetupMass, SetupMass_loc, &qf_setup);
100*b6972d74SZach Atkins       CeedQFunctionAddInput(qf_setup, "x", dim * dim, CEED_EVAL_GRAD);
101*b6972d74SZach Atkins       CeedQFunctionAddInput(qf_setup, "weight", 1, CEED_EVAL_WEIGHT);
102*b6972d74SZach Atkins       CeedQFunctionAddOutput(qf_setup, "rho", 1, CEED_EVAL_NONE);
103*b6972d74SZach Atkins 
104*b6972d74SZach Atkins       CeedOperatorCreateAtPoints((*ctx)->ceed, qf_setup, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_setup);
105*b6972d74SZach Atkins       CeedOperatorSetField(op_setup, "x", elem_restr_x_mesh, basis_x, CEED_VECTOR_ACTIVE);
106*b6972d74SZach Atkins       CeedOperatorSetField(op_setup, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, CEED_VECTOR_NONE);
107*b6972d74SZach Atkins       CeedOperatorSetField(op_setup, "rho", elem_restr_q_data_points, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
108*b6972d74SZach Atkins       CeedOperatorAtPointsSetPoints(op_setup, elem_restr_x_points, x_ref_points);
109*b6972d74SZach Atkins 
110*b6972d74SZach Atkins       CeedOperatorApply(op_setup, x_coord, q_data_points, CEED_REQUEST_IMMEDIATE);
111*b6972d74SZach Atkins 
112*b6972d74SZach Atkins       // Cleanup
113*b6972d74SZach Atkins       CeedVectorDestroy(&x_coord);
114*b6972d74SZach Atkins       CeedQFunctionDestroy(&qf_setup);
115*b6972d74SZach Atkins       CeedOperatorDestroy(&op_setup);
116*b6972d74SZach Atkins     }
117*b6972d74SZach Atkins 
118*b6972d74SZach Atkins     // -- Cleanup
119*b6972d74SZach Atkins     PetscCall(ISDestroy(&is_points));
120*b6972d74SZach Atkins     PetscCall(VecDestroy(&X_ref));
121*b6972d74SZach Atkins   }
122*b6972d74SZach Atkins 
123*b6972d74SZach Atkins   PetscCall(DMSetApplicationContext(dm_mesh, (void *)(*ctx)));
124*b6972d74SZach Atkins 
125*b6972d74SZach Atkins   // Create operators
126*b6972d74SZach Atkins   // Mesh to points interpolation operator
127*b6972d74SZach Atkins   {
128*b6972d74SZach Atkins     CeedQFunction qf_mesh_to_points;
129*b6972d74SZach Atkins 
130*b6972d74SZach Atkins     // -- Create operator
131*b6972d74SZach Atkins     CeedQFunctionCreateIdentity((*ctx)->ceed, num_comp, CEED_EVAL_INTERP, CEED_EVAL_NONE, &qf_mesh_to_points);
132*b6972d74SZach Atkins 
133*b6972d74SZach Atkins     CeedOperatorCreateAtPoints((*ctx)->ceed, qf_mesh_to_points, NULL, NULL, &(*ctx)->op_mesh_to_points);
134*b6972d74SZach Atkins     CeedOperatorSetField((*ctx)->op_mesh_to_points, "input", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE);
135*b6972d74SZach Atkins     CeedOperatorSetField((*ctx)->op_mesh_to_points, "output", elem_restr_u_points, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
136*b6972d74SZach Atkins     CeedOperatorAtPointsSetPoints((*ctx)->op_mesh_to_points, elem_restr_x_points, x_ref_points);
137*b6972d74SZach Atkins 
138*b6972d74SZach Atkins     // -- Cleanup
139*b6972d74SZach Atkins     CeedQFunctionDestroy(&qf_mesh_to_points);
140*b6972d74SZach Atkins   }
141*b6972d74SZach Atkins 
142*b6972d74SZach Atkins   // RHS operator
143*b6972d74SZach Atkins   {
144*b6972d74SZach Atkins     CeedQFunction        qf_pts_to_mesh;
145*b6972d74SZach Atkins     CeedQFunctionContext qf_ctx;
146*b6972d74SZach Atkins 
147*b6972d74SZach Atkins     // -- Mass QFunction
148*b6972d74SZach Atkins     CeedQFunctionCreateInterior((*ctx)->ceed, 1, Mass, Mass_loc, &qf_pts_to_mesh);
149*b6972d74SZach Atkins     CeedQFunctionAddInput(qf_pts_to_mesh, "q data", 1, CEED_EVAL_NONE);
150*b6972d74SZach Atkins     CeedQFunctionAddInput(qf_pts_to_mesh, "u", num_comp, CEED_EVAL_NONE);
151*b6972d74SZach Atkins     CeedQFunctionAddOutput(qf_pts_to_mesh, "v", num_comp, CEED_EVAL_INTERP);
152*b6972d74SZach Atkins 
153*b6972d74SZach Atkins     // -- QFunction context
154*b6972d74SZach Atkins     CeedQFunctionContextCreate((*ctx)->ceed, &qf_ctx);
155*b6972d74SZach Atkins     CeedQFunctionContextSetData(qf_ctx, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof(num_comp), &num_comp);
156*b6972d74SZach Atkins     CeedQFunctionSetContext(qf_pts_to_mesh, qf_ctx);
157*b6972d74SZach Atkins 
158*b6972d74SZach Atkins     // -- Mass Operator
159*b6972d74SZach Atkins     CeedOperatorCreateAtPoints((*ctx)->ceed, qf_pts_to_mesh, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &(*ctx)->op_points_to_mesh);
160*b6972d74SZach Atkins     CeedOperatorSetField((*ctx)->op_points_to_mesh, "q data", elem_restr_q_data_points, CEED_BASIS_NONE, q_data_points);
161*b6972d74SZach Atkins     CeedOperatorSetField((*ctx)->op_points_to_mesh, "u", elem_restr_u_points, CEED_BASIS_NONE, CEED_VECTOR_ACTIVE);
162*b6972d74SZach Atkins     CeedOperatorSetField((*ctx)->op_points_to_mesh, "v", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE);
163*b6972d74SZach Atkins     CeedOperatorAtPointsSetPoints((*ctx)->op_points_to_mesh, elem_restr_x_points, x_ref_points);
164*b6972d74SZach Atkins 
165*b6972d74SZach Atkins     // -- Cleanup
166*b6972d74SZach Atkins     CeedQFunctionContextDestroy(&qf_ctx);
167*b6972d74SZach Atkins     CeedQFunctionDestroy(&qf_pts_to_mesh);
168*b6972d74SZach Atkins   }
169*b6972d74SZach Atkins 
170*b6972d74SZach Atkins   // Mass operator
171*b6972d74SZach Atkins   {
172*b6972d74SZach Atkins     CeedQFunction        qf_mass;
173*b6972d74SZach Atkins     CeedQFunctionContext ctx_mass;
174*b6972d74SZach Atkins 
175*b6972d74SZach Atkins     // -- Mass QFunction
176*b6972d74SZach Atkins     CeedQFunctionCreateInterior((*ctx)->ceed, 1, Mass, Mass_loc, &qf_mass);
177*b6972d74SZach Atkins     CeedQFunctionAddInput(qf_mass, "q data", 1, CEED_EVAL_NONE);
178*b6972d74SZach Atkins     CeedQFunctionAddInput(qf_mass, "u", num_comp, CEED_EVAL_INTERP);
179*b6972d74SZach Atkins     CeedQFunctionAddOutput(qf_mass, "v", num_comp, CEED_EVAL_INTERP);
180*b6972d74SZach Atkins 
181*b6972d74SZach Atkins     // -- QFunction context
182*b6972d74SZach Atkins     CeedQFunctionContextCreate((*ctx)->ceed, &ctx_mass);
183*b6972d74SZach Atkins     CeedQFunctionContextSetData(ctx_mass, CEED_MEM_HOST, CEED_COPY_VALUES, sizeof(num_comp), &num_comp);
184*b6972d74SZach Atkins     CeedQFunctionSetContext(qf_mass, ctx_mass);
185*b6972d74SZach Atkins 
186*b6972d74SZach Atkins     // -- Mass Operator
187*b6972d74SZach Atkins     CeedOperatorCreateAtPoints((*ctx)->ceed, qf_mass, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &(*ctx)->op_mass);
188*b6972d74SZach Atkins     CeedOperatorSetField((*ctx)->op_mass, "q data", elem_restr_q_data_points, CEED_BASIS_NONE, q_data_points);
189*b6972d74SZach Atkins     CeedOperatorSetField((*ctx)->op_mass, "u", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE);
190*b6972d74SZach Atkins     CeedOperatorSetField((*ctx)->op_mass, "v", elem_restr_u_mesh, basis_u, CEED_VECTOR_ACTIVE);
191*b6972d74SZach Atkins     CeedOperatorAtPointsSetPoints((*ctx)->op_mass, elem_restr_x_points, x_ref_points);
192*b6972d74SZach Atkins 
193*b6972d74SZach Atkins     // -- Cleanup
194*b6972d74SZach Atkins     CeedQFunctionContextDestroy(&ctx_mass);
195*b6972d74SZach Atkins     CeedQFunctionDestroy(&qf_mass);
196*b6972d74SZach Atkins   }
197*b6972d74SZach Atkins 
198*b6972d74SZach Atkins   // Cleanup
199*b6972d74SZach Atkins   CeedElemRestrictionDestroy(&elem_restr_u_mesh);
200*b6972d74SZach Atkins   CeedElemRestrictionDestroy(&elem_restr_x_mesh);
201*b6972d74SZach Atkins   CeedElemRestrictionDestroy(&elem_restr_u_points);
202*b6972d74SZach Atkins   CeedElemRestrictionDestroy(&elem_restr_x_points);
203*b6972d74SZach Atkins   CeedElemRestrictionDestroy(&elem_restr_q_data_points);
204*b6972d74SZach Atkins   CeedBasisDestroy(&basis_u);
205*b6972d74SZach Atkins   CeedBasisDestroy(&basis_x);
206*b6972d74SZach Atkins   CeedVectorDestroy(&x_ref_points);
207*b6972d74SZach Atkins   CeedVectorDestroy(&q_data_points);
208*b6972d74SZach Atkins 
209*b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
210*b6972d74SZach Atkins }
211*b6972d74SZach Atkins 
212*b6972d74SZach Atkins PetscErrorCode DMSwarmCeedContextDestroy(DMSwarmCeedContext *ctx) {
213*b6972d74SZach Atkins   PetscFunctionBeginUser;
214*b6972d74SZach Atkins   CeedDestroy(&(*ctx)->ceed);
215*b6972d74SZach Atkins   CeedVectorDestroy(&(*ctx)->u_mesh);
216*b6972d74SZach Atkins   CeedVectorDestroy(&(*ctx)->v_mesh);
217*b6972d74SZach Atkins   CeedVectorDestroy(&(*ctx)->u_points);
218*b6972d74SZach Atkins   CeedOperatorDestroy(&(*ctx)->op_mesh_to_points);
219*b6972d74SZach Atkins   CeedOperatorDestroy(&(*ctx)->op_points_to_mesh);
220*b6972d74SZach Atkins   CeedOperatorDestroy(&(*ctx)->op_mass);
221*b6972d74SZach Atkins   PetscCall(PetscFree(*ctx));
222*b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
223*b6972d74SZach Atkins }
224*b6972d74SZach Atkins 
225*b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
226*b6972d74SZach Atkins // PETSc-libCEED memory space utilities
227*b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
228*b6972d74SZach Atkins PetscErrorCode DMSwarmPICFieldP2C(DM dm_swarm, const char *field, CeedVector x_ceed) {
229*b6972d74SZach Atkins   PetscScalar *x;
230*b6972d74SZach Atkins 
231*b6972d74SZach Atkins   PetscFunctionBeginUser;
232*b6972d74SZach Atkins   PetscCall(DMSwarmGetField(dm_swarm, field, NULL, NULL, (void **)&x));
233*b6972d74SZach Atkins   CeedVectorSetArray(x_ceed, CEED_MEM_HOST, CEED_USE_POINTER, (CeedScalar *)x);
234*b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
235*b6972d74SZach Atkins }
236*b6972d74SZach Atkins 
237*b6972d74SZach Atkins PetscErrorCode DMSwarmPICFieldC2P(DM dm_swarm, const char *field, CeedVector x_ceed) {
238*b6972d74SZach Atkins   PetscScalar *x;
239*b6972d74SZach Atkins 
240*b6972d74SZach Atkins   PetscFunctionBeginUser;
241*b6972d74SZach Atkins   CeedVectorTakeArray(x_ceed, CEED_MEM_HOST, (CeedScalar **)&x);
242*b6972d74SZach Atkins   PetscCall(DMSwarmRestoreField(dm_swarm, field, NULL, NULL, (void **)&x));
243*b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
244*b6972d74SZach Atkins }
245*b6972d74SZach Atkins 
246*b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
247*b6972d74SZach Atkins // Swarm point location utility
248*b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
249*b6972d74SZach Atkins PetscErrorCode DMSwarmInitalizePointLocations(DM dm_swarm, PointSwarmType point_swarm_type, PetscInt num_points, PetscInt num_points_per_cell) {
250*b6972d74SZach Atkins   PetscFunctionBeginUser;
251*b6972d74SZach Atkins 
252*b6972d74SZach Atkins   switch (point_swarm_type) {
253*b6972d74SZach Atkins     case SWARM_GAUSS:
254*b6972d74SZach Atkins     case SWARM_UNIFORM: {
255*b6972d74SZach Atkins       // -- Set gauss or uniform point locations in each cell
256*b6972d74SZach Atkins       PetscInt    num_points_per_cell_1d = round(cbrt(num_points_per_cell * 1.0)), dim = 3;
257*b6972d74SZach Atkins       PetscScalar point_coords[num_points_per_cell * 3];
258*b6972d74SZach Atkins       CeedScalar  points_1d[num_points_per_cell_1d], weights_1d[num_points_per_cell_1d];
259*b6972d74SZach Atkins 
260*b6972d74SZach Atkins       if (point_swarm_type == SWARM_GAUSS) {
261*b6972d74SZach Atkins         PetscCall(CeedGaussQuadrature(num_points_per_cell_1d, points_1d, weights_1d));
262*b6972d74SZach Atkins       } else {
263*b6972d74SZach Atkins         for (PetscInt i = 0; i < num_points_per_cell_1d; i++) points_1d[i] = 2.0 * (PetscReal)(i + 1) / (PetscReal)(num_points_per_cell_1d + 1) - 1;
264*b6972d74SZach Atkins       }
265*b6972d74SZach Atkins       for (PetscInt i = 0; i < num_points_per_cell_1d; i++) {
266*b6972d74SZach Atkins         for (PetscInt j = 0; j < num_points_per_cell_1d; j++) {
267*b6972d74SZach Atkins           for (PetscInt k = 0; k < num_points_per_cell_1d; k++) {
268*b6972d74SZach Atkins             PetscInt p = (i * num_points_per_cell_1d + j) * num_points_per_cell_1d + k;
269*b6972d74SZach Atkins 
270*b6972d74SZach Atkins             point_coords[p * dim + 0] = points_1d[i];
271*b6972d74SZach Atkins             point_coords[p * dim + 1] = points_1d[j];
272*b6972d74SZach Atkins             point_coords[p * dim + 2] = points_1d[k];
273*b6972d74SZach Atkins           }
274*b6972d74SZach Atkins         }
275*b6972d74SZach Atkins       }
276*b6972d74SZach Atkins       PetscCall(DMSwarmSetPointCoordinatesCellwise(dm_swarm, num_points_per_cell_1d * num_points_per_cell_1d * num_points_per_cell_1d, point_coords));
277*b6972d74SZach Atkins     } break;
278*b6972d74SZach Atkins     case SWARM_CELL_RANDOM: {
279*b6972d74SZach Atkins       // -- Set points randomly in each cell
280*b6972d74SZach Atkins       PetscInt     dim = 3, num_cells_total = 1, num_cells[] = {1, 1, 1};
281*b6972d74SZach Atkins       PetscScalar *point_coords;
282*b6972d74SZach Atkins       PetscRandom  rng;
283*b6972d74SZach Atkins 
284*b6972d74SZach Atkins       PetscOptionsBegin(PetscObjectComm((PetscObject)dm_swarm), NULL, "libCEED example using PETSc with DMSwarm", NULL);
285*b6972d74SZach Atkins 
286*b6972d74SZach Atkins       PetscCall(PetscOptionsInt("-dm_plex_dim", "Background mesh dimension", NULL, dim, &dim, NULL));
287*b6972d74SZach Atkins       PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of cells", NULL, num_cells, &dim, NULL));
288*b6972d74SZach Atkins 
289*b6972d74SZach Atkins       PetscOptionsEnd();
290*b6972d74SZach Atkins 
291*b6972d74SZach Atkins       PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm_swarm), &rng));
292*b6972d74SZach Atkins 
293*b6972d74SZach Atkins       num_cells_total = num_cells[0] * num_cells[1] * num_cells[2];
294*b6972d74SZach Atkins       PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords));
295*b6972d74SZach Atkins       for (PetscInt c = 0; c < num_cells_total; c++) {
296*b6972d74SZach Atkins         PetscInt cell_index[3] = {c % num_cells[0], (c / num_cells[0]) % num_cells[1], (c / num_cells[0] / num_cells[1]) % num_cells[2]};
297*b6972d74SZach Atkins 
298*b6972d74SZach Atkins         for (PetscInt p = 0; p < num_points_per_cell; p++) {
299*b6972d74SZach Atkins           PetscInt    point_index = c * num_points_per_cell + p;
300*b6972d74SZach Atkins           PetscScalar random_value;
301*b6972d74SZach Atkins 
302*b6972d74SZach Atkins           for (PetscInt i = 0; i < dim; i++) {
303*b6972d74SZach Atkins             PetscCall(PetscRandomGetValue(rng, &random_value));
304*b6972d74SZach Atkins             point_coords[point_index * dim + i] = -1.0 + cell_index[i] * 2.0 / (num_cells[i] + 1.0) + random_value;
305*b6972d74SZach Atkins           }
306*b6972d74SZach Atkins         }
307*b6972d74SZach Atkins       }
308*b6972d74SZach Atkins       PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords));
309*b6972d74SZach Atkins       PetscCall(PetscRandomDestroy(&rng));
310*b6972d74SZach Atkins     } break;
311*b6972d74SZach Atkins     case SWARM_SINUSOIDAL: {
312*b6972d74SZach Atkins       // -- Set points distributed per sinusoidal functions
313*b6972d74SZach Atkins       PetscInt     dim = 3;
314*b6972d74SZach Atkins       PetscScalar *point_coords;
315*b6972d74SZach Atkins       DM           dm_mesh;
316*b6972d74SZach Atkins 
317*b6972d74SZach Atkins       PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
318*b6972d74SZach Atkins       PetscCall(DMGetDimension(dm_mesh, &dim));
319*b6972d74SZach Atkins 
320*b6972d74SZach Atkins       PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords));
321*b6972d74SZach Atkins       for (PetscInt p = 0; p < num_points; p++) {
322*b6972d74SZach Atkins         point_coords[p * dim + 0] = -PetscCosReal((PetscReal)(p + 1) / (PetscReal)(num_points + 1) * PETSC_PI);
323*b6972d74SZach Atkins         if (dim > 1) point_coords[p * dim + 1] = -PetscSinReal((PetscReal)(p + 1) / (PetscReal)(num_points + 1) * PETSC_PI);
324*b6972d74SZach Atkins         if (dim > 2) point_coords[p * dim + 2] = PetscSinReal((PetscReal)(p + 1) / (PetscReal)(num_points + 1) * PETSC_PI);
325*b6972d74SZach Atkins       }
326*b6972d74SZach Atkins       PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords));
327*b6972d74SZach Atkins     } break;
328*b6972d74SZach Atkins   }
329*b6972d74SZach Atkins   PetscCall(DMSwarmMigrate(dm_swarm, PETSC_TRUE));
330*b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
331*b6972d74SZach Atkins }
332*b6972d74SZach Atkins 
333*b6972d74SZach Atkins /*@C
334*b6972d74SZach Atkins   DMSwarmCreateReferenceCoordinates - Compute the cell reference coordinates for local DMSwarm points.
335*b6972d74SZach Atkins 
336*b6972d74SZach Atkins   Collective
337*b6972d74SZach Atkins 
338*b6972d74SZach Atkins   Input Parameter:
339*b6972d74SZach Atkins . dm_swarm  - the `DMSwarm`
340*b6972d74SZach Atkins 
341*b6972d74SZach Atkins   Output Parameters:
342*b6972d74SZach Atkins + is_points    - The IS object for indexing into points per cell
343*b6972d74SZach Atkins - X_points_ref - Vec holding the cell reference coordinates for local DMSwarm points
344*b6972d74SZach Atkins 
345*b6972d74SZach Atkins The index set contains ranges of indices for each local cell. This range contains the indices of every point in the cell.
346*b6972d74SZach Atkins 
347*b6972d74SZach Atkins ```
348*b6972d74SZach Atkins total_num_cells
349*b6972d74SZach Atkins cell_0_start_index
350*b6972d74SZach Atkins cell_1_start_index
351*b6972d74SZach Atkins cell_2_start_index
352*b6972d74SZach Atkins ...
353*b6972d74SZach Atkins cell_n_start_index
354*b6972d74SZach Atkins cell_n_stop_index
355*b6972d74SZach Atkins cell_0_point_0
356*b6972d74SZach Atkins cell_0_point_0
357*b6972d74SZach Atkins ...
358*b6972d74SZach Atkins cell_n_point_m
359*b6972d74SZach Atkins ```
360*b6972d74SZach Atkins 
361*b6972d74SZach Atkins   Level: beginner
362*b6972d74SZach Atkins 
363*b6972d74SZach Atkins .seealso: `DMSwarm`
364*b6972d74SZach Atkins @*/
365*b6972d74SZach Atkins PetscErrorCode DMSwarmCreateReferenceCoordinates(DM dm_swarm, IS *is_points, Vec *X_points_ref) {
366*b6972d74SZach Atkins   PetscInt           cell_start, cell_end, num_cells_local, dim, num_points_local, *cell_points, points_offset;
367*b6972d74SZach Atkins   PetscScalar       *coords_points_ref;
368*b6972d74SZach Atkins   const PetscScalar *coords_points_true;
369*b6972d74SZach Atkins   DM                 dm_mesh;
370*b6972d74SZach Atkins 
371*b6972d74SZach Atkins   PetscFunctionBeginUser;
372*b6972d74SZach Atkins   PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
373*b6972d74SZach Atkins 
374*b6972d74SZach Atkins   // Create vector to hold reference coordinates
375*b6972d74SZach Atkins   {
376*b6972d74SZach Atkins     Vec X_points_true;
377*b6972d74SZach Atkins 
378*b6972d74SZach Atkins     PetscCall(DMSwarmCreateLocalVectorFromField(dm_swarm, DMSwarmPICField_coor, &X_points_true));
379*b6972d74SZach Atkins     PetscCall(VecDuplicate(X_points_true, X_points_ref));
380*b6972d74SZach Atkins     PetscCall(DMSwarmDestroyLocalVectorFromField(dm_swarm, DMSwarmPICField_coor, &X_points_true));
381*b6972d74SZach Atkins   }
382*b6972d74SZach Atkins 
383*b6972d74SZach Atkins   // Allocate index set array
384*b6972d74SZach Atkins   PetscCall(DMPlexGetHeightStratum(dm_mesh, 0, &cell_start, &cell_end));
385*b6972d74SZach Atkins   num_cells_local = cell_end - cell_start;
386*b6972d74SZach Atkins   points_offset   = num_cells_local + 2;
387*b6972d74SZach Atkins   PetscCall(VecGetLocalSize(*X_points_ref, &num_points_local));
388*b6972d74SZach Atkins   PetscCall(DMGetDimension(dm_mesh, &dim));
389*b6972d74SZach Atkins   num_points_local /= dim;
390*b6972d74SZach Atkins   PetscCall(PetscMalloc1(num_points_local + num_cells_local + 2, &cell_points));
391*b6972d74SZach Atkins   cell_points[0] = num_cells_local;
392*b6972d74SZach Atkins 
393*b6972d74SZach Atkins   // Get reference coordinates for each swarm point wrt the elements in the background mesh
394*b6972d74SZach Atkins   PetscCall(DMSwarmSortGetAccess(dm_swarm));
395*b6972d74SZach Atkins   PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points_true));
396*b6972d74SZach Atkins   PetscCall(VecGetArray(*X_points_ref, &coords_points_ref));
397*b6972d74SZach Atkins   for (PetscInt cell = cell_start, num_points_processed = 0; cell < cell_end; cell++) {
398*b6972d74SZach Atkins     PetscInt *points_in_cell, num_points_in_cell, local_cell = cell - cell_start;
399*b6972d74SZach Atkins     PetscReal v[3], J[9], invJ[9], detJ, v0_ref[3] = {-1.0, -1.0, -1.0};
400*b6972d74SZach Atkins 
401*b6972d74SZach Atkins     PetscCall(DMSwarmSortGetPointsPerCell(dm_swarm, cell, &num_points_in_cell, &points_in_cell));
402*b6972d74SZach Atkins     // -- Reference coordinates for swarm points in background mesh element
403*b6972d74SZach Atkins     PetscCall(DMPlexComputeCellGeometryFEM(dm_mesh, cell, NULL, v, J, invJ, &detJ));
404*b6972d74SZach Atkins     cell_points[local_cell + 1] = num_points_processed + points_offset;
405*b6972d74SZach Atkins     for (PetscInt p = 0; p < num_points_in_cell; p++) {
406*b6972d74SZach Atkins       const CeedInt point = points_in_cell[p];
407*b6972d74SZach Atkins 
408*b6972d74SZach Atkins       cell_points[points_offset + (num_points_processed++)] = point;
409*b6972d74SZach Atkins       CoordinatesRealToRef(dim, dim, v0_ref, v, invJ, &coords_points_true[point * dim], &coords_points_ref[point * dim]);
410*b6972d74SZach Atkins     }
411*b6972d74SZach Atkins 
412*b6972d74SZach Atkins     // -- Cleanup
413*b6972d74SZach Atkins     PetscCall(PetscFree(points_in_cell));
414*b6972d74SZach Atkins   }
415*b6972d74SZach Atkins   cell_points[points_offset - 1] = num_points_local + points_offset;
416*b6972d74SZach Atkins 
417*b6972d74SZach Atkins   // Cleanup
418*b6972d74SZach Atkins   PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points_true));
419*b6972d74SZach Atkins   PetscCall(VecRestoreArray(*X_points_ref, &coords_points_ref));
420*b6972d74SZach Atkins   PetscCall(DMSwarmSortRestoreAccess(dm_swarm));
421*b6972d74SZach Atkins 
422*b6972d74SZach Atkins   // Create index set
423*b6972d74SZach Atkins   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, num_points_local + points_offset, cell_points, PETSC_OWN_POINTER, is_points));
424*b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
425*b6972d74SZach Atkins }
426*b6972d74SZach Atkins 
427*b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
428*b6972d74SZach Atkins // RHS for Swarm to Mesh projection
429*b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
430*b6972d74SZach Atkins PetscErrorCode DMSwarmCreateProjectionRHS(DM dm_swarm, const char *field, Vec B_mesh) {
431*b6972d74SZach Atkins   PetscMemType       B_mem_type;
432*b6972d74SZach Atkins   DM                 dm_mesh;
433*b6972d74SZach Atkins   Vec                B_mesh_loc;
434*b6972d74SZach Atkins   DMSwarmCeedContext swarm_ceed_context;
435*b6972d74SZach Atkins 
436*b6972d74SZach Atkins   PetscFunctionBeginUser;
437*b6972d74SZach Atkins   // Get mesh DM
438*b6972d74SZach Atkins   PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
439*b6972d74SZach Atkins   PetscCall(DMGetApplicationContext(dm_mesh, (void *)&swarm_ceed_context));
440*b6972d74SZach Atkins 
441*b6972d74SZach Atkins   // Get mesh values
442*b6972d74SZach Atkins   PetscCall(DMGetLocalVector(dm_mesh, &B_mesh_loc));
443*b6972d74SZach Atkins   PetscCall(VecZeroEntries(B_mesh_loc));
444*b6972d74SZach Atkins   PetscCall(VecP2C(B_mesh_loc, &B_mem_type, swarm_ceed_context->v_mesh));
445*b6972d74SZach Atkins 
446*b6972d74SZach Atkins   // Get swarm access
447*b6972d74SZach Atkins   PetscCall(DMSwarmSortGetAccess(dm_swarm));
448*b6972d74SZach Atkins   PetscCall(DMSwarmPICFieldP2C(dm_swarm, field, swarm_ceed_context->u_points));
449*b6972d74SZach Atkins 
450*b6972d74SZach Atkins   // Interpolate field from swarm points to mesh
451*b6972d74SZach Atkins   CeedOperatorApply(swarm_ceed_context->op_points_to_mesh, swarm_ceed_context->u_points, swarm_ceed_context->v_mesh, CEED_REQUEST_IMMEDIATE);
452*b6972d74SZach Atkins 
453*b6972d74SZach Atkins   // Restore PETSc Vecs and Local to Global
454*b6972d74SZach Atkins   PetscCall(VecC2P(swarm_ceed_context->v_mesh, B_mem_type, B_mesh_loc));
455*b6972d74SZach Atkins   PetscCall(VecZeroEntries(B_mesh));
456*b6972d74SZach Atkins   PetscCall(DMLocalToGlobal(dm_mesh, B_mesh_loc, ADD_VALUES, B_mesh));
457*b6972d74SZach Atkins 
458*b6972d74SZach Atkins   // Cleanup
459*b6972d74SZach Atkins   PetscCall(DMSwarmPICFieldC2P(dm_swarm, field, swarm_ceed_context->u_points));
460*b6972d74SZach Atkins   PetscCall(DMSwarmSortRestoreAccess(dm_swarm));
461*b6972d74SZach Atkins   PetscCall(DMRestoreLocalVector(dm_mesh, &B_mesh_loc));
462*b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
463*b6972d74SZach Atkins }
464*b6972d74SZach Atkins 
465*b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
466*b6972d74SZach Atkins // Swarm "mass matrix"
467*b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
468*b6972d74SZach Atkins PetscErrorCode MatMult_SwarmMass(Mat A, Vec U_mesh, Vec V_mesh) {
469*b6972d74SZach Atkins   PetscMemType       U_mem_type, V_mem_type;
470*b6972d74SZach Atkins   DM                 dm_mesh;
471*b6972d74SZach Atkins   Vec                U_mesh_loc, V_mesh_loc;
472*b6972d74SZach Atkins   DMSwarmCeedContext swarm_ceed_context;
473*b6972d74SZach Atkins 
474*b6972d74SZach Atkins   PetscFunctionBeginUser;
475*b6972d74SZach Atkins   // Get mesh DM
476*b6972d74SZach Atkins   PetscCall(MatGetDM(A, &dm_mesh));
477*b6972d74SZach Atkins   PetscCall(DMGetApplicationContext(dm_mesh, (void *)&swarm_ceed_context));
478*b6972d74SZach Atkins 
479*b6972d74SZach Atkins   // Global to Local and get PETSc Vec access
480*b6972d74SZach Atkins   PetscCall(DMGetLocalVector(dm_mesh, &U_mesh_loc));
481*b6972d74SZach Atkins   PetscCall(VecZeroEntries(U_mesh_loc));
482*b6972d74SZach Atkins   PetscCall(DMGlobalToLocal(dm_mesh, U_mesh, INSERT_VALUES, U_mesh_loc));
483*b6972d74SZach Atkins   PetscCall(VecReadP2C(U_mesh_loc, &U_mem_type, swarm_ceed_context->u_mesh));
484*b6972d74SZach Atkins   PetscCall(DMGetLocalVector(dm_mesh, &V_mesh_loc));
485*b6972d74SZach Atkins   PetscCall(VecZeroEntries(V_mesh_loc));
486*b6972d74SZach Atkins   PetscCall(VecP2C(V_mesh_loc, &V_mem_type, swarm_ceed_context->v_mesh));
487*b6972d74SZach Atkins 
488*b6972d74SZach Atkins   // Apply swarm mass operator
489*b6972d74SZach Atkins   CeedOperatorApply(swarm_ceed_context->op_mass, swarm_ceed_context->u_mesh, swarm_ceed_context->v_mesh, CEED_REQUEST_IMMEDIATE);
490*b6972d74SZach Atkins 
491*b6972d74SZach Atkins   // Restore PETSc Vecs and Local to Global
492*b6972d74SZach Atkins   PetscCall(VecReadC2P(swarm_ceed_context->u_mesh, U_mem_type, U_mesh_loc));
493*b6972d74SZach Atkins   PetscCall(VecC2P(swarm_ceed_context->v_mesh, V_mem_type, V_mesh_loc));
494*b6972d74SZach Atkins   PetscCall(VecZeroEntries(V_mesh));
495*b6972d74SZach Atkins   PetscCall(DMLocalToGlobal(dm_mesh, V_mesh_loc, ADD_VALUES, V_mesh));
496*b6972d74SZach Atkins 
497*b6972d74SZach Atkins   // Cleanup
498*b6972d74SZach Atkins   PetscCall(DMRestoreLocalVector(dm_mesh, &U_mesh_loc));
499*b6972d74SZach Atkins   PetscCall(DMRestoreLocalVector(dm_mesh, &V_mesh_loc));
500*b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
501*b6972d74SZach Atkins }
502*b6972d74SZach Atkins 
503*b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
504*b6972d74SZach Atkins // Swarm to mesh projection
505*b6972d74SZach Atkins // ------------------------------------------------------------------------------------------------
506*b6972d74SZach Atkins PetscErrorCode DMSwarmProjectFromSwarmToCells(DM dm_swarm, const char *field, Vec U_mesh) {
507*b6972d74SZach Atkins   PetscBool          test_mode;
508*b6972d74SZach Atkins   Vec                B_mesh;
509*b6972d74SZach Atkins   Mat                M;
510*b6972d74SZach Atkins   KSP                ksp;
511*b6972d74SZach Atkins   DM                 dm_mesh;
512*b6972d74SZach Atkins   DMSwarmCeedContext swarm_ceed_context;
513*b6972d74SZach Atkins   MPI_Comm           comm;
514*b6972d74SZach Atkins 
515*b6972d74SZach Atkins   PetscFunctionBeginUser;
516*b6972d74SZach Atkins 
517*b6972d74SZach Atkins   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swarm-to-Mesh Projection Options", NULL);
518*b6972d74SZach Atkins   PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL));
519*b6972d74SZach Atkins   PetscOptionsEnd();
520*b6972d74SZach Atkins 
521*b6972d74SZach Atkins   comm = PetscObjectComm((PetscObject)dm_swarm);
522*b6972d74SZach Atkins   PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
523*b6972d74SZach Atkins   PetscCall(DMGetApplicationContext(dm_mesh, (void *)&swarm_ceed_context));
524*b6972d74SZach Atkins   PetscCall(VecDuplicate(U_mesh, &B_mesh));
525*b6972d74SZach Atkins 
526*b6972d74SZach Atkins   // Setup "mass matrix"
527*b6972d74SZach Atkins   {
528*b6972d74SZach Atkins     PetscInt l_size, g_size;
529*b6972d74SZach Atkins 
530*b6972d74SZach Atkins     PetscCall(VecGetLocalSize(U_mesh, &l_size));
531*b6972d74SZach Atkins     PetscCall(VecGetSize(U_mesh, &g_size));
532*b6972d74SZach Atkins     PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, swarm_ceed_context, &M));
533*b6972d74SZach Atkins     PetscCall(MatSetDM(M, dm_mesh));
534*b6972d74SZach Atkins     PetscCall(MatShellSetOperation(M, MATOP_MULT, (void (*)(void))MatMult_SwarmMass));
535*b6972d74SZach Atkins   }
536*b6972d74SZach Atkins 
537*b6972d74SZach Atkins   // Setup KSP
538*b6972d74SZach Atkins   {
539*b6972d74SZach Atkins     PC pc;
540*b6972d74SZach Atkins 
541*b6972d74SZach Atkins     PetscCall(KSPCreate(comm, &ksp));
542*b6972d74SZach Atkins     PetscCall(KSPGetPC(ksp, &pc));
543*b6972d74SZach Atkins     PetscCall(PCSetType(pc, PCJACOBI));
544*b6972d74SZach Atkins     PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
545*b6972d74SZach Atkins     PetscCall(KSPSetType(ksp, KSPCG));
546*b6972d74SZach Atkins     PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
547*b6972d74SZach Atkins     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
548*b6972d74SZach Atkins     PetscCall(KSPSetOperators(ksp, M, M));
549*b6972d74SZach Atkins     PetscCall(KSPSetFromOptions(ksp));
550*b6972d74SZach Atkins     PetscCall(PetscObjectSetName((PetscObject)ksp, "Swarm-to-Mesh Projection"));
551*b6972d74SZach Atkins     PetscCall(KSPViewFromOptions(ksp, NULL, "-ksp_projection_view"));
552*b6972d74SZach Atkins   }
553*b6972d74SZach Atkins 
554*b6972d74SZach Atkins   // Setup RHS
555*b6972d74SZach Atkins   PetscCall(DMSwarmCreateProjectionRHS(dm_swarm, field, B_mesh));
556*b6972d74SZach Atkins 
557*b6972d74SZach Atkins   // Solve
558*b6972d74SZach Atkins   PetscCall(VecZeroEntries(U_mesh));
559*b6972d74SZach Atkins   PetscCall(KSPSolve(ksp, B_mesh, U_mesh));
560*b6972d74SZach Atkins 
561*b6972d74SZach Atkins   // KSP summary
562*b6972d74SZach Atkins   {
563*b6972d74SZach Atkins     KSPType            ksp_type;
564*b6972d74SZach Atkins     KSPConvergedReason reason;
565*b6972d74SZach Atkins     PetscReal          rnorm;
566*b6972d74SZach Atkins     PetscInt           its;
567*b6972d74SZach Atkins     PetscCall(KSPGetType(ksp, &ksp_type));
568*b6972d74SZach Atkins     PetscCall(KSPGetConvergedReason(ksp, &reason));
569*b6972d74SZach Atkins     PetscCall(KSPGetIterationNumber(ksp, &its));
570*b6972d74SZach Atkins     PetscCall(KSPGetResidualNorm(ksp, &rnorm));
571*b6972d74SZach Atkins 
572*b6972d74SZach Atkins     if (!test_mode || reason < 0 || rnorm > 1e-8) {
573*b6972d74SZach Atkins       PetscCall(PetscPrintf(comm,
574*b6972d74SZach Atkins                             "Swarm-to-Mesh Projection KSP Solve:\n"
575*b6972d74SZach Atkins                             "  KSP type: %s\n"
576*b6972d74SZach Atkins                             "  KSP convergence: %s\n"
577*b6972d74SZach Atkins                             "  Total KSP iterations: %" PetscInt_FMT "\n"
578*b6972d74SZach Atkins                             "  Final rnorm: %e\n",
579*b6972d74SZach Atkins                             ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
580*b6972d74SZach Atkins     }
581*b6972d74SZach Atkins   }
582*b6972d74SZach Atkins 
583*b6972d74SZach Atkins   // Optional viewing
584*b6972d74SZach Atkins   PetscCall(KSPViewFromOptions(ksp, NULL, "-ksp_view"));
585*b6972d74SZach Atkins 
586*b6972d74SZach Atkins   // Cleanup
587*b6972d74SZach Atkins   PetscCall(VecDestroy(&B_mesh));
588*b6972d74SZach Atkins   PetscCall(MatDestroy(&M));
589*b6972d74SZach Atkins   PetscCall(KSPDestroy(&ksp));
590*b6972d74SZach Atkins 
591*b6972d74SZach Atkins   PetscFunctionReturn(PETSC_SUCCESS);
592*b6972d74SZach Atkins }
593