xref: /libCEED/examples/petsc/bps.c (revision cb32e2e7f026784d97a57f1901677e9727def907)
1819eb1b3Sjeremylt // Copyright (c) 2017, Lawrence Livermore National Security, LLC. Produced at
2819eb1b3Sjeremylt // the Lawrence Livermore National Laboratory. LLNL-CODE-734707. All Rights
3819eb1b3Sjeremylt // reserved. See files LICENSE and NOTICE for details.
4819eb1b3Sjeremylt //
5819eb1b3Sjeremylt // This file is part of CEED, a collection of benchmarks, miniapps, software
6819eb1b3Sjeremylt // libraries and APIs for efficient high-order finite element and spectral
7819eb1b3Sjeremylt // element discretizations for exascale applications. For more information and
8819eb1b3Sjeremylt // source code availability see http://github.com/ceed.
9819eb1b3Sjeremylt //
10819eb1b3Sjeremylt // The CEED research is supported by the Exascale Computing Project 17-SC-20-SC,
11819eb1b3Sjeremylt // a collaborative effort of two U.S. Department of Energy organizations (Office
12819eb1b3Sjeremylt // of Science and the National Nuclear Security Administration) responsible for
13819eb1b3Sjeremylt // the planning and preparation of a capable exascale ecosystem, including
14819eb1b3Sjeremylt // software, applications, hardware, advanced system engineering and early
15819eb1b3Sjeremylt // testbed platforms, in support of the nation's exascale computing imperative.
16819eb1b3Sjeremylt 
170c59ef15SJeremy L Thompson //                        libCEED + PETSc Example: CEED BPs
180c59ef15SJeremy L Thompson //
190c59ef15SJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve the
200c59ef15SJeremy L Thompson // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
210c59ef15SJeremy L Thompson //
22*cb32e2e7SValeria Barra // The code uses higher level communication protocols in DMPlex.
230c59ef15SJeremy L Thompson //
240c59ef15SJeremy L Thompson // Build with:
250c59ef15SJeremy L Thompson //
260c59ef15SJeremy L Thompson //     make bps [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
270c59ef15SJeremy L Thompson //
280c59ef15SJeremy L Thompson // Sample runs:
290c59ef15SJeremy L Thompson //
30*cb32e2e7SValeria Barra //     bps -problem bp1 -degree 3
31*cb32e2e7SValeria Barra //     bps -problem bp2 -ceed /cpu/self -degree 3
32*cb32e2e7SValeria Barra //     bps -problem bp3 -ceed /gpu/occa -degree 3
33*cb32e2e7SValeria Barra //     bps -problem bp4 -ceed /cpu/occa -degree 3
34*cb32e2e7SValeria Barra //     bps -problem bp5 -ceed /omp/occa -degree 3
35*cb32e2e7SValeria Barra //     bps -problem bp6 -ceed /ocl/occa -degree 3
360c59ef15SJeremy L Thompson //
37*cb32e2e7SValeria Barra //TESTARGS -ceed {ceed_resource} -test -problem bp5 -degree 3
380c59ef15SJeremy L Thompson 
390c59ef15SJeremy L Thompson /// @file
40*cb32e2e7SValeria Barra /// CEED BPs example using PETSc with DMPlex
41*cb32e2e7SValeria Barra /// See bpsraw.c for a "raw" implementation using a structured grid.
42*cb32e2e7SValeria Barra const char help[] = "Solve CEED BPs using PETSc with DMPlex\n";
430c59ef15SJeremy L Thompson 
440c59ef15SJeremy L Thompson #include <stdbool.h>
450c59ef15SJeremy L Thompson #include <string.h>
464d537eeaSYohann #include <petscksp.h>
47*cb32e2e7SValeria Barra #include <petscdmplex.h>
484d537eeaSYohann #include <ceed.h>
49*cb32e2e7SValeria Barra #include "setup.h"
500c59ef15SJeremy L Thompson 
510c59ef15SJeremy L Thompson int main(int argc, char **argv) {
520c59ef15SJeremy L Thompson   PetscInt ierr;
530c59ef15SJeremy L Thompson   MPI_Comm comm;
54*cb32e2e7SValeria Barra   char ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self",
55*cb32e2e7SValeria Barra        filename[PETSC_MAX_PATH_LEN];
56819eb1b3Sjeremylt   double my_rt_start, my_rt, rt_min, rt_max;
57*cb32e2e7SValeria Barra   PetscInt degree = 3, qextra, lsize, gsize, dim = 3, melem[3] = {3, 3, 3},
58*cb32e2e7SValeria Barra            ncompu = 1, xlsize;
590c59ef15SJeremy L Thompson   PetscScalar *r;
60*cb32e2e7SValeria Barra   PetscBool test_mode, benchmark_mode, read_mesh, write_solution;
61819eb1b3Sjeremylt   Vec X, Xloc, rhs, rhsloc;
62*cb32e2e7SValeria Barra   Mat matO;
63819eb1b3Sjeremylt   KSP ksp;
64*cb32e2e7SValeria Barra   DM  dm;
65*cb32e2e7SValeria Barra   UserO userO;
660c59ef15SJeremy L Thompson   Ceed ceed;
67*cb32e2e7SValeria Barra   CeedData ceeddata;
68*cb32e2e7SValeria Barra   CeedQFunction qf_error;
69*cb32e2e7SValeria Barra   CeedOperator op_error;
70*cb32e2e7SValeria Barra   CeedVector rhsceed, target;
710c59ef15SJeremy L Thompson   bpType bpChoice;
720c59ef15SJeremy L Thompson 
730c59ef15SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
740c59ef15SJeremy L Thompson   if (ierr) return ierr;
750c59ef15SJeremy L Thompson   comm = PETSC_COMM_WORLD;
76*cb32e2e7SValeria Barra 
77*cb32e2e7SValeria Barra   // Read CL options
780c59ef15SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
790c59ef15SJeremy L Thompson   bpChoice = CEED_BP1;
800c59ef15SJeremy L Thompson   ierr = PetscOptionsEnum("-problem",
810c59ef15SJeremy L Thompson                           "CEED benchmark problem to solve", NULL,
820c59ef15SJeremy L Thompson                           bpTypes, (PetscEnum)bpChoice, (PetscEnum *)&bpChoice,
830c59ef15SJeremy L Thompson                           NULL); CHKERRQ(ierr);
8434d77899SValeria Barra   ncompu = bpOptions[bpChoice].ncompu;
850c59ef15SJeremy L Thompson   test_mode = PETSC_FALSE;
860c59ef15SJeremy L Thompson   ierr = PetscOptionsBool("-test",
870c59ef15SJeremy L Thompson                           "Testing mode (do not print unless error is large)",
880c59ef15SJeremy L Thompson                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
890c59ef15SJeremy L Thompson   benchmark_mode = PETSC_FALSE;
900c59ef15SJeremy L Thompson   ierr = PetscOptionsBool("-benchmark",
910c59ef15SJeremy L Thompson                           "Benchmarking mode (prints benchmark statistics)",
920c59ef15SJeremy L Thompson                           NULL, benchmark_mode, &benchmark_mode, NULL);
930c59ef15SJeremy L Thompson   CHKERRQ(ierr);
946bd9afcaSjeremylt   write_solution = PETSC_FALSE;
956bd9afcaSjeremylt   ierr = PetscOptionsBool("-write_solution",
966bd9afcaSjeremylt                           "Write solution for visualization",
976bd9afcaSjeremylt                           NULL, write_solution, &write_solution, NULL);
986bd9afcaSjeremylt   CHKERRQ(ierr);
99*cb32e2e7SValeria Barra   degree = test_mode ? 3 : 2;
1000c59ef15SJeremy L Thompson   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
1010c59ef15SJeremy L Thompson                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
1020c59ef15SJeremy L Thompson   qextra = bpOptions[bpChoice].qextra;
1030c59ef15SJeremy L Thompson   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
1040c59ef15SJeremy L Thompson                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
1050c59ef15SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
1060c59ef15SJeremy L Thompson                             NULL, ceedresource, ceedresource,
1070c59ef15SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
108*cb32e2e7SValeria Barra   read_mesh = PETSC_FALSE;
109*cb32e2e7SValeria Barra   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
110*cb32e2e7SValeria Barra                             filename, filename, sizeof(filename), &read_mesh);
1110c59ef15SJeremy L Thompson   CHKERRQ(ierr);
112*cb32e2e7SValeria Barra   if (!read_mesh) {
113*cb32e2e7SValeria Barra     PetscInt tmp = dim;
114*cb32e2e7SValeria Barra     ierr = PetscOptionsIntArray("-cells","Number of cells per dimension", NULL,
115*cb32e2e7SValeria Barra                                 melem, &tmp, NULL); CHKERRQ(ierr);
116*cb32e2e7SValeria Barra   }
117*cb32e2e7SValeria Barra   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
118*cb32e2e7SValeria Barra 
119*cb32e2e7SValeria Barra   // Setup DM
120*cb32e2e7SValeria Barra   if (read_mesh) {
121*cb32e2e7SValeria Barra     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm);
122*cb32e2e7SValeria Barra     CHKERRQ(ierr);
123*cb32e2e7SValeria Barra   } else {
124*cb32e2e7SValeria Barra     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL,
125*cb32e2e7SValeria Barra                                NULL, NULL, PETSC_TRUE, &dm); CHKERRQ(ierr);
126*cb32e2e7SValeria Barra   }
127*cb32e2e7SValeria Barra 
128*cb32e2e7SValeria Barra   {
129*cb32e2e7SValeria Barra     DM dmDist = NULL;
130*cb32e2e7SValeria Barra     PetscPartitioner part;
131*cb32e2e7SValeria Barra 
132*cb32e2e7SValeria Barra     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
133*cb32e2e7SValeria Barra     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
134*cb32e2e7SValeria Barra     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
135*cb32e2e7SValeria Barra     if (dmDist) {
136*cb32e2e7SValeria Barra       ierr = DMDestroy(&dm); CHKERRQ(ierr);
137*cb32e2e7SValeria Barra       dm  = dmDist;
138*cb32e2e7SValeria Barra     }
139*cb32e2e7SValeria Barra   }
140*cb32e2e7SValeria Barra 
141*cb32e2e7SValeria Barra   // Create DM
142*cb32e2e7SValeria Barra   ierr = SetupDMByDegree(dm, degree, ncompu, bpChoice);
143*cb32e2e7SValeria Barra   CHKERRQ(ierr);
144*cb32e2e7SValeria Barra 
145*cb32e2e7SValeria Barra   // Create vectors
146*cb32e2e7SValeria Barra   ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr);
147*cb32e2e7SValeria Barra   ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr);
148*cb32e2e7SValeria Barra   ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
149*cb32e2e7SValeria Barra   ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr);
150*cb32e2e7SValeria Barra   ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr);
151*cb32e2e7SValeria Barra   ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr);
152*cb32e2e7SValeria Barra 
153*cb32e2e7SValeria Barra   // Operator
154*cb32e2e7SValeria Barra   ierr = PetscMalloc1(1, &userO); CHKERRQ(ierr);
155*cb32e2e7SValeria Barra   ierr = MatCreateShell(comm, lsize, lsize, gsize, gsize,
156*cb32e2e7SValeria Barra                         userO, &matO); CHKERRQ(ierr);
157*cb32e2e7SValeria Barra   ierr = MatShellSetOperation(matO, MATOP_MULT,
158*cb32e2e7SValeria Barra                               (void(*)(void))MatMult_Ceed);
159*cb32e2e7SValeria Barra   CHKERRQ(ierr);
1600c59ef15SJeremy L Thompson 
1612d03409cSjeremylt   // Set up libCEED
1622d03409cSjeremylt   CeedInit(ceedresource, &ceed);
1632d03409cSjeremylt 
164819eb1b3Sjeremylt   // Print summary
1650c59ef15SJeremy L Thompson   if (!test_mode) {
166*cb32e2e7SValeria Barra     PetscInt P = degree + 1, Q = P + qextra;
1672d03409cSjeremylt     const char *usedresource;
1682d03409cSjeremylt     CeedGetResource(ceed, &usedresource);
1693ce86546SJeremy L Thompson     ierr = PetscPrintf(comm,
1703ce86546SJeremy L Thompson                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n"
1713ce86546SJeremy L Thompson                        "  libCEED:\n"
1723ce86546SJeremy L Thompson                        "    libCEED Backend                    : %s\n"
1733ce86546SJeremy L Thompson                        "  Mesh:\n"
1743ce86546SJeremy L Thompson                        "    Number of 1D Basis Nodes (p)       : %d\n"
1753ce86546SJeremy L Thompson                        "    Number of 1D Quadrature Points (q) : %d\n"
176b52ddc88Sjeremylt                        "    Global nodes                       : %D\n"
177*cb32e2e7SValeria Barra                        "    Owned nodes                        : %D\n",
178*cb32e2e7SValeria Barra                        bpChoice+1, usedresource, P, Q, gsize/ncompu,
179*cb32e2e7SValeria Barra                        lsize/ncompu); CHKERRQ(ierr);
1800c59ef15SJeremy L Thompson   }
1810c59ef15SJeremy L Thompson 
182*cb32e2e7SValeria Barra   // Create RHS vector
183*cb32e2e7SValeria Barra   ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr);
184*cb32e2e7SValeria Barra   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
185*cb32e2e7SValeria Barra   ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr);
186*cb32e2e7SValeria Barra   CeedVectorCreate(ceed, xlsize, &rhsceed);
187*cb32e2e7SValeria Barra   CeedVectorSetArray(rhsceed, CEED_MEM_HOST, CEED_USE_POINTER, r);
1880c59ef15SJeremy L Thompson 
189*cb32e2e7SValeria Barra   ierr = PetscMalloc1(1, &ceeddata); CHKERRQ(ierr);
190*cb32e2e7SValeria Barra   ierr = SetupLibceedByDegree(dm, ceed, degree, dim, qextra,
191*cb32e2e7SValeria Barra                               ncompu, gsize, xlsize, bpChoice, ceeddata,
192*cb32e2e7SValeria Barra                               true, rhsceed, &target); CHKERRQ(ierr);
1930c59ef15SJeremy L Thompson 
194*cb32e2e7SValeria Barra   // Gather RHS
195*cb32e2e7SValeria Barra   ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
196*cb32e2e7SValeria Barra   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
197*cb32e2e7SValeria Barra   ierr = DMLocalToGlobalBegin(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
198*cb32e2e7SValeria Barra   ierr = DMLocalToGlobalEnd(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
199*cb32e2e7SValeria Barra   CeedVectorDestroy(&rhsceed);
2000c59ef15SJeremy L Thompson 
201*cb32e2e7SValeria Barra   // Create the error Q-function
2020c59ef15SJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpChoice].error,
2030c59ef15SJeremy L Thompson                               bpOptions[bpChoice].errorfname, &qf_error);
20434d77899SValeria Barra   CeedQFunctionAddInput(qf_error, "u", ncompu, CEED_EVAL_INTERP);
20534d77899SValeria Barra   CeedQFunctionAddInput(qf_error, "true_soln", ncompu, CEED_EVAL_NONE);
20634d77899SValeria Barra   CeedQFunctionAddOutput(qf_error, "error", ncompu, CEED_EVAL_NONE);
2070c59ef15SJeremy L Thompson 
2080c59ef15SJeremy L Thompson   // Create the error operator
2090c59ef15SJeremy L Thompson   CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
210*cb32e2e7SValeria Barra   CeedOperatorSetField(op_error, "u", ceeddata->Erestrictu,
211*cb32e2e7SValeria Barra                        CEED_TRANSPOSE, ceeddata->basisu,
212*cb32e2e7SValeria Barra                        CEED_VECTOR_ACTIVE);
213*cb32e2e7SValeria Barra   CeedOperatorSetField(op_error, "true_soln", ceeddata->Erestrictui,
214*cb32e2e7SValeria Barra                        CEED_NOTRANSPOSE, CEED_BASIS_COLLOCATED, target);
215*cb32e2e7SValeria Barra   CeedOperatorSetField(op_error, "error", ceeddata->Erestrictui,
216*cb32e2e7SValeria Barra                        CEED_NOTRANSPOSE, CEED_BASIS_COLLOCATED,
217*cb32e2e7SValeria Barra                        CEED_VECTOR_ACTIVE);
2180c59ef15SJeremy L Thompson 
2190c59ef15SJeremy L Thompson   // Set up Mat
220*cb32e2e7SValeria Barra   userO->comm = comm;
221*cb32e2e7SValeria Barra   userO->dm = dm;
222*cb32e2e7SValeria Barra   userO->Xloc = Xloc;
223*cb32e2e7SValeria Barra   ierr = VecDuplicate(Xloc, &userO->Yloc); CHKERRQ(ierr);
224*cb32e2e7SValeria Barra   userO->xceed = ceeddata->xceed;
225*cb32e2e7SValeria Barra   userO->yceed = ceeddata->yceed;
226*cb32e2e7SValeria Barra   userO->op = ceeddata->op_apply;
227*cb32e2e7SValeria Barra   userO->ceed = ceed;
2280c59ef15SJeremy L Thompson 
2290c59ef15SJeremy L Thompson   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
2300c59ef15SJeremy L Thompson   {
2310c59ef15SJeremy L Thompson     PC pc;
2320c59ef15SJeremy L Thompson     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
2330c59ef15SJeremy L Thompson     if (bpChoice == CEED_BP1 || bpChoice == CEED_BP2) {
2340c59ef15SJeremy L Thompson       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
2350c59ef15SJeremy L Thompson       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
2360c59ef15SJeremy L Thompson     } else {
2370c59ef15SJeremy L Thompson       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
2380c59ef15SJeremy L Thompson     }
2390c59ef15SJeremy L Thompson     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
2400c59ef15SJeremy L Thompson     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
2410c59ef15SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
2420c59ef15SJeremy L Thompson                             PETSC_DEFAULT); CHKERRQ(ierr);
2430c59ef15SJeremy L Thompson   }
2440c59ef15SJeremy L Thompson   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
245*cb32e2e7SValeria Barra   ierr = KSPSetOperators(ksp, matO, matO); CHKERRQ(ierr);
246*cb32e2e7SValeria Barra 
2470c59ef15SJeremy L Thompson   // First run, if benchmarking
2480c59ef15SJeremy L Thompson   if (benchmark_mode) {
2490c59ef15SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
2500c59ef15SJeremy L Thompson     CHKERRQ(ierr);
2510c59ef15SJeremy L Thompson     my_rt_start = MPI_Wtime();
2520c59ef15SJeremy L Thompson     ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
2530c59ef15SJeremy L Thompson     my_rt = MPI_Wtime() - my_rt_start;
254b5009ec9Sjeremylt     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
255b5009ec9Sjeremylt     CHKERRQ(ierr);
2560c59ef15SJeremy L Thompson     // Set maxits based on first iteration timing
2570c59ef15SJeremy L Thompson     if (my_rt > 0.02) {
2580c59ef15SJeremy L Thompson       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
2590c59ef15SJeremy L Thompson       CHKERRQ(ierr);
2600c59ef15SJeremy L Thompson     } else {
2610c59ef15SJeremy L Thompson       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
2620c59ef15SJeremy L Thompson       CHKERRQ(ierr);
2630c59ef15SJeremy L Thompson     }
2640c59ef15SJeremy L Thompson   }
265*cb32e2e7SValeria Barra 
2660c59ef15SJeremy L Thompson   // Timed solve
267*cb32e2e7SValeria Barra   ierr = VecZeroEntries(X); CHKERRQ(ierr);
2683b3d4a15SJeremy L Thompson   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
2690c59ef15SJeremy L Thompson   my_rt_start = MPI_Wtime();
2700c59ef15SJeremy L Thompson   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
2710c59ef15SJeremy L Thompson   my_rt = MPI_Wtime() - my_rt_start;
272*cb32e2e7SValeria Barra 
273*cb32e2e7SValeria Barra   // Output results
2740c59ef15SJeremy L Thompson   {
2750c59ef15SJeremy L Thompson     KSPType ksptype;
2760c59ef15SJeremy L Thompson     KSPConvergedReason reason;
2770c59ef15SJeremy L Thompson     PetscReal rnorm;
2780c59ef15SJeremy L Thompson     PetscInt its;
2790c59ef15SJeremy L Thompson     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
2800c59ef15SJeremy L Thompson     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
2810c59ef15SJeremy L Thompson     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
2820c59ef15SJeremy L Thompson     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
2830c59ef15SJeremy L Thompson     if (!test_mode || reason < 0 || rnorm > 1e-8) {
2843ce86546SJeremy L Thompson       ierr = PetscPrintf(comm,
2853ce86546SJeremy L Thompson                          "  KSP:\n"
2863ce86546SJeremy L Thompson                          "    KSP Type                           : %s\n"
2873ce86546SJeremy L Thompson                          "    KSP Convergence                    : %s\n"
2883ce86546SJeremy L Thompson                          "    Total KSP Iterations               : %D\n"
2893ce86546SJeremy L Thompson                          "    Final rnorm                        : %e\n",
2903ce86546SJeremy L Thompson                          ksptype, KSPConvergedReasons[reason], its,
2913ce86546SJeremy L Thompson                          (double)rnorm); CHKERRQ(ierr);
2920c59ef15SJeremy L Thompson     }
293*cb32e2e7SValeria Barra     if (!test_mode) {
294*cb32e2e7SValeria Barra       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
295*cb32e2e7SValeria Barra     }
296*cb32e2e7SValeria Barra     {
297*cb32e2e7SValeria Barra       PetscReal maxerror;
298*cb32e2e7SValeria Barra       ierr = ComputeErrorMax(userO, op_error, X, target, &maxerror);
299*cb32e2e7SValeria Barra       CHKERRQ(ierr);
300*cb32e2e7SValeria Barra       PetscReal tol = 5e-2;
301*cb32e2e7SValeria Barra       if (!test_mode || maxerror > tol) {
302b5009ec9Sjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
303b5009ec9Sjeremylt         CHKERRQ(ierr);
304b5009ec9Sjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
305b5009ec9Sjeremylt         CHKERRQ(ierr);
3060c59ef15SJeremy L Thompson         ierr = PetscPrintf(comm,
307*cb32e2e7SValeria Barra                            "    Pointwise Error (max)              : %e\n"
308*cb32e2e7SValeria Barra                            "    CG Solve Time                      : %g (%g) sec\n",
309*cb32e2e7SValeria Barra                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
3100c59ef15SJeremy L Thompson       }
3110c59ef15SJeremy L Thompson     }
312*cb32e2e7SValeria Barra     if (benchmark_mode && (!test_mode)) {
3133ce86546SJeremy L Thompson       ierr = PetscPrintf(comm,
314*cb32e2e7SValeria Barra                          "    DoFs/Sec in CG                     : %g (%g) million\n",
315*cb32e2e7SValeria Barra                          1e-6*gsize*its/rt_max,
316*cb32e2e7SValeria Barra                          1e-6*gsize*its/rt_min); CHKERRQ(ierr);
3170c59ef15SJeremy L Thompson     }
3180c59ef15SJeremy L Thompson   }
3190c59ef15SJeremy L Thompson 
3206bd9afcaSjeremylt   if (write_solution) {
3216bd9afcaSjeremylt     PetscViewer vtkviewersoln;
3226bd9afcaSjeremylt 
3236bd9afcaSjeremylt     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
3246bd9afcaSjeremylt     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
3256bd9afcaSjeremylt     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
3266bd9afcaSjeremylt     ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr);
3276bd9afcaSjeremylt     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
3286bd9afcaSjeremylt   }
3296bd9afcaSjeremylt 
330*cb32e2e7SValeria Barra   // Cleanup
331*cb32e2e7SValeria Barra   ierr = VecDestroy(&X); CHKERRQ(ierr);
332*cb32e2e7SValeria Barra   ierr = VecDestroy(&Xloc); CHKERRQ(ierr);
333*cb32e2e7SValeria Barra   ierr = VecDestroy(&userO->Yloc); CHKERRQ(ierr);
334*cb32e2e7SValeria Barra   ierr = MatDestroy(&matO); CHKERRQ(ierr);
335*cb32e2e7SValeria Barra   ierr = PetscFree(userO); CHKERRQ(ierr);
336*cb32e2e7SValeria Barra   ierr = CeedDataDestroy(0, ceeddata); CHKERRQ(ierr);
337*cb32e2e7SValeria Barra   ierr = DMDestroy(&dm); CHKERRQ(ierr);
338*cb32e2e7SValeria Barra 
3390c59ef15SJeremy L Thompson   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
3400c59ef15SJeremy L Thompson   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
3410c59ef15SJeremy L Thompson   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
3420c59ef15SJeremy L Thompson   CeedVectorDestroy(&target);
3430c59ef15SJeremy L Thompson   CeedQFunctionDestroy(&qf_error);
344*cb32e2e7SValeria Barra   CeedOperatorDestroy(&op_error);
3450c59ef15SJeremy L Thompson   CeedDestroy(&ceed);
3460c59ef15SJeremy L Thompson   return PetscFinalize();
3470c59ef15SJeremy L Thompson }
348