xref: /libCEED/examples/petsc/bps.c (revision 153bcb04e7732932dbed1d2a53d9a5bc67c3a7e1)
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 //
22cb32e2e7SValeria 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 //
3032d2ee49SValeria Barra //     ./bps -problem bp1 -degree 3
3132d2ee49SValeria Barra //     ./bps -problem bp2 -ceed /cpu/self -degree 3
3232d2ee49SValeria Barra //     ./bps -problem bp3 -ceed /gpu/occa -degree 3
3332d2ee49SValeria Barra //     ./bps -problem bp4 -ceed /cpu/occa -degree 3
3432d2ee49SValeria Barra //     ./bps -problem bp5 -ceed /omp/occa -degree 3
3532d2ee49SValeria Barra //     ./bps -problem bp6 -ceed /ocl/occa -degree 3
360c59ef15SJeremy L Thompson //
37cb32e2e7SValeria Barra //TESTARGS -ceed {ceed_resource} -test -problem bp5 -degree 3
380c59ef15SJeremy L Thompson 
390c59ef15SJeremy L Thompson /// @file
40cb32e2e7SValeria Barra /// CEED BPs example using PETSc with DMPlex
41cb32e2e7SValeria Barra /// See bpsraw.c for a "raw" implementation using a structured grid.
42cb32e2e7SValeria 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>
47cb32e2e7SValeria Barra #include <petscdmplex.h>
484d537eeaSYohann #include <ceed.h>
49cb32e2e7SValeria 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;
54cb0b5415Sjeremylt   char filename[PETSC_MAX_PATH_LEN],
55cb0b5415Sjeremylt        ceedresource[PETSC_MAX_PATH_LEN] = "/cpu/self";
56819eb1b3Sjeremylt   double my_rt_start, my_rt, rt_min, rt_max;
57cb32e2e7SValeria Barra   PetscInt degree = 3, qextra, lsize, gsize, dim = 3, melem[3] = {3, 3, 3},
58*153bcb04Svaleriabarra            ncompu = 1, xlsize, localnodes, lelem;
590c59ef15SJeremy L Thompson   PetscScalar *r;
60*153bcb04Svaleriabarra   PetscBool test_mode, benchmark_mode, read_mesh, write_solution,
61*153bcb04Svaleriabarra             userlnodes = PETSC_FALSE;
6209a940d7Sjeremylt   PetscLogStage solvestage;
63819eb1b3Sjeremylt   Vec X, Xloc, rhs, rhsloc;
64cb32e2e7SValeria Barra   Mat matO;
65819eb1b3Sjeremylt   KSP ksp;
66cb32e2e7SValeria Barra   DM  dm;
67cb32e2e7SValeria Barra   UserO userO;
680c59ef15SJeremy L Thompson   Ceed ceed;
69cb32e2e7SValeria Barra   CeedData ceeddata;
70c70bd2a0Sjeremylt   CeedQFunction qferror;
71c70bd2a0Sjeremylt   CeedOperator operror;
72cb32e2e7SValeria Barra   CeedVector rhsceed, target;
739396343dSjeremylt   CeedMemType memtyperequested;
74199551f5Sjeremylt   bpType bpchoice;
750c59ef15SJeremy L Thompson 
769396343dSjeremylt   // Check PETSc CUDA support
779396343dSjeremylt   PetscBool petschavecuda, setmemtyperequest = PETSC_FALSE;
789396343dSjeremylt   // *INDENT-OFF*
799396343dSjeremylt   #ifdef PETSC_HAVE_CUDA
809396343dSjeremylt   petschavecuda = PETSC_TRUE;
819396343dSjeremylt   #else
829396343dSjeremylt   petschavecuda = PETSC_FALSE;
839396343dSjeremylt   #endif
849396343dSjeremylt   // *INDENT-ON*
859396343dSjeremylt 
860c59ef15SJeremy L Thompson   ierr = PetscInitialize(&argc, &argv, NULL, help);
870c59ef15SJeremy L Thompson   if (ierr) return ierr;
880c59ef15SJeremy L Thompson   comm = PETSC_COMM_WORLD;
89cb32e2e7SValeria Barra 
9032d2ee49SValeria Barra   // Read command line options
910c59ef15SJeremy L Thompson   ierr = PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); CHKERRQ(ierr);
92199551f5Sjeremylt   bpchoice = CEED_BP1;
930c59ef15SJeremy L Thompson   ierr = PetscOptionsEnum("-problem",
940c59ef15SJeremy L Thompson                           "CEED benchmark problem to solve", NULL,
95199551f5Sjeremylt                           bpTypes, (PetscEnum)bpchoice, (PetscEnum *)&bpchoice,
960c59ef15SJeremy L Thompson                           NULL); CHKERRQ(ierr);
97199551f5Sjeremylt   ncompu = bpOptions[bpchoice].ncompu;
980c59ef15SJeremy L Thompson   test_mode = PETSC_FALSE;
990c59ef15SJeremy L Thompson   ierr = PetscOptionsBool("-test",
1000c59ef15SJeremy L Thompson                           "Testing mode (do not print unless error is large)",
1010c59ef15SJeremy L Thompson                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
1020c59ef15SJeremy L Thompson   benchmark_mode = PETSC_FALSE;
1030c59ef15SJeremy L Thompson   ierr = PetscOptionsBool("-benchmark",
1040c59ef15SJeremy L Thompson                           "Benchmarking mode (prints benchmark statistics)",
1050c59ef15SJeremy L Thompson                           NULL, benchmark_mode, &benchmark_mode, NULL);
1060c59ef15SJeremy L Thompson   CHKERRQ(ierr);
1076bd9afcaSjeremylt   write_solution = PETSC_FALSE;
1086bd9afcaSjeremylt   ierr = PetscOptionsBool("-write_solution",
1096bd9afcaSjeremylt                           "Write solution for visualization",
1106bd9afcaSjeremylt                           NULL, write_solution, &write_solution, NULL);
1116bd9afcaSjeremylt   CHKERRQ(ierr);
112cb32e2e7SValeria Barra   degree = test_mode ? 3 : 2;
1130c59ef15SJeremy L Thompson   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
1140c59ef15SJeremy L Thompson                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
115199551f5Sjeremylt   qextra = bpOptions[bpchoice].qextra;
1160c59ef15SJeremy L Thompson   ierr = PetscOptionsInt("-qextra", "Number of extra quadrature points",
1170c59ef15SJeremy L Thompson                          NULL, qextra, &qextra, NULL); CHKERRQ(ierr);
1180c59ef15SJeremy L Thompson   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
1190c59ef15SJeremy L Thompson                             NULL, ceedresource, ceedresource,
1200c59ef15SJeremy L Thompson                             sizeof(ceedresource), NULL); CHKERRQ(ierr);
121cb32e2e7SValeria Barra   read_mesh = PETSC_FALSE;
122cb32e2e7SValeria Barra   ierr = PetscOptionsString("-mesh", "Read mesh from file", NULL,
123cb32e2e7SValeria Barra                             filename, filename, sizeof(filename), &read_mesh);
1240c59ef15SJeremy L Thompson   CHKERRQ(ierr);
125cb32e2e7SValeria Barra   if (!read_mesh) {
126cb32e2e7SValeria Barra     PetscInt tmp = dim;
127cb32e2e7SValeria Barra     ierr = PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL,
128cb32e2e7SValeria Barra                                 melem, &tmp, NULL); CHKERRQ(ierr);
129cb32e2e7SValeria Barra   }
1309396343dSjeremylt   memtyperequested = petschavecuda ? CEED_MEM_DEVICE : CEED_MEM_HOST;
1319396343dSjeremylt   ierr = PetscOptionsEnum("-memtype",
1329396343dSjeremylt                           "CEED MemType requested", NULL,
1339396343dSjeremylt                           memTypes, (PetscEnum)memtyperequested,
1349396343dSjeremylt                           (PetscEnum *)&memtyperequested, &setmemtyperequest);
1359396343dSjeremylt   CHKERRQ(ierr);
136*153bcb04Svaleriabarra   localnodes = 1000;
137*153bcb04Svaleriabarra   ierr = PetscOptionsInt("-local_nodes",
138*153bcb04Svaleriabarra                          "Target number of locally owned nodes per process",
139*153bcb04Svaleriabarra                          NULL, localnodes, &localnodes, &userlnodes);
140*153bcb04Svaleriabarra   CHKERRQ(ierr);
1419396343dSjeremylt 
142cb32e2e7SValeria Barra   ierr = PetscOptionsEnd(); CHKERRQ(ierr);
143cb32e2e7SValeria Barra 
144cb32e2e7SValeria Barra   // Setup DM
145cb32e2e7SValeria Barra   if (read_mesh) {
146cb32e2e7SValeria Barra     ierr = DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, PETSC_TRUE, &dm);
147cb32e2e7SValeria Barra     CHKERRQ(ierr);
148cb32e2e7SValeria Barra   } else {
149*153bcb04Svaleriabarra     if (userlnodes) {
150*153bcb04Svaleriabarra       // Find a nicely composite number of elements no less than lnodes
151*153bcb04Svaleriabarra       for (lelem = PetscMax(1, localnodes / (degree*degree*degree)); ;
152*153bcb04Svaleriabarra            lelem++) {
153*153bcb04Svaleriabarra         Split3(lelem, melem, true);
154*153bcb04Svaleriabarra         if (Max3(melem) / Min3(melem) <= 2) break;
155*153bcb04Svaleriabarra       }
156*153bcb04Svaleriabarra     } else {
157*153bcb04Svaleriabarra       lelem = melem[0]*melem[1]*melem[2];
158*153bcb04Svaleriabarra     }
159cb32e2e7SValeria Barra     ierr = DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, melem, NULL,
160cb32e2e7SValeria Barra                                NULL, NULL, PETSC_TRUE, &dm); CHKERRQ(ierr);
161cb32e2e7SValeria Barra   }
162cb32e2e7SValeria Barra 
163cb32e2e7SValeria Barra   {
164cb32e2e7SValeria Barra     DM dmDist = NULL;
165cb32e2e7SValeria Barra     PetscPartitioner part;
166cb32e2e7SValeria Barra 
167cb32e2e7SValeria Barra     ierr = DMPlexGetPartitioner(dm, &part); CHKERRQ(ierr);
168cb32e2e7SValeria Barra     ierr = PetscPartitionerSetFromOptions(part); CHKERRQ(ierr);
169cb32e2e7SValeria Barra     ierr = DMPlexDistribute(dm, 0, NULL, &dmDist); CHKERRQ(ierr);
170cb32e2e7SValeria Barra     if (dmDist) {
171cb32e2e7SValeria Barra       ierr = DMDestroy(&dm); CHKERRQ(ierr);
172cb32e2e7SValeria Barra       dm  = dmDist;
173cb32e2e7SValeria Barra     }
174cb32e2e7SValeria Barra   }
175cb32e2e7SValeria Barra 
1769396343dSjeremylt   // Set up libCEED
1779396343dSjeremylt   CeedInit(ceedresource, &ceed);
1789396343dSjeremylt   CeedMemType memtypebackend;
1799396343dSjeremylt   CeedGetPreferredMemType(ceed, &memtypebackend);
1809396343dSjeremylt 
1819396343dSjeremylt   // Check memtype compatibility
1829396343dSjeremylt   if (!setmemtyperequest)
1839396343dSjeremylt     memtyperequested = memtypebackend;
1849396343dSjeremylt   else if (!petschavecuda && memtyperequested == CEED_MEM_DEVICE)
1859396343dSjeremylt     SETERRQ1(PETSC_COMM_WORLD, PETSC_ERR_SUP_SYS,
1869396343dSjeremylt              "PETSc was not built with CUDA. "
1879396343dSjeremylt              "Requested MemType CEED_MEM_DEVICE is not supported.", NULL);
1889396343dSjeremylt 
189cb32e2e7SValeria Barra   // Create DM
190199551f5Sjeremylt   ierr = SetupDMByDegree(dm, degree, ncompu, bpchoice);
191cb32e2e7SValeria Barra   CHKERRQ(ierr);
192cb32e2e7SValeria Barra 
193cb32e2e7SValeria Barra   // Create vectors
1949396343dSjeremylt   if (memtyperequested == CEED_MEM_DEVICE) {
1959396343dSjeremylt     ierr = DMSetVecType(dm, VECCUDA); CHKERRQ(ierr);
1969396343dSjeremylt   }
197cb32e2e7SValeria Barra   ierr = DMCreateGlobalVector(dm, &X); CHKERRQ(ierr);
198cb32e2e7SValeria Barra   ierr = VecGetLocalSize(X, &lsize); CHKERRQ(ierr);
199cb32e2e7SValeria Barra   ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
200cb32e2e7SValeria Barra   ierr = DMCreateLocalVector(dm, &Xloc); CHKERRQ(ierr);
201cb32e2e7SValeria Barra   ierr = VecGetSize(Xloc, &xlsize); CHKERRQ(ierr);
202cb32e2e7SValeria Barra   ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr);
203cb32e2e7SValeria Barra 
204cb32e2e7SValeria Barra   // Operator
205cb32e2e7SValeria Barra   ierr = PetscMalloc1(1, &userO); CHKERRQ(ierr);
206cb32e2e7SValeria Barra   ierr = MatCreateShell(comm, lsize, lsize, gsize, gsize,
207cb32e2e7SValeria Barra                         userO, &matO); CHKERRQ(ierr);
208cb32e2e7SValeria Barra   ierr = MatShellSetOperation(matO, MATOP_MULT,
209ce74dcefSjeremylt                               (void(*)(void))MatMult_Ceed); CHKERRQ(ierr);
210ce74dcefSjeremylt   ierr = MatShellSetOperation(matO, MATOP_GET_DIAGONAL,
211ce74dcefSjeremylt                               (void(*)(void))MatGetDiag); CHKERRQ(ierr);
2129396343dSjeremylt   if (memtyperequested == CEED_MEM_DEVICE) {
2139396343dSjeremylt     ierr = MatShellSetVecType(matO, VECCUDA); CHKERRQ(ierr);
2149396343dSjeremylt   }
2152d03409cSjeremylt 
216819eb1b3Sjeremylt   // Print summary
2170c59ef15SJeremy L Thompson   if (!test_mode) {
218cb32e2e7SValeria Barra     PetscInt P = degree + 1, Q = P + qextra;
2199396343dSjeremylt 
2202d03409cSjeremylt     const char *usedresource;
2212d03409cSjeremylt     CeedGetResource(ceed, &usedresource);
2229396343dSjeremylt 
2239396343dSjeremylt     VecType vectype;
2249396343dSjeremylt     ierr = VecGetType(X, &vectype); CHKERRQ(ierr);
2259396343dSjeremylt 
2263ce86546SJeremy L Thompson     ierr = PetscPrintf(comm,
2273ce86546SJeremy L Thompson                        "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n"
2289396343dSjeremylt                        "  PETSc:\n"
2299396343dSjeremylt                        "    PETSc Vec Type                     : %s\n"
2303ce86546SJeremy L Thompson                        "  libCEED:\n"
2313ce86546SJeremy L Thompson                        "    libCEED Backend                    : %s\n"
2329396343dSjeremylt                        "    libCEED Backend MemType            : %s\n"
2339396343dSjeremylt                        "    libCEED User Requested MemType     : %s\n"
2343ce86546SJeremy L Thompson                        "  Mesh:\n"
2353ce86546SJeremy L Thompson                        "    Number of 1D Basis Nodes (p)       : %d\n"
2363ce86546SJeremy L Thompson                        "    Number of 1D Quadrature Points (q) : %d\n"
237b52ddc88Sjeremylt                        "    Global nodes                       : %D\n"
238db419314Sjeremylt                        "    Owned nodes                        : %D\n"
239*153bcb04Svaleriabarra                        "    DoF per node                       : %D\n"
240*153bcb04Svaleriabarra                        "    Owned elements                     : %D\n",
2419396343dSjeremylt                        bpchoice+1, vectype, usedresource,
2429396343dSjeremylt                        CeedMemTypes[memtypebackend],
2439396343dSjeremylt                        (setmemtyperequest) ?
2449396343dSjeremylt                        CeedMemTypes[memtyperequested] : "none",
245*153bcb04Svaleriabarra                        P, Q, gsize/ncompu, lsize/ncompu, ncompu, lelem);
246*153bcb04Svaleriabarra     CHKERRQ(ierr);
2470c59ef15SJeremy L Thompson   }
2480c59ef15SJeremy L Thompson 
249cb32e2e7SValeria Barra   // Create RHS vector
250cb32e2e7SValeria Barra   ierr = VecDuplicate(Xloc, &rhsloc); CHKERRQ(ierr);
251cb32e2e7SValeria Barra   ierr = VecZeroEntries(rhsloc); CHKERRQ(ierr);
2529396343dSjeremylt   if (memtyperequested == CEED_MEM_HOST) {
253cb32e2e7SValeria Barra     ierr = VecGetArray(rhsloc, &r); CHKERRQ(ierr);
2549396343dSjeremylt   } else {
2559396343dSjeremylt     ierr = VecCUDAGetArray(rhsloc, &r); CHKERRQ(ierr);
2569396343dSjeremylt   }
257cb32e2e7SValeria Barra   CeedVectorCreate(ceed, xlsize, &rhsceed);
2589396343dSjeremylt   CeedVectorSetArray(rhsceed, memtyperequested, CEED_USE_POINTER, r);
2590c59ef15SJeremy L Thompson 
260cb32e2e7SValeria Barra   ierr = PetscMalloc1(1, &ceeddata); CHKERRQ(ierr);
261cb32e2e7SValeria Barra   ierr = SetupLibceedByDegree(dm, ceed, degree, dim, qextra,
262199551f5Sjeremylt                               ncompu, gsize, xlsize, bpchoice, ceeddata,
263cb32e2e7SValeria Barra                               true, rhsceed, &target); CHKERRQ(ierr);
2640c59ef15SJeremy L Thompson 
265cb32e2e7SValeria Barra   // Gather RHS
2669396343dSjeremylt   CeedVectorSyncArray(rhsceed, memtyperequested);
2679396343dSjeremylt   if (memtyperequested == CEED_MEM_HOST) {
268cb32e2e7SValeria Barra     ierr = VecRestoreArray(rhsloc, &r); CHKERRQ(ierr);
2699396343dSjeremylt   } else {
2709396343dSjeremylt     ierr = VecCUDARestoreArray(rhsloc, &r); CHKERRQ(ierr);
2719396343dSjeremylt   }
272cb32e2e7SValeria Barra   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
273cb32e2e7SValeria Barra   ierr = DMLocalToGlobalBegin(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
274cb32e2e7SValeria Barra   ierr = DMLocalToGlobalEnd(dm, rhsloc, ADD_VALUES, rhs); CHKERRQ(ierr);
275cb32e2e7SValeria Barra   CeedVectorDestroy(&rhsceed);
2760c59ef15SJeremy L Thompson 
277cb32e2e7SValeria Barra   // Create the error Q-function
278199551f5Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bpOptions[bpchoice].error,
279199551f5Sjeremylt                               bpOptions[bpchoice].errorfname, &qferror);
280c70bd2a0Sjeremylt   CeedQFunctionAddInput(qferror, "u", ncompu, CEED_EVAL_INTERP);
281c70bd2a0Sjeremylt   CeedQFunctionAddInput(qferror, "true_soln", ncompu, CEED_EVAL_NONE);
282c70bd2a0Sjeremylt   CeedQFunctionAddOutput(qferror, "error", ncompu, CEED_EVAL_NONE);
2830c59ef15SJeremy L Thompson 
2840c59ef15SJeremy L Thompson   // Create the error operator
285c70bd2a0Sjeremylt   CeedOperatorCreate(ceed, qferror, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
286c70bd2a0Sjeremylt                      &operror);
287c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "u", ceeddata->Erestrictu,
288a8d32208Sjeremylt                        ceeddata->basisu, CEED_VECTOR_ACTIVE);
289c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "true_soln", ceeddata->Erestrictui,
290a8d32208Sjeremylt                        CEED_BASIS_COLLOCATED, target);
291c70bd2a0Sjeremylt   CeedOperatorSetField(operror, "error", ceeddata->Erestrictui,
292a8d32208Sjeremylt                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
2930c59ef15SJeremy L Thompson 
2940c59ef15SJeremy L Thompson   // Set up Mat
295cb32e2e7SValeria Barra   userO->comm = comm;
296cb32e2e7SValeria Barra   userO->dm = dm;
297cb32e2e7SValeria Barra   userO->Xloc = Xloc;
298cb32e2e7SValeria Barra   ierr = VecDuplicate(Xloc, &userO->Yloc); CHKERRQ(ierr);
299cb32e2e7SValeria Barra   userO->xceed = ceeddata->xceed;
300cb32e2e7SValeria Barra   userO->yceed = ceeddata->yceed;
301c70bd2a0Sjeremylt   userO->op = ceeddata->opapply;
302cb32e2e7SValeria Barra   userO->ceed = ceed;
3039396343dSjeremylt   userO->memtype = memtyperequested;
3049396343dSjeremylt   if (memtyperequested == CEED_MEM_HOST) {
3059396343dSjeremylt     userO->VecGetArray = VecGetArray;
3069396343dSjeremylt     userO->VecGetArrayRead = VecGetArrayRead;
3079396343dSjeremylt     userO->VecRestoreArray = VecRestoreArray;
3089396343dSjeremylt     userO->VecRestoreArrayRead = VecRestoreArrayRead;
3099396343dSjeremylt   } else {
3109396343dSjeremylt     userO->VecGetArray = VecCUDAGetArray;
3119396343dSjeremylt     userO->VecGetArrayRead = VecCUDAGetArrayRead;
3129396343dSjeremylt     userO->VecRestoreArray = VecCUDARestoreArray;
3139396343dSjeremylt     userO->VecRestoreArrayRead = VecCUDARestoreArrayRead;
3149396343dSjeremylt   }
3150c59ef15SJeremy L Thompson 
3160c59ef15SJeremy L Thompson   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
3170c59ef15SJeremy L Thompson   {
3180c59ef15SJeremy L Thompson     PC pc;
3190c59ef15SJeremy L Thompson     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
320199551f5Sjeremylt     if (bpchoice == CEED_BP1 || bpchoice == CEED_BP2) {
3210c59ef15SJeremy L Thompson       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
3220c59ef15SJeremy L Thompson       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
3230c59ef15SJeremy L Thompson     } else {
3240c59ef15SJeremy L Thompson       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
3250c59ef15SJeremy L Thompson     }
3260c59ef15SJeremy L Thompson     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
3270c59ef15SJeremy L Thompson     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
3280c59ef15SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
3290c59ef15SJeremy L Thompson                             PETSC_DEFAULT); CHKERRQ(ierr);
3300c59ef15SJeremy L Thompson   }
331cb32e2e7SValeria Barra   ierr = KSPSetOperators(ksp, matO, matO); CHKERRQ(ierr);
332cb32e2e7SValeria Barra 
3330c59ef15SJeremy L Thompson   // First run, if benchmarking
3340c59ef15SJeremy L Thompson   if (benchmark_mode) {
3350c59ef15SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
3360c59ef15SJeremy L Thompson     CHKERRQ(ierr);
3370c59ef15SJeremy L Thompson     my_rt_start = MPI_Wtime();
3380c59ef15SJeremy L Thompson     ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
3390c59ef15SJeremy L Thompson     my_rt = MPI_Wtime() - my_rt_start;
340b5009ec9Sjeremylt     ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
341b5009ec9Sjeremylt     CHKERRQ(ierr);
3420c59ef15SJeremy L Thompson     // Set maxits based on first iteration timing
3430c59ef15SJeremy L Thompson     if (my_rt > 0.02) {
3440c59ef15SJeremy L Thompson       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5);
3450c59ef15SJeremy L Thompson       CHKERRQ(ierr);
3460c59ef15SJeremy L Thompson     } else {
3470c59ef15SJeremy L Thompson       ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20);
3480c59ef15SJeremy L Thompson       CHKERRQ(ierr);
3490c59ef15SJeremy L Thompson     }
3500c59ef15SJeremy L Thompson   }
351debcf919SJed Brown   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
352cb32e2e7SValeria Barra 
3530c59ef15SJeremy L Thompson   // Timed solve
354cb32e2e7SValeria Barra   ierr = VecZeroEntries(X); CHKERRQ(ierr);
3553b3d4a15SJeremy L Thompson   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
35609a940d7Sjeremylt 
35709a940d7Sjeremylt   // -- Performance logging
35809a940d7Sjeremylt   ierr = PetscLogStageRegister("Solve Stage", &solvestage); CHKERRQ(ierr);
35909a940d7Sjeremylt   ierr = PetscLogStagePush(solvestage); CHKERRQ(ierr);
36009a940d7Sjeremylt 
36109a940d7Sjeremylt   // -- Solve
3620c59ef15SJeremy L Thompson   my_rt_start = MPI_Wtime();
3630c59ef15SJeremy L Thompson   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
3640c59ef15SJeremy L Thompson   my_rt = MPI_Wtime() - my_rt_start;
365cb32e2e7SValeria Barra 
36609a940d7Sjeremylt   // -- Performance logging
36709a940d7Sjeremylt   ierr = PetscLogStagePop();
36809a940d7Sjeremylt 
369cb32e2e7SValeria Barra   // Output results
3700c59ef15SJeremy L Thompson   {
3710c59ef15SJeremy L Thompson     KSPType ksptype;
3720c59ef15SJeremy L Thompson     KSPConvergedReason reason;
3730c59ef15SJeremy L Thompson     PetscReal rnorm;
3740c59ef15SJeremy L Thompson     PetscInt its;
3750c59ef15SJeremy L Thompson     ierr = KSPGetType(ksp, &ksptype); CHKERRQ(ierr);
3760c59ef15SJeremy L Thompson     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
3770c59ef15SJeremy L Thompson     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
3780c59ef15SJeremy L Thompson     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
3790c59ef15SJeremy L Thompson     if (!test_mode || reason < 0 || rnorm > 1e-8) {
3803ce86546SJeremy L Thompson       ierr = PetscPrintf(comm,
3813ce86546SJeremy L Thompson                          "  KSP:\n"
3823ce86546SJeremy L Thompson                          "    KSP Type                           : %s\n"
3833ce86546SJeremy L Thompson                          "    KSP Convergence                    : %s\n"
3843ce86546SJeremy L Thompson                          "    Total KSP Iterations               : %D\n"
3853ce86546SJeremy L Thompson                          "    Final rnorm                        : %e\n",
3863ce86546SJeremy L Thompson                          ksptype, KSPConvergedReasons[reason], its,
3873ce86546SJeremy L Thompson                          (double)rnorm); CHKERRQ(ierr);
3880c59ef15SJeremy L Thompson     }
389cb32e2e7SValeria Barra     if (!test_mode) {
390cb32e2e7SValeria Barra       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
391cb32e2e7SValeria Barra     }
392cb32e2e7SValeria Barra     {
393cb32e2e7SValeria Barra       PetscReal maxerror;
394c70bd2a0Sjeremylt       ierr = ComputeErrorMax(userO, operror, X, target, &maxerror);
395cb32e2e7SValeria Barra       CHKERRQ(ierr);
396cb32e2e7SValeria Barra       PetscReal tol = 5e-2;
397cb32e2e7SValeria Barra       if (!test_mode || maxerror > tol) {
398b5009ec9Sjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
399b5009ec9Sjeremylt         CHKERRQ(ierr);
400b5009ec9Sjeremylt         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
401b5009ec9Sjeremylt         CHKERRQ(ierr);
4020c59ef15SJeremy L Thompson         ierr = PetscPrintf(comm,
403cb32e2e7SValeria Barra                            "    Pointwise Error (max)              : %e\n"
404cb32e2e7SValeria Barra                            "    CG Solve Time                      : %g (%g) sec\n",
405cb32e2e7SValeria Barra                            (double)maxerror, rt_max, rt_min); CHKERRQ(ierr);
4060c59ef15SJeremy L Thompson       }
4070c59ef15SJeremy L Thompson     }
408cb32e2e7SValeria Barra     if (benchmark_mode && (!test_mode)) {
4093ce86546SJeremy L Thompson       ierr = PetscPrintf(comm,
410cb32e2e7SValeria Barra                          "    DoFs/Sec in CG                     : %g (%g) million\n",
411cb32e2e7SValeria Barra                          1e-6*gsize*its/rt_max,
412cb32e2e7SValeria Barra                          1e-6*gsize*its/rt_min); CHKERRQ(ierr);
4130c59ef15SJeremy L Thompson     }
4140c59ef15SJeremy L Thompson   }
4150c59ef15SJeremy L Thompson 
4166bd9afcaSjeremylt   if (write_solution) {
4176bd9afcaSjeremylt     PetscViewer vtkviewersoln;
4186bd9afcaSjeremylt 
4196bd9afcaSjeremylt     ierr = PetscViewerCreate(comm, &vtkviewersoln); CHKERRQ(ierr);
4206bd9afcaSjeremylt     ierr = PetscViewerSetType(vtkviewersoln, PETSCVIEWERVTK); CHKERRQ(ierr);
4216bd9afcaSjeremylt     ierr = PetscViewerFileSetName(vtkviewersoln, "solution.vtk"); CHKERRQ(ierr);
4226bd9afcaSjeremylt     ierr = VecView(X, vtkviewersoln); CHKERRQ(ierr);
4236bd9afcaSjeremylt     ierr = PetscViewerDestroy(&vtkviewersoln); CHKERRQ(ierr);
4246bd9afcaSjeremylt   }
4256bd9afcaSjeremylt 
426cb32e2e7SValeria Barra   // Cleanup
427cb32e2e7SValeria Barra   ierr = VecDestroy(&X); CHKERRQ(ierr);
428cb32e2e7SValeria Barra   ierr = VecDestroy(&Xloc); CHKERRQ(ierr);
429cb32e2e7SValeria Barra   ierr = VecDestroy(&userO->Yloc); CHKERRQ(ierr);
430cb32e2e7SValeria Barra   ierr = MatDestroy(&matO); CHKERRQ(ierr);
431cb32e2e7SValeria Barra   ierr = PetscFree(userO); CHKERRQ(ierr);
432cb32e2e7SValeria Barra   ierr = CeedDataDestroy(0, ceeddata); CHKERRQ(ierr);
433cb32e2e7SValeria Barra   ierr = DMDestroy(&dm); CHKERRQ(ierr);
434cb32e2e7SValeria Barra 
4350c59ef15SJeremy L Thompson   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
4360c59ef15SJeremy L Thompson   ierr = VecDestroy(&rhsloc); CHKERRQ(ierr);
4370c59ef15SJeremy L Thompson   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
4380c59ef15SJeremy L Thompson   CeedVectorDestroy(&target);
439c70bd2a0Sjeremylt   CeedQFunctionDestroy(&qferror);
440c70bd2a0Sjeremylt   CeedOperatorDestroy(&operror);
4410c59ef15SJeremy L Thompson   CeedDestroy(&ceed);
4420c59ef15SJeremy L Thompson   return PetscFinalize();
4430c59ef15SJeremy L Thompson }
444