xref: /libCEED/examples/petsc/dmswarm.c (revision 9dbcead6862c6d965f40afcd1255316d96af311a)
1*9dbcead6SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2*9dbcead6SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3*9dbcead6SJeremy L Thompson //
4*9dbcead6SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5*9dbcead6SJeremy L Thompson //
6*9dbcead6SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7*9dbcead6SJeremy L Thompson 
8*9dbcead6SJeremy L Thompson //                        libCEED + PETSc DMSwarm Example
9*9dbcead6SJeremy L Thompson //
10*9dbcead6SJeremy L Thompson // This example demonstrates a simple usage of libCEED with DMSwarm.
11*9dbcead6SJeremy L Thompson // This example combines elements of PETSc src/impls/dm/swam/tutorials/ex1.c and src/impls/dm/swarm/tests/ex6.c
12*9dbcead6SJeremy L Thompson //
13*9dbcead6SJeremy L Thompson // Build with:
14*9dbcead6SJeremy L Thompson //
15*9dbcead6SJeremy L Thompson //     make dmswarm [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
16*9dbcead6SJeremy L Thompson //
17*9dbcead6SJeremy L Thompson // Sample runs:
18*9dbcead6SJeremy L Thompson //
19*9dbcead6SJeremy L Thompson //  ./dmswarm -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_plex_box_lower -1.0,-1.0,-1.0 -dm_plex_simplex 0 -num_comp 2 -swarm gauss
20*9dbcead6SJeremy L Thompson //
21*9dbcead6SJeremy L Thompson //TESTARGS(name="Uniform swarm, CG projection") -ceed {ceed_resource} -test -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_plex_box_lower -1.0,-1.0,-1.0 -dm_plex_simplex 0 -dm_plex_hash_location true -num_comp 2 -swarm uniform -solution_order 3 -points_per_cell 125
22*9dbcead6SJeremy L Thompson //TESTARGS(name="Gauss swarm, lumped projection") -ceed {ceed_resource} -test -dm_plex_dim 3 -dm_plex_box_faces 3,3,3 -dm_plex_box_lower -1.0,-1.0,-1.0 -dm_plex_simplex 0 -dm_plex_hash_location true -num_comp 2 -swarm gauss -ksp_type preonly -pc_type jacobi -pc_jacobi_type rowsum -tolerance 9e-2
23*9dbcead6SJeremy L Thompson 
24*9dbcead6SJeremy L Thompson /// @file
25*9dbcead6SJeremy L Thompson /// libCEED example using PETSc with DMSwarm
26*9dbcead6SJeremy L Thompson const char help[] = "libCEED example using PETSc with DMSwarm\n";
27*9dbcead6SJeremy L Thompson 
28*9dbcead6SJeremy L Thompson #include <ceed.h>
29*9dbcead6SJeremy L Thompson #include <math.h>
30*9dbcead6SJeremy L Thompson #include <petscdmplex.h>
31*9dbcead6SJeremy L Thompson #include <petscdmswarm.h>
32*9dbcead6SJeremy L Thompson #include <petscds.h>
33*9dbcead6SJeremy L Thompson #include <petscfe.h>
34*9dbcead6SJeremy L Thompson #include <petscksp.h>
35*9dbcead6SJeremy L Thompson #include <petsc/private/petscfeimpl.h> /* For interpolation */
36*9dbcead6SJeremy L Thompson 
37*9dbcead6SJeremy L Thompson #include "include/petscutils.h"
38*9dbcead6SJeremy L Thompson 
39*9dbcead6SJeremy L Thompson const char DMSwarmPICField_u[] = "u";
40*9dbcead6SJeremy L Thompson 
41*9dbcead6SJeremy L Thompson // libCEED context data
42*9dbcead6SJeremy L Thompson typedef struct DMSwarmCeedContext_ *DMSwarmCeedContext;
43*9dbcead6SJeremy L Thompson struct DMSwarmCeedContext_ {
44*9dbcead6SJeremy L Thompson   Ceed                ceed;
45*9dbcead6SJeremy L Thompson   CeedVector          u_mesh_loc, v_mesh_loc, u_mesh_elem, u_points_loc, u_points_elem, x_ref_points_loc, x_ref_points_elem;
46*9dbcead6SJeremy L Thompson   CeedElemRestriction restriction_u_mesh, restriction_x_points, restriction_u_points;
47*9dbcead6SJeremy L Thompson   CeedBasis           basis_u;
48*9dbcead6SJeremy L Thompson };
49*9dbcead6SJeremy L Thompson 
50*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmCeedContextCreate(DM dm_swarm, const char *ceed_resource, DMSwarmCeedContext *ctx);
51*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmCeedContextDestroy(DMSwarmCeedContext *ctx);
52*9dbcead6SJeremy L Thompson 
53*9dbcead6SJeremy L Thompson // Target functions
54*9dbcead6SJeremy L Thompson typedef enum { TARGET_TANH = 0, TARGET_POLYNOMIAL = 1, TARGET_SPHERE = 2 } TargetType;
55*9dbcead6SJeremy L Thompson static const char *const target_types[] = {"tanh", "polynomial", "sphere", "TargetType", "TARGET", 0};
56*9dbcead6SJeremy L Thompson 
57*9dbcead6SJeremy L Thompson typedef PetscErrorCode (*TargetFunc)(PetscInt dim, const PetscScalar x[]);
58*9dbcead6SJeremy L Thompson typedef PetscErrorCode (*TargetFuncProj)(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx);
59*9dbcead6SJeremy L Thompson 
60*9dbcead6SJeremy L Thompson PetscScalar    EvalU_Tanh(PetscInt dim, const PetscScalar x[]);
61*9dbcead6SJeremy L Thompson PetscScalar    EvalU_Poly(PetscInt dim, const PetscScalar x[]);
62*9dbcead6SJeremy L Thompson PetscScalar    EvalU_Sphere(PetscInt dim, const PetscScalar x[]);
63*9dbcead6SJeremy L Thompson PetscErrorCode EvalU_Tanh_proj(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt num_comp, PetscScalar *u, void *ctx);
64*9dbcead6SJeremy L Thompson PetscErrorCode EvalU_Poly_proj(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt num_comp, PetscScalar *u, void *ctx);
65*9dbcead6SJeremy L Thompson PetscErrorCode EvalU_Sphere_proj(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt num_comp, PetscScalar *u, void *ctx);
66*9dbcead6SJeremy L Thompson 
67*9dbcead6SJeremy L Thompson // Swarm point distribution
68*9dbcead6SJeremy L Thompson typedef enum { SWARM_GAUSS = 0, SWARM_UNIFORM = 1, SWARM_CELL_RANDOM = 2, SWARM_SINUSOIDAL = 3 } PointSwarmType;
69*9dbcead6SJeremy L Thompson static const char *const point_swarm_types[] = {"gauss", "uniform", "cell_random", "sinusoidal", "PointSwarmType", "SWARM", 0};
70*9dbcead6SJeremy L Thompson 
71*9dbcead6SJeremy L Thompson // Swarm helper function
72*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmInitalizePointLocations(DM dm_swarm, PointSwarmType point_swarm_type, PetscInt num_points, PetscInt num_points_per_cell);
73*9dbcead6SJeremy L Thompson 
74*9dbcead6SJeremy L Thompson // Memory utilities
75*9dbcead6SJeremy L Thompson PetscErrorCode VecP2C(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed);
76*9dbcead6SJeremy L Thompson PetscErrorCode VecC2P(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc);
77*9dbcead6SJeremy L Thompson PetscErrorCode VecReadP2C(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed);
78*9dbcead6SJeremy L Thompson PetscErrorCode VecReadC2P(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc);
79*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmPICFieldP2C(DM dm_swarm, const char *field, CeedVector x_ceed);
80*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmPICFieldC2P(DM dm_swarm, const char *field, CeedVector x_ceed);
81*9dbcead6SJeremy L Thompson 
82*9dbcead6SJeremy L Thompson // Swarm to mesh and mesh to swarm
83*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmCreateReferenceCoordinates(DM dm_swarm, IS *is_points, Vec *ref_coords);
84*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmInterpolateFromCellToSwarm_Petsc(DM dm_swarm, const char *field, Vec U_mesh);
85*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmInterpolateFromCellToSwarm_Ceed(DM dm_swarm, const char *field, Vec U_mesh);
86*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmCheckSwarmValues(DM dm_swarm, const char *field, PetscScalar tolerance, TargetFuncProj TrueSolution);
87*9dbcead6SJeremy L Thompson 
88*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmCreateProjectionRHS(DM dm_swarm, const char *field, Vec B_mesh);
89*9dbcead6SJeremy L Thompson PetscErrorCode MatMult_SwarmMass(Mat A, Vec U_mesh, Vec V_mesh);
90*9dbcead6SJeremy L Thompson 
91*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
92*9dbcead6SJeremy L Thompson // main driver
93*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
94*9dbcead6SJeremy L Thompson int main(int argc, char **argv) {
95*9dbcead6SJeremy L Thompson   MPI_Comm           comm;
96*9dbcead6SJeremy L Thompson   char               ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self";
97*9dbcead6SJeremy L Thompson   PetscBool          test_mode = PETSC_FALSE, view_petsc_swarm = PETSC_FALSE, view_ceed_swarm = PETSC_FALSE;
98*9dbcead6SJeremy L Thompson   PetscInt           dim = 3, num_comp = 1, num_points = 1728, num_points_per_cell = 64, mesh_order = 1, solution_order = 2, q_extra = 3;
99*9dbcead6SJeremy L Thompson   PetscScalar        tolerance = 1E-3;
100*9dbcead6SJeremy L Thompson   DM                 dm_mesh, dm_swarm;
101*9dbcead6SJeremy L Thompson   Vec                U_mesh;
102*9dbcead6SJeremy L Thompson   DMSwarmCeedContext swarm_ceed_context;
103*9dbcead6SJeremy L Thompson   PointSwarmType     point_swarm_type     = SWARM_UNIFORM;
104*9dbcead6SJeremy L Thompson   TargetType         target_type          = TARGET_TANH;
105*9dbcead6SJeremy L Thompson   TargetFuncProj     target_function_proj = EvalU_Tanh_proj;
106*9dbcead6SJeremy L Thompson 
107*9dbcead6SJeremy L Thompson   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
108*9dbcead6SJeremy L Thompson   comm = PETSC_COMM_WORLD;
109*9dbcead6SJeremy L Thompson 
110*9dbcead6SJeremy L Thompson   // Read command line options
111*9dbcead6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "libCEED example using PETSc with DMSwarm", NULL);
112*9dbcead6SJeremy L Thompson 
113*9dbcead6SJeremy L Thompson   PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL));
114*9dbcead6SJeremy L Thompson   PetscCall(
115*9dbcead6SJeremy L Thompson       PetscOptionsBool("-u_petsc_swarm_view", "View XDMF of swarm values interpolated by PETSc", NULL, view_petsc_swarm, &view_petsc_swarm, NULL));
116*9dbcead6SJeremy L Thompson   PetscCall(
117*9dbcead6SJeremy L Thompson       PetscOptionsBool("-u_ceed_swarm_view", "View XDMF of swarm values interpolated by libCEED", NULL, view_ceed_swarm, &view_ceed_swarm, NULL));
118*9dbcead6SJeremy L Thompson   PetscCall(PetscOptionsEnum("-target", "Target field function", NULL, target_types, (PetscEnum)target_type, (PetscEnum *)&target_type, NULL));
119*9dbcead6SJeremy L Thompson   PetscCall(PetscOptionsInt("-solution_order", "Order of mesh solution space", NULL, solution_order, &solution_order, NULL));
120*9dbcead6SJeremy L Thompson   PetscCall(PetscOptionsInt("-mesh_order", "Order of mesh coordinate space", NULL, mesh_order, &mesh_order, NULL));
121*9dbcead6SJeremy L Thompson   PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra, NULL));
122*9dbcead6SJeremy L Thompson   PetscCall(PetscOptionsInt("-num_comp", "Number of components in solution", NULL, num_comp, &num_comp, NULL));
123*9dbcead6SJeremy L Thompson   PetscCall(PetscOptionsEnum("-swarm", "Swarm points distribution", NULL, point_swarm_types, (PetscEnum)point_swarm_type,
124*9dbcead6SJeremy L Thompson                              (PetscEnum *)&point_swarm_type, NULL));
125*9dbcead6SJeremy L Thompson   {
126*9dbcead6SJeremy L Thompson     PetscBool user_set_num_points_per_cell = PETSC_FALSE;
127*9dbcead6SJeremy L Thompson     PetscInt  dim = 3, num_cells_total = 1;
128*9dbcead6SJeremy L Thompson     PetscInt  num_cells[] = {1, 1, 1};
129*9dbcead6SJeremy L Thompson 
130*9dbcead6SJeremy L Thompson     PetscCall(PetscOptionsInt("-points_per_cell", "Total number of swarm points in each cell", NULL, num_points_per_cell, &num_points_per_cell,
131*9dbcead6SJeremy L Thompson                               &user_set_num_points_per_cell));
132*9dbcead6SJeremy L Thompson     PetscCall(PetscOptionsInt("-dm_plex_dim", "Background mesh dimension", NULL, dim, &dim, NULL));
133*9dbcead6SJeremy L Thompson     PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of cells", NULL, num_cells, &dim, NULL));
134*9dbcead6SJeremy L Thompson 
135*9dbcead6SJeremy L Thompson     num_cells_total = num_cells[0] * num_cells[1] * num_cells[2];
136*9dbcead6SJeremy L Thompson     PetscCheck(!user_set_num_points_per_cell || point_swarm_type != SWARM_SINUSOIDAL, comm, PETSC_ERR_USER,
137*9dbcead6SJeremy L Thompson                "Cannot specify points per cell with sinusoidal points locations");
138*9dbcead6SJeremy L Thompson     if (!user_set_num_points_per_cell) {
139*9dbcead6SJeremy L Thompson       PetscCall(PetscOptionsInt("-points", "Total number of swarm points", NULL, num_points, &num_points, NULL));
140*9dbcead6SJeremy L Thompson       num_points_per_cell = PetscCeilInt(num_points, num_cells_total);
141*9dbcead6SJeremy L Thompson     }
142*9dbcead6SJeremy L Thompson     if (point_swarm_type != SWARM_SINUSOIDAL) {
143*9dbcead6SJeremy L Thompson       PetscInt num_points_per_cell_1d = round(cbrt(num_points_per_cell * 1.0));
144*9dbcead6SJeremy L Thompson 
145*9dbcead6SJeremy L Thompson       num_points_per_cell = 1;
146*9dbcead6SJeremy L Thompson       for (PetscInt i = 0; i < dim; i++) num_points_per_cell *= num_points_per_cell_1d;
147*9dbcead6SJeremy L Thompson     }
148*9dbcead6SJeremy L Thompson     num_points = num_points_per_cell * num_cells_total;
149*9dbcead6SJeremy L Thompson   }
150*9dbcead6SJeremy L Thompson   PetscCall(PetscOptionsScalar("-tolerance", "Tolerance for swarm point values and projection relative L2 error", NULL, tolerance, &tolerance, NULL));
151*9dbcead6SJeremy L Thompson   PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resource, sizeof(ceed_resource), NULL));
152*9dbcead6SJeremy L Thompson 
153*9dbcead6SJeremy L Thompson   PetscOptionsEnd();
154*9dbcead6SJeremy L Thompson 
155*9dbcead6SJeremy L Thompson   // Create background mesh
156*9dbcead6SJeremy L Thompson   {
157*9dbcead6SJeremy L Thompson     PetscCall(DMCreate(comm, &dm_mesh));
158*9dbcead6SJeremy L Thompson     PetscCall(DMSetType(dm_mesh, DMPLEX));
159*9dbcead6SJeremy L Thompson     PetscCall(DMSetFromOptions(dm_mesh));
160*9dbcead6SJeremy L Thompson 
161*9dbcead6SJeremy L Thompson     // -- Check for tensor product mesh
162*9dbcead6SJeremy L Thompson     {
163*9dbcead6SJeremy L Thompson       PetscBool is_simplex;
164*9dbcead6SJeremy L Thompson 
165*9dbcead6SJeremy L Thompson       PetscCall(DMPlexIsSimplex(dm_mesh, &is_simplex));
166*9dbcead6SJeremy L Thompson       PetscCheck(!is_simplex, comm, PETSC_ERR_USER, "Only tensor-product background meshes supported");
167*9dbcead6SJeremy L Thompson     }
168*9dbcead6SJeremy L Thompson 
169*9dbcead6SJeremy L Thompson     // -- Mesh FE space
170*9dbcead6SJeremy L Thompson     PetscCall(DMGetDimension(dm_mesh, &dim));
171*9dbcead6SJeremy L Thompson     {
172*9dbcead6SJeremy L Thompson       PetscFE fe;
173*9dbcead6SJeremy L Thompson 
174*9dbcead6SJeremy L Thompson       PetscCall(DMGetDimension(dm_mesh, &dim));
175*9dbcead6SJeremy L Thompson       PetscCall(PetscFECreateLagrange(comm, dim, num_comp, PETSC_FALSE, solution_order, solution_order + q_extra, &fe));
176*9dbcead6SJeremy L Thompson       PetscCall(DMAddField(dm_mesh, NULL, (PetscObject)fe));
177*9dbcead6SJeremy L Thompson       PetscCall(PetscFEDestroy(&fe));
178*9dbcead6SJeremy L Thompson     }
179*9dbcead6SJeremy L Thompson     PetscCall(DMCreateDS(dm_mesh));
180*9dbcead6SJeremy L Thompson 
181*9dbcead6SJeremy L Thompson     // -- Coordinate FE space
182*9dbcead6SJeremy L Thompson     {
183*9dbcead6SJeremy L Thompson       PetscFE fe_coord;
184*9dbcead6SJeremy L Thompson 
185*9dbcead6SJeremy L Thompson       PetscCall(PetscFECreateLagrange(comm, dim, dim, PETSC_FALSE, mesh_order, solution_order + q_extra, &fe_coord));
186*9dbcead6SJeremy L Thompson       PetscCall(DMProjectCoordinates(dm_mesh, fe_coord));
187*9dbcead6SJeremy L Thompson       PetscCall(PetscFEDestroy(&fe_coord));
188*9dbcead6SJeremy L Thompson     }
189*9dbcead6SJeremy L Thompson 
190*9dbcead6SJeremy L Thompson     // -- Set tensor permutation
191*9dbcead6SJeremy L Thompson     {
192*9dbcead6SJeremy L Thompson       DM dm_coord;
193*9dbcead6SJeremy L Thompson 
194*9dbcead6SJeremy L Thompson       PetscCall(DMGetCoordinateDM(dm_mesh, &dm_coord));
195*9dbcead6SJeremy L Thompson       PetscCall(DMPlexSetClosurePermutationTensor(dm_mesh, PETSC_DETERMINE, NULL));
196*9dbcead6SJeremy L Thompson       PetscCall(DMPlexSetClosurePermutationTensor(dm_coord, PETSC_DETERMINE, NULL));
197*9dbcead6SJeremy L Thompson     }
198*9dbcead6SJeremy L Thompson 
199*9dbcead6SJeremy L Thompson     // -- Final background mesh
200*9dbcead6SJeremy L Thompson     PetscCall(PetscObjectSetName((PetscObject)dm_mesh, "Background Mesh"));
201*9dbcead6SJeremy L Thompson     PetscCall(DMViewFromOptions(dm_mesh, NULL, "-dm_mesh_view"));
202*9dbcead6SJeremy L Thompson   }
203*9dbcead6SJeremy L Thompson 
204*9dbcead6SJeremy L Thompson   // Create particle swarm
205*9dbcead6SJeremy L Thompson   {
206*9dbcead6SJeremy L Thompson     PetscCall(DMCreate(comm, &dm_swarm));
207*9dbcead6SJeremy L Thompson     PetscCall(DMSetType(dm_swarm, DMSWARM));
208*9dbcead6SJeremy L Thompson     PetscCall(DMSetDimension(dm_swarm, dim));
209*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmSetType(dm_swarm, DMSWARM_PIC));
210*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmSetCellDM(dm_swarm, dm_mesh));
211*9dbcead6SJeremy L Thompson 
212*9dbcead6SJeremy L Thompson     // -- Swarm field
213*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmRegisterPetscDatatypeField(dm_swarm, DMSwarmPICField_u, num_comp, PETSC_SCALAR));
214*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmFinalizeFieldRegister(dm_swarm));
215*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmSetLocalSizes(dm_swarm, num_points, 0));
216*9dbcead6SJeremy L Thompson     PetscCall(DMSetFromOptions(dm_swarm));
217*9dbcead6SJeremy L Thompson 
218*9dbcead6SJeremy L Thompson     // -- Set swarm point locations
219*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmInitalizePointLocations(dm_swarm, point_swarm_type, num_points, num_points_per_cell));
220*9dbcead6SJeremy L Thompson 
221*9dbcead6SJeremy L Thompson     // -- Final particle swarm
222*9dbcead6SJeremy L Thompson     PetscCall(PetscObjectSetName((PetscObject)dm_swarm, "Particle Swarm"));
223*9dbcead6SJeremy L Thompson     PetscCall(DMViewFromOptions(dm_swarm, NULL, "-dm_swarm_view"));
224*9dbcead6SJeremy L Thompson   }
225*9dbcead6SJeremy L Thompson 
226*9dbcead6SJeremy L Thompson   // Set field values on background mesh
227*9dbcead6SJeremy L Thompson   PetscCall(DMCreateGlobalVector(dm_mesh, &U_mesh));
228*9dbcead6SJeremy L Thompson   switch (target_type) {
229*9dbcead6SJeremy L Thompson     case TARGET_TANH:
230*9dbcead6SJeremy L Thompson       target_function_proj = EvalU_Tanh_proj;
231*9dbcead6SJeremy L Thompson       break;
232*9dbcead6SJeremy L Thompson     case TARGET_POLYNOMIAL:
233*9dbcead6SJeremy L Thompson       target_function_proj = EvalU_Poly_proj;
234*9dbcead6SJeremy L Thompson       break;
235*9dbcead6SJeremy L Thompson     case TARGET_SPHERE:
236*9dbcead6SJeremy L Thompson       target_function_proj = EvalU_Sphere_proj;
237*9dbcead6SJeremy L Thompson       break;
238*9dbcead6SJeremy L Thompson   }
239*9dbcead6SJeremy L Thompson   {
240*9dbcead6SJeremy L Thompson     TargetFuncProj mesh_solution[1] = {target_function_proj};
241*9dbcead6SJeremy L Thompson 
242*9dbcead6SJeremy L Thompson     PetscCall(DMProjectFunction(dm_mesh, 0.0, mesh_solution, NULL, INSERT_VALUES, U_mesh));
243*9dbcead6SJeremy L Thompson   }
244*9dbcead6SJeremy L Thompson 
245*9dbcead6SJeremy L Thompson   // Visualize background mesh
246*9dbcead6SJeremy L Thompson   PetscCall(PetscObjectSetName((PetscObject)U_mesh, "U on Background Mesh"));
247*9dbcead6SJeremy L Thompson   PetscCall(VecViewFromOptions(U_mesh, NULL, "-u_mesh_view"));
248*9dbcead6SJeremy L Thompson 
249*9dbcead6SJeremy L Thompson   // libCEED objects for swarm and background mesh
250*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmCeedContextCreate(dm_swarm, ceed_resource, &swarm_ceed_context));
251*9dbcead6SJeremy L Thompson 
252*9dbcead6SJeremy L Thompson   // Interpolate from mesh to points via PETSc
253*9dbcead6SJeremy L Thompson   {
254*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmInterpolateFromCellToSwarm_Petsc(dm_swarm, DMSwarmPICField_u, U_mesh));
255*9dbcead6SJeremy L Thompson     if (view_petsc_swarm) PetscCall(DMSwarmViewXDMF(dm_swarm, "swarm_petsc.xmf"));
256*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmCheckSwarmValues(dm_swarm, DMSwarmPICField_u, tolerance, target_function_proj));
257*9dbcead6SJeremy L Thompson   }
258*9dbcead6SJeremy L Thompson 
259*9dbcead6SJeremy L Thompson   // Interpolate from mesh to points via libCEED
260*9dbcead6SJeremy L Thompson   {
261*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmInterpolateFromCellToSwarm_Ceed(dm_swarm, DMSwarmPICField_u, U_mesh));
262*9dbcead6SJeremy L Thompson     if (view_ceed_swarm) PetscCall(DMSwarmViewXDMF(dm_swarm, "swarm_ceed.xmf"));
263*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmCheckSwarmValues(dm_swarm, DMSwarmPICField_u, tolerance, target_function_proj));
264*9dbcead6SJeremy L Thompson   }
265*9dbcead6SJeremy L Thompson 
266*9dbcead6SJeremy L Thompson   // Project from points to mesh via libCEED
267*9dbcead6SJeremy L Thompson   {
268*9dbcead6SJeremy L Thompson     Vec B_mesh, U_projected;
269*9dbcead6SJeremy L Thompson     Mat M;
270*9dbcead6SJeremy L Thompson     KSP ksp;
271*9dbcead6SJeremy L Thompson 
272*9dbcead6SJeremy L Thompson     PetscCall(VecDuplicate(U_mesh, &B_mesh));
273*9dbcead6SJeremy L Thompson     PetscCall(VecDuplicate(U_mesh, &U_projected));
274*9dbcead6SJeremy L Thompson 
275*9dbcead6SJeremy L Thompson     // -- Setup "mass matrix"
276*9dbcead6SJeremy L Thompson     {
277*9dbcead6SJeremy L Thompson       PetscInt l_size, g_size;
278*9dbcead6SJeremy L Thompson 
279*9dbcead6SJeremy L Thompson       PetscCall(VecGetLocalSize(U_mesh, &l_size));
280*9dbcead6SJeremy L Thompson       PetscCall(VecGetSize(U_mesh, &g_size));
281*9dbcead6SJeremy L Thompson       PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, swarm_ceed_context, &M));
282*9dbcead6SJeremy L Thompson       PetscCall(MatSetDM(M, dm_mesh));
283*9dbcead6SJeremy L Thompson       PetscCall(MatShellSetOperation(M, MATOP_MULT, (void (*)(void))MatMult_SwarmMass));
284*9dbcead6SJeremy L Thompson     }
285*9dbcead6SJeremy L Thompson 
286*9dbcead6SJeremy L Thompson     // -- Setup KSP
287*9dbcead6SJeremy L Thompson     {
288*9dbcead6SJeremy L Thompson       PC pc;
289*9dbcead6SJeremy L Thompson 
290*9dbcead6SJeremy L Thompson       PetscCall(KSPCreate(comm, &ksp));
291*9dbcead6SJeremy L Thompson       PetscCall(KSPGetPC(ksp, &pc));
292*9dbcead6SJeremy L Thompson       PetscCall(PCSetType(pc, PCJACOBI));
293*9dbcead6SJeremy L Thompson       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
294*9dbcead6SJeremy L Thompson       PetscCall(KSPSetType(ksp, KSPCG));
295*9dbcead6SJeremy L Thompson       PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
296*9dbcead6SJeremy L Thompson       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
297*9dbcead6SJeremy L Thompson       PetscCall(KSPSetOperators(ksp, M, M));
298*9dbcead6SJeremy L Thompson       PetscCall(KSPSetFromOptions(ksp));
299*9dbcead6SJeremy L Thompson       PetscCall(PetscObjectSetName((PetscObject)ksp, "Swarm-to-Mesh Projection"));
300*9dbcead6SJeremy L Thompson       PetscCall(KSPViewFromOptions(ksp, NULL, "-ksp_projection_view"));
301*9dbcead6SJeremy L Thompson     }
302*9dbcead6SJeremy L Thompson 
303*9dbcead6SJeremy L Thompson     // -- Setup RHS
304*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmCreateProjectionRHS(dm_swarm, DMSwarmPICField_u, B_mesh));
305*9dbcead6SJeremy L Thompson 
306*9dbcead6SJeremy L Thompson     // -- Solve
307*9dbcead6SJeremy L Thompson     PetscCall(VecZeroEntries(U_projected));
308*9dbcead6SJeremy L Thompson     PetscCall(KSPSolve(ksp, B_mesh, U_projected));
309*9dbcead6SJeremy L Thompson 
310*9dbcead6SJeremy L Thompson     // -- KSP summary
311*9dbcead6SJeremy L Thompson     {
312*9dbcead6SJeremy L Thompson       KSPType            ksp_type;
313*9dbcead6SJeremy L Thompson       KSPConvergedReason reason;
314*9dbcead6SJeremy L Thompson       PetscReal          rnorm;
315*9dbcead6SJeremy L Thompson       PetscInt           its;
316*9dbcead6SJeremy L Thompson       PetscCall(KSPGetType(ksp, &ksp_type));
317*9dbcead6SJeremy L Thompson       PetscCall(KSPGetConvergedReason(ksp, &reason));
318*9dbcead6SJeremy L Thompson       PetscCall(KSPGetIterationNumber(ksp, &its));
319*9dbcead6SJeremy L Thompson       PetscCall(KSPGetResidualNorm(ksp, &rnorm));
320*9dbcead6SJeremy L Thompson 
321*9dbcead6SJeremy L Thompson       if (!test_mode || reason < 0 || rnorm > 1e-8) {
322*9dbcead6SJeremy L Thompson         PetscCall(PetscPrintf(comm,
323*9dbcead6SJeremy L Thompson                               "Swarm-to-Mesh Projection KSP Solve:\n"
324*9dbcead6SJeremy L Thompson                               "  KSP type: %s\n"
325*9dbcead6SJeremy L Thompson                               "  KSP convergence: %s\n"
326*9dbcead6SJeremy L Thompson                               "  Total KSP iterations: %" PetscInt_FMT "\n"
327*9dbcead6SJeremy L Thompson                               "  Final rnorm: %e\n",
328*9dbcead6SJeremy L Thompson                               ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
329*9dbcead6SJeremy L Thompson       }
330*9dbcead6SJeremy L Thompson     }
331*9dbcead6SJeremy L Thompson 
332*9dbcead6SJeremy L Thompson     // -- Check error
333*9dbcead6SJeremy L Thompson     PetscCall(KSPViewFromOptions(ksp, NULL, "-ksp_view"));
334*9dbcead6SJeremy L Thompson     PetscCall(PetscObjectSetName((PetscObject)U_mesh, "U projected to Background Mesh"));
335*9dbcead6SJeremy L Thompson     PetscCall(VecViewFromOptions(U_projected, NULL, "-u_projected_view"));
336*9dbcead6SJeremy L Thompson     PetscCall(VecAXPY(U_projected, -1.0, U_mesh));
337*9dbcead6SJeremy L Thompson     PetscCall(PetscObjectSetName((PetscObject)U_projected, "U Projection Error"));
338*9dbcead6SJeremy L Thompson     PetscCall(VecViewFromOptions(U_projected, NULL, "-u_error_view"));
339*9dbcead6SJeremy L Thompson     {
340*9dbcead6SJeremy L Thompson       PetscScalar error, norm_u_mesh;
341*9dbcead6SJeremy L Thompson 
342*9dbcead6SJeremy L Thompson       PetscCall(VecNorm(U_projected, NORM_2, &error));
343*9dbcead6SJeremy L Thompson       PetscCall(VecNorm(U_mesh, NORM_2, &norm_u_mesh));
344*9dbcead6SJeremy L Thompson       PetscCheck(error / norm_u_mesh < tolerance, comm, PETSC_ERR_USER, "Projection error too high: %e\n", error / norm_u_mesh);
345*9dbcead6SJeremy L Thompson       if (!test_mode) PetscCall(PetscPrintf(comm, "  Projection error: %e\n", error / norm_u_mesh));
346*9dbcead6SJeremy L Thompson     }
347*9dbcead6SJeremy L Thompson 
348*9dbcead6SJeremy L Thompson     // -- Cleanup
349*9dbcead6SJeremy L Thompson     PetscCall(VecDestroy(&B_mesh));
350*9dbcead6SJeremy L Thompson     PetscCall(VecDestroy(&U_projected));
351*9dbcead6SJeremy L Thompson     PetscCall(MatDestroy(&M));
352*9dbcead6SJeremy L Thompson     PetscCall(KSPDestroy(&ksp));
353*9dbcead6SJeremy L Thompson   }
354*9dbcead6SJeremy L Thompson 
355*9dbcead6SJeremy L Thompson   // Cleanup
356*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmCeedContextDestroy(&swarm_ceed_context));
357*9dbcead6SJeremy L Thompson   PetscCall(DMDestroy(&dm_swarm));
358*9dbcead6SJeremy L Thompson   PetscCall(DMDestroy(&dm_mesh));
359*9dbcead6SJeremy L Thompson   PetscCall(VecDestroy(&U_mesh));
360*9dbcead6SJeremy L Thompson   return PetscFinalize();
361*9dbcead6SJeremy L Thompson }
362*9dbcead6SJeremy L Thompson 
363*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
364*9dbcead6SJeremy L Thompson // Context utilities
365*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
366*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmCeedContextCreate(DM dm_swarm, const char *ceed_resource, DMSwarmCeedContext *ctx) {
367*9dbcead6SJeremy L Thompson   DM dm_mesh;
368*9dbcead6SJeremy L Thompson 
369*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
370*9dbcead6SJeremy L Thompson   PetscCall(PetscNew(ctx));
371*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
372*9dbcead6SJeremy L Thompson 
373*9dbcead6SJeremy L Thompson   CeedInit(ceed_resource, &(*ctx)->ceed);
374*9dbcead6SJeremy L Thompson   // Background mesh objects
375*9dbcead6SJeremy L Thompson   {
376*9dbcead6SJeremy L Thompson     CeedInt elem_size, num_comp;
377*9dbcead6SJeremy L Thompson     BPData  bp_data = {.q_mode = CEED_GAUSS};
378*9dbcead6SJeremy L Thompson 
379*9dbcead6SJeremy L Thompson     PetscCall(CreateBasisFromPlex((*ctx)->ceed, dm_mesh, NULL, 0, 0, 0, bp_data, &(*ctx)->basis_u));
380*9dbcead6SJeremy L Thompson     PetscCall(CreateRestrictionFromPlex((*ctx)->ceed, dm_mesh, 0, NULL, 0, &(*ctx)->restriction_u_mesh));
381*9dbcead6SJeremy L Thompson 
382*9dbcead6SJeremy L Thompson     // -- U vector
383*9dbcead6SJeremy L Thompson     CeedElemRestrictionCreateVector((*ctx)->restriction_u_mesh, &(*ctx)->u_mesh_loc, NULL);
384*9dbcead6SJeremy L Thompson     CeedElemRestrictionCreateVector((*ctx)->restriction_u_mesh, &(*ctx)->v_mesh_loc, NULL);
385*9dbcead6SJeremy L Thompson     CeedElemRestrictionGetElementSize((*ctx)->restriction_u_mesh, &elem_size);
386*9dbcead6SJeremy L Thompson     CeedElemRestrictionGetNumComponents((*ctx)->restriction_u_mesh, &num_comp);
387*9dbcead6SJeremy L Thompson     CeedVectorCreate((*ctx)->ceed, elem_size * num_comp, &(*ctx)->u_mesh_elem);
388*9dbcead6SJeremy L Thompson   }
389*9dbcead6SJeremy L Thompson   // Swarm objects
390*9dbcead6SJeremy L Thompson   {
391*9dbcead6SJeremy L Thompson     PetscInt        dim;
392*9dbcead6SJeremy L Thompson     const PetscInt *cell_points;
393*9dbcead6SJeremy L Thompson     IS              is_points;
394*9dbcead6SJeremy L Thompson     Vec             X_ref;
395*9dbcead6SJeremy L Thompson     CeedInt         num_elem, num_comp, max_points_in_cell;
396*9dbcead6SJeremy L Thompson 
397*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmCreateReferenceCoordinates(dm_swarm, &is_points, &X_ref));
398*9dbcead6SJeremy L Thompson     PetscCall(DMGetDimension(dm_mesh, &dim));
399*9dbcead6SJeremy L Thompson     CeedElemRestrictionGetNumElements((*ctx)->restriction_u_mesh, &num_elem);
400*9dbcead6SJeremy L Thompson     CeedElemRestrictionGetNumComponents((*ctx)->restriction_u_mesh, &num_comp);
401*9dbcead6SJeremy L Thompson 
402*9dbcead6SJeremy L Thompson     PetscCall(ISGetIndices(is_points, &cell_points));
403*9dbcead6SJeremy L Thompson     PetscInt num_points = cell_points[num_elem + 1] - num_elem - 2;
404*9dbcead6SJeremy L Thompson     CeedInt  offsets[num_elem + 1 + num_points];
405*9dbcead6SJeremy L Thompson 
406*9dbcead6SJeremy L Thompson     for (PetscInt i = 0; i < num_elem + 1; i++) offsets[i] = cell_points[i + 1] - 1;
407*9dbcead6SJeremy L Thompson     for (PetscInt i = num_elem + 1; i < num_points + num_elem + 1; i++) offsets[i] = cell_points[i + 1];
408*9dbcead6SJeremy L Thompson     PetscCall(ISRestoreIndices(is_points, &cell_points));
409*9dbcead6SJeremy L Thompson 
410*9dbcead6SJeremy L Thompson     // -- Points restrictions
411*9dbcead6SJeremy L Thompson     CeedElemRestrictionCreateAtPoints((*ctx)->ceed, num_elem, num_points, num_comp, num_points * num_comp, CEED_MEM_HOST, CEED_COPY_VALUES, offsets,
412*9dbcead6SJeremy L Thompson                                       &(*ctx)->restriction_u_points);
413*9dbcead6SJeremy L Thompson     CeedElemRestrictionCreateAtPoints((*ctx)->ceed, num_elem, num_points, dim, num_points * dim, CEED_MEM_HOST, CEED_COPY_VALUES, offsets,
414*9dbcead6SJeremy L Thompson                                       &(*ctx)->restriction_x_points);
415*9dbcead6SJeremy L Thompson 
416*9dbcead6SJeremy L Thompson     // -- U vector
417*9dbcead6SJeremy L Thompson     CeedElemRestrictionGetMaxPointsInElement((*ctx)->restriction_u_points, &max_points_in_cell);
418*9dbcead6SJeremy L Thompson     CeedElemRestrictionCreateVector((*ctx)->restriction_u_points, &(*ctx)->u_points_loc, NULL);
419*9dbcead6SJeremy L Thompson     CeedVectorCreate((*ctx)->ceed, max_points_in_cell * num_comp, &(*ctx)->u_points_elem);
420*9dbcead6SJeremy L Thompson 
421*9dbcead6SJeremy L Thompson     // -- Ref coordinates
422*9dbcead6SJeremy L Thompson     {
423*9dbcead6SJeremy L Thompson       PetscMemType       X_mem_type;
424*9dbcead6SJeremy L Thompson       const PetscScalar *x;
425*9dbcead6SJeremy L Thompson 
426*9dbcead6SJeremy L Thompson       CeedVectorCreate((*ctx)->ceed, num_points * dim, &(*ctx)->x_ref_points_loc);
427*9dbcead6SJeremy L Thompson       CeedVectorCreate((*ctx)->ceed, max_points_in_cell * dim, &(*ctx)->x_ref_points_elem);
428*9dbcead6SJeremy L Thompson 
429*9dbcead6SJeremy L Thompson       PetscCall(VecGetArrayReadAndMemType(X_ref, (const PetscScalar **)&x, &X_mem_type));
430*9dbcead6SJeremy L Thompson       CeedVectorSetArray((*ctx)->x_ref_points_loc, MemTypeP2C(X_mem_type), CEED_COPY_VALUES, (CeedScalar *)x);
431*9dbcead6SJeremy L Thompson       PetscCall(VecRestoreArrayReadAndMemType(X_ref, (const PetscScalar **)&x));
432*9dbcead6SJeremy L Thompson     }
433*9dbcead6SJeremy L Thompson 
434*9dbcead6SJeremy L Thompson     // -- Cleanup
435*9dbcead6SJeremy L Thompson     PetscCall(ISDestroy(&is_points));
436*9dbcead6SJeremy L Thompson     PetscCall(VecDestroy(&X_ref));
437*9dbcead6SJeremy L Thompson   }
438*9dbcead6SJeremy L Thompson 
439*9dbcead6SJeremy L Thompson   PetscCall(DMSetApplicationContext(dm_mesh, (void *)(*ctx)));
440*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
441*9dbcead6SJeremy L Thompson }
442*9dbcead6SJeremy L Thompson 
443*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmCeedContextDestroy(DMSwarmCeedContext *ctx) {
444*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
445*9dbcead6SJeremy L Thompson   CeedDestroy(&(*ctx)->ceed);
446*9dbcead6SJeremy L Thompson   CeedVectorDestroy(&(*ctx)->u_mesh_loc);
447*9dbcead6SJeremy L Thompson   CeedVectorDestroy(&(*ctx)->v_mesh_loc);
448*9dbcead6SJeremy L Thompson   CeedVectorDestroy(&(*ctx)->u_mesh_elem);
449*9dbcead6SJeremy L Thompson   CeedVectorDestroy(&(*ctx)->u_points_loc);
450*9dbcead6SJeremy L Thompson   CeedVectorDestroy(&(*ctx)->u_points_elem);
451*9dbcead6SJeremy L Thompson   CeedVectorDestroy(&(*ctx)->x_ref_points_loc);
452*9dbcead6SJeremy L Thompson   CeedVectorDestroy(&(*ctx)->x_ref_points_elem);
453*9dbcead6SJeremy L Thompson   CeedElemRestrictionDestroy(&(*ctx)->restriction_u_mesh);
454*9dbcead6SJeremy L Thompson   CeedElemRestrictionDestroy(&(*ctx)->restriction_x_points);
455*9dbcead6SJeremy L Thompson   CeedElemRestrictionDestroy(&(*ctx)->restriction_u_points);
456*9dbcead6SJeremy L Thompson   CeedBasisDestroy(&(*ctx)->basis_u);
457*9dbcead6SJeremy L Thompson   PetscCall(PetscFree(*ctx));
458*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
459*9dbcead6SJeremy L Thompson }
460*9dbcead6SJeremy L Thompson 
461*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
462*9dbcead6SJeremy L Thompson // PETSc-libCEED memory space utilities
463*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
464*9dbcead6SJeremy L Thompson PetscErrorCode VecP2C(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed) {
465*9dbcead6SJeremy L Thompson   PetscScalar *x;
466*9dbcead6SJeremy L Thompson 
467*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
468*9dbcead6SJeremy L Thompson   PetscCall(VecGetArrayAndMemType(X_petsc, &x, mem_type));
469*9dbcead6SJeremy L Thompson   CeedVectorSetArray(x_ceed, MemTypeP2C(*mem_type), CEED_USE_POINTER, x);
470*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
471*9dbcead6SJeremy L Thompson }
472*9dbcead6SJeremy L Thompson 
473*9dbcead6SJeremy L Thompson PetscErrorCode VecC2P(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc) {
474*9dbcead6SJeremy L Thompson   PetscScalar *x;
475*9dbcead6SJeremy L Thompson 
476*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
477*9dbcead6SJeremy L Thompson   CeedVectorTakeArray(x_ceed, MemTypeP2C(mem_type), &x);
478*9dbcead6SJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(X_petsc, &x));
479*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
480*9dbcead6SJeremy L Thompson }
481*9dbcead6SJeremy L Thompson 
482*9dbcead6SJeremy L Thompson PetscErrorCode VecReadP2C(Vec X_petsc, PetscMemType *mem_type, CeedVector x_ceed) {
483*9dbcead6SJeremy L Thompson   PetscScalar *x;
484*9dbcead6SJeremy L Thompson 
485*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
486*9dbcead6SJeremy L Thompson   PetscCall(VecGetArrayReadAndMemType(X_petsc, (const PetscScalar **)&x, mem_type));
487*9dbcead6SJeremy L Thompson   CeedVectorSetArray(x_ceed, MemTypeP2C(*mem_type), CEED_USE_POINTER, x);
488*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
489*9dbcead6SJeremy L Thompson }
490*9dbcead6SJeremy L Thompson 
491*9dbcead6SJeremy L Thompson PetscErrorCode VecReadC2P(CeedVector x_ceed, PetscMemType mem_type, Vec X_petsc) {
492*9dbcead6SJeremy L Thompson   PetscScalar *x;
493*9dbcead6SJeremy L Thompson 
494*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
495*9dbcead6SJeremy L Thompson   CeedVectorTakeArray(x_ceed, MemTypeP2C(mem_type), &x);
496*9dbcead6SJeremy L Thompson   PetscCall(VecRestoreArrayReadAndMemType(X_petsc, (const PetscScalar **)&x));
497*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
498*9dbcead6SJeremy L Thompson }
499*9dbcead6SJeremy L Thompson 
500*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmPICFieldP2C(DM dm_swarm, const char *field, CeedVector x_ceed) {
501*9dbcead6SJeremy L Thompson   PetscScalar *x;
502*9dbcead6SJeremy L Thompson 
503*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
504*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmGetField(dm_swarm, field, NULL, NULL, (void **)&x));
505*9dbcead6SJeremy L Thompson   CeedVectorSetArray(x_ceed, CEED_MEM_HOST, CEED_USE_POINTER, (CeedScalar *)x);
506*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
507*9dbcead6SJeremy L Thompson }
508*9dbcead6SJeremy L Thompson 
509*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmPICFieldC2P(DM dm_swarm, const char *field, CeedVector x_ceed) {
510*9dbcead6SJeremy L Thompson   PetscScalar *x;
511*9dbcead6SJeremy L Thompson 
512*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
513*9dbcead6SJeremy L Thompson   CeedVectorTakeArray(x_ceed, CEED_MEM_HOST, (CeedScalar **)&x);
514*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmRestoreField(dm_swarm, field, NULL, NULL, (void **)&x));
515*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
516*9dbcead6SJeremy L Thompson }
517*9dbcead6SJeremy L Thompson 
518*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
519*9dbcead6SJeremy L Thompson // Solution functions
520*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
521*9dbcead6SJeremy L Thompson PetscScalar EvalU_Poly(PetscInt dim, const PetscScalar x[]) {
522*9dbcead6SJeremy L Thompson   PetscScalar       result = 0.0;
523*9dbcead6SJeremy L Thompson   const PetscScalar p[5]   = {3, 1, 4, 1, 5};
524*9dbcead6SJeremy L Thompson 
525*9dbcead6SJeremy L Thompson   for (PetscInt d = 0; d < dim; d++) {
526*9dbcead6SJeremy L Thompson     PetscScalar result_1d = 1.0;
527*9dbcead6SJeremy L Thompson 
528*9dbcead6SJeremy L Thompson     for (PetscInt i = 4; i >= 0; i--) result_1d = result_1d * x[d] + p[i];
529*9dbcead6SJeremy L Thompson     result += result_1d;
530*9dbcead6SJeremy L Thompson   }
531*9dbcead6SJeremy L Thompson   return result * 1E-3;
532*9dbcead6SJeremy L Thompson }
533*9dbcead6SJeremy L Thompson 
534*9dbcead6SJeremy L Thompson PetscScalar EvalU_Tanh(PetscInt dim, const PetscScalar x[]) {
535*9dbcead6SJeremy L Thompson   PetscScalar result = 1.0, center = 0.1;
536*9dbcead6SJeremy L Thompson 
537*9dbcead6SJeremy L Thompson   for (PetscInt d = 0; d < dim; d++) {
538*9dbcead6SJeremy L Thompson     result *= tanh(x[d] - center);
539*9dbcead6SJeremy L Thompson     center += 0.1;
540*9dbcead6SJeremy L Thompson   }
541*9dbcead6SJeremy L Thompson   return result;
542*9dbcead6SJeremy L Thompson }
543*9dbcead6SJeremy L Thompson 
544*9dbcead6SJeremy L Thompson PetscScalar EvalU_Sphere(PetscInt dim, const PetscScalar x[]) {
545*9dbcead6SJeremy L Thompson   PetscScalar distance = sqrt(x[0] * x[0] + x[1] * x[1] + x[2] * x[2]);
546*9dbcead6SJeremy L Thompson 
547*9dbcead6SJeremy L Thompson   return distance < 1.0 ? 1.0 : 0.1;
548*9dbcead6SJeremy L Thompson }
549*9dbcead6SJeremy L Thompson 
550*9dbcead6SJeremy L Thompson PetscErrorCode EvalU_Poly_proj(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt num_comp, PetscScalar *u, void *ctx) {
551*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
552*9dbcead6SJeremy L Thompson 
553*9dbcead6SJeremy L Thompson   const PetscScalar f_x = EvalU_Poly(dim, x);
554*9dbcead6SJeremy L Thompson 
555*9dbcead6SJeremy L Thompson   for (PetscInt c = 0; c < num_comp; c++) u[c] = (c + 1.0) * f_x;
556*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
557*9dbcead6SJeremy L Thompson }
558*9dbcead6SJeremy L Thompson 
559*9dbcead6SJeremy L Thompson PetscErrorCode EvalU_Tanh_proj(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt num_comp, PetscScalar *u, void *ctx) {
560*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
561*9dbcead6SJeremy L Thompson 
562*9dbcead6SJeremy L Thompson   const PetscScalar f_x = EvalU_Tanh(dim, x);
563*9dbcead6SJeremy L Thompson 
564*9dbcead6SJeremy L Thompson   for (PetscInt c = 0; c < num_comp; c++) u[c] = (c + 1.0) * f_x;
565*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
566*9dbcead6SJeremy L Thompson }
567*9dbcead6SJeremy L Thompson 
568*9dbcead6SJeremy L Thompson PetscErrorCode EvalU_Sphere_proj(PetscInt dim, PetscReal t, const PetscReal x[], PetscInt num_comp, PetscScalar *u, void *ctx) {
569*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
570*9dbcead6SJeremy L Thompson 
571*9dbcead6SJeremy L Thompson   const PetscScalar f_x = EvalU_Sphere(dim, x);
572*9dbcead6SJeremy L Thompson 
573*9dbcead6SJeremy L Thompson   for (PetscInt c = 0; c < num_comp; c++) u[c] = (c + 1.0) * f_x;
574*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
575*9dbcead6SJeremy L Thompson }
576*9dbcead6SJeremy L Thompson 
577*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
578*9dbcead6SJeremy L Thompson // Swarm point location utility
579*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
580*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmInitalizePointLocations(DM dm_swarm, PointSwarmType point_swarm_type, PetscInt num_points, PetscInt num_points_per_cell) {
581*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
582*9dbcead6SJeremy L Thompson 
583*9dbcead6SJeremy L Thompson   switch (point_swarm_type) {
584*9dbcead6SJeremy L Thompson     case SWARM_GAUSS:
585*9dbcead6SJeremy L Thompson     case SWARM_UNIFORM: {
586*9dbcead6SJeremy L Thompson       // -- Set gauss or uniform point locations in each cell
587*9dbcead6SJeremy L Thompson       PetscInt    num_points_per_cell_1d = round(cbrt(num_points_per_cell * 1.0)), dim = 3;
588*9dbcead6SJeremy L Thompson       PetscScalar point_coords[num_points_per_cell * 3];
589*9dbcead6SJeremy L Thompson       CeedScalar  points_1d[num_points_per_cell_1d], weights_1d[num_points_per_cell_1d];
590*9dbcead6SJeremy L Thompson 
591*9dbcead6SJeremy L Thompson       if (point_swarm_type == SWARM_GAUSS) {
592*9dbcead6SJeremy L Thompson         PetscCall(CeedGaussQuadrature(num_points_per_cell_1d, points_1d, weights_1d));
593*9dbcead6SJeremy L Thompson       } else {
594*9dbcead6SJeremy L Thompson         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;
595*9dbcead6SJeremy L Thompson       }
596*9dbcead6SJeremy L Thompson       for (PetscInt i = 0; i < num_points_per_cell_1d; i++) {
597*9dbcead6SJeremy L Thompson         for (PetscInt j = 0; j < num_points_per_cell_1d; j++) {
598*9dbcead6SJeremy L Thompson           for (PetscInt k = 0; k < num_points_per_cell_1d; k++) {
599*9dbcead6SJeremy L Thompson             PetscInt p = (i * num_points_per_cell_1d + j) * num_points_per_cell_1d + k;
600*9dbcead6SJeremy L Thompson 
601*9dbcead6SJeremy L Thompson             point_coords[p * dim + 0] = points_1d[i];
602*9dbcead6SJeremy L Thompson             point_coords[p * dim + 1] = points_1d[j];
603*9dbcead6SJeremy L Thompson             point_coords[p * dim + 2] = points_1d[k];
604*9dbcead6SJeremy L Thompson           }
605*9dbcead6SJeremy L Thompson         }
606*9dbcead6SJeremy L Thompson       }
607*9dbcead6SJeremy L Thompson       PetscCall(DMSwarmSetPointCoordinatesCellwise(dm_swarm, num_points_per_cell_1d * num_points_per_cell_1d * num_points_per_cell_1d, point_coords));
608*9dbcead6SJeremy L Thompson     } break;
609*9dbcead6SJeremy L Thompson     case SWARM_CELL_RANDOM: {
610*9dbcead6SJeremy L Thompson       // -- Set points randomly in each cell
611*9dbcead6SJeremy L Thompson       PetscInt     dim = 3, num_cells_total = 1, num_cells[] = {1, 1, 1};
612*9dbcead6SJeremy L Thompson       PetscScalar *point_coords;
613*9dbcead6SJeremy L Thompson       PetscRandom  rng;
614*9dbcead6SJeremy L Thompson 
615*9dbcead6SJeremy L Thompson       PetscOptionsBegin(PetscObjectComm((PetscObject)dm_swarm), NULL, "libCEED example using PETSc with DMSwarm", NULL);
616*9dbcead6SJeremy L Thompson 
617*9dbcead6SJeremy L Thompson       PetscCall(PetscOptionsInt("-dm_plex_dim", "Background mesh dimension", NULL, dim, &dim, NULL));
618*9dbcead6SJeremy L Thompson       PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of cells", NULL, num_cells, &dim, NULL));
619*9dbcead6SJeremy L Thompson 
620*9dbcead6SJeremy L Thompson       PetscOptionsEnd();
621*9dbcead6SJeremy L Thompson 
622*9dbcead6SJeremy L Thompson       PetscCall(PetscRandomCreate(PetscObjectComm((PetscObject)dm_swarm), &rng));
623*9dbcead6SJeremy L Thompson 
624*9dbcead6SJeremy L Thompson       num_cells_total = num_cells[0] * num_cells[1] * num_cells[2];
625*9dbcead6SJeremy L Thompson       PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords));
626*9dbcead6SJeremy L Thompson       for (PetscInt c = 0; c < num_cells_total; c++) {
627*9dbcead6SJeremy L Thompson         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]};
628*9dbcead6SJeremy L Thompson 
629*9dbcead6SJeremy L Thompson         for (PetscInt p = 0; p < num_points_per_cell; p++) {
630*9dbcead6SJeremy L Thompson           PetscInt    point_index = c * num_points_per_cell + p;
631*9dbcead6SJeremy L Thompson           PetscScalar random_value;
632*9dbcead6SJeremy L Thompson 
633*9dbcead6SJeremy L Thompson           for (PetscInt i = 0; i < dim; i++) {
634*9dbcead6SJeremy L Thompson             PetscCall(PetscRandomGetValue(rng, &random_value));
635*9dbcead6SJeremy L Thompson             point_coords[point_index * dim + i] = -1.0 + cell_index[i] * 2.0 / (num_cells[i] + 1.0) + random_value;
636*9dbcead6SJeremy L Thompson           }
637*9dbcead6SJeremy L Thompson         }
638*9dbcead6SJeremy L Thompson       }
639*9dbcead6SJeremy L Thompson       PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords));
640*9dbcead6SJeremy L Thompson       PetscCall(PetscRandomDestroy(&rng));
641*9dbcead6SJeremy L Thompson     } break;
642*9dbcead6SJeremy L Thompson     case SWARM_SINUSOIDAL: {
643*9dbcead6SJeremy L Thompson       // -- Set points distributed per sinusoidal functions
644*9dbcead6SJeremy L Thompson       PetscInt     dim = 3;
645*9dbcead6SJeremy L Thompson       PetscScalar *point_coords;
646*9dbcead6SJeremy L Thompson       DM           dm_mesh;
647*9dbcead6SJeremy L Thompson 
648*9dbcead6SJeremy L Thompson       PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
649*9dbcead6SJeremy L Thompson       PetscCall(DMGetDimension(dm_mesh, &dim));
650*9dbcead6SJeremy L Thompson 
651*9dbcead6SJeremy L Thompson       PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords));
652*9dbcead6SJeremy L Thompson       for (PetscInt p = 0; p < num_points; p++) {
653*9dbcead6SJeremy L Thompson         point_coords[p * dim + 0] = -PetscCosReal((PetscReal)(p + 1) / (PetscReal)(num_points + 1) * PETSC_PI);
654*9dbcead6SJeremy L Thompson         if (dim > 1) point_coords[p * dim + 1] = -PetscSinReal((PetscReal)(p + 1) / (PetscReal)(num_points + 1) * PETSC_PI);
655*9dbcead6SJeremy L Thompson         if (dim > 2) point_coords[p * dim + 2] = PetscSinReal((PetscReal)(p + 1) / (PetscReal)(num_points + 1) * PETSC_PI);
656*9dbcead6SJeremy L Thompson       }
657*9dbcead6SJeremy L Thompson       PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&point_coords));
658*9dbcead6SJeremy L Thompson     } break;
659*9dbcead6SJeremy L Thompson   }
660*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmMigrate(dm_swarm, PETSC_TRUE));
661*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
662*9dbcead6SJeremy L Thompson }
663*9dbcead6SJeremy L Thompson 
664*9dbcead6SJeremy L Thompson /*@C
665*9dbcead6SJeremy L Thompson   DMSwarmCreateReferenceCoordinates - Compute the cell reference coordinates for local DMSwarm points.
666*9dbcead6SJeremy L Thompson 
667*9dbcead6SJeremy L Thompson   Collective
668*9dbcead6SJeremy L Thompson 
669*9dbcead6SJeremy L Thompson   Input Parameter:
670*9dbcead6SJeremy L Thompson . dm_swarm  - the `DMSwarm`
671*9dbcead6SJeremy L Thompson 
672*9dbcead6SJeremy L Thompson   Output Parameters:
673*9dbcead6SJeremy L Thompson + is_points    - The IS object for indexing into points per cell
674*9dbcead6SJeremy L Thompson - X_points_ref - Vec holding the cell reference coordinates for local DMSwarm points
675*9dbcead6SJeremy L Thompson 
676*9dbcead6SJeremy L Thompson The index set contains ranges of indices for each local cell. This range contains the indices of every point in the cell.
677*9dbcead6SJeremy L Thompson 
678*9dbcead6SJeremy L Thompson ```
679*9dbcead6SJeremy L Thompson total_num_cells
680*9dbcead6SJeremy L Thompson cell_0_start_index
681*9dbcead6SJeremy L Thompson cell_1_start_index
682*9dbcead6SJeremy L Thompson cell_2_start_index
683*9dbcead6SJeremy L Thompson ...
684*9dbcead6SJeremy L Thompson cell_n_start_index
685*9dbcead6SJeremy L Thompson cell_n_stop_index
686*9dbcead6SJeremy L Thompson cell_0_point_0
687*9dbcead6SJeremy L Thompson cell_0_point_0
688*9dbcead6SJeremy L Thompson ...
689*9dbcead6SJeremy L Thompson cell_n_point_m
690*9dbcead6SJeremy L Thompson ```
691*9dbcead6SJeremy L Thompson 
692*9dbcead6SJeremy L Thompson   Level: beginner
693*9dbcead6SJeremy L Thompson 
694*9dbcead6SJeremy L Thompson .seealso: `DMSwarm`
695*9dbcead6SJeremy L Thompson @*/
696*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmCreateReferenceCoordinates(DM dm_swarm, IS *is_points, Vec *X_points_ref) {
697*9dbcead6SJeremy L Thompson   PetscInt           cell_start, cell_end, num_cells_local, dim, num_points_local, *cell_points, points_offset;
698*9dbcead6SJeremy L Thompson   PetscScalar       *coords_points_ref;
699*9dbcead6SJeremy L Thompson   const PetscScalar *coords_points_true;
700*9dbcead6SJeremy L Thompson   DM                 dm_mesh;
701*9dbcead6SJeremy L Thompson 
702*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
703*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
704*9dbcead6SJeremy L Thompson 
705*9dbcead6SJeremy L Thompson   // Create vector to hold reference coordinates
706*9dbcead6SJeremy L Thompson   {
707*9dbcead6SJeremy L Thompson     Vec X_points_true;
708*9dbcead6SJeremy L Thompson 
709*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmCreateLocalVectorFromField(dm_swarm, DMSwarmPICField_coor, &X_points_true));
710*9dbcead6SJeremy L Thompson     PetscCall(VecDuplicate(X_points_true, X_points_ref));
711*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmDestroyLocalVectorFromField(dm_swarm, DMSwarmPICField_coor, &X_points_true));
712*9dbcead6SJeremy L Thompson   }
713*9dbcead6SJeremy L Thompson 
714*9dbcead6SJeremy L Thompson   // Allocate index set array
715*9dbcead6SJeremy L Thompson   PetscCall(DMPlexGetHeightStratum(dm_mesh, 0, &cell_start, &cell_end));
716*9dbcead6SJeremy L Thompson   num_cells_local = cell_end - cell_start;
717*9dbcead6SJeremy L Thompson   points_offset   = num_cells_local + 2;
718*9dbcead6SJeremy L Thompson   PetscCall(VecGetLocalSize(*X_points_ref, &num_points_local));
719*9dbcead6SJeremy L Thompson   PetscCall(DMGetDimension(dm_mesh, &dim));
720*9dbcead6SJeremy L Thompson   num_points_local /= dim;
721*9dbcead6SJeremy L Thompson   PetscCall(PetscMalloc1(num_points_local + num_cells_local + 2, &cell_points));
722*9dbcead6SJeremy L Thompson   cell_points[0] = num_cells_local;
723*9dbcead6SJeremy L Thompson 
724*9dbcead6SJeremy L Thompson   // Get reference coordinates for each swarm point wrt the elements in the background mesh
725*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmSortGetAccess(dm_swarm));
726*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points_true));
727*9dbcead6SJeremy L Thompson   PetscCall(VecGetArray(*X_points_ref, &coords_points_ref));
728*9dbcead6SJeremy L Thompson   for (PetscInt cell = cell_start, num_points_processed = 0; cell < cell_end; cell++) {
729*9dbcead6SJeremy L Thompson     PetscInt *points_in_cell, num_points_in_cell, local_cell = cell - cell_start;
730*9dbcead6SJeremy L Thompson     PetscReal v[3], J[9], invJ[9], detJ, v0_ref[3] = {-1.0, -1.0, -1.0};
731*9dbcead6SJeremy L Thompson 
732*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmSortGetPointsPerCell(dm_swarm, cell, &num_points_in_cell, &points_in_cell));
733*9dbcead6SJeremy L Thompson     // -- Reference coordinates for swarm points in background mesh element
734*9dbcead6SJeremy L Thompson     PetscCall(DMPlexComputeCellGeometryFEM(dm_mesh, cell, NULL, v, J, invJ, &detJ));
735*9dbcead6SJeremy L Thompson     cell_points[local_cell + 1] = num_points_processed + points_offset;
736*9dbcead6SJeremy L Thompson     for (PetscInt p = 0; p < num_points_in_cell; p++) {
737*9dbcead6SJeremy L Thompson       const CeedInt point = points_in_cell[p];
738*9dbcead6SJeremy L Thompson 
739*9dbcead6SJeremy L Thompson       cell_points[points_offset + (num_points_processed++)] = point;
740*9dbcead6SJeremy L Thompson       CoordinatesRealToRef(dim, dim, v0_ref, v, invJ, &coords_points_true[point * dim], &coords_points_ref[point * dim]);
741*9dbcead6SJeremy L Thompson     }
742*9dbcead6SJeremy L Thompson 
743*9dbcead6SJeremy L Thompson     // -- Cleanup
744*9dbcead6SJeremy L Thompson     PetscCall(PetscFree(points_in_cell));
745*9dbcead6SJeremy L Thompson   }
746*9dbcead6SJeremy L Thompson   cell_points[points_offset - 1] = num_points_local + points_offset;
747*9dbcead6SJeremy L Thompson 
748*9dbcead6SJeremy L Thompson   // Cleanup
749*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points_true));
750*9dbcead6SJeremy L Thompson   PetscCall(VecRestoreArray(*X_points_ref, &coords_points_ref));
751*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmSortRestoreAccess(dm_swarm));
752*9dbcead6SJeremy L Thompson 
753*9dbcead6SJeremy L Thompson   // Create index set
754*9dbcead6SJeremy L Thompson   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, num_points_local + points_offset, cell_points, PETSC_OWN_POINTER, is_points));
755*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
756*9dbcead6SJeremy L Thompson }
757*9dbcead6SJeremy L Thompson 
758*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
759*9dbcead6SJeremy L Thompson // Projection via PETSc
760*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
761*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmInterpolateFromCellToSwarm_Petsc(DM dm_swarm, const char *field, Vec U_mesh) {
762*9dbcead6SJeremy L Thompson   PetscInt           dim, num_comp, cell_start, cell_end;
763*9dbcead6SJeremy L Thompson   PetscScalar       *u_points;
764*9dbcead6SJeremy L Thompson   const PetscScalar *coords_points;
765*9dbcead6SJeremy L Thompson   const PetscReal    v0_ref[3] = {-1.0, -1.0, -1.0};
766*9dbcead6SJeremy L Thompson   DM                 dm_mesh;
767*9dbcead6SJeremy L Thompson   PetscSection       section_u_mesh_loc;
768*9dbcead6SJeremy L Thompson   PetscDS            ds;
769*9dbcead6SJeremy L Thompson   PetscFE            fe;
770*9dbcead6SJeremy L Thompson   PetscFEGeom        fe_geometry;
771*9dbcead6SJeremy L Thompson   PetscQuadrature    quadrature;
772*9dbcead6SJeremy L Thompson   Vec                U_loc;
773*9dbcead6SJeremy L Thompson 
774*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
775*9dbcead6SJeremy L Thompson   // Get mesh DM
776*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
777*9dbcead6SJeremy L Thompson   PetscCall(DMGetDimension(dm_mesh, &dim));
778*9dbcead6SJeremy L Thompson   {
779*9dbcead6SJeremy L Thompson     PetscSection section_u_mesh_loc_closure_permutation;
780*9dbcead6SJeremy L Thompson 
781*9dbcead6SJeremy L Thompson     PetscCall(DMGetLocalSection(dm_mesh, &section_u_mesh_loc_closure_permutation));
782*9dbcead6SJeremy L Thompson     PetscCall(PetscSectionClone(section_u_mesh_loc_closure_permutation, &section_u_mesh_loc));
783*9dbcead6SJeremy L Thompson     PetscCall(PetscSectionResetClosurePermutation(section_u_mesh_loc));
784*9dbcead6SJeremy L Thompson   }
785*9dbcead6SJeremy L Thompson 
786*9dbcead6SJeremy L Thompson   // Get local mesh values
787*9dbcead6SJeremy L Thompson   PetscCall(DMGetLocalVector(dm_mesh, &U_loc));
788*9dbcead6SJeremy L Thompson   PetscCall(VecZeroEntries(U_loc));
789*9dbcead6SJeremy L Thompson   PetscCall(DMGlobalToLocal(dm_mesh, U_mesh, INSERT_VALUES, U_loc));
790*9dbcead6SJeremy L Thompson 
791*9dbcead6SJeremy L Thompson   // Get local swarm data
792*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmSortGetAccess(dm_swarm));
793*9dbcead6SJeremy L Thompson   PetscCall(DMPlexGetHeightStratum(dm_mesh, 0, &cell_start, &cell_end));
794*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points));
795*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmGetField(dm_swarm, field, &num_comp, NULL, (void **)&u_points));
796*9dbcead6SJeremy L Thompson 
797*9dbcead6SJeremy L Thompson   // Interpolate values to each swarm point, one element in the background mesh at a time
798*9dbcead6SJeremy L Thompson   PetscCall(DMGetDS(dm_mesh, &ds));
799*9dbcead6SJeremy L Thompson   PetscCall(PetscDSGetDiscretization(ds, 0, (PetscObject *)&fe));
800*9dbcead6SJeremy L Thompson   for (PetscInt cell = cell_start; cell < cell_end; cell++) {
801*9dbcead6SJeremy L Thompson     PetscTabulation tabulation;
802*9dbcead6SJeremy L Thompson     PetscScalar    *u_cell = NULL, *coords_points_cell_true, *coords_points_cell_ref;
803*9dbcead6SJeremy L Thompson     PetscReal       v[dim], J[dim * dim], invJ[dim * dim], detJ;
804*9dbcead6SJeremy L Thompson     PetscInt       *points_cell;
805*9dbcead6SJeremy L Thompson     PetscInt        num_points_in_cell;
806*9dbcead6SJeremy L Thompson 
807*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmSortGetPointsPerCell(dm_swarm, cell, &num_points_in_cell, &points_cell));
808*9dbcead6SJeremy L Thompson     PetscCall(DMGetWorkArray(dm_mesh, num_points_in_cell * dim, MPIU_REAL, &coords_points_cell_true));
809*9dbcead6SJeremy L Thompson     PetscCall(DMGetWorkArray(dm_mesh, num_points_in_cell * dim, MPIU_REAL, &coords_points_cell_ref));
810*9dbcead6SJeremy L Thompson     // -- Reference coordinates for swarm points in background mesh element
811*9dbcead6SJeremy L Thompson     for (PetscInt p = 0; p < num_points_in_cell; p++) {
812*9dbcead6SJeremy L Thompson       for (PetscInt d = 0; d < dim; d++) coords_points_cell_true[p * dim + d] = coords_points[points_cell[p] * dim + d];
813*9dbcead6SJeremy L Thompson     }
814*9dbcead6SJeremy L Thompson     PetscCall(DMPlexComputeCellGeometryFEM(dm_mesh, cell, NULL, v, J, invJ, &detJ));
815*9dbcead6SJeremy L Thompson     for (PetscInt p = 0; p < num_points_in_cell; p++) {
816*9dbcead6SJeremy L Thompson       CoordinatesRealToRef(dim, dim, v0_ref, v, invJ, &coords_points_cell_true[p * dim], &coords_points_cell_ref[p * dim]);
817*9dbcead6SJeremy L Thompson     }
818*9dbcead6SJeremy L Thompson     // -- Interpolate values from current element in background mesh to swarm points
819*9dbcead6SJeremy L Thompson     PetscCall(PetscFECreateTabulation(fe, 1, num_points_in_cell, coords_points_cell_ref, 1, &tabulation));
820*9dbcead6SJeremy L Thompson     PetscCall(DMPlexVecGetClosure(dm_mesh, section_u_mesh_loc, U_loc, cell, NULL, &u_cell));
821*9dbcead6SJeremy L Thompson     PetscCall(PetscFEGetQuadrature(fe, &quadrature));
822*9dbcead6SJeremy L Thompson     PetscCall(PetscFECreateCellGeometry(fe, quadrature, &fe_geometry));
823*9dbcead6SJeremy L Thompson     for (PetscInt p = 0; p < num_points_in_cell; p++) {
824*9dbcead6SJeremy L Thompson       PetscCall(PetscFEInterpolateAtPoints_Static(fe, tabulation, u_cell, &fe_geometry, p, &u_points[points_cell[p] * num_comp]));
825*9dbcead6SJeremy L Thompson     }
826*9dbcead6SJeremy L Thompson 
827*9dbcead6SJeremy L Thompson     // -- Cleanup
828*9dbcead6SJeremy L Thompson     PetscCall(PetscFEDestroyCellGeometry(fe, &fe_geometry));
829*9dbcead6SJeremy L Thompson     PetscCall(DMPlexVecRestoreClosure(dm_mesh, section_u_mesh_loc, U_loc, cell, NULL, &u_cell));
830*9dbcead6SJeremy L Thompson     PetscCall(DMRestoreWorkArray(dm_mesh, num_points_in_cell * dim, MPIU_REAL, &coords_points_cell_true));
831*9dbcead6SJeremy L Thompson     PetscCall(DMRestoreWorkArray(dm_mesh, num_points_in_cell * dim, MPIU_REAL, &coords_points_cell_ref));
832*9dbcead6SJeremy L Thompson     PetscCall(PetscTabulationDestroy(&tabulation));
833*9dbcead6SJeremy L Thompson     PetscCall(PetscFree(points_cell));
834*9dbcead6SJeremy L Thompson   }
835*9dbcead6SJeremy L Thompson 
836*9dbcead6SJeremy L Thompson   // Cleanup
837*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points));
838*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmRestoreField(dm_swarm, field, NULL, NULL, (void **)&u_points));
839*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmSortRestoreAccess(dm_swarm));
840*9dbcead6SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm_mesh, &U_loc));
841*9dbcead6SJeremy L Thompson   PetscCall(PetscSectionDestroy(&section_u_mesh_loc));
842*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
843*9dbcead6SJeremy L Thompson }
844*9dbcead6SJeremy L Thompson 
845*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
846*9dbcead6SJeremy L Thompson // Projection via libCEED
847*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
848*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmInterpolateFromCellToSwarm_Ceed(DM dm_swarm, const char *field, Vec U_mesh) {
849*9dbcead6SJeremy L Thompson   PetscInt           num_elem;
850*9dbcead6SJeremy L Thompson   PetscMemType       U_mem_type;
851*9dbcead6SJeremy L Thompson   DM                 dm_mesh;
852*9dbcead6SJeremy L Thompson   Vec                U_mesh_loc;
853*9dbcead6SJeremy L Thompson   DMSwarmCeedContext swarm_ceed_context;
854*9dbcead6SJeremy L Thompson 
855*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
856*9dbcead6SJeremy L Thompson   // Get mesh DM
857*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
858*9dbcead6SJeremy L Thompson   PetscCall(DMGetApplicationContext(dm_mesh, (void *)&swarm_ceed_context));
859*9dbcead6SJeremy L Thompson 
860*9dbcead6SJeremy L Thompson   // Get mesh values
861*9dbcead6SJeremy L Thompson   PetscCall(DMGetLocalVector(dm_mesh, &U_mesh_loc));
862*9dbcead6SJeremy L Thompson   PetscCall(VecZeroEntries(U_mesh_loc));
863*9dbcead6SJeremy L Thompson   PetscCall(DMGlobalToLocal(dm_mesh, U_mesh, INSERT_VALUES, U_mesh_loc));
864*9dbcead6SJeremy L Thompson   PetscCall(VecReadP2C(U_mesh_loc, &U_mem_type, swarm_ceed_context->u_mesh_loc));
865*9dbcead6SJeremy L Thompson 
866*9dbcead6SJeremy L Thompson   // Get swarm access
867*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmSortGetAccess(dm_swarm));
868*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmPICFieldP2C(dm_swarm, field, swarm_ceed_context->u_points_loc));
869*9dbcead6SJeremy L Thompson 
870*9dbcead6SJeremy L Thompson   // Interpolate values to each swarm point, one element in the background mesh at a time
871*9dbcead6SJeremy L Thompson   CeedElemRestrictionGetNumElements(swarm_ceed_context->restriction_u_mesh, &num_elem);
872*9dbcead6SJeremy L Thompson   for (PetscInt e = 0; e < num_elem; e++) {
873*9dbcead6SJeremy L Thompson     PetscInt num_points_in_elem;
874*9dbcead6SJeremy L Thompson 
875*9dbcead6SJeremy L Thompson     CeedElemRestrictionGetNumPointsInElement(swarm_ceed_context->restriction_u_points, e, &num_points_in_elem);
876*9dbcead6SJeremy L Thompson 
877*9dbcead6SJeremy L Thompson     // -- Reference coordinates for swarm points in background mesh element
878*9dbcead6SJeremy L Thompson     CeedElemRestrictionApplyAtPointsInElement(swarm_ceed_context->restriction_x_points, e, CEED_NOTRANSPOSE, swarm_ceed_context->x_ref_points_loc,
879*9dbcead6SJeremy L Thompson                                               swarm_ceed_context->x_ref_points_elem, CEED_REQUEST_IMMEDIATE);
880*9dbcead6SJeremy L Thompson 
881*9dbcead6SJeremy L Thompson     // -- Interpolate values from current element in background mesh to swarm points
882*9dbcead6SJeremy L Thompson     // Note: This will only work for CPU backends at this time, as only CPU backends support ApplyBlock and ApplyAtPoints
883*9dbcead6SJeremy L Thompson     CeedElemRestrictionApplyBlock(swarm_ceed_context->restriction_u_mesh, e, CEED_NOTRANSPOSE, swarm_ceed_context->u_mesh_loc,
884*9dbcead6SJeremy L Thompson                                   swarm_ceed_context->u_mesh_elem, CEED_REQUEST_IMMEDIATE);
885*9dbcead6SJeremy L Thompson     CeedBasisApplyAtPoints(swarm_ceed_context->basis_u, num_points_in_elem, CEED_NOTRANSPOSE, CEED_EVAL_INTERP, swarm_ceed_context->x_ref_points_elem,
886*9dbcead6SJeremy L Thompson                            swarm_ceed_context->u_mesh_elem, swarm_ceed_context->u_points_elem);
887*9dbcead6SJeremy L Thompson 
888*9dbcead6SJeremy L Thompson     // -- Insert result back into local vector
889*9dbcead6SJeremy L Thompson     CeedElemRestrictionApplyAtPointsInElement(swarm_ceed_context->restriction_u_points, e, CEED_TRANSPOSE, swarm_ceed_context->u_points_elem,
890*9dbcead6SJeremy L Thompson                                               swarm_ceed_context->u_points_loc, CEED_REQUEST_IMMEDIATE);
891*9dbcead6SJeremy L Thompson   }
892*9dbcead6SJeremy L Thompson 
893*9dbcead6SJeremy L Thompson   // Cleanup
894*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmPICFieldC2P(dm_swarm, field, swarm_ceed_context->u_points_loc));
895*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmSortRestoreAccess(dm_swarm));
896*9dbcead6SJeremy L Thompson   PetscCall(VecReadC2P(swarm_ceed_context->u_mesh_loc, U_mem_type, U_mesh_loc));
897*9dbcead6SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm_mesh, &U_mesh_loc));
898*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
899*9dbcead6SJeremy L Thompson }
900*9dbcead6SJeremy L Thompson 
901*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
902*9dbcead6SJeremy L Thompson // Error checking utility
903*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
904*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmCheckSwarmValues(DM dm_swarm, const char *field, PetscScalar tolerance, TargetFuncProj TrueSolution) {
905*9dbcead6SJeremy L Thompson   PetscBool          within_tolerance = PETSC_TRUE;
906*9dbcead6SJeremy L Thompson   PetscInt           dim, num_comp, cell_start, cell_end;
907*9dbcead6SJeremy L Thompson   const PetscScalar *u_points, *coords_points;
908*9dbcead6SJeremy L Thompson   DM                 dm_mesh;
909*9dbcead6SJeremy L Thompson 
910*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
911*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmSortGetAccess(dm_swarm));
912*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
913*9dbcead6SJeremy L Thompson   PetscCall(DMGetDimension(dm_mesh, &dim));
914*9dbcead6SJeremy L Thompson   PetscCall(DMPlexGetHeightStratum(dm_mesh, 0, &cell_start, &cell_end));
915*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmGetField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points));
916*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmGetField(dm_swarm, field, &num_comp, NULL, (void **)&u_points));
917*9dbcead6SJeremy L Thompson 
918*9dbcead6SJeremy L Thompson   // Interpolate values to each swarm point, one element in the background mesh at a time
919*9dbcead6SJeremy L Thompson   for (PetscInt cell = cell_start; cell < cell_end; cell++) {
920*9dbcead6SJeremy L Thompson     PetscInt *points;
921*9dbcead6SJeremy L Thompson     PetscInt  num_points_in_cell;
922*9dbcead6SJeremy L Thompson 
923*9dbcead6SJeremy L Thompson     PetscCall(DMSwarmSortGetPointsPerCell(dm_swarm, cell, &num_points_in_cell, &points));
924*9dbcead6SJeremy L Thompson     // -- Reference coordinates for swarm points in background mesh element
925*9dbcead6SJeremy L Thompson     for (PetscInt p = 0; p < num_points_in_cell; p++) {
926*9dbcead6SJeremy L Thompson       PetscScalar x[dim], u_true[num_comp];
927*9dbcead6SJeremy L Thompson 
928*9dbcead6SJeremy L Thompson       for (PetscInt d = 0; d < dim; d++) x[d] = coords_points[points[p] * dim + d];
929*9dbcead6SJeremy L Thompson       PetscCall(TrueSolution(dim, 0.0, x, num_comp, u_true, NULL));
930*9dbcead6SJeremy L Thompson       for (PetscInt i = 0; i < num_comp; i++) {
931*9dbcead6SJeremy L Thompson         if (PetscAbs(u_points[points[p] * num_comp + i] - u_true[i]) > tolerance) {
932*9dbcead6SJeremy L Thompson           within_tolerance = PETSC_FALSE;
933*9dbcead6SJeremy L Thompson           PetscCall(PetscPrintf(PetscObjectComm((PetscObject)dm_swarm),
934*9dbcead6SJeremy L Thompson                                 "Incorrect interpolated value, cell %" PetscInt_FMT " point %" PetscInt_FMT " component %" PetscInt_FMT
935*9dbcead6SJeremy L Thompson                                 ", found %f expected %f\n",
936*9dbcead6SJeremy L Thompson                                 cell, p, i, u_points[points[p] * num_comp + i], u_true[i]));
937*9dbcead6SJeremy L Thompson         }
938*9dbcead6SJeremy L Thompson       }
939*9dbcead6SJeremy L Thompson     }
940*9dbcead6SJeremy L Thompson 
941*9dbcead6SJeremy L Thompson     // -- Cleanup
942*9dbcead6SJeremy L Thompson     PetscCall(PetscFree(points));
943*9dbcead6SJeremy L Thompson   }
944*9dbcead6SJeremy L Thompson 
945*9dbcead6SJeremy L Thompson   // Cleanup
946*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmRestoreField(dm_swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&coords_points));
947*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmRestoreField(dm_swarm, field, NULL, NULL, (void **)&u_points));
948*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmSortRestoreAccess(dm_swarm));
949*9dbcead6SJeremy L Thompson   PetscCheck(within_tolerance, PetscObjectComm((PetscObject)dm_swarm), PETSC_ERR_USER, "Interpolation to swarm points not within tolerance");
950*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
951*9dbcead6SJeremy L Thompson }
952*9dbcead6SJeremy L Thompson 
953*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
954*9dbcead6SJeremy L Thompson // RHS for Swarm to Mesh projection
955*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
956*9dbcead6SJeremy L Thompson PetscErrorCode DMSwarmCreateProjectionRHS(DM dm_swarm, const char *field, Vec B_mesh) {
957*9dbcead6SJeremy L Thompson   PetscInt           num_elem;
958*9dbcead6SJeremy L Thompson   PetscMemType       B_mem_type;
959*9dbcead6SJeremy L Thompson   DM                 dm_mesh;
960*9dbcead6SJeremy L Thompson   Vec                B_mesh_loc;
961*9dbcead6SJeremy L Thompson   DMSwarmCeedContext swarm_ceed_context;
962*9dbcead6SJeremy L Thompson 
963*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
964*9dbcead6SJeremy L Thompson   // Get mesh DM
965*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmGetCellDM(dm_swarm, &dm_mesh));
966*9dbcead6SJeremy L Thompson   PetscCall(DMGetApplicationContext(dm_mesh, (void *)&swarm_ceed_context));
967*9dbcead6SJeremy L Thompson 
968*9dbcead6SJeremy L Thompson   // Get mesh values
969*9dbcead6SJeremy L Thompson   PetscCall(DMGetLocalVector(dm_mesh, &B_mesh_loc));
970*9dbcead6SJeremy L Thompson   PetscCall(VecZeroEntries(B_mesh_loc));
971*9dbcead6SJeremy L Thompson   PetscCall(VecP2C(B_mesh_loc, &B_mem_type, swarm_ceed_context->v_mesh_loc));
972*9dbcead6SJeremy L Thompson 
973*9dbcead6SJeremy L Thompson   // Get swarm access
974*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmSortGetAccess(dm_swarm));
975*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmPICFieldP2C(dm_swarm, field, swarm_ceed_context->u_points_loc));
976*9dbcead6SJeremy L Thompson 
977*9dbcead6SJeremy L Thompson   // Interpolate values to each swarm point, one element in the background mesh at a time
978*9dbcead6SJeremy L Thompson   CeedElemRestrictionGetNumElements(swarm_ceed_context->restriction_u_mesh, &num_elem);
979*9dbcead6SJeremy L Thompson   for (PetscInt e = 0; e < num_elem; e++) {
980*9dbcead6SJeremy L Thompson     PetscInt num_points_in_elem;
981*9dbcead6SJeremy L Thompson 
982*9dbcead6SJeremy L Thompson     CeedElemRestrictionGetNumPointsInElement(swarm_ceed_context->restriction_u_points, e, &num_points_in_elem);
983*9dbcead6SJeremy L Thompson 
984*9dbcead6SJeremy L Thompson     // -- Reference coordinates for swarm points in background mesh element
985*9dbcead6SJeremy L Thompson     CeedElemRestrictionApplyAtPointsInElement(swarm_ceed_context->restriction_x_points, e, CEED_NOTRANSPOSE, swarm_ceed_context->x_ref_points_loc,
986*9dbcead6SJeremy L Thompson                                               swarm_ceed_context->x_ref_points_elem, CEED_REQUEST_IMMEDIATE);
987*9dbcead6SJeremy L Thompson 
988*9dbcead6SJeremy L Thompson     // -- Interpolate values from current element in background mesh to swarm points
989*9dbcead6SJeremy L Thompson     // Note: This will only work for CPU backends at this time, as only CPU backends support ApplyBlock and ApplyAtPoints
990*9dbcead6SJeremy L Thompson     CeedElemRestrictionApplyAtPointsInElement(swarm_ceed_context->restriction_u_points, e, CEED_NOTRANSPOSE, swarm_ceed_context->u_points_loc,
991*9dbcead6SJeremy L Thompson                                               swarm_ceed_context->u_points_elem, CEED_REQUEST_IMMEDIATE);
992*9dbcead6SJeremy L Thompson     CeedBasisApplyAtPoints(swarm_ceed_context->basis_u, num_points_in_elem, CEED_TRANSPOSE, CEED_EVAL_INTERP, swarm_ceed_context->x_ref_points_elem,
993*9dbcead6SJeremy L Thompson                            swarm_ceed_context->u_points_elem, swarm_ceed_context->u_mesh_elem);
994*9dbcead6SJeremy L Thompson 
995*9dbcead6SJeremy L Thompson     // -- Insert result back into local vector
996*9dbcead6SJeremy L Thompson     CeedElemRestrictionApplyBlock(swarm_ceed_context->restriction_u_mesh, e, CEED_TRANSPOSE, swarm_ceed_context->u_mesh_elem,
997*9dbcead6SJeremy L Thompson                                   swarm_ceed_context->v_mesh_loc, CEED_REQUEST_IMMEDIATE);
998*9dbcead6SJeremy L Thompson   }
999*9dbcead6SJeremy L Thompson 
1000*9dbcead6SJeremy L Thompson   // Restore PETSc Vecs and Local to Global
1001*9dbcead6SJeremy L Thompson   PetscCall(VecC2P(swarm_ceed_context->v_mesh_loc, B_mem_type, B_mesh_loc));
1002*9dbcead6SJeremy L Thompson   PetscCall(VecZeroEntries(B_mesh));
1003*9dbcead6SJeremy L Thompson   PetscCall(DMLocalToGlobal(dm_mesh, B_mesh_loc, ADD_VALUES, B_mesh));
1004*9dbcead6SJeremy L Thompson 
1005*9dbcead6SJeremy L Thompson   // Cleanup
1006*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmPICFieldC2P(dm_swarm, field, swarm_ceed_context->u_points_loc));
1007*9dbcead6SJeremy L Thompson   PetscCall(DMSwarmSortRestoreAccess(dm_swarm));
1008*9dbcead6SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm_mesh, &B_mesh_loc));
1009*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
1010*9dbcead6SJeremy L Thompson }
1011*9dbcead6SJeremy L Thompson 
1012*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
1013*9dbcead6SJeremy L Thompson // Swarm "mass matrix"
1014*9dbcead6SJeremy L Thompson // ------------------------------------------------------------------------------------------------
1015*9dbcead6SJeremy L Thompson PetscErrorCode MatMult_SwarmMass(Mat A, Vec U_mesh, Vec V_mesh) {
1016*9dbcead6SJeremy L Thompson   PetscInt           num_elem;
1017*9dbcead6SJeremy L Thompson   PetscMemType       U_mem_type, V_mem_type;
1018*9dbcead6SJeremy L Thompson   DM                 dm_mesh;
1019*9dbcead6SJeremy L Thompson   Vec                U_mesh_loc, V_mesh_loc;
1020*9dbcead6SJeremy L Thompson   DMSwarmCeedContext swarm_ceed_context;
1021*9dbcead6SJeremy L Thompson 
1022*9dbcead6SJeremy L Thompson   PetscFunctionBeginUser;
1023*9dbcead6SJeremy L Thompson   // Get mesh DM
1024*9dbcead6SJeremy L Thompson   PetscCall(MatGetDM(A, &dm_mesh));
1025*9dbcead6SJeremy L Thompson   PetscCall(DMGetApplicationContext(dm_mesh, (void *)&swarm_ceed_context));
1026*9dbcead6SJeremy L Thompson 
1027*9dbcead6SJeremy L Thompson   // Global to Local and get PETSc Vec access
1028*9dbcead6SJeremy L Thompson   PetscCall(DMGetLocalVector(dm_mesh, &U_mesh_loc));
1029*9dbcead6SJeremy L Thompson   PetscCall(VecZeroEntries(U_mesh_loc));
1030*9dbcead6SJeremy L Thompson   PetscCall(DMGlobalToLocal(dm_mesh, U_mesh, INSERT_VALUES, U_mesh_loc));
1031*9dbcead6SJeremy L Thompson   PetscCall(VecReadP2C(U_mesh_loc, &U_mem_type, swarm_ceed_context->u_mesh_loc));
1032*9dbcead6SJeremy L Thompson   PetscCall(DMGetLocalVector(dm_mesh, &V_mesh_loc));
1033*9dbcead6SJeremy L Thompson   PetscCall(VecZeroEntries(V_mesh_loc));
1034*9dbcead6SJeremy L Thompson   PetscCall(VecP2C(V_mesh_loc, &V_mem_type, swarm_ceed_context->v_mesh_loc));
1035*9dbcead6SJeremy L Thompson 
1036*9dbcead6SJeremy L Thompson   // Interpolate values to each swarm point, one element in the background mesh at a time
1037*9dbcead6SJeremy L Thompson   CeedElemRestrictionGetNumElements(swarm_ceed_context->restriction_u_mesh, &num_elem);
1038*9dbcead6SJeremy L Thompson   for (PetscInt e = 0; e < num_elem; e++) {
1039*9dbcead6SJeremy L Thompson     PetscInt num_points_in_elem;
1040*9dbcead6SJeremy L Thompson 
1041*9dbcead6SJeremy L Thompson     CeedElemRestrictionGetNumPointsInElement(swarm_ceed_context->restriction_u_points, e, &num_points_in_elem);
1042*9dbcead6SJeremy L Thompson 
1043*9dbcead6SJeremy L Thompson     // -- Reference coordinates for swarm points in background mesh element
1044*9dbcead6SJeremy L Thompson     CeedElemRestrictionApplyAtPointsInElement(swarm_ceed_context->restriction_x_points, e, CEED_NOTRANSPOSE, swarm_ceed_context->x_ref_points_loc,
1045*9dbcead6SJeremy L Thompson                                               swarm_ceed_context->x_ref_points_elem, CEED_REQUEST_IMMEDIATE);
1046*9dbcead6SJeremy L Thompson 
1047*9dbcead6SJeremy L Thompson     // -- Interpolate values from current element in background mesh to swarm points
1048*9dbcead6SJeremy L Thompson     // Note: This will only work for CPU backends at this time, as only CPU backends support ApplyBlock and ApplyAtPoints
1049*9dbcead6SJeremy L Thompson     CeedElemRestrictionApplyBlock(swarm_ceed_context->restriction_u_mesh, e, CEED_NOTRANSPOSE, swarm_ceed_context->u_mesh_loc,
1050*9dbcead6SJeremy L Thompson                                   swarm_ceed_context->u_mesh_elem, CEED_REQUEST_IMMEDIATE);
1051*9dbcead6SJeremy L Thompson     CeedBasisApplyAtPoints(swarm_ceed_context->basis_u, num_points_in_elem, CEED_NOTRANSPOSE, CEED_EVAL_INTERP, swarm_ceed_context->x_ref_points_elem,
1052*9dbcead6SJeremy L Thompson                            swarm_ceed_context->u_mesh_elem, swarm_ceed_context->u_points_elem);
1053*9dbcead6SJeremy L Thompson 
1054*9dbcead6SJeremy L Thompson     // -- Interpolate transpose back from swarm points to mesh
1055*9dbcead6SJeremy L Thompson     CeedBasisApplyAtPoints(swarm_ceed_context->basis_u, num_points_in_elem, CEED_TRANSPOSE, CEED_EVAL_INTERP, swarm_ceed_context->x_ref_points_elem,
1056*9dbcead6SJeremy L Thompson                            swarm_ceed_context->u_points_elem, swarm_ceed_context->u_mesh_elem);
1057*9dbcead6SJeremy L Thompson     CeedElemRestrictionApplyBlock(swarm_ceed_context->restriction_u_mesh, e, CEED_TRANSPOSE, swarm_ceed_context->u_mesh_elem,
1058*9dbcead6SJeremy L Thompson                                   swarm_ceed_context->v_mesh_loc, CEED_REQUEST_IMMEDIATE);
1059*9dbcead6SJeremy L Thompson   }
1060*9dbcead6SJeremy L Thompson 
1061*9dbcead6SJeremy L Thompson   // Restore PETSc Vecs and Local to Global
1062*9dbcead6SJeremy L Thompson   PetscCall(VecReadC2P(swarm_ceed_context->u_mesh_loc, U_mem_type, U_mesh_loc));
1063*9dbcead6SJeremy L Thompson   PetscCall(VecC2P(swarm_ceed_context->v_mesh_loc, V_mem_type, V_mesh_loc));
1064*9dbcead6SJeremy L Thompson   PetscCall(VecZeroEntries(V_mesh));
1065*9dbcead6SJeremy L Thompson   PetscCall(DMLocalToGlobal(dm_mesh, V_mesh_loc, ADD_VALUES, V_mesh));
1066*9dbcead6SJeremy L Thompson 
1067*9dbcead6SJeremy L Thompson   // Cleanup
1068*9dbcead6SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm_mesh, &U_mesh_loc));
1069*9dbcead6SJeremy L Thompson   PetscCall(DMRestoreLocalVector(dm_mesh, &V_mesh_loc));
1070*9dbcead6SJeremy L Thompson   PetscFunctionReturn(PETSC_SUCCESS);
1071*9dbcead6SJeremy L Thompson }
1072