xref: /libCEED/examples/petsc/multigrid.c (revision 53b04fa6da4b584296fd192b66fdd6cbd83ea08e)
16c5df90dSjeremylt // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
26c5df90dSjeremylt // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
36c5df90dSjeremylt // reserved. See files LICENSE and NOTICE for details.
46c5df90dSjeremylt //
56c5df90dSjeremylt // This file is part of CEED, a collection of benchmarks, miniapps, software
66c5df90dSjeremylt // libraries and APIs for efficient high-order finite element and spectral
76c5df90dSjeremylt // element discretizations for exascale applications. For more information and
86c5df90dSjeremylt // source code availability see http://github.com/ceed.
96c5df90dSjeremylt //
106c5df90dSjeremylt // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
116c5df90dSjeremylt // a collaborative effort of two U.S. Department of Energy organizations (Office
126c5df90dSjeremylt // of Science and the National Nuclear Security Administration) responsible for
136c5df90dSjeremylt // the planning and preparation of a capable exascale ecosystem, including
146c5df90dSjeremylt // software, applications, hardware, advanced system engineering and early
156c5df90dSjeremylt // testbed platforms, in support of the nation's exascale computing imperative.
166c5df90dSjeremylt 
176c5df90dSjeremylt //                        libCEED + PETSc Example: CEED BPs 3-6 with Multigrid
186c5df90dSjeremylt //
196c5df90dSjeremylt // This example demonstrates a simple usage of libCEED with PETSc to solve the
206c5df90dSjeremylt // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
216c5df90dSjeremylt //
226c5df90dSjeremylt // The code uses higher level communication protocols in DMPlex.
236c5df90dSjeremylt //
246c5df90dSjeremylt // Build with:
256c5df90dSjeremylt //
266c5df90dSjeremylt //     make multigrid [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
276c5df90dSjeremylt //
286c5df90dSjeremylt // Sample runs:
296c5df90dSjeremylt //
306c5df90dSjeremylt //     multigrid -problem bp3
3128688798Sjeremylt //     multigrid -problem bp4
3228688798Sjeremylt //     multigrid -problem bp5 -ceed /cpu/self
336c5df90dSjeremylt //     multigrid -problem bp6 -ceed /gpu/cuda
346c5df90dSjeremylt //
356c5df90dSjeremylt //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3
366c5df90dSjeremylt 
376c5df90dSjeremylt /// @file
386c5df90dSjeremylt /// CEED BPs 1-6 multigrid example using PETSc
396c5df90dSjeremylt const char help[] = "Solve CEED BPs using p-multigrid with PETSc and DMPlex\n";
406c5df90dSjeremylt 
41636cccdbSjeremylt #include <stdbool.h>
42636cccdbSjeremylt #include <string.h>
43636cccdbSjeremylt #include <ceed.h>
44636cccdbSjeremylt #include <petsc.h>
45636cccdbSjeremylt #include <petscdmplex.h>
46636cccdbSjeremylt #include <petscksp.h>
47636cccdbSjeremylt #include <petscsys.h>
48636cccdbSjeremylt 
49e83e87a5Sjeremylt #include "bps.h"
50636cccdbSjeremylt #include "include/bpsproblemdata.h"
51636cccdbSjeremylt #include "include/petscmacros.h"
52636cccdbSjeremylt #include "include/petscutils.h"
53636cccdbSjeremylt #include "include/matops.h"
54636cccdbSjeremylt #include "include/structs.h"
55636cccdbSjeremylt #include "include/libceedsetup.h"
56636cccdbSjeremylt 
57636cccdbSjeremylt #if PETSC_VERSION_LT(3,12,0)
58636cccdbSjeremylt #ifdef PETSC_HAVE_CUDA
59636cccdbSjeremylt #include <petsccuda.h>
60636cccdbSjeremylt // Note: With PETSc prior to version 3.12.0, providing the source path to
61636cccdbSjeremylt //       include 'cublas_v2.h' will be needed to use 'petsccuda.h'.
62636cccdbSjeremylt #endif
63636cccdbSjeremylt #endif
646c5df90dSjeremylt 
656c5df90dSjeremylt int main(int argc, char **argv) {
666c5df90dSjeremylt   PetscInt ierr;
676c5df90dSjeremylt   MPI_Comm comm;
68cb0b5415Sjeremylt   char filename[PETSC_MAX_PATH_LEN],
699b072555Sjeremylt        ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self";
706c5df90dSjeremylt   double my_rt_start, my_rt, rt_min, rt_max;
719b072555Sjeremylt   PetscInt degree = 3, q_extra, *l_size, *xl_size, *g_size, dim = 3, fine_level,
729b072555Sjeremylt            mesh_elem[3] = {3, 3, 3}, num_comp_u = 1, num_levels = degree, *level_degrees;
736c5df90dSjeremylt   PetscScalar *r;
74cfa59c5bSRey   PetscScalar eps = 1.0;
756c5df90dSjeremylt   PetscBool test_mode, benchmark_mode, read_mesh, write_solution;
769b072555Sjeremylt   PetscLogStage solve_stage;
779b072555Sjeremylt   DM  *dm, dm_orig;
786c5df90dSjeremylt   KSP ksp;
796c5df90dSjeremylt   PC pc;
809b072555Sjeremylt   Mat *mat_O, *mat_pr, mat_coarse;
819b072555Sjeremylt   Vec *X, *X_loc, *mult, rhs, rhs_loc;
829b072555Sjeremylt   PetscMemType mem_type;
839b072555Sjeremylt   UserO *user_O;
849b072555Sjeremylt   UserProlongRestr *user_pr;
856c5df90dSjeremylt   Ceed ceed;
869b072555Sjeremylt   CeedData *ceed_data;
879b072555Sjeremylt   CeedVector rhs_ceed, target;
889b072555Sjeremylt   CeedQFunction qf_error, qf_restrict, qf_prolong;
899b072555Sjeremylt   CeedOperator op_error;
909b072555Sjeremylt   BPType bp_choice;
919b072555Sjeremylt   CoarsenType coarsen;
926c5df90dSjeremylt 
936c5df90dSjeremylt   ierr = PetscInitialize(&argc, &argv, NULL, help);
946c5df90dSjeremylt   if (ierr) return ierr;
956c5df90dSjeremylt   comm = PETSC_COMM_WORLD;
966c5df90dSjeremylt 
976c5df90dSjeremylt   // Parse command line options
986c5df90dSjeremylt   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
999b072555Sjeremylt   bp_choice = CEED_BP3;
1006c5df90dSjeremylt   ierr = PetscOptionsEnum("-problem",
1016c5df90dSjeremylt                           "CEED benchmark problem to solve", NULL,
1029b072555Sjeremylt                           bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice,
1036c5df90dSjeremylt                           NULL); CHKERRQ(ierr);
1049b072555Sjeremylt   num_comp_u = bp_options[bp_choice].num_comp_u;
1056c5df90dSjeremylt   test_mode = PETSC_FALSE;
1066c5df90dSjeremylt   ierr = PetscOptionsBool("-test",
1076c5df90dSjeremylt                           "Testing mode (do not print unless error is large)",
1086c5df90dSjeremylt                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
1096c5df90dSjeremylt   benchmark_mode = PETSC_FALSE;
1106c5df90dSjeremylt   ierr = PetscOptionsBool("-benchmark",
1116c5df90dSjeremylt                           "Benchmarking mode (prints benchmark statistics)",
1126c5df90dSjeremylt                           NULL, benchmark_mode, &benchmark_mode, NULL);
1136c5df90dSjeremylt   CHKERRQ(ierr);
1146c5df90dSjeremylt   write_solution = PETSC_FALSE;
1156c5df90dSjeremylt   ierr = PetscOptionsBool("-write_solution",
1166c5df90dSjeremylt                           "Write solution for visualization",
1176c5df90dSjeremylt                           NULL, write_solution, &write_solution, NULL);
1186c5df90dSjeremylt   CHKERRQ(ierr);
119cfa59c5bSRey   ierr = PetscOptionsScalar("-eps",
120cfa59c5bSRey                             "Epsilon parameter for Kershaw mesh transformation",
121cfa59c5bSRey                             NULL, eps, &eps, NULL);
122cfa59c5bSRey   if (eps > 1 || eps <= 0) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
123cfa59c5bSRey                                       "-eps %D must be (0,1]", eps);
1246c5df90dSjeremylt   degree = test_mode ? 3 : 2;
1256c5df90dSjeremylt   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
1266c5df90dSjeremylt                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
1276c5df90dSjeremylt   if (degree < 1) SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE,
1286c5df90dSjeremylt                              "-degree %D must be at least 1", degree);
1299b072555Sjeremylt   q_extra = bp_options[bp_choice].q_extra;
1309b072555Sjeremylt   ierr = PetscOptionsInt("-q_extra", "Number of extra quadrature points",
1319b072555Sjeremylt                          NULL, q_extra, &q_extra, NULL); CHKERRQ(ierr);
1326c5df90dSjeremylt   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
1339b072555Sjeremylt                             NULL, ceed_resource, ceed_resource,
1349b072555Sjeremylt                             sizeof(ceed_resource), NULL); CHKERRQ(ierr);
1356c5df90dSjeremylt   coarsen = COARSEN_UNIFORM;
1366c5df90dSjeremylt   ierr = PetscOptionsEnum("-coarsen",
1376c5df90dSjeremylt                           "Coarsening strategy to use", NULL,
1389b072555Sjeremylt                           coarsen_types, (PetscEnum)coarsen,
1396c5df90dSjeremylt                           (PetscEnum *)&coarsen, NULL); CHKERRQ(ierr);
140cb32e2e7SValeria Barra   read_mesh = PETSC_FALSE;
1416c5df90dSjeremylt   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
1426c5df90dSjeremylt                             filename, filename, sizeof(filename), &read_mesh);
1436c5df90dSjeremylt   CHKERRQ(ierr);
1446c5df90dSjeremylt   if (!read_mesh) {
1456c5df90dSjeremylt     PetscInt tmp = dim;
1466c5df90dSjeremylt     ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL,
1479b072555Sjeremylt                                 mesh_elem, &tmp, NULL); CHKERRQ(ierr);
1486c5df90dSjeremylt   }
1496c5df90dSjeremylt   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
1506c5df90dSjeremylt 
1519396343dSjeremylt   // Set up libCEED
1529b072555Sjeremylt   CeedInit(ceed_resource, &ceed);
1539b072555Sjeremylt   CeedMemType mem_type_backend;
1549b072555Sjeremylt   CeedGetPreferredMemType(ceed, &mem_type_backend);
1559396343dSjeremylt 
1566c5df90dSjeremylt   // Setup DM
1576c5df90dSjeremylt   if (read_mesh) {
1587ed3e4cdSJeremy L Thompson     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, NULL, PETSC_TRUE,
1597ed3e4cdSJeremy L Thompson                                 &dm_orig);
1606c5df90dSjeremylt     CHKERRQ(ierr);
1616c5df90dSjeremylt   } else {
1629b072555Sjeremylt     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, mesh_elem, NULL,
1639b072555Sjeremylt                                NULL, NULL, PETSC_TRUE, &dm_orig); CHKERRQ(ierr);
1646c5df90dSjeremylt   }
1656c5df90dSjeremylt 
1666c5df90dSjeremylt   {
1679b072555Sjeremylt     DM dm_dist = NULL;
1686c5df90dSjeremylt     PetscPartitioner part;
1696c5df90dSjeremylt 
1709b072555Sjeremylt     ierr = DMPlexGetPartitioner(dm_orig, &part); CHKERRQ(ierr);
1716c5df90dSjeremylt     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
1729b072555Sjeremylt     ierr = DMPlexDistribute(dm_orig, 0, NULL, &dm_dist); CHKERRQ(ierr);
1739b072555Sjeremylt     if (dm_dist) {
1749b072555Sjeremylt       ierr = DMDestroy(&dm_orig); CHKERRQ(ierr);
1759b072555Sjeremylt       dm_orig = dm_dist;
1766c5df90dSjeremylt     }
1776c5df90dSjeremylt   }
1786c5df90dSjeremylt 
1799b072555Sjeremylt   // Apply Kershaw mesh transformation
1809b072555Sjeremylt   ierr = Kershaw(dm_orig, eps); CHKERRQ(ierr);
181cfa59c5bSRey 
1829b072555Sjeremylt   VecType vec_type;
1839b072555Sjeremylt   switch (mem_type_backend) {
1849b072555Sjeremylt   case CEED_MEM_HOST: vec_type = VECSTANDARD; break;
185b68a8d79SJed Brown   case CEED_MEM_DEVICE: {
186b68a8d79SJed Brown     const char *resolved;
187b68a8d79SJed Brown     CeedGetResource(ceed, &resolved);
1889b072555Sjeremylt     if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
189ca2d516cSJed Brown     else if (strstr(resolved, "/gpu/hip/occa"))
1909b072555Sjeremylt       vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678
1919b072555Sjeremylt     else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
1929b072555Sjeremylt     else vec_type = VECSTANDARD;
193b68a8d79SJed Brown   }
194b68a8d79SJed Brown   }
1959b072555Sjeremylt   ierr = DMSetVecType(dm_orig, vec_type); CHKERRQ(ierr);
1969b072555Sjeremylt   ierr = DMSetFromOptions(dm_orig); CHKERRQ(ierr);
197b68a8d79SJed Brown 
1986c5df90dSjeremylt   // Allocate arrays for PETSc objects for each level
1996c5df90dSjeremylt   switch (coarsen) {
2006c5df90dSjeremylt   case COARSEN_UNIFORM:
2019b072555Sjeremylt     num_levels = degree;
2026c5df90dSjeremylt     break;
203dc7d240cSValeria Barra   case COARSEN_LOGARITHMIC:
2049b072555Sjeremylt     num_levels = ceil(log(degree)/log(2)) + 1;
2056c5df90dSjeremylt     break;
2066c5df90dSjeremylt   }
2079b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &level_degrees); CHKERRQ(ierr);
2089b072555Sjeremylt   fine_level = num_levels - 1;
20961608365Sjeremylt 
2106c5df90dSjeremylt   switch (coarsen) {
2116c5df90dSjeremylt   case COARSEN_UNIFORM:
2129b072555Sjeremylt     for (int i=0; i<num_levels; i++) level_degrees[i] = i + 1;
2136c5df90dSjeremylt     break;
214dc7d240cSValeria Barra   case COARSEN_LOGARITHMIC:
2159b072555Sjeremylt     for (int i=0; i<num_levels - 1; i++) level_degrees[i] = pow(2,i);
2169b072555Sjeremylt     level_degrees[fine_level] = degree;
2176c5df90dSjeremylt     break;
2186c5df90dSjeremylt   }
2199b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &dm); CHKERRQ(ierr);
2209b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &X); CHKERRQ(ierr);
2219b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &X_loc); CHKERRQ(ierr);
2229b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &mult); CHKERRQ(ierr);
2239b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &user_O); CHKERRQ(ierr);
2249b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &user_pr); CHKERRQ(ierr);
2259b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &mat_O); CHKERRQ(ierr);
2269b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &mat_pr); CHKERRQ(ierr);
2279b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &l_size); CHKERRQ(ierr);
2289b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &xl_size); CHKERRQ(ierr);
2299b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &g_size); CHKERRQ(ierr);
2306c5df90dSjeremylt 
2316c5df90dSjeremylt   // Setup DM and Operator Mat Shells for each level
2329b072555Sjeremylt   for (CeedInt i=0; i<num_levels; i++) {
2336c5df90dSjeremylt     // Create DM
2349b072555Sjeremylt     ierr = DMClone(dm_orig, &dm[i]); CHKERRQ(ierr);
2359b072555Sjeremylt     ierr = DMGetVecType(dm_orig, &vec_type); CHKERRQ(ierr);
2369b072555Sjeremylt     ierr = DMSetVecType(dm[i], vec_type); CHKERRQ(ierr);
237e83e87a5Sjeremylt     PetscInt dim;
238e83e87a5Sjeremylt     ierr = DMGetDimension(dm[i], &dim); CHKERRQ(ierr);
2399b072555Sjeremylt     ierr = SetupDMByDegree(dm[i], level_degrees[i], num_comp_u, dim,
2409b072555Sjeremylt                            bp_options[bp_choice].enforce_bc, bp_options[bp_choice].bc_func);
2416c5df90dSjeremylt     CHKERRQ(ierr);
2426c5df90dSjeremylt 
2436c5df90dSjeremylt     // Create vectors
2446c5df90dSjeremylt     ierr = DMCreateGlobalVector(dm[i], &X[i]); CHKERRQ(ierr);
2459b072555Sjeremylt     ierr = VecGetLocalSize(X[i], &l_size[i]); CHKERRQ(ierr);
2469b072555Sjeremylt     ierr = VecGetSize(X[i], &g_size[i]); CHKERRQ(ierr);
2479b072555Sjeremylt     ierr = DMCreateLocalVector(dm[i], &X_loc[i]); CHKERRQ(ierr);
2489b072555Sjeremylt     ierr = VecGetSize(X_loc[i], &xl_size[i]); CHKERRQ(ierr);
2496c5df90dSjeremylt 
2506c5df90dSjeremylt     // Operator
2519b072555Sjeremylt     ierr = PetscMalloc1(1, &user_O[i]); CHKERRQ(ierr);
2529b072555Sjeremylt     ierr = MatCreateShell(comm, l_size[i], l_size[i], g_size[i], g_size[i],
2539b072555Sjeremylt                           user_O[i], &mat_O[i]); CHKERRQ(ierr);
2549b072555Sjeremylt     ierr = MatShellSetOperation(mat_O[i], MATOP_MULT,
255ce74dcefSjeremylt                                 (void(*)(void))MatMult_Ceed); CHKERRQ(ierr);
2569b072555Sjeremylt     ierr = MatShellSetOperation(mat_O[i], MATOP_GET_DIAGONAL,
257ce74dcefSjeremylt                                 (void(*)(void))MatGetDiag); CHKERRQ(ierr);
2589b072555Sjeremylt     ierr = MatShellSetVecType(mat_O[i], vec_type); CHKERRQ(ierr);
2596c5df90dSjeremylt 
2606c5df90dSjeremylt     // Level transfers
2616c5df90dSjeremylt     if (i > 0) {
2626c5df90dSjeremylt       // Interp
2639b072555Sjeremylt       ierr = PetscMalloc1(1, &user_pr[i]); CHKERRQ(ierr);
2649b072555Sjeremylt       ierr = MatCreateShell(comm, l_size[i], l_size[i-1], g_size[i], g_size[i-1],
2659b072555Sjeremylt                             user_pr[i], &mat_pr[i]); CHKERRQ(ierr);
2669b072555Sjeremylt       ierr = MatShellSetOperation(mat_pr[i], MATOP_MULT,
267a97643b0Sjeremylt                                   (void(*)(void))MatMult_Prolong);
2686c5df90dSjeremylt       CHKERRQ(ierr);
2699b072555Sjeremylt       ierr = MatShellSetOperation(mat_pr[i], MATOP_MULT_TRANSPOSE,
2706c5df90dSjeremylt                                   (void(*)(void))MatMult_Restrict);
2716c5df90dSjeremylt       CHKERRQ(ierr);
2729b072555Sjeremylt       ierr = MatShellSetVecType(mat_pr[i], vec_type); CHKERRQ(ierr);
2736c5df90dSjeremylt     }
2746c5df90dSjeremylt   }
2759b072555Sjeremylt   ierr = VecDuplicate(X[fine_level], &rhs); CHKERRQ(ierr);
2766c5df90dSjeremylt 
2776c5df90dSjeremylt   // Print global grid information
2786c5df90dSjeremylt   if (!test_mode) {
2799b072555Sjeremylt     PetscInt P = degree + 1, Q = P + q_extra;
2809396343dSjeremylt 
2819b072555Sjeremylt     const char *used_resource;
2829b072555Sjeremylt     CeedGetResource(ceed, &used_resource);
2839396343dSjeremylt 
2849b072555Sjeremylt     ierr = VecGetType(X[0], &vec_type); CHKERRQ(ierr);
2859396343dSjeremylt 
2866c5df90dSjeremylt     ierr = PetscPrintf(comm,
2876c5df90dSjeremylt                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc + PCMG --\n"
2889396343dSjeremylt                        "  PETSc:\n"
2899396343dSjeremylt                        "    PETSc Vec Type                     : %s\n"
2906c5df90dSjeremylt                        "  libCEED:\n"
2916c5df90dSjeremylt                        "    libCEED Backend                    : %s\n"
2929396343dSjeremylt                        "    libCEED Backend MemType            : %s\n"
2936c5df90dSjeremylt                        "  Mesh:\n"
2946c5df90dSjeremylt                        "    Number of 1D Basis Nodes (p)       : %d\n"
2956c5df90dSjeremylt                        "    Number of 1D Quadrature Points (q) : %d\n"
2966c5df90dSjeremylt                        "    Global Nodes                       : %D\n"
2976c5df90dSjeremylt                        "    Owned Nodes                        : %D\n"
298db419314Sjeremylt                        "    DoF per node                       : %D\n"
2996c5df90dSjeremylt                        "  Multigrid:\n"
3006c5df90dSjeremylt                        "    Number of Levels                   : %d\n",
3019b072555Sjeremylt                        bp_choice+1, vec_type, used_resource,
3029b072555Sjeremylt                        CeedMemTypes[mem_type_backend],
3039b072555Sjeremylt                        P, Q, g_size[fine_level]/num_comp_u, l_size[fine_level]/num_comp_u,
3049b072555Sjeremylt                        num_comp_u, num_levels); CHKERRQ(ierr);
3056c5df90dSjeremylt   }
3066c5df90dSjeremylt 
3076c5df90dSjeremylt   // Create RHS vector
3089b072555Sjeremylt   ierr = VecDuplicate(X_loc[fine_level], &rhs_loc); CHKERRQ(ierr);
3099b072555Sjeremylt   ierr = VecZeroEntries(rhs_loc); CHKERRQ(ierr);
3109b072555Sjeremylt   ierr = VecGetArrayAndMemType(rhs_loc, &r, &mem_type); CHKERRQ(ierr);
3119b072555Sjeremylt   CeedVectorCreate(ceed, xl_size[fine_level], &rhs_ceed);
3129b072555Sjeremylt   CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r);
3136c5df90dSjeremylt 
3146c5df90dSjeremylt   // Set up libCEED operators on each level
3159b072555Sjeremylt   ierr = PetscMalloc1(num_levels, &ceed_data); CHKERRQ(ierr);
3169b072555Sjeremylt   for (int i=0; i<num_levels; i++) {
3176c5df90dSjeremylt     // Print level information
3189b072555Sjeremylt     if (!test_mode && (i == 0 || i == fine_level)) {
3196c5df90dSjeremylt       ierr = PetscPrintf(comm,"    Level %D (%s):\n"
3206c5df90dSjeremylt                          "      Number of 1D Basis Nodes (p)     : %d\n"
3216c5df90dSjeremylt                          "      Global Nodes                     : %D\n"
3226c5df90dSjeremylt                          "      Owned Nodes                      : %D\n",
3239b072555Sjeremylt                          i, (i? "fine" : "coarse"), level_degrees[i] + 1,
3249b072555Sjeremylt                          g_size[i]/num_comp_u, l_size[i]/num_comp_u); CHKERRQ(ierr);
3256c5df90dSjeremylt     }
3269b072555Sjeremylt     ierr = PetscMalloc1(1, &ceed_data[i]); CHKERRQ(ierr);
3279b072555Sjeremylt     ierr = SetupLibceedByDegree(dm[i], ceed, level_degrees[i], dim, q_extra,
3289b072555Sjeremylt                                 dim, num_comp_u, g_size[i], xl_size[i], bp_options[bp_choice],
3299b072555Sjeremylt                                 ceed_data[i], i==(fine_level), rhs_ceed, &target);
330f9342789SJeremy L Thompson     CHKERRQ(ierr);
3316c5df90dSjeremylt   }
3326c5df90dSjeremylt 
3336c5df90dSjeremylt   // Gather RHS
3349b072555Sjeremylt   CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL);
3359b072555Sjeremylt   ierr = VecRestoreArrayAndMemType(rhs_loc, &r); CHKERRQ(ierr);
3366c5df90dSjeremylt   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
3379b072555Sjeremylt   ierr = DMLocalToGlobal(dm[fine_level], rhs_loc, ADD_VALUES, rhs); CHKERRQ(ierr);
3389b072555Sjeremylt   CeedVectorDestroy(&rhs_ceed);
3396c5df90dSjeremylt 
34043eb8658SJeremy L Thompson   // Create the restriction/interpolation QFunction
3419b072555Sjeremylt   CeedQFunctionCreateIdentity(ceed, num_comp_u, CEED_EVAL_NONE, CEED_EVAL_INTERP,
3429b072555Sjeremylt                               &qf_restrict);
3439b072555Sjeremylt   CeedQFunctionCreateIdentity(ceed, num_comp_u, CEED_EVAL_INTERP, CEED_EVAL_NONE,
3449b072555Sjeremylt                               &qf_prolong);
3456c5df90dSjeremylt 
3466c5df90dSjeremylt   // Set up libCEED level transfer operators
3479b072555Sjeremylt   ierr = CeedLevelTransferSetup(ceed, num_levels, num_comp_u, ceed_data,
3489b072555Sjeremylt                                 level_degrees,
3499b072555Sjeremylt                                 qf_restrict, qf_prolong); CHKERRQ(ierr);
3506c5df90dSjeremylt 
35143eb8658SJeremy L Thompson   // Create the error QFunction
3529b072555Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error,
3539b072555Sjeremylt                               bp_options[bp_choice].error_loc, &qf_error);
3549b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
3559b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
3569b072555Sjeremylt   CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_NONE);
3576c5df90dSjeremylt 
3586c5df90dSjeremylt   // Create the error operator
3599b072555Sjeremylt   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
3609b072555Sjeremylt                      &op_error);
3619b072555Sjeremylt   CeedOperatorSetField(op_error, "u", ceed_data[fine_level]->elem_restr_u,
3629b072555Sjeremylt                        ceed_data[fine_level]->basis_u, CEED_VECTOR_ACTIVE);
3639b072555Sjeremylt   CeedOperatorSetField(op_error, "true_soln",
3649b072555Sjeremylt                        ceed_data[fine_level]->elem_restr_u_i,
365a8d32208Sjeremylt                        CEED_BASIS_COLLOCATED, target);
3669b072555Sjeremylt   CeedOperatorSetField(op_error, "error", ceed_data[fine_level]->elem_restr_u_i,
367a8d32208Sjeremylt                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
3686c5df90dSjeremylt 
3696c5df90dSjeremylt   // Calculate multiplicity
3709b072555Sjeremylt   for (int i=0; i<num_levels; i++) {
3716c5df90dSjeremylt     PetscScalar *x;
3726c5df90dSjeremylt 
3736c5df90dSjeremylt     // CEED vector
3749b072555Sjeremylt     ierr = VecZeroEntries(X_loc[i]); CHKERRQ(ierr);
3759b072555Sjeremylt     ierr = VecGetArray(X_loc[i], &x); CHKERRQ(ierr);
3769b072555Sjeremylt     CeedVectorSetArray(ceed_data[i]->x_ceed, CEED_MEM_HOST, CEED_USE_POINTER, x);
3776c5df90dSjeremylt 
3786c5df90dSjeremylt     // Multiplicity
3799b072555Sjeremylt     CeedElemRestrictionGetMultiplicity(ceed_data[i]->elem_restr_u,
3809b072555Sjeremylt                                        ceed_data[i]->x_ceed);
3819b072555Sjeremylt     CeedVectorSyncArray(ceed_data[i]->x_ceed, CEED_MEM_HOST);
3826c5df90dSjeremylt 
3836c5df90dSjeremylt     // Restore vector
3849b072555Sjeremylt     ierr = VecRestoreArray(X_loc[i], &x); CHKERRQ(ierr);
3856c5df90dSjeremylt 
3866c5df90dSjeremylt     // Creat mult vector
3879b072555Sjeremylt     ierr = VecDuplicate(X_loc[i], &mult[i]); CHKERRQ(ierr);
3886c5df90dSjeremylt 
3896c5df90dSjeremylt     // Local-to-global
3906c5df90dSjeremylt     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
3919b072555Sjeremylt     ierr = DMLocalToGlobal(dm[i], X_loc[i], ADD_VALUES, X[i]);
3926c5df90dSjeremylt     CHKERRQ(ierr);
3939b072555Sjeremylt     ierr = VecZeroEntries(X_loc[i]); CHKERRQ(ierr);
3946c5df90dSjeremylt 
3956c5df90dSjeremylt     // Global-to-local
396483f8b0dSjeremylt     ierr = DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]);
3976c5df90dSjeremylt     CHKERRQ(ierr);
3986c5df90dSjeremylt     ierr = VecZeroEntries(X[i]); CHKERRQ(ierr);
3996c5df90dSjeremylt 
4006c5df90dSjeremylt     // Multiplicity scaling
4016c5df90dSjeremylt     ierr = VecReciprocal(mult[i]);
4026c5df90dSjeremylt   }
4036c5df90dSjeremylt 
4046c5df90dSjeremylt   // Set up Mat
4059b072555Sjeremylt   for (int i=0; i<num_levels; i++) {
406226c3a8fSjeremylt     // User Operator
4079b072555Sjeremylt     user_O[i]->comm = comm;
4089b072555Sjeremylt     user_O[i]->dm = dm[i];
4099b072555Sjeremylt     user_O[i]->X_loc = X_loc[i];
4109b072555Sjeremylt     ierr = VecDuplicate(X_loc[i], &user_O[i]->Y_loc); CHKERRQ(ierr);
4119b072555Sjeremylt     user_O[i]->x_ceed = ceed_data[i]->x_ceed;
4129b072555Sjeremylt     user_O[i]->y_ceed = ceed_data[i]->y_ceed;
4139b072555Sjeremylt     user_O[i]->op = ceed_data[i]->op_apply;
4149b072555Sjeremylt     user_O[i]->ceed = ceed;
4156c5df90dSjeremylt 
4166c5df90dSjeremylt     if (i > 0) {
417a97643b0Sjeremylt       // Prolongation/Restriction Operator
4189b072555Sjeremylt       user_pr[i]->comm = comm;
4199b072555Sjeremylt       user_pr[i]->dmf = dm[i];
4209b072555Sjeremylt       user_pr[i]->dmc = dm[i-1];
4219b072555Sjeremylt       user_pr[i]->loc_vec_c = X_loc[i-1];
4229b072555Sjeremylt       user_pr[i]->loc_vec_f = user_O[i]->Y_loc;
4239b072555Sjeremylt       user_pr[i]->mult_vec = mult[i];
4249b072555Sjeremylt       user_pr[i]->ceed_vec_c = user_O[i-1]->x_ceed;
4259b072555Sjeremylt       user_pr[i]->ceed_vec_f = user_O[i]->y_ceed;
4269b072555Sjeremylt       user_pr[i]->op_prolong = ceed_data[i]->op_prolong;
4279b072555Sjeremylt       user_pr[i]->op_restrict = ceed_data[i]->op_restrict;
4289b072555Sjeremylt       user_pr[i]->ceed = ceed;
4296c5df90dSjeremylt     }
4306c5df90dSjeremylt   }
4316c5df90dSjeremylt 
432*53b04fa6SJed Brown   // Assemble coarse grid Jacobian for AMG (or other sparse matrix) solve
433*53b04fa6SJed Brown   ierr = DMCreateMatrix(dm[0], &mat_coarse); CHKERRQ(ierr);
434*53b04fa6SJed Brown 
435*53b04fa6SJed Brown   if (PETSC_VERSION_GE(3,16,0)) { // Assemble matrix analytically
436*53b04fa6SJed Brown     CeedInt num_entries, *rows, *cols;
437*53b04fa6SJed Brown     CeedVector coo_values;
438*53b04fa6SJed Brown     CeedOperatorLinearAssembleSymbolic(user_O[0]->op, &num_entries, &rows, &cols);
439*53b04fa6SJed Brown     ISLocalToGlobalMapping ltog_row, ltog_col;
440*53b04fa6SJed Brown     ierr = MatGetLocalToGlobalMapping(mat_coarse, &ltog_row, &ltog_col);
441*53b04fa6SJed Brown     CHKERRQ(ierr);
442*53b04fa6SJed Brown     ierr = ISLocalToGlobalMappingApply(ltog_row, num_entries, rows, rows);
443*53b04fa6SJed Brown     CHKERRQ(ierr);
444*53b04fa6SJed Brown     ierr = ISLocalToGlobalMappingApply(ltog_col, num_entries, cols, cols);
445*53b04fa6SJed Brown     CHKERRQ(ierr);
446*53b04fa6SJed Brown     ierr = MatSetPreallocationCOO(mat_coarse, num_entries, rows, cols);
447*53b04fa6SJed Brown     CHKERRQ(ierr);
448*53b04fa6SJed Brown     free(rows);
449*53b04fa6SJed Brown     free(cols);
450*53b04fa6SJed Brown     CeedVectorCreate(ceed, num_entries, &coo_values);
451*53b04fa6SJed Brown     CeedOperatorLinearAssemble(user_O[0]->op, coo_values);
452*53b04fa6SJed Brown     const CeedScalar *values;
453*53b04fa6SJed Brown     CeedVectorGetArrayRead(coo_values, CEED_MEM_HOST, &values);
454*53b04fa6SJed Brown     ierr = MatSetValuesCOO(mat_coarse, values, ADD_VALUES); CHKERRQ(ierr);
455*53b04fa6SJed Brown     CeedVectorRestoreArrayRead(coo_values, &values);
456*53b04fa6SJed Brown     CeedVectorDestroy(&coo_values);
457*53b04fa6SJed Brown   } else { // Assemble matrix using coloring. Our interfaces match SNES so it's convenient
458*53b04fa6SJed Brown     SNES snes_dummy;
459*53b04fa6SJed Brown 
46015ce0ef0Sjeremylt     // Setup dummy SNES for AMG coarse solve
4619b072555Sjeremylt     ierr = SNESCreate(comm, &snes_dummy); CHKERRQ(ierr);
4629b072555Sjeremylt     ierr = SNESSetDM(snes_dummy, dm[0]); CHKERRQ(ierr);
4639b072555Sjeremylt     ierr = SNESSetSolution(snes_dummy, X[0]); CHKERRQ(ierr);
4649b072555Sjeremylt     ierr = SNESSetFunction(snes_dummy, X[0], FormResidual_Ceed,
4659b072555Sjeremylt                            user_O[0]); CHKERRQ(ierr);
466*53b04fa6SJed Brown     ierr = SNESSetJacobian(snes_dummy, mat_coarse, mat_coarse, NULL,
467*53b04fa6SJed Brown                            NULL); CHKERRQ(ierr);
4689b072555Sjeremylt     ierr = SNESComputeJacobianDefaultColor(snes_dummy, X[0], mat_O[0],
4699b072555Sjeremylt                                            mat_coarse, NULL); CHKERRQ(ierr);
470*53b04fa6SJed Brown     ierr = SNESDestroy(&snes_dummy); CHKERRQ(ierr);
471*53b04fa6SJed Brown   }
47215ce0ef0Sjeremylt 
4736c5df90dSjeremylt   // Set up KSP
4746c5df90dSjeremylt   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
4756c5df90dSjeremylt   {
4766c5df90dSjeremylt     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
4776c5df90dSjeremylt     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
4786c5df90dSjeremylt     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
4796c5df90dSjeremylt                             PETSC_DEFAULT); CHKERRQ(ierr);
4806c5df90dSjeremylt   }
4816c5df90dSjeremylt   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
4829b072555Sjeremylt   ierr = KSPSetOperators(ksp, mat_O[fine_level], mat_O[fine_level]);
4836c5df90dSjeremylt   CHKERRQ(ierr);
4846c5df90dSjeremylt 
4856c5df90dSjeremylt   // Set up PCMG
4866c5df90dSjeremylt   ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
4879b072555Sjeremylt   PCMGCycleType pcmg_cycle_type = PC_MG_CYCLE_V;
4886c5df90dSjeremylt   {
4896c5df90dSjeremylt     ierr = PCSetType(pc, PCMG); CHKERRQ(ierr);
4906c5df90dSjeremylt 
4916c5df90dSjeremylt     // PCMG levels
4929b072555Sjeremylt     ierr = PCMGSetLevels(pc, num_levels, NULL); CHKERRQ(ierr);
4939b072555Sjeremylt     for (int i=0; i<num_levels; i++) {
4946c5df90dSjeremylt       // Smoother
4956c5df90dSjeremylt       KSP smoother;
4966c5df90dSjeremylt       PC smoother_pc;
4976c5df90dSjeremylt       ierr = PCMGGetSmoother(pc, i, &smoother); CHKERRQ(ierr);
4986c5df90dSjeremylt       ierr = KSPSetType(smoother, KSPCHEBYSHEV); CHKERRQ(ierr);
4996c5df90dSjeremylt       ierr = KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1); CHKERRQ(ierr);
5006c5df90dSjeremylt       ierr = KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE); CHKERRQ(ierr);
5019b072555Sjeremylt       ierr = KSPSetOperators(smoother, mat_O[i], mat_O[i]); CHKERRQ(ierr);
5026c5df90dSjeremylt       ierr = KSPGetPC(smoother, &smoother_pc); CHKERRQ(ierr);
5036c5df90dSjeremylt       ierr = PCSetType(smoother_pc, PCJACOBI); CHKERRQ(ierr);
5046c5df90dSjeremylt       ierr = PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL); CHKERRQ(ierr);
5056c5df90dSjeremylt 
5066c5df90dSjeremylt       // Work vector
5079b072555Sjeremylt       if (i < num_levels - 1) {
5086c5df90dSjeremylt         ierr = PCMGSetX(pc, i, X[i]); CHKERRQ(ierr);
5096c5df90dSjeremylt       }
5106c5df90dSjeremylt 
5116c5df90dSjeremylt       // Level transfers
5126c5df90dSjeremylt       if (i > 0) {
5136c5df90dSjeremylt         // Interpolation
5149b072555Sjeremylt         ierr = PCMGSetInterpolation(pc, i, mat_pr[i]); CHKERRQ(ierr);
5156c5df90dSjeremylt       }
5166c5df90dSjeremylt 
5176c5df90dSjeremylt       // Coarse solve
5186c5df90dSjeremylt       KSP coarse;
5196c5df90dSjeremylt       PC coarse_pc;
5206c5df90dSjeremylt       ierr = PCMGGetCoarseSolve(pc, &coarse); CHKERRQ(ierr);
52115ce0ef0Sjeremylt       ierr = KSPSetType(coarse, KSPPREONLY); CHKERRQ(ierr);
5229b072555Sjeremylt       ierr = KSPSetOperators(coarse, mat_coarse, mat_coarse); CHKERRQ(ierr);
52315ce0ef0Sjeremylt 
5246c5df90dSjeremylt       ierr = KSPGetPC(coarse, &coarse_pc); CHKERRQ(ierr);
52515ce0ef0Sjeremylt       ierr = PCSetType(coarse_pc, PCGAMG); CHKERRQ(ierr);
52615ce0ef0Sjeremylt 
52715ce0ef0Sjeremylt       ierr = KSPSetOptionsPrefix(coarse, "coarse_"); CHKERRQ(ierr);
52815ce0ef0Sjeremylt       ierr = PCSetOptionsPrefix(coarse_pc, "coarse_"); CHKERRQ(ierr);
52915ce0ef0Sjeremylt       ierr = KSPSetFromOptions(coarse); CHKERRQ(ierr);
53015ce0ef0Sjeremylt       ierr = PCSetFromOptions(coarse_pc); CHKERRQ(ierr);
5316c5df90dSjeremylt     }
5326c5df90dSjeremylt 
5336c5df90dSjeremylt     // PCMG options
5346c5df90dSjeremylt     ierr = PCMGSetType(pc, PC_MG_MULTIPLICATIVE); CHKERRQ(ierr);
5356c5df90dSjeremylt     ierr = PCMGSetNumberSmooth(pc, 3); CHKERRQ(ierr);
5369b072555Sjeremylt     ierr = PCMGSetCycleType(pc, pcmg_cycle_type); CHKERRQ(ierr);
5376c5df90dSjeremylt   }
5386c5df90dSjeremylt 
5396c5df90dSjeremylt   // First run, if benchmarking
5406c5df90dSjeremylt   if (benchmark_mode) {
5416c5df90dSjeremylt     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
5426c5df90dSjeremylt     CHKERRQ(ierr);
5439b072555Sjeremylt     ierr = VecZeroEntries(X[fine_level]); CHKERRQ(ierr);
5446c5df90dSjeremylt     my_rt_start = MPI_Wtime();
5459b072555Sjeremylt     ierr = KSPSolve(ksp, rhs, X[fine_level]); CHKERRQ(ierr);
5466c5df90dSjeremylt     my_rt = MPI_Wtime() - my_rt_start;
5476c5df90dSjeremylt     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
5486c5df90dSjeremylt     CHKERRQ(ierr);
5496c5df90dSjeremylt     // Set maxits based on first iteration timing
5506c5df90dSjeremylt     if (my_rt > 0.02) {
5516c5df90dSjeremylt       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
5526c5df90dSjeremylt       CHKERRQ(ierr);
5536c5df90dSjeremylt     } else {
5546c5df90dSjeremylt       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
5556c5df90dSjeremylt       CHKERRQ(ierr);
5566c5df90dSjeremylt     }
5576c5df90dSjeremylt   }
5586c5df90dSjeremylt 
5596c5df90dSjeremylt   // Timed solve
5609b072555Sjeremylt   ierr = VecZeroEntries(X[fine_level]); CHKERRQ(ierr);
5616c5df90dSjeremylt   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
56209a940d7Sjeremylt 
56309a940d7Sjeremylt   // -- Performance logging
5649b072555Sjeremylt   ierr = PetscLogStageRegister("Solve Stage", &solve_stage); CHKERRQ(ierr);
5659b072555Sjeremylt   ierr = PetscLogStagePush(solve_stage); CHKERRQ(ierr);
56609a940d7Sjeremylt 
56709a940d7Sjeremylt   // -- Solve
5686c5df90dSjeremylt   my_rt_start = MPI_Wtime();
5699b072555Sjeremylt   ierr = KSPSolve(ksp, rhs, X[fine_level]); CHKERRQ(ierr);
5706c5df90dSjeremylt   my_rt = MPI_Wtime() - my_rt_start;
5716c5df90dSjeremylt 
57209a940d7Sjeremylt 
57309a940d7Sjeremylt   // -- Performance logging
57409a940d7Sjeremylt   ierr = PetscLogStagePop();
57509a940d7Sjeremylt 
5766c5df90dSjeremylt   // Output results
5776c5df90dSjeremylt   {
5789b072555Sjeremylt     KSPType ksp_type;
5799b072555Sjeremylt     PCMGType pcmg_type;
5806c5df90dSjeremylt     KSPConvergedReason reason;
5816c5df90dSjeremylt     PetscReal rnorm;
5826c5df90dSjeremylt     PetscInt its;
5839b072555Sjeremylt     ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr);
5846c5df90dSjeremylt     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
5856c5df90dSjeremylt     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
5866c5df90dSjeremylt     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
5879b072555Sjeremylt     ierr = PCMGGetType(pc, &pcmg_type); CHKERRQ(ierr);
5886c5df90dSjeremylt     if (!test_mode || reason < 0 || rnorm > 1e-8) {
5896c5df90dSjeremylt       ierr = PetscPrintf(comm,
5906c5df90dSjeremylt                          "  KSP:\n"
5916c5df90dSjeremylt                          "    KSP Type                           : %s\n"
5926c5df90dSjeremylt                          "    KSP Convergence                    : %s\n"
5936c5df90dSjeremylt                          "    Total KSP Iterations               : %D\n"
5946c5df90dSjeremylt                          "    Final rnorm                        : %e\n",
5959b072555Sjeremylt                          ksp_type, KSPConvergedReasons[reason], its,
5966c5df90dSjeremylt                          (double)rnorm); CHKERRQ(ierr);
5976c5df90dSjeremylt       ierr = PetscPrintf(comm,
5986c5df90dSjeremylt                          "  PCMG:\n"
5996c5df90dSjeremylt                          "    PCMG Type                          : %s\n"
6006c5df90dSjeremylt                          "    PCMG Cycle Type                    : %s\n",
6019b072555Sjeremylt                          PCMGTypes[pcmg_type],
6029b072555Sjeremylt                          PCMGCycleTypes[pcmg_cycle_type]); CHKERRQ(ierr);
6036c5df90dSjeremylt     }
6046c5df90dSjeremylt     if (!test_mode) {
6056c5df90dSjeremylt       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
6066c5df90dSjeremylt     }
6076c5df90dSjeremylt     {
6089b072555Sjeremylt       PetscReal max_error;
6099b072555Sjeremylt       ierr = ComputeErrorMax(user_O[fine_level], op_error, X[fine_level], target,
6109b072555Sjeremylt                              &max_error); CHKERRQ(ierr);
6116c5df90dSjeremylt       PetscReal tol = 5e-2;
6129b072555Sjeremylt       if (!test_mode || max_error > tol) {
6136c5df90dSjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
6146c5df90dSjeremylt         CHKERRQ(ierr);
6156c5df90dSjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
6166c5df90dSjeremylt         CHKERRQ(ierr);
6176c5df90dSjeremylt         ierr = PetscPrintf(comm,
6186c5df90dSjeremylt                            "    Pointwise Error (max)              : %e\n"
6196c5df90dSjeremylt                            "    CG Solve Time                      : %g (%g) sec\n",
6209b072555Sjeremylt                            (double)max_error, rt_max, rt_min); CHKERRQ(ierr);
6216c5df90dSjeremylt       }
6226c5df90dSjeremylt     }
6236c5df90dSjeremylt     if (benchmark_mode && (!test_mode)) {
6246c5df90dSjeremylt       ierr = PetscPrintf(comm,
6256c5df90dSjeremylt                          "    DoFs/Sec in CG                     : %g (%g) million\n",
6269b072555Sjeremylt                          1e-6*g_size[fine_level]*its/rt_max,
6279b072555Sjeremylt                          1e-6*g_size[fine_level]*its/rt_min);
6286c5df90dSjeremylt       CHKERRQ(ierr);
6296c5df90dSjeremylt     }
6306c5df90dSjeremylt   }
6316c5df90dSjeremylt 
6326c5df90dSjeremylt   if (write_solution) {
6339b072555Sjeremylt     PetscViewer vtk_viewer_soln;
6346c5df90dSjeremylt 
6359b072555Sjeremylt     ierr = PetscViewerCreate(comm, &vtk_viewer_soln); CHKERRQ(ierr);
6369b072555Sjeremylt     ierr = PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK); CHKERRQ(ierr);
6379b072555Sjeremylt     ierr = PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"); CHKERRQ(ierr);
6389b072555Sjeremylt     ierr = VecView(X[fine_level], vtk_viewer_soln); CHKERRQ(ierr);
6399b072555Sjeremylt     ierr = PetscViewerDestroy(&vtk_viewer_soln); CHKERRQ(ierr);
6406c5df90dSjeremylt   }
6416c5df90dSjeremylt 
6426c5df90dSjeremylt   // Cleanup
6439b072555Sjeremylt   for (int i=0; i<num_levels; i++) {
6446c5df90dSjeremylt     ierr = VecDestroy(&X[i]); CHKERRQ(ierr);
6459b072555Sjeremylt     ierr = VecDestroy(&X_loc[i]); CHKERRQ(ierr);
6466c5df90dSjeremylt     ierr = VecDestroy(&mult[i]); CHKERRQ(ierr);
6479b072555Sjeremylt     ierr = VecDestroy(&user_O[i]->Y_loc); CHKERRQ(ierr);
6489b072555Sjeremylt     ierr = MatDestroy(&mat_O[i]); CHKERRQ(ierr);
6499b072555Sjeremylt     ierr = PetscFree(user_O[i]); CHKERRQ(ierr);
6506c5df90dSjeremylt     if (i > 0) {
6519b072555Sjeremylt       ierr = MatDestroy(&mat_pr[i]); CHKERRQ(ierr);
6529b072555Sjeremylt       ierr = PetscFree(user_pr[i]); CHKERRQ(ierr);
6536c5df90dSjeremylt     }
6549b072555Sjeremylt     ierr = CeedDataDestroy(i, ceed_data[i]); CHKERRQ(ierr);
6556c5df90dSjeremylt     ierr = DMDestroy(&dm[i]); CHKERRQ(ierr);
6566c5df90dSjeremylt   }
6579b072555Sjeremylt   ierr = PetscFree(level_degrees); CHKERRQ(ierr);
6586c5df90dSjeremylt   ierr = PetscFree(dm); CHKERRQ(ierr);
6596c5df90dSjeremylt   ierr = PetscFree(X); CHKERRQ(ierr);
6609b072555Sjeremylt   ierr = PetscFree(X_loc); CHKERRQ(ierr);
6616c5df90dSjeremylt   ierr = PetscFree(mult); CHKERRQ(ierr);
6629b072555Sjeremylt   ierr = PetscFree(mat_O); CHKERRQ(ierr);
6639b072555Sjeremylt   ierr = PetscFree(mat_pr); CHKERRQ(ierr);
6649b072555Sjeremylt   ierr = PetscFree(ceed_data); CHKERRQ(ierr);
6659b072555Sjeremylt   ierr = PetscFree(user_O); CHKERRQ(ierr);
6669b072555Sjeremylt   ierr = PetscFree(user_pr); CHKERRQ(ierr);
6679b072555Sjeremylt   ierr = PetscFree(l_size); CHKERRQ(ierr);
6689b072555Sjeremylt   ierr = PetscFree(xl_size); CHKERRQ(ierr);
6699b072555Sjeremylt   ierr = PetscFree(g_size); CHKERRQ(ierr);
6706c5df90dSjeremylt   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
6719b072555Sjeremylt   ierr = VecDestroy(&rhs_loc); CHKERRQ(ierr);
6729b072555Sjeremylt   ierr = MatDestroy(&mat_coarse); CHKERRQ(ierr);
6736c5df90dSjeremylt   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
6749b072555Sjeremylt   ierr = DMDestroy(&dm_orig); CHKERRQ(ierr);
6756c5df90dSjeremylt   CeedVectorDestroy(&target);
6769b072555Sjeremylt   CeedQFunctionDestroy(&qf_error);
6779b072555Sjeremylt   CeedQFunctionDestroy(&qf_restrict);
6789b072555Sjeremylt   CeedQFunctionDestroy(&qf_prolong);
6799b072555Sjeremylt   CeedOperatorDestroy(&op_error);
6806c5df90dSjeremylt   CeedDestroy(&ceed);
6816c5df90dSjeremylt   return PetscFinalize();
6826c5df90dSjeremylt }
683