xref: /libCEED/examples/petsc/bpssphere.c (revision 356036fa84f714fa73ef64c9a80ce2028dde816f)
13d8e8822SJeremy L Thompson // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
23d8e8822SJeremy L Thompson // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3ed264d09SValeria Barra //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5ed264d09SValeria Barra //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7ed264d09SValeria Barra 
8ed264d09SValeria Barra //                        libCEED + PETSc Example: CEED BPs
9ed264d09SValeria Barra //
10ea61e9acSJeremy L Thompson // This example demonstrates a simple usage of libCEED with PETSc to solve the CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps, on
11ea61e9acSJeremy L Thompson // a closed surface, such as the one of a discrete sphere.
12ed264d09SValeria Barra //
13ed264d09SValeria Barra // The code uses higher level communication protocols in DMPlex.
14ed264d09SValeria Barra //
15ed264d09SValeria Barra // Build with:
16ed264d09SValeria Barra //
17ed264d09SValeria Barra //     make bpssphere [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
18ed264d09SValeria Barra //
19ed264d09SValeria Barra // Sample runs:
20ed264d09SValeria Barra //
21ed264d09SValeria Barra //     bpssphere -problem bp1 -degree 3
2228688798Sjeremylt //     bpssphere -problem bp2 -degree 3
2328688798Sjeremylt //     bpssphere -problem bp3 -degree 3
2428688798Sjeremylt //     bpssphere -problem bp4 -degree 3
2528688798Sjeremylt //     bpssphere -problem bp5 -degree 3 -ceed /cpu/self
2628688798Sjeremylt //     bpssphere -problem bp6 -degree 3 -ceed /gpu/cuda
27ed264d09SValeria Barra //
28587be3cdSvaleriabarra //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3 -dm_refine 2
29ed264d09SValeria Barra 
30ed264d09SValeria Barra /// @file
31ed264d09SValeria Barra /// CEED BPs example using PETSc with DMPlex
32ea61e9acSJeremy L Thompson /// See bps.c for a "raw" implementation using a structured grid and bpsdmplex.c for an implementation using an unstructured grid.
33ed264d09SValeria Barra static const char help[] = "Solve CEED BPs on a sphere using DMPlex in PETSc\n";
34ed264d09SValeria Barra 
352b730f8bSJeremy L Thompson #include "bpssphere.h"
362b730f8bSJeremy L Thompson 
37636cccdbSjeremylt #include <ceed.h>
38636cccdbSjeremylt #include <petscdmplex.h>
39636cccdbSjeremylt #include <petscksp.h>
402b730f8bSJeremy L Thompson #include <stdbool.h>
412b730f8bSJeremy L Thompson #include <string.h>
42636cccdbSjeremylt 
432b730f8bSJeremy L Thompson #include "include/libceedsetup.h"
442b730f8bSJeremy L Thompson #include "include/matops.h"
45636cccdbSjeremylt #include "include/petscutils.h"
46b8962995SJeremy L Thompson #include "include/petscversion.h"
472b730f8bSJeremy L Thompson #include "include/sphereproblemdata.h"
48636cccdbSjeremylt 
49ed264d09SValeria Barra int main(int argc, char **argv) {
50ed264d09SValeria Barra   MPI_Comm             comm;
512b730f8bSJeremy L Thompson   char                 ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self", filename[PETSC_MAX_PATH_LEN];
52ed264d09SValeria Barra   double               my_rt_start, my_rt, rt_min, rt_max;
532b730f8bSJeremy L Thompson   PetscInt             degree = 3, q_extra, l_size, g_size, topo_dim = 2, num_comp_x = 3, num_comp_u = 1, xl_size;
54ed264d09SValeria Barra   PetscScalar         *r;
5509a940d7Sjeremylt   PetscBool            test_mode, benchmark_mode, read_mesh, write_solution, simplex;
569b072555Sjeremylt   PetscLogStage        solve_stage;
579b072555Sjeremylt   Vec                  X, X_loc, rhs, rhs_loc;
589b072555Sjeremylt   Mat                  mat_O;
59ed264d09SValeria Barra   KSP                  ksp;
60ed264d09SValeria Barra   DM                   dm;
616c88e6a2Srezgarshakeri   OperatorApplyContext op_apply_ctx, op_error_ctx;
62ed264d09SValeria Barra   Ceed                 ceed;
639b072555Sjeremylt   CeedData             ceed_data;
649b072555Sjeremylt   CeedQFunction        qf_error;
659b072555Sjeremylt   CeedOperator         op_error;
669b072555Sjeremylt   CeedVector           rhs_ceed, target;
679b072555Sjeremylt   BPType               bp_choice;
689b072555Sjeremylt   VecType              vec_type;
699b072555Sjeremylt   PetscMemType         mem_type;
70ed264d09SValeria Barra 
712b730f8bSJeremy L Thompson   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
72ed264d09SValeria Barra   comm = PETSC_COMM_WORLD;
73ed264d09SValeria Barra 
74ed264d09SValeria Barra   // Read command line options
7567490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
769b072555Sjeremylt   bp_choice = CEED_BP1;
772b730f8bSJeremy L Thompson   PetscCall(PetscOptionsEnum("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice, NULL));
789b072555Sjeremylt   num_comp_u = bp_options[bp_choice].num_comp_u;
79ed264d09SValeria Barra   test_mode  = PETSC_FALSE;
802b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL));
81ed264d09SValeria Barra   benchmark_mode = PETSC_FALSE;
822b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-benchmark", "Benchmarking mode (prints benchmark statistics)", NULL, benchmark_mode, &benchmark_mode, NULL));
83ed264d09SValeria Barra   write_solution = PETSC_FALSE;
842b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, write_solution, &write_solution, NULL));
85ed264d09SValeria Barra   degree = test_mode ? 3 : 2;
862b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, degree, &degree, NULL));
879b072555Sjeremylt   q_extra = bp_options[bp_choice].q_extra;
882b730f8bSJeremy L Thompson   PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra, NULL));
892b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resource, sizeof(ceed_resource), NULL));
90ed264d09SValeria Barra   read_mesh = PETSC_FALSE;
912b730f8bSJeremy L Thompson   PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &read_mesh));
92ed264d09SValeria Barra   simplex = PETSC_FALSE;
932b730f8bSJeremy L Thompson   PetscCall(PetscOptionsBool("-simplex", "Use simplices, or tensor product cells", NULL, simplex, &simplex, NULL));
9467490bc6SJeremy L Thompson   PetscOptionsEnd();
95ed264d09SValeria Barra 
96ed264d09SValeria Barra   // Setup DM
97ed264d09SValeria Barra   if (read_mesh) {
982b730f8bSJeremy L Thompson     PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, NULL, PETSC_TRUE, &dm));
99ed264d09SValeria Barra   } else {
100ea61e9acSJeremy L Thompson     // Create the mesh as a 0-refined sphere.
101ea61e9acSJeremy L Thompson     // This will create a cubic surface, not a box, and will snap to the unit sphere upon refinement.
1022b730f8bSJeremy L Thompson     PetscCall(DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topo_dim, simplex, 1., &dm));
103ed264d09SValeria Barra     // Set the object name
1042b730f8bSJeremy L Thompson     PetscCall(PetscObjectSetName((PetscObject)dm, "Sphere"));
105ed264d09SValeria Barra     // Refine DMPlex with uniform refinement using runtime option -dm_refine
1062b730f8bSJeremy L Thompson     PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE));
1073fc8a154SJed Brown   }
1082b730f8bSJeremy L Thompson   PetscCall(DMSetFromOptions(dm));
109ed264d09SValeria Barra   // View DMPlex via runtime option
1102b730f8bSJeremy L Thompson   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
111ed264d09SValeria Barra 
112ed264d09SValeria Barra   // Create DM
1132b730f8bSJeremy L Thompson   PetscCall(SetupDMByDegree(dm, degree, q_extra, num_comp_u, topo_dim, false));
114ed264d09SValeria Barra 
115ed264d09SValeria Barra   // Create vectors
1162b730f8bSJeremy L Thompson   PetscCall(DMCreateGlobalVector(dm, &X));
1172b730f8bSJeremy L Thompson   PetscCall(VecGetLocalSize(X, &l_size));
1182b730f8bSJeremy L Thompson   PetscCall(VecGetSize(X, &g_size));
1192b730f8bSJeremy L Thompson   PetscCall(DMCreateLocalVector(dm, &X_loc));
1202b730f8bSJeremy L Thompson   PetscCall(VecGetSize(X_loc, &xl_size));
1212b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(X, &rhs));
122ed264d09SValeria Barra 
123ed264d09SValeria Barra   // Operator
1242b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(1, &op_apply_ctx));
1252b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(1, &op_error_ctx));
1262b730f8bSJeremy L Thompson   PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, op_apply_ctx, &mat_O));
1272b730f8bSJeremy L Thompson   PetscCall(MatShellSetOperation(mat_O, MATOP_MULT, (void (*)(void))MatMult_Ceed));
128ed264d09SValeria Barra 
129ed264d09SValeria Barra   // Set up libCEED
1309b072555Sjeremylt   CeedInit(ceed_resource, &ceed);
1319b072555Sjeremylt   CeedMemType mem_type_backend;
1329b072555Sjeremylt   CeedGetPreferredMemType(ceed, &mem_type_backend);
133e83e87a5Sjeremylt 
1342b730f8bSJeremy L Thompson   PetscCall(DMGetVecType(dm, &vec_type));
1359b072555Sjeremylt   if (!vec_type) {  // Not yet set by user -dm_vec_type
1369b072555Sjeremylt     switch (mem_type_backend) {
1372b730f8bSJeremy L Thompson       case CEED_MEM_HOST:
1382b730f8bSJeremy L Thompson         vec_type = VECSTANDARD;
1392b730f8bSJeremy L Thompson         break;
140e83e87a5Sjeremylt       case CEED_MEM_DEVICE: {
141e83e87a5Sjeremylt         const char *resolved;
142e83e87a5Sjeremylt         CeedGetResource(ceed, &resolved);
1439b072555Sjeremylt         if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
1442b730f8bSJeremy L Thompson         else if (strstr(resolved, "/gpu/hip/occa")) vec_type = VECSTANDARD;  // https://github.com/CEED/libCEED/issues/678
1459b072555Sjeremylt         else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
1469b072555Sjeremylt         else vec_type = VECSTANDARD;
147e83e87a5Sjeremylt       }
148e83e87a5Sjeremylt     }
1492b730f8bSJeremy L Thompson     PetscCall(DMSetVecType(dm, vec_type));
150e83e87a5Sjeremylt   }
151ed264d09SValeria Barra 
152ed264d09SValeria Barra   // Print summary
153ed264d09SValeria Barra   if (!test_mode) {
1549b072555Sjeremylt     PetscInt    P = degree + 1, Q = P + q_extra;
1559b072555Sjeremylt     const char *used_resource;
1569b072555Sjeremylt     CeedGetResource(ceed, &used_resource);
1572b730f8bSJeremy L Thompson     PetscCall(PetscPrintf(comm,
1582b730f8bSJeremy L Thompson                           "\n-- CEED Benchmark Problem %" CeedInt_FMT " on the Sphere -- libCEED + PETSc --\n"
159ed264d09SValeria Barra                           "  libCEED:\n"
160ed264d09SValeria Barra                           "    libCEED Backend                         : %s\n"
161e83e87a5Sjeremylt                           "    libCEED Backend MemType                 : %s\n"
162ed264d09SValeria Barra                           "  Mesh:\n"
163751eb813Srezgarshakeri                           "    Solution Order (P)                      : %" CeedInt_FMT "\n"
164751eb813Srezgarshakeri                           "    Quadrature  Order (Q)                   : %" CeedInt_FMT "\n"
16551ad7d5bSrezgarshakeri                           "    Additional quadrature points (q_extra)  : %" CeedInt_FMT "\n"
16608140895SJed Brown                           "    Global nodes                            : %" PetscInt_FMT "\n",
1672b730f8bSJeremy L Thompson                           bp_choice + 1, ceed_resource, CeedMemTypes[mem_type_backend], P, Q, q_extra, g_size / num_comp_u));
168ed264d09SValeria Barra   }
169ed264d09SValeria Barra 
170ed264d09SValeria Barra   // Create RHS vector
1712b730f8bSJeremy L Thompson   PetscCall(VecDuplicate(X_loc, &rhs_loc));
1722b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(rhs_loc));
1732b730f8bSJeremy L Thompson   PetscCall(VecGetArrayAndMemType(rhs_loc, &r, &mem_type));
1749b072555Sjeremylt   CeedVectorCreate(ceed, xl_size, &rhs_ceed);
1759b072555Sjeremylt   CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r);
176ed264d09SValeria Barra 
177ed264d09SValeria Barra   // Setup libCEED's objects
1782b730f8bSJeremy L Thompson   PetscCall(PetscMalloc1(1, &ceed_data));
1792b730f8bSJeremy L Thompson   PetscCall(SetupLibceedByDegree(dm, ceed, degree, topo_dim, q_extra, num_comp_x, num_comp_u, g_size, xl_size, bp_options[bp_choice], ceed_data, true,
1802b730f8bSJeremy L Thompson                                  rhs_ceed, &target));
181ed264d09SValeria Barra 
182ed264d09SValeria Barra   // Gather RHS
1839b072555Sjeremylt   CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL);
1842b730f8bSJeremy L Thompson   PetscCall(VecRestoreArrayAndMemType(rhs_loc, &r));
1852b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(rhs));
1862b730f8bSJeremy L Thompson   PetscCall(DMLocalToGlobal(dm, rhs_loc, ADD_VALUES, rhs));
1879b072555Sjeremylt   CeedVectorDestroy(&rhs_ceed);
188ed264d09SValeria Barra 
189ed264d09SValeria Barra   // Create the error Q-function
1902b730f8bSJeremy L Thompson   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error, bp_options[bp_choice].error_loc, &qf_error);
1919b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
1929b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
1932b730f8bSJeremy L Thompson   CeedQFunctionAddInput(qf_error, "qdata", ceed_data->q_data_size, CEED_EVAL_NONE);
19438f32c05Srezgarshakeri   CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_INTERP);
195ed264d09SValeria Barra 
196ed264d09SValeria Barra   // Create the error operator
1979b072555Sjeremylt   CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error);
1982b730f8bSJeremy L Thompson   CeedOperatorSetField(op_error, "u", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE);
199*356036faSJeremy L Thompson   CeedOperatorSetField(op_error, "true_soln", ceed_data->elem_restr_u_i, CEED_BASIS_NONE, target);
200*356036faSJeremy L Thompson   CeedOperatorSetField(op_error, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data);
2012b730f8bSJeremy L Thompson   CeedOperatorSetField(op_error, "error", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE);
202ed264d09SValeria Barra 
2036c88e6a2Srezgarshakeri   // Set up apply operator context
2042b730f8bSJeremy L Thompson   PetscCall(SetupApplyOperatorCtx(comm, dm, ceed, ceed_data, X_loc, op_apply_ctx));
205ed264d09SValeria Barra 
206ed264d09SValeria Barra   // Setup solver
2072b730f8bSJeremy L Thompson   PetscCall(KSPCreate(comm, &ksp));
208ed264d09SValeria Barra   {
209ed264d09SValeria Barra     PC pc;
2102b730f8bSJeremy L Thompson     PetscCall(KSPGetPC(ksp, &pc));
2119b072555Sjeremylt     if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) {
2122b730f8bSJeremy L Thompson       PetscCall(PCSetType(pc, PCJACOBI));
2132b730f8bSJeremy L Thompson       PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
214ed264d09SValeria Barra     } else {
2152b730f8bSJeremy L Thompson       PetscCall(PCSetType(pc, PCNONE));
216ed264d09SValeria Barra       MatNullSpace nullspace;
217ed264d09SValeria Barra 
2182b730f8bSJeremy L Thompson       PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, 0, &nullspace));
2192b730f8bSJeremy L Thompson       PetscCall(MatSetNullSpace(mat_O, nullspace));
2202b730f8bSJeremy L Thompson       PetscCall(MatNullSpaceDestroy(&nullspace));
221ed264d09SValeria Barra     }
2222b730f8bSJeremy L Thompson     PetscCall(KSPSetType(ksp, KSPCG));
2232b730f8bSJeremy L Thompson     PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
2242b730f8bSJeremy L Thompson     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
225ed264d09SValeria Barra   }
2262b730f8bSJeremy L Thompson   PetscCall(KSPSetFromOptions(ksp));
2272b730f8bSJeremy L Thompson   PetscCall(KSPSetOperators(ksp, mat_O, mat_O));
228ed264d09SValeria Barra 
229ed264d09SValeria Barra   // First run, if benchmarking
230ed264d09SValeria Barra   if (benchmark_mode) {
2312b730f8bSJeremy L Thompson     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1));
232ed264d09SValeria Barra     my_rt_start = MPI_Wtime();
2332b730f8bSJeremy L Thompson     PetscCall(KSPSolve(ksp, rhs, X));
234ed264d09SValeria Barra     my_rt = MPI_Wtime() - my_rt_start;
2352b730f8bSJeremy L Thompson     PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm));
236ed264d09SValeria Barra     // Set maxits based on first iteration timing
237ed264d09SValeria Barra     if (my_rt > 0.02) {
2382b730f8bSJeremy L Thompson       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5));
239ed264d09SValeria Barra     } else {
2402b730f8bSJeremy L Thompson       PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20));
241ed264d09SValeria Barra     }
242ed264d09SValeria Barra   }
243ed264d09SValeria Barra 
244ed264d09SValeria Barra   // Timed solve
2452b730f8bSJeremy L Thompson   PetscCall(VecZeroEntries(X));
2462b730f8bSJeremy L Thompson   PetscCall(PetscBarrier((PetscObject)ksp));
24709a940d7Sjeremylt 
24809a940d7Sjeremylt   // -- Performance logging
2492b730f8bSJeremy L Thompson   PetscCall(PetscLogStageRegister("Solve Stage", &solve_stage));
2502b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePush(solve_stage));
25109a940d7Sjeremylt 
25209a940d7Sjeremylt   // -- Solve
253ed264d09SValeria Barra   my_rt_start = MPI_Wtime();
2542b730f8bSJeremy L Thompson   PetscCall(KSPSolve(ksp, rhs, X));
255ed264d09SValeria Barra   my_rt = MPI_Wtime() - my_rt_start;
256ed264d09SValeria Barra 
25709a940d7Sjeremylt   // -- Performance logging
2582b730f8bSJeremy L Thompson   PetscCall(PetscLogStagePop());
25909a940d7Sjeremylt 
260ed264d09SValeria Barra   // Output results
261ed264d09SValeria Barra   {
2629b072555Sjeremylt     KSPType            ksp_type;
263ed264d09SValeria Barra     KSPConvergedReason reason;
264ed264d09SValeria Barra     PetscReal          rnorm;
265ed264d09SValeria Barra     PetscInt           its;
2662b730f8bSJeremy L Thompson     PetscCall(KSPGetType(ksp, &ksp_type));
2672b730f8bSJeremy L Thompson     PetscCall(KSPGetConvergedReason(ksp, &reason));
2682b730f8bSJeremy L Thompson     PetscCall(KSPGetIterationNumber(ksp, &its));
2692b730f8bSJeremy L Thompson     PetscCall(KSPGetResidualNorm(ksp, &rnorm));
270ed264d09SValeria Barra     if (!test_mode || reason < 0 || rnorm > 1e-8) {
2712b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm,
272ed264d09SValeria Barra                             "  KSP:\n"
273ed264d09SValeria Barra                             "    KSP Type                                : %s\n"
274ed264d09SValeria Barra                             "    KSP Convergence                         : %s\n"
275a9b2c5ddSrezgarshakeri                             "    Total KSP Iterations                    : %" PetscInt_FMT "\n"
276ed264d09SValeria Barra                             "    Final rnorm                             : %e\n",
2772b730f8bSJeremy L Thompson                             ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
278ed264d09SValeria Barra     }
279ed264d09SValeria Barra     if (!test_mode) {
2802b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm, "  Performance:\n"));
281ed264d09SValeria Barra     }
282ed264d09SValeria Barra     {
2836c88e6a2Srezgarshakeri       // Set up error operator context
2842b730f8bSJeremy L Thompson       PetscCall(SetupErrorOperatorCtx(comm, dm, ceed, ceed_data, X_loc, op_error, op_error_ctx));
28538f32c05Srezgarshakeri       PetscScalar l2_error;
2862b730f8bSJeremy L Thompson       PetscCall(ComputeL2Error(X, &l2_error, op_error_ctx));
287587be3cdSvaleriabarra       PetscReal tol = 5e-4;
28838f32c05Srezgarshakeri       if (!test_mode || l2_error > tol) {
2892b730f8bSJeremy L Thompson         PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm));
2902b730f8bSJeremy L Thompson         PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm));
2912b730f8bSJeremy L Thompson         PetscCall(PetscPrintf(comm,
29238f32c05Srezgarshakeri                               "    L2 Error                                : %e\n"
293ed264d09SValeria Barra                               "    CG Solve Time                           : %g (%g) sec\n",
2942b730f8bSJeremy L Thompson                               (double)l2_error, rt_max, rt_min));
295ed264d09SValeria Barra       }
296ed264d09SValeria Barra     }
297ed264d09SValeria Barra     if (benchmark_mode && (!test_mode)) {
2982b730f8bSJeremy L Thompson       PetscCall(PetscPrintf(comm, "    DoFs/Sec in CG                            : %g (%g) million\n", 1e-6 * g_size * its / rt_max,
2992b730f8bSJeremy L Thompson                             1e-6 * g_size * its / rt_min));
300ed264d09SValeria Barra     }
301ed264d09SValeria Barra   }
302ed264d09SValeria Barra 
303ed264d09SValeria Barra   // Output solution
304ed264d09SValeria Barra   if (write_solution) {
3059b072555Sjeremylt     PetscViewer vtk_viewer_soln;
306ed264d09SValeria Barra 
3072b730f8bSJeremy L Thompson     PetscCall(PetscViewerCreate(comm, &vtk_viewer_soln));
3082b730f8bSJeremy L Thompson     PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK));
3092b730f8bSJeremy L Thompson     PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"));
3102b730f8bSJeremy L Thompson     PetscCall(VecView(X, vtk_viewer_soln));
3112b730f8bSJeremy L Thompson     PetscCall(PetscViewerDestroy(&vtk_viewer_soln));
312ed264d09SValeria Barra   }
313ed264d09SValeria Barra 
314ed264d09SValeria Barra   // Cleanup
3152b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&X));
3162b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&X_loc));
3172b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&op_apply_ctx->Y_loc));
3182b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&op_error_ctx->Y_loc));
3192b730f8bSJeremy L Thompson   PetscCall(MatDestroy(&mat_O));
3202b730f8bSJeremy L Thompson   PetscCall(PetscFree(op_apply_ctx));
3212b730f8bSJeremy L Thompson   PetscCall(PetscFree(op_error_ctx));
3222b730f8bSJeremy L Thompson   PetscCall(CeedDataDestroy(0, ceed_data));
3232b730f8bSJeremy L Thompson   PetscCall(DMDestroy(&dm));
324ed264d09SValeria Barra 
3252b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&rhs));
3262b730f8bSJeremy L Thompson   PetscCall(VecDestroy(&rhs_loc));
3272b730f8bSJeremy L Thompson   PetscCall(KSPDestroy(&ksp));
328ed264d09SValeria Barra   CeedVectorDestroy(&target);
3299b072555Sjeremylt   CeedQFunctionDestroy(&qf_error);
3309b072555Sjeremylt   CeedOperatorDestroy(&op_error);
331ed264d09SValeria Barra   CeedDestroy(&ceed);
332ed264d09SValeria Barra   return PetscFinalize();
333ed264d09SValeria Barra }
334