xref: /libCEED/examples/petsc/bpsswarm.c (revision 78f7fce354c760f980d0580f63d29ea51c63cedc)
1*78f7fce3SZach Atkins // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2*78f7fce3SZach Atkins // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3*78f7fce3SZach Atkins //
4*78f7fce3SZach Atkins // SPDX-License-Identifier: BSD-2-Clause
5*78f7fce3SZach Atkins //
6*78f7fce3SZach Atkins // This file is part of CEED:  http://github.com/ceed
7*78f7fce3SZach Atkins 
8*78f7fce3SZach Atkins //                        libCEED + PETSc Example: CEED BPs
9*78f7fce3SZach Atkins //
10*78f7fce3SZach Atkins // This example demonstrates a simple usage of libCEED with PETSc to solve the CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps, on
11*78f7fce3SZach Atkins // a particle swarm.
12*78f7fce3SZach Atkins //
13*78f7fce3SZach Atkins // The code uses higher level communication protocols in DMPlex and DMSwarm.
14*78f7fce3SZach Atkins //
15*78f7fce3SZach Atkins // Build with:
16*78f7fce3SZach Atkins //
17*78f7fce3SZach Atkins //     make bpsswarm [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
18*78f7fce3SZach Atkins //
19*78f7fce3SZach Atkins // Sample runs:
20*78f7fce3SZach Atkins //
21*78f7fce3SZach Atkins //     bpssphere -problem bp1 -degree 3
22*78f7fce3SZach Atkins //     bpssphere -problem bp2 -degree 3
23*78f7fce3SZach Atkins //     bpssphere -problem bp3 -degree 3
24*78f7fce3SZach Atkins //
25*78f7fce3SZach Atkins //TESTARGS -ceed {ceed_resource} -test -problem bp3 -dm_plex_dim 3 -dm_plex_box_faces 10,10,10 -dm_plex_simplex 0 -swarm uniform -points_per_cell 500
26*78f7fce3SZach Atkins 
27*78f7fce3SZach Atkins /// @file
28*78f7fce3SZach Atkins /// CEED BPs example using PETSc with DMPlex
29*78f7fce3SZach Atkins /// See bpsraw.c for a "raw" implementation using a structured grid and bps.c for an implementation using an unstructured grid.
30*78f7fce3SZach Atkins static const char help[]              = "Solve CEED BPs on a particle swarm using DMPlex and DMSwarm in PETSc\n";
31*78f7fce3SZach Atkins const char        DMSwarmPICField_u[] = "u";
32*78f7fce3SZach Atkins 
33*78f7fce3SZach Atkins #include "bps.h"
34*78f7fce3SZach Atkins 
35*78f7fce3SZach Atkins #include <ceed.h>
36*78f7fce3SZach Atkins #include <petscdmplex.h>
37*78f7fce3SZach Atkins #include <petscksp.h>
38*78f7fce3SZach Atkins #include <stdbool.h>
39*78f7fce3SZach Atkins #include <string.h>
40*78f7fce3SZach Atkins 
41*78f7fce3SZach Atkins #include "include/bpsproblemdata.h"
42*78f7fce3SZach Atkins #include "include/libceedsetup.h"
43*78f7fce3SZach Atkins #include "include/matops.h"
44*78f7fce3SZach Atkins #include "include/petscutils.h"
45*78f7fce3SZach Atkins #include "include/petscversion.h"
46*78f7fce3SZach Atkins #include "include/swarmutils.h"
47*78f7fce3SZach Atkins 
48*78f7fce3SZach Atkins int main(int argc, char **argv) {
49*78f7fce3SZach Atkins   MPI_Comm             comm;
50*78f7fce3SZach Atkins   char                 ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self", filename[PETSC_MAX_PATH_LEN];
51*78f7fce3SZach Atkins   double               my_rt_start, my_rt, rt_min, rt_max;
52*78f7fce3SZach Atkins   PetscInt             comm_size, degree = 3, q_extra, l_size, g_size, dim = 3, num_comp_u = 1, xl_size, num_points = 1728, num_points_per_cell = 64;
53*78f7fce3SZach Atkins   PetscBool            test_mode, benchmark_mode, read_mesh, write_solution, write_true_solution_swarm;
54*78f7fce3SZach Atkins   PetscLogStage        solve_stage;
55*78f7fce3SZach Atkins   Vec                  X, X_loc, rhs;
56*78f7fce3SZach Atkins   Mat                  mat_O;
57*78f7fce3SZach Atkins   KSP                  ksp;
58*78f7fce3SZach Atkins   DM                   dm_mesh, dm_swarm;
59*78f7fce3SZach Atkins   OperatorApplyContext op_apply_ctx, op_error_ctx;
60*78f7fce3SZach Atkins   Ceed                 ceed;
61*78f7fce3SZach Atkins   CeedData             ceed_data;
62*78f7fce3SZach Atkins   CeedOperator         op_error;
63*78f7fce3SZach Atkins   BPType               bp_choice;
64*78f7fce3SZach Atkins   VecType              vec_type;
65*78f7fce3SZach Atkins   PointSwarmType       point_swarm_type = SWARM_GAUSS;
66*78f7fce3SZach Atkins   PetscMPIInt          ranks_per_node;
67*78f7fce3SZach Atkins   char                 hostname[PETSC_MAX_PATH_LEN];
68*78f7fce3SZach Atkins 
69*78f7fce3SZach Atkins   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
70*78f7fce3SZach Atkins   comm = PETSC_COMM_WORLD;
71*78f7fce3SZach Atkins   PetscCall(MPI_Comm_size(comm, &comm_size));
72*78f7fce3SZach Atkins #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
73*78f7fce3SZach Atkins   {
74*78f7fce3SZach Atkins     MPI_Comm splitcomm;
75*78f7fce3SZach Atkins     PetscCall(MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL, &splitcomm));
76*78f7fce3SZach Atkins     PetscCall(MPI_Comm_size(splitcomm, &ranks_per_node));
77*78f7fce3SZach Atkins     PetscCall(MPI_Comm_free(&splitcomm));
78*78f7fce3SZach Atkins   }
79*78f7fce3SZach Atkins #else
80*78f7fce3SZach Atkins   ranks_per_node = -1;  // Unknown
81*78f7fce3SZach Atkins #endif
82*78f7fce3SZach Atkins 
83*78f7fce3SZach Atkins   // Read command line options
84*78f7fce3SZach Atkins   PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
85*78f7fce3SZach Atkins   bp_choice = CEED_BP1;
86*78f7fce3SZach Atkins   PetscCall(PetscOptionsEnum("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice, NULL));
87*78f7fce3SZach Atkins   num_comp_u = bp_options[bp_choice].num_comp_u;
88*78f7fce3SZach Atkins   test_mode  = PETSC_FALSE;
89*78f7fce3SZach Atkins   PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL));
90*78f7fce3SZach Atkins   benchmark_mode = PETSC_FALSE;
91*78f7fce3SZach Atkins   PetscCall(PetscOptionsBool("-benchmark", "Benchmarking mode (prints benchmark statistics)", NULL, benchmark_mode, &benchmark_mode, NULL));
92*78f7fce3SZach Atkins   write_solution = PETSC_FALSE;
93*78f7fce3SZach Atkins   PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, write_solution, &write_solution, NULL));
94*78f7fce3SZach Atkins   write_true_solution_swarm = PETSC_FALSE;
95*78f7fce3SZach Atkins   PetscCall(PetscOptionsBool("-write_true_solution_swarm", "Write true solution at swarm points for visualization", NULL, write_true_solution_swarm,
96*78f7fce3SZach Atkins                              &write_true_solution_swarm, NULL));
97*78f7fce3SZach Atkins   degree = test_mode ? 3 : 2;
98*78f7fce3SZach Atkins   PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, degree, &degree, NULL));
99*78f7fce3SZach Atkins   q_extra = bp_options[bp_choice].q_extra;
100*78f7fce3SZach Atkins   PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra, NULL));
101*78f7fce3SZach Atkins   PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resource, sizeof(ceed_resource), NULL));
102*78f7fce3SZach Atkins   PetscCall(PetscGetHostName(hostname, sizeof hostname));
103*78f7fce3SZach Atkins   PetscCall(PetscOptionsString("-hostname", "Hostname for output", NULL, hostname, hostname, sizeof(hostname), NULL));
104*78f7fce3SZach Atkins   read_mesh = PETSC_FALSE;
105*78f7fce3SZach Atkins   PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &read_mesh));
106*78f7fce3SZach Atkins   PetscCall(PetscOptionsEnum("-swarm", "Swarm points distribution", NULL, point_swarm_types, (PetscEnum)point_swarm_type,
107*78f7fce3SZach Atkins                              (PetscEnum *)&point_swarm_type, NULL));
108*78f7fce3SZach Atkins   {
109*78f7fce3SZach Atkins     PetscBool user_set_num_points_per_cell = PETSC_FALSE;
110*78f7fce3SZach Atkins     PetscInt  num_cells_total = 1, tmp = dim;
111*78f7fce3SZach Atkins     PetscInt  num_cells[] = {1, 1, 1};
112*78f7fce3SZach Atkins 
113*78f7fce3SZach Atkins     PetscCall(PetscOptionsInt("-points_per_cell", "Total number of swarm points in each cell", NULL, num_points_per_cell, &num_points_per_cell,
114*78f7fce3SZach Atkins                               &user_set_num_points_per_cell));
115*78f7fce3SZach Atkins     PetscCall(PetscOptionsInt("-dm_plex_dim", "Background mesh dimension", NULL, dim, &dim, NULL));
116*78f7fce3SZach Atkins     PetscCall(PetscOptionsIntArray("-dm_plex_box_faces", "Number of cells", NULL, num_cells, &tmp, NULL));
117*78f7fce3SZach Atkins 
118*78f7fce3SZach Atkins     PetscCheck(tmp == dim, comm, PETSC_ERR_USER, "Number of values for -dm_plex_box_faces must match dimension");
119*78f7fce3SZach Atkins 
120*78f7fce3SZach Atkins     num_cells_total = num_cells[0] * num_cells[1] * num_cells[2];
121*78f7fce3SZach Atkins     PetscCheck(!user_set_num_points_per_cell || point_swarm_type != SWARM_SINUSOIDAL, comm, PETSC_ERR_USER,
122*78f7fce3SZach Atkins                "Cannot specify points per cell with sinusoidal points locations");
123*78f7fce3SZach Atkins     if (!user_set_num_points_per_cell) {
124*78f7fce3SZach Atkins       PetscCall(PetscOptionsInt("-points", "Total number of swarm points", NULL, num_points, &num_points, NULL));
125*78f7fce3SZach Atkins       num_points_per_cell = PetscCeilInt(num_points, num_cells_total);
126*78f7fce3SZach Atkins     }
127*78f7fce3SZach Atkins     if (point_swarm_type != SWARM_SINUSOIDAL) {
128*78f7fce3SZach Atkins       PetscInt num_points_per_cell_1d = round(cbrt(num_points_per_cell * 1.0));
129*78f7fce3SZach Atkins 
130*78f7fce3SZach Atkins       num_points_per_cell = 1;
131*78f7fce3SZach Atkins       for (PetscInt i = 0; i < dim; i++) num_points_per_cell *= num_points_per_cell_1d;
132*78f7fce3SZach Atkins     }
133*78f7fce3SZach Atkins     num_points = num_points_per_cell * num_cells_total;
134*78f7fce3SZach Atkins   }
135*78f7fce3SZach Atkins   {
136*78f7fce3SZach Atkins     PetscBool flg;
137*78f7fce3SZach Atkins     PetscInt  p = ranks_per_node;
138*78f7fce3SZach Atkins     PetscCall(PetscOptionsInt("-p", "Number of MPI ranks per node", NULL, p, &p, &flg));
139*78f7fce3SZach Atkins     if (flg) ranks_per_node = p;
140*78f7fce3SZach Atkins   }
141*78f7fce3SZach Atkins   PetscOptionsEnd();
142*78f7fce3SZach Atkins 
143*78f7fce3SZach Atkins   // Setup DM
144*78f7fce3SZach Atkins   if (read_mesh) {
145*78f7fce3SZach Atkins     PetscCall(DMPlexCreateFromFile(comm, filename, NULL, PETSC_TRUE, &dm_mesh));
146*78f7fce3SZach Atkins   } else {
147*78f7fce3SZach Atkins     PetscCall(DMCreate(comm, &dm_mesh));
148*78f7fce3SZach Atkins     PetscCall(DMSetType(dm_mesh, DMPLEX));
149*78f7fce3SZach Atkins     PetscCall(DMSetFromOptions(dm_mesh));
150*78f7fce3SZach Atkins 
151*78f7fce3SZach Atkins     // -- Check for tensor product mesh
152*78f7fce3SZach Atkins     {
153*78f7fce3SZach Atkins       PetscBool is_simplex;
154*78f7fce3SZach Atkins 
155*78f7fce3SZach Atkins       PetscCall(DMPlexIsSimplex(dm_mesh, &is_simplex));
156*78f7fce3SZach Atkins       PetscCheck(!is_simplex, comm, PETSC_ERR_USER, "Only tensor-product background meshes supported");
157*78f7fce3SZach Atkins     }
158*78f7fce3SZach Atkins   }
159*78f7fce3SZach Atkins   PetscCall(DMGetDimension(dm_mesh, &dim));
160*78f7fce3SZach Atkins   PetscCall(SetupDMByDegree(dm_mesh, degree, q_extra, num_comp_u, dim, bp_options[bp_choice].enforce_bc));
161*78f7fce3SZach Atkins 
162*78f7fce3SZach Atkins   // View mesh
163*78f7fce3SZach Atkins   PetscCall(DMSetOptionsPrefix(dm_mesh, "final_"));
164*78f7fce3SZach Atkins   PetscCall(DMViewFromOptions(dm_mesh, NULL, "-dm_view"));
165*78f7fce3SZach Atkins 
166*78f7fce3SZach Atkins   // Create particle swarm
167*78f7fce3SZach Atkins   PetscCall(DMCreate(comm, &dm_swarm));
168*78f7fce3SZach Atkins   PetscCall(DMSetType(dm_swarm, DMSWARM));
169*78f7fce3SZach Atkins   PetscCall(DMSetDimension(dm_swarm, dim));
170*78f7fce3SZach Atkins   PetscCall(DMSwarmSetType(dm_swarm, DMSWARM_PIC));
171*78f7fce3SZach Atkins   PetscCall(DMSwarmSetCellDM(dm_swarm, dm_mesh));
172*78f7fce3SZach Atkins 
173*78f7fce3SZach Atkins   // -- Swarm field
174*78f7fce3SZach Atkins   PetscCall(DMSwarmRegisterPetscDatatypeField(dm_swarm, DMSwarmPICField_u, num_comp_u, PETSC_SCALAR));
175*78f7fce3SZach Atkins   PetscCall(DMSwarmFinalizeFieldRegister(dm_swarm));
176*78f7fce3SZach Atkins   PetscCall(DMSwarmSetLocalSizes(dm_swarm, num_points, 0));
177*78f7fce3SZach Atkins   PetscCall(DMSetFromOptions(dm_swarm));
178*78f7fce3SZach Atkins 
179*78f7fce3SZach Atkins   // -- Set swarm point locations
180*78f7fce3SZach Atkins   PetscCall(DMSwarmInitalizePointLocations(dm_swarm, point_swarm_type, num_points, num_points_per_cell));
181*78f7fce3SZach Atkins   PetscCall(DMSwarmVectorDefineField(dm_swarm, DMSwarmPICField_u));
182*78f7fce3SZach Atkins 
183*78f7fce3SZach Atkins   // -- Final particle swarm
184*78f7fce3SZach Atkins   PetscCall(PetscObjectSetName((PetscObject)dm_swarm, "Particle Swarm"));
185*78f7fce3SZach Atkins   PetscCall(DMViewFromOptions(dm_swarm, NULL, "-dm_swarm_view"));
186*78f7fce3SZach Atkins 
187*78f7fce3SZach Atkins   // Create vectors
188*78f7fce3SZach Atkins   PetscCall(DMCreateGlobalVector(dm_mesh, &X));
189*78f7fce3SZach Atkins   PetscCall(VecGetLocalSize(X, &l_size));
190*78f7fce3SZach Atkins   PetscCall(VecGetSize(X, &g_size));
191*78f7fce3SZach Atkins   PetscCall(DMCreateLocalVector(dm_mesh, &X_loc));
192*78f7fce3SZach Atkins   PetscCall(VecGetSize(X_loc, &xl_size));
193*78f7fce3SZach Atkins   PetscCall(VecDuplicate(X, &rhs));
194*78f7fce3SZach Atkins 
195*78f7fce3SZach Atkins   // Operator
196*78f7fce3SZach Atkins   PetscCall(PetscMalloc1(1, &op_apply_ctx));
197*78f7fce3SZach Atkins   PetscCall(PetscMalloc1(1, &op_error_ctx));
198*78f7fce3SZach Atkins   PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, op_apply_ctx, &mat_O));
199*78f7fce3SZach Atkins   PetscCall(MatSetDM(mat_O, dm_mesh));
200*78f7fce3SZach Atkins   PetscCall(MatShellSetOperation(mat_O, MATOP_MULT, (void (*)(void))MatMult_Ceed));
201*78f7fce3SZach Atkins 
202*78f7fce3SZach Atkins   // Set up libCEED
203*78f7fce3SZach Atkins   CeedInit(ceed_resource, &ceed);
204*78f7fce3SZach Atkins   CeedMemType mem_type_backend;
205*78f7fce3SZach Atkins   CeedGetPreferredMemType(ceed, &mem_type_backend);
206*78f7fce3SZach Atkins 
207*78f7fce3SZach Atkins   PetscCall(DMGetVecType(dm_mesh, &vec_type));
208*78f7fce3SZach Atkins   if (!vec_type) {  // Not yet set by user -dm_vec_type
209*78f7fce3SZach Atkins     switch (mem_type_backend) {
210*78f7fce3SZach Atkins       case CEED_MEM_HOST:
211*78f7fce3SZach Atkins         vec_type = VECSTANDARD;
212*78f7fce3SZach Atkins         break;
213*78f7fce3SZach Atkins       case CEED_MEM_DEVICE: {
214*78f7fce3SZach Atkins         const char *resolved;
215*78f7fce3SZach Atkins         CeedGetResource(ceed, &resolved);
216*78f7fce3SZach Atkins         if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
217*78f7fce3SZach Atkins         else if (strstr(resolved, "/gpu/hip/occa")) vec_type = VECSTANDARD;  // https://github.com/CEED/libCEED/issues/678
218*78f7fce3SZach Atkins         else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
219*78f7fce3SZach Atkins         else vec_type = VECSTANDARD;
220*78f7fce3SZach Atkins       }
221*78f7fce3SZach Atkins     }
222*78f7fce3SZach Atkins     PetscCall(DMSetVecType(dm_mesh, vec_type));
223*78f7fce3SZach Atkins   }
224*78f7fce3SZach Atkins 
225*78f7fce3SZach Atkins   // Print summary
226*78f7fce3SZach Atkins   if (!test_mode) {
227*78f7fce3SZach Atkins     PetscInt P = degree + 1, Q = P + q_extra;
228*78f7fce3SZach Atkins 
229*78f7fce3SZach Atkins     const char *used_resource;
230*78f7fce3SZach Atkins     CeedGetResource(ceed, &used_resource);
231*78f7fce3SZach Atkins 
232*78f7fce3SZach Atkins     VecType vec_type;
233*78f7fce3SZach Atkins     PetscCall(VecGetType(X, &vec_type));
234*78f7fce3SZach Atkins 
235*78f7fce3SZach Atkins     PetscInt c_start, c_end, num_cells_local;
236*78f7fce3SZach Atkins     PetscCall(DMPlexGetHeightStratum(dm_mesh, 0, &c_start, &c_end));
237*78f7fce3SZach Atkins     num_cells_local = c_end - c_start;
238*78f7fce3SZach Atkins     DMPolytopeType cell_type;
239*78f7fce3SZach Atkins     PetscCall(DMPlexGetCellType(dm_mesh, c_start, &cell_type));
240*78f7fce3SZach Atkins     PetscMPIInt comm_size;
241*78f7fce3SZach Atkins     PetscCall(MPI_Comm_size(comm, &comm_size));
242*78f7fce3SZach Atkins 
243*78f7fce3SZach Atkins     PetscInt num_points_local;
244*78f7fce3SZach Atkins     PetscCall(DMSwarmGetLocalSize(dm_swarm, &num_points_local));
245*78f7fce3SZach Atkins 
246*78f7fce3SZach Atkins     PetscCall(PetscPrintf(comm,
247*78f7fce3SZach Atkins                           "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc --\n"
248*78f7fce3SZach Atkins                           "  MPI:\n"
249*78f7fce3SZach Atkins                           "    Hostname                                : %s\n"
250*78f7fce3SZach Atkins                           "    Total ranks                             : %d\n"
251*78f7fce3SZach Atkins                           "    Ranks per compute node                  : %d\n"
252*78f7fce3SZach Atkins                           "  PETSc:\n"
253*78f7fce3SZach Atkins                           "    PETSc Vec Type                          : %s\n"
254*78f7fce3SZach Atkins                           "  libCEED:\n"
255*78f7fce3SZach Atkins                           "    libCEED Backend                         : %s\n"
256*78f7fce3SZach Atkins                           "    libCEED Backend MemType                 : %s\n"
257*78f7fce3SZach Atkins                           "  Mesh:\n"
258*78f7fce3SZach Atkins                           "    Solution Order (P)                      : %" CeedInt_FMT "\n"
259*78f7fce3SZach Atkins                           "    Quadrature  Order (Q)                   : %" CeedInt_FMT "\n"
260*78f7fce3SZach Atkins                           "    Additional quadrature points (q_extra)  : %" CeedInt_FMT "\n"
261*78f7fce3SZach Atkins                           "    Global nodes                            : %" PetscInt_FMT "\n"
262*78f7fce3SZach Atkins                           "    Local Elements                          : %" PetscInt_FMT "\n"
263*78f7fce3SZach Atkins                           "    Owned nodes                             : %" PetscInt_FMT "\n"
264*78f7fce3SZach Atkins                           "    DoF per node                            : %" PetscInt_FMT "\n"
265*78f7fce3SZach Atkins                           "  Swarm:\n"
266*78f7fce3SZach Atkins                           "    Global points                           : %" PetscInt_FMT "\n"
267*78f7fce3SZach Atkins                           "    Local points                            : %" PetscInt_FMT "\n"
268*78f7fce3SZach Atkins                           "    Avg points per cell                     : %" PetscInt_FMT "\n"
269*78f7fce3SZach Atkins                           "    Point distribution                      : %s\n",
270*78f7fce3SZach Atkins                           bp_choice + 1, hostname, comm_size, ranks_per_node, vec_type, used_resource, CeedMemTypes[mem_type_backend], P, Q, q_extra,
271*78f7fce3SZach Atkins                           g_size / num_comp_u, num_cells_local, l_size / num_comp_u, num_comp_u, num_points, num_points_local,
272*78f7fce3SZach Atkins                           num_cells_local > 0 ? num_points_local / num_cells_local : 0, point_swarm_types[point_swarm_type]));
273*78f7fce3SZach Atkins   }
274*78f7fce3SZach Atkins 
275*78f7fce3SZach Atkins   // Setup libCEED's objects
276*78f7fce3SZach Atkins   Vec target;
277*78f7fce3SZach Atkins 
278*78f7fce3SZach Atkins   PetscCall(DMCreateLocalVector(dm_swarm, &target));
279*78f7fce3SZach Atkins   PetscCall(PetscMalloc1(1, &ceed_data));
280*78f7fce3SZach Atkins   PetscCall(SetupProblemSwarm(dm_swarm, ceed, bp_options[bp_choice], ceed_data, true, rhs, target));
281*78f7fce3SZach Atkins   PetscCall(SetupErrorOperator(dm_mesh, ceed, bp_options[bp_choice], dim, dim, num_comp_u, &op_error));
282*78f7fce3SZach Atkins 
283*78f7fce3SZach Atkins   // Set up apply operator context
284*78f7fce3SZach Atkins   PetscCall(SetupApplyOperatorCtx(comm, dm_mesh, ceed, ceed_data, X_loc, op_apply_ctx));
285*78f7fce3SZach Atkins 
286*78f7fce3SZach Atkins   // Setup solver
287*78f7fce3SZach Atkins   PetscCall(KSPCreate(comm, &ksp));
288*78f7fce3SZach Atkins   {
289*78f7fce3SZach Atkins     PC pc;
290*78f7fce3SZach Atkins     PetscCall(KSPGetPC(ksp, &pc));
291*78f7fce3SZach Atkins     if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) {
292*78f7fce3SZach Atkins       PetscCall(PCSetType(pc, PCJACOBI));
293*78f7fce3SZach Atkins       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
294*78f7fce3SZach Atkins     } else {
295*78f7fce3SZach Atkins       PetscCall(PCSetType(pc, PCNONE));
296*78f7fce3SZach Atkins       MatNullSpace nullspace;
297*78f7fce3SZach Atkins 
298*78f7fce3SZach Atkins       PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, 0, &nullspace));
299*78f7fce3SZach Atkins       PetscCall(MatSetNullSpace(mat_O, nullspace));
300*78f7fce3SZach Atkins       PetscCall(MatNullSpaceDestroy(&nullspace));
301*78f7fce3SZach Atkins     }
302*78f7fce3SZach Atkins     PetscCall(KSPSetType(ksp, KSPCG));
303*78f7fce3SZach Atkins     PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
304*78f7fce3SZach Atkins     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
305*78f7fce3SZach Atkins   }
306*78f7fce3SZach Atkins   PetscCall(KSPSetFromOptions(ksp));
307*78f7fce3SZach Atkins   PetscCall(KSPSetOperators(ksp, mat_O, mat_O));
308*78f7fce3SZach Atkins 
309*78f7fce3SZach Atkins   // First run, if benchmarking
310*78f7fce3SZach Atkins   if (benchmark_mode) {
311*78f7fce3SZach Atkins     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1));
312*78f7fce3SZach Atkins     my_rt_start = MPI_Wtime();
313*78f7fce3SZach Atkins     PetscCall(KSPSolve(ksp, rhs, X));
314*78f7fce3SZach Atkins     my_rt = MPI_Wtime() - my_rt_start;
315*78f7fce3SZach Atkins     PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm));
316*78f7fce3SZach Atkins     // Set maxits based on first iteration timing
317*78f7fce3SZach Atkins     if (my_rt > 0.02) {
318*78f7fce3SZach Atkins       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5));
319*78f7fce3SZach Atkins     } else {
320*78f7fce3SZach Atkins       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20));
321*78f7fce3SZach Atkins     }
322*78f7fce3SZach Atkins   }
323*78f7fce3SZach Atkins 
324*78f7fce3SZach Atkins   // Timed solve
325*78f7fce3SZach Atkins   PetscCall(VecZeroEntries(X));
326*78f7fce3SZach Atkins   PetscCall(PetscBarrier((PetscObject)ksp));
327*78f7fce3SZach Atkins 
328*78f7fce3SZach Atkins   // -- Performance logging
329*78f7fce3SZach Atkins   PetscCall(PetscLogStageRegister("Solve Stage", &solve_stage));
330*78f7fce3SZach Atkins   PetscCall(PetscLogStagePush(solve_stage));
331*78f7fce3SZach Atkins 
332*78f7fce3SZach Atkins   // -- Solve
333*78f7fce3SZach Atkins   my_rt_start = MPI_Wtime();
334*78f7fce3SZach Atkins   PetscCall(KSPSolve(ksp, rhs, X));
335*78f7fce3SZach Atkins   my_rt = MPI_Wtime() - my_rt_start;
336*78f7fce3SZach Atkins 
337*78f7fce3SZach Atkins   // -- Performance logging
338*78f7fce3SZach Atkins   PetscCall(PetscLogStagePop());
339*78f7fce3SZach Atkins 
340*78f7fce3SZach Atkins   // Output results
341*78f7fce3SZach Atkins   {
342*78f7fce3SZach Atkins     KSPType            ksp_type;
343*78f7fce3SZach Atkins     KSPConvergedReason reason;
344*78f7fce3SZach Atkins     PetscReal          rnorm;
345*78f7fce3SZach Atkins     PetscInt           its;
346*78f7fce3SZach Atkins     PetscCall(KSPGetType(ksp, &ksp_type));
347*78f7fce3SZach Atkins     PetscCall(KSPGetConvergedReason(ksp, &reason));
348*78f7fce3SZach Atkins     PetscCall(KSPGetIterationNumber(ksp, &its));
349*78f7fce3SZach Atkins     PetscCall(KSPGetResidualNorm(ksp, &rnorm));
350*78f7fce3SZach Atkins     if (!test_mode || reason < 0 || rnorm > 1e-8) {
351*78f7fce3SZach Atkins       PetscCall(PetscPrintf(comm,
352*78f7fce3SZach Atkins                             "  KSP:\n"
353*78f7fce3SZach Atkins                             "    KSP Type                                : %s\n"
354*78f7fce3SZach Atkins                             "    KSP Convergence                         : %s\n"
355*78f7fce3SZach Atkins                             "    Total KSP Iterations                    : %" PetscInt_FMT "\n"
356*78f7fce3SZach Atkins                             "    Final rnorm                             : %e\n",
357*78f7fce3SZach Atkins                             ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
358*78f7fce3SZach Atkins     }
359*78f7fce3SZach Atkins     if (!test_mode) {
360*78f7fce3SZach Atkins       PetscCall(PetscPrintf(comm, "  Performance:\n"));
361*78f7fce3SZach Atkins     }
362*78f7fce3SZach Atkins 
363*78f7fce3SZach Atkins     // View true solution at particles
364*78f7fce3SZach Atkins     if (write_true_solution_swarm) {
365*78f7fce3SZach Atkins       Vec u_swarm, u_swarm_old;
366*78f7fce3SZach Atkins       PetscCall(DMSwarmSortGetAccess(dm_swarm));
367*78f7fce3SZach Atkins       PetscCall(DMSwarmCreateLocalVectorFromField(dm_swarm, DMSwarmPICField_u, &u_swarm));
368*78f7fce3SZach Atkins       PetscCall(VecDuplicate(u_swarm, &u_swarm_old));
369*78f7fce3SZach Atkins       PetscCall(VecCopy(u_swarm, u_swarm_old));
370*78f7fce3SZach Atkins       PetscCall(VecCopy(target, u_swarm));
371*78f7fce3SZach Atkins       PetscCall(DMSwarmDestroyLocalVectorFromField(dm_swarm, DMSwarmPICField_u, &u_swarm));
372*78f7fce3SZach Atkins       PetscCall(DMSwarmSortRestoreAccess(dm_swarm));
373*78f7fce3SZach Atkins 
374*78f7fce3SZach Atkins       PetscCall(DMSwarmViewXDMF(dm_swarm, "swarm.xmf"));
375*78f7fce3SZach Atkins 
376*78f7fce3SZach Atkins       PetscCall(DMSwarmSortGetAccess(dm_swarm));
377*78f7fce3SZach Atkins       PetscCall(DMSwarmCreateLocalVectorFromField(dm_swarm, DMSwarmPICField_u, &u_swarm));
378*78f7fce3SZach Atkins       PetscCall(VecCopy(u_swarm_old, u_swarm));
379*78f7fce3SZach Atkins       PetscCall(DMSwarmDestroyLocalVectorFromField(dm_swarm, DMSwarmPICField_u, &u_swarm));
380*78f7fce3SZach Atkins       PetscCall(DMSwarmSortRestoreAccess(dm_swarm));
381*78f7fce3SZach Atkins       PetscCall(VecDestroy(&u_swarm_old));
382*78f7fce3SZach Atkins     }
383*78f7fce3SZach Atkins 
384*78f7fce3SZach Atkins     // View solution at mesh points
385*78f7fce3SZach Atkins     PetscCall(VecViewFromOptions(X, NULL, "-solution_view"));
386*78f7fce3SZach Atkins 
387*78f7fce3SZach Atkins     // Compute L2 Error
388*78f7fce3SZach Atkins     {
389*78f7fce3SZach Atkins       // Set up error operator context
390*78f7fce3SZach Atkins       PetscCall(SetupErrorOperatorCtx(comm, dm_mesh, ceed, ceed_data, X_loc, op_error, op_error_ctx));
391*78f7fce3SZach Atkins       PetscScalar l2_error;
392*78f7fce3SZach Atkins       PetscCall(ComputeL2Error(X, &l2_error, op_error_ctx));
393*78f7fce3SZach Atkins 
394*78f7fce3SZach Atkins       PetscReal tol = 5e-4;
395*78f7fce3SZach Atkins       if (!test_mode || l2_error > tol) {
396*78f7fce3SZach Atkins         PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm));
397*78f7fce3SZach Atkins         PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm));
398*78f7fce3SZach Atkins         PetscCall(PetscPrintf(comm,
399*78f7fce3SZach Atkins                               "    L2 Error                                : %e\n"
400*78f7fce3SZach Atkins                               "    CG Solve Time                           : %g (%g) sec\n",
401*78f7fce3SZach Atkins                               (double)l2_error, rt_max, rt_min));
402*78f7fce3SZach Atkins       }
403*78f7fce3SZach Atkins     }
404*78f7fce3SZach Atkins     if (benchmark_mode && (!test_mode)) {
405*78f7fce3SZach Atkins       PetscCall(PetscPrintf(comm, "    DoFs/Sec in CG                            : %g (%g) million\n", 1e-6 * g_size * its / rt_max,
406*78f7fce3SZach Atkins                             1e-6 * g_size * its / rt_min));
407*78f7fce3SZach Atkins     }
408*78f7fce3SZach Atkins   }
409*78f7fce3SZach Atkins 
410*78f7fce3SZach Atkins   // Output solution
411*78f7fce3SZach Atkins   if (write_solution) {
412*78f7fce3SZach Atkins     PetscViewer vtk_viewer_soln;
413*78f7fce3SZach Atkins 
414*78f7fce3SZach Atkins     PetscCall(PetscViewerCreate(comm, &vtk_viewer_soln));
415*78f7fce3SZach Atkins     PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK));
416*78f7fce3SZach Atkins     PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"));
417*78f7fce3SZach Atkins     PetscCall(VecView(X, vtk_viewer_soln));
418*78f7fce3SZach Atkins     PetscCall(PetscViewerDestroy(&vtk_viewer_soln));
419*78f7fce3SZach Atkins   }
420*78f7fce3SZach Atkins 
421*78f7fce3SZach Atkins   // Cleanup
422*78f7fce3SZach Atkins   PetscCall(VecDestroy(&X));
423*78f7fce3SZach Atkins   PetscCall(VecDestroy(&X_loc));
424*78f7fce3SZach Atkins   PetscCall(VecDestroy(&op_apply_ctx->Y_loc));
425*78f7fce3SZach Atkins   PetscCall(VecDestroy(&op_error_ctx->Y_loc));
426*78f7fce3SZach Atkins   PetscCall(MatDestroy(&mat_O));
427*78f7fce3SZach Atkins   PetscCall(PetscFree(op_apply_ctx));
428*78f7fce3SZach Atkins   PetscCall(PetscFree(op_error_ctx));
429*78f7fce3SZach Atkins   PetscCall(CeedDataDestroy(0, ceed_data));
430*78f7fce3SZach Atkins   PetscCall(DMDestroy(&dm_mesh));
431*78f7fce3SZach Atkins   PetscCall(DMDestroy(&dm_swarm));
432*78f7fce3SZach Atkins 
433*78f7fce3SZach Atkins   PetscCall(VecDestroy(&rhs));
434*78f7fce3SZach Atkins   PetscCall(VecDestroy(&target));
435*78f7fce3SZach Atkins   PetscCall(KSPDestroy(&ksp));
436*78f7fce3SZach Atkins   CeedOperatorDestroy(&op_error);
437*78f7fce3SZach Atkins   CeedDestroy(&ceed);
438*78f7fce3SZach Atkins   return PetscFinalize();
439*78f7fce3SZach Atkins }
440