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 // 10ed264d09SValeria Barra // This example demonstrates a simple usage of libCEED with PETSc to solve the 11ed264d09SValeria Barra // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps, 12ed264d09SValeria Barra // on a closed surface, such as the one of a discrete sphere. 13ed264d09SValeria Barra // 14ed264d09SValeria Barra // The code uses higher level communication protocols in DMPlex. 15ed264d09SValeria Barra // 16ed264d09SValeria Barra // Build with: 17ed264d09SValeria Barra // 18ed264d09SValeria Barra // make bpssphere [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 19ed264d09SValeria Barra // 20ed264d09SValeria Barra // Sample runs: 21ed264d09SValeria Barra // 22ed264d09SValeria Barra // bpssphere -problem bp1 -degree 3 2328688798Sjeremylt // bpssphere -problem bp2 -degree 3 2428688798Sjeremylt // bpssphere -problem bp3 -degree 3 2528688798Sjeremylt // bpssphere -problem bp4 -degree 3 2628688798Sjeremylt // bpssphere -problem bp5 -degree 3 -ceed /cpu/self 2728688798Sjeremylt // bpssphere -problem bp6 -degree 3 -ceed /gpu/cuda 28ed264d09SValeria Barra // 29587be3cdSvaleriabarra //TESTARGS -ceed {ceed_resource} -test -problem bp3 -degree 3 -dm_refine 2 30ed264d09SValeria Barra 31ed264d09SValeria Barra /// @file 32ed264d09SValeria Barra /// CEED BPs example using PETSc with DMPlex 33ed264d09SValeria Barra /// See bps.c for a "raw" implementation using a structured grid. 34ed264d09SValeria Barra /// and bpsdmplex.c for an implementation using an unstructured grid. 35ed264d09SValeria Barra static const char help[] = "Solve CEED BPs on a sphere using DMPlex in PETSc\n"; 36ed264d09SValeria Barra 37*2b730f8bSJeremy L Thompson #include "bpssphere.h" 38*2b730f8bSJeremy L Thompson 39636cccdbSjeremylt #include <ceed.h> 40636cccdbSjeremylt #include <petsc.h> 41636cccdbSjeremylt #include <petscdmplex.h> 42636cccdbSjeremylt #include <petscksp.h> 43*2b730f8bSJeremy L Thompson #include <stdbool.h> 44*2b730f8bSJeremy L Thompson #include <string.h> 45636cccdbSjeremylt 46*2b730f8bSJeremy L Thompson #include "include/libceedsetup.h" 47*2b730f8bSJeremy L Thompson #include "include/matops.h" 48636cccdbSjeremylt #include "include/petscutils.h" 49b8962995SJeremy L Thompson #include "include/petscversion.h" 50*2b730f8bSJeremy L Thompson #include "include/sphereproblemdata.h" 51636cccdbSjeremylt 52636cccdbSjeremylt #if PETSC_VERSION_LT(3, 12, 0) 53636cccdbSjeremylt #ifdef PETSC_HAVE_CUDA 54636cccdbSjeremylt #include <petsccuda.h> 55636cccdbSjeremylt // Note: With PETSc prior to version 3.12.0, providing the source path to 56636cccdbSjeremylt // include 'cublas_v2.h' will be needed to use 'petsccuda.h'. 57636cccdbSjeremylt #endif 58636cccdbSjeremylt #endif 59ed264d09SValeria Barra 60ed264d09SValeria Barra int main(int argc, char **argv) { 61ed264d09SValeria Barra MPI_Comm comm; 62*2b730f8bSJeremy L Thompson char ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self", filename[PETSC_MAX_PATH_LEN]; 63ed264d09SValeria Barra double my_rt_start, my_rt, rt_min, rt_max; 64*2b730f8bSJeremy L Thompson PetscInt degree = 3, q_extra, l_size, g_size, topo_dim = 2, num_comp_x = 3, num_comp_u = 1, xl_size; 65ed264d09SValeria Barra PetscScalar *r; 6609a940d7Sjeremylt PetscBool test_mode, benchmark_mode, read_mesh, write_solution, simplex; 679b072555Sjeremylt PetscLogStage solve_stage; 689b072555Sjeremylt Vec X, X_loc, rhs, rhs_loc; 699b072555Sjeremylt Mat mat_O; 70ed264d09SValeria Barra KSP ksp; 71ed264d09SValeria Barra DM dm; 726c88e6a2Srezgarshakeri OperatorApplyContext op_apply_ctx, op_error_ctx; 73ed264d09SValeria Barra Ceed ceed; 749b072555Sjeremylt CeedData ceed_data; 759b072555Sjeremylt CeedQFunction qf_error; 769b072555Sjeremylt CeedOperator op_error; 779b072555Sjeremylt CeedVector rhs_ceed, target; 789b072555Sjeremylt BPType bp_choice; 799b072555Sjeremylt VecType vec_type; 809b072555Sjeremylt PetscMemType mem_type; 81ed264d09SValeria Barra 82*2b730f8bSJeremy L Thompson PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 83ed264d09SValeria Barra comm = PETSC_COMM_WORLD; 84ed264d09SValeria Barra 85ed264d09SValeria Barra // Read command line options 8667490bc6SJeremy L Thompson PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); 879b072555Sjeremylt bp_choice = CEED_BP1; 88*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsEnum("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice, NULL)); 899b072555Sjeremylt num_comp_u = bp_options[bp_choice].num_comp_u; 90ed264d09SValeria Barra test_mode = PETSC_FALSE; 91*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL)); 92ed264d09SValeria Barra benchmark_mode = PETSC_FALSE; 93*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-benchmark", "Benchmarking mode (prints benchmark statistics)", NULL, benchmark_mode, &benchmark_mode, NULL)); 94ed264d09SValeria Barra write_solution = PETSC_FALSE; 95*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, write_solution, &write_solution, NULL)); 96ed264d09SValeria Barra degree = test_mode ? 3 : 2; 97*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, degree, °ree, NULL)); 989b072555Sjeremylt q_extra = bp_options[bp_choice].q_extra; 99*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra, NULL)); 100*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resource, sizeof(ceed_resource), NULL)); 101ed264d09SValeria Barra read_mesh = PETSC_FALSE; 102*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &read_mesh)); 103ed264d09SValeria Barra simplex = PETSC_FALSE; 104*2b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-simplex", "Use simplices, or tensor product cells", NULL, simplex, &simplex, NULL)); 10567490bc6SJeremy L Thompson PetscOptionsEnd(); 106ed264d09SValeria Barra 107ed264d09SValeria Barra // Setup DM 108ed264d09SValeria Barra if (read_mesh) { 109*2b730f8bSJeremy L Thompson PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, NULL, PETSC_TRUE, &dm)); 110ed264d09SValeria Barra } else { 1113fc8a154SJed Brown // Create the mesh as a 0-refined sphere. This will create a cubic surface, 1123fc8a154SJed Brown // not a box, and will snap to the unit sphere upon refinement. 113*2b730f8bSJeremy L Thompson PetscCall(DMPlexCreateSphereMesh(PETSC_COMM_WORLD, topo_dim, simplex, 1., &dm)); 114ed264d09SValeria Barra // Set the object name 115*2b730f8bSJeremy L Thompson PetscCall(PetscObjectSetName((PetscObject)dm, "Sphere")); 116ed264d09SValeria Barra // Refine DMPlex with uniform refinement using runtime option -dm_refine 117*2b730f8bSJeremy L Thompson PetscCall(DMPlexSetRefinementUniform(dm, PETSC_TRUE)); 1183fc8a154SJed Brown } 119*2b730f8bSJeremy L Thompson PetscCall(DMSetFromOptions(dm)); 120ed264d09SValeria Barra // View DMPlex via runtime option 121*2b730f8bSJeremy L Thompson PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 122ed264d09SValeria Barra 123ed264d09SValeria Barra // Create DM 124*2b730f8bSJeremy L Thompson PetscCall(SetupDMByDegree(dm, degree, q_extra, num_comp_u, topo_dim, false)); 125ed264d09SValeria Barra 126ed264d09SValeria Barra // Create vectors 127*2b730f8bSJeremy L Thompson PetscCall(DMCreateGlobalVector(dm, &X)); 128*2b730f8bSJeremy L Thompson PetscCall(VecGetLocalSize(X, &l_size)); 129*2b730f8bSJeremy L Thompson PetscCall(VecGetSize(X, &g_size)); 130*2b730f8bSJeremy L Thompson PetscCall(DMCreateLocalVector(dm, &X_loc)); 131*2b730f8bSJeremy L Thompson PetscCall(VecGetSize(X_loc, &xl_size)); 132*2b730f8bSJeremy L Thompson PetscCall(VecDuplicate(X, &rhs)); 133ed264d09SValeria Barra 134ed264d09SValeria Barra // Operator 135*2b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &op_apply_ctx)); 136*2b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &op_error_ctx)); 137*2b730f8bSJeremy L Thompson PetscCall(MatCreateShell(comm, l_size, l_size, g_size, g_size, op_apply_ctx, &mat_O)); 138*2b730f8bSJeremy L Thompson PetscCall(MatShellSetOperation(mat_O, MATOP_MULT, (void (*)(void))MatMult_Ceed)); 139ed264d09SValeria Barra 140ed264d09SValeria Barra // Set up libCEED 1419b072555Sjeremylt CeedInit(ceed_resource, &ceed); 1429b072555Sjeremylt CeedMemType mem_type_backend; 1439b072555Sjeremylt CeedGetPreferredMemType(ceed, &mem_type_backend); 144e83e87a5Sjeremylt 145*2b730f8bSJeremy L Thompson PetscCall(DMGetVecType(dm, &vec_type)); 1469b072555Sjeremylt if (!vec_type) { // Not yet set by user -dm_vec_type 1479b072555Sjeremylt switch (mem_type_backend) { 148*2b730f8bSJeremy L Thompson case CEED_MEM_HOST: 149*2b730f8bSJeremy L Thompson vec_type = VECSTANDARD; 150*2b730f8bSJeremy L Thompson break; 151e83e87a5Sjeremylt case CEED_MEM_DEVICE: { 152e83e87a5Sjeremylt const char *resolved; 153e83e87a5Sjeremylt CeedGetResource(ceed, &resolved); 1549b072555Sjeremylt if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA; 155*2b730f8bSJeremy L Thompson else if (strstr(resolved, "/gpu/hip/occa")) vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678 1569b072555Sjeremylt else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP; 1579b072555Sjeremylt else vec_type = VECSTANDARD; 158e83e87a5Sjeremylt } 159e83e87a5Sjeremylt } 160*2b730f8bSJeremy L Thompson PetscCall(DMSetVecType(dm, vec_type)); 161e83e87a5Sjeremylt } 162ed264d09SValeria Barra 163ed264d09SValeria Barra // Print summary 164ed264d09SValeria Barra if (!test_mode) { 1659b072555Sjeremylt PetscInt P = degree + 1, Q = P + q_extra; 1669b072555Sjeremylt const char *used_resource; 1679b072555Sjeremylt CeedGetResource(ceed, &used_resource); 168*2b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 169*2b730f8bSJeremy L Thompson "\n-- CEED Benchmark Problem %" CeedInt_FMT " on the Sphere -- libCEED + PETSc --\n" 170ed264d09SValeria Barra " libCEED:\n" 171ed264d09SValeria Barra " libCEED Backend : %s\n" 172e83e87a5Sjeremylt " libCEED Backend MemType : %s\n" 173ed264d09SValeria Barra " Mesh:\n" 174751eb813Srezgarshakeri " Solution Order (P) : %" CeedInt_FMT "\n" 175751eb813Srezgarshakeri " Quadrature Order (Q) : %" CeedInt_FMT "\n" 17651ad7d5bSrezgarshakeri " Additional quadrature points (q_extra) : %" CeedInt_FMT "\n" 17708140895SJed Brown " Global nodes : %" PetscInt_FMT "\n", 178*2b730f8bSJeremy L Thompson bp_choice + 1, ceed_resource, CeedMemTypes[mem_type_backend], P, Q, q_extra, g_size / num_comp_u)); 179ed264d09SValeria Barra } 180ed264d09SValeria Barra 181ed264d09SValeria Barra // Create RHS vector 182*2b730f8bSJeremy L Thompson PetscCall(VecDuplicate(X_loc, &rhs_loc)); 183*2b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(rhs_loc)); 184*2b730f8bSJeremy L Thompson PetscCall(VecGetArrayAndMemType(rhs_loc, &r, &mem_type)); 1859b072555Sjeremylt CeedVectorCreate(ceed, xl_size, &rhs_ceed); 1869b072555Sjeremylt CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r); 187ed264d09SValeria Barra 188ed264d09SValeria Barra // Setup libCEED's objects 189*2b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &ceed_data)); 190*2b730f8bSJeremy 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, 191*2b730f8bSJeremy L Thompson rhs_ceed, &target)); 192ed264d09SValeria Barra 193ed264d09SValeria Barra // Gather RHS 1949b072555Sjeremylt CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL); 195*2b730f8bSJeremy L Thompson PetscCall(VecRestoreArrayAndMemType(rhs_loc, &r)); 196*2b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(rhs)); 197*2b730f8bSJeremy L Thompson PetscCall(DMLocalToGlobal(dm, rhs_loc, ADD_VALUES, rhs)); 1989b072555Sjeremylt CeedVectorDestroy(&rhs_ceed); 199ed264d09SValeria Barra 200ed264d09SValeria Barra // Create the error Q-function 201*2b730f8bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error, bp_options[bp_choice].error_loc, &qf_error); 2029b072555Sjeremylt CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP); 2039b072555Sjeremylt CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE); 204*2b730f8bSJeremy L Thompson CeedQFunctionAddInput(qf_error, "qdata", ceed_data->q_data_size, CEED_EVAL_NONE); 20538f32c05Srezgarshakeri CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_INTERP); 206ed264d09SValeria Barra 207ed264d09SValeria Barra // Create the error operator 2089b072555Sjeremylt CeedOperatorCreate(ceed, qf_error, NULL, NULL, &op_error); 209*2b730f8bSJeremy L Thompson CeedOperatorSetField(op_error, "u", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE); 210*2b730f8bSJeremy L Thompson CeedOperatorSetField(op_error, "true_soln", ceed_data->elem_restr_u_i, CEED_BASIS_COLLOCATED, target); 211*2b730f8bSJeremy L Thompson CeedOperatorSetField(op_error, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_COLLOCATED, ceed_data->q_data); 212*2b730f8bSJeremy L Thompson CeedOperatorSetField(op_error, "error", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE); 213ed264d09SValeria Barra 2146c88e6a2Srezgarshakeri // Set up apply operator context 215*2b730f8bSJeremy L Thompson PetscCall(SetupApplyOperatorCtx(comm, dm, ceed, ceed_data, X_loc, op_apply_ctx)); 216ed264d09SValeria Barra 217ed264d09SValeria Barra // Setup solver 218*2b730f8bSJeremy L Thompson PetscCall(KSPCreate(comm, &ksp)); 219ed264d09SValeria Barra { 220ed264d09SValeria Barra PC pc; 221*2b730f8bSJeremy L Thompson PetscCall(KSPGetPC(ksp, &pc)); 2229b072555Sjeremylt if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) { 223*2b730f8bSJeremy L Thompson PetscCall(PCSetType(pc, PCJACOBI)); 224*2b730f8bSJeremy L Thompson PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM)); 225ed264d09SValeria Barra } else { 226*2b730f8bSJeremy L Thompson PetscCall(PCSetType(pc, PCNONE)); 227ed264d09SValeria Barra MatNullSpace nullspace; 228ed264d09SValeria Barra 229*2b730f8bSJeremy L Thompson PetscCall(MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, 0, &nullspace)); 230*2b730f8bSJeremy L Thompson PetscCall(MatSetNullSpace(mat_O, nullspace)); 231*2b730f8bSJeremy L Thompson PetscCall(MatNullSpaceDestroy(&nullspace)); 232ed264d09SValeria Barra } 233*2b730f8bSJeremy L Thompson PetscCall(KSPSetType(ksp, KSPCG)); 234*2b730f8bSJeremy L Thompson PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL)); 235*2b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 236ed264d09SValeria Barra } 237*2b730f8bSJeremy L Thompson PetscCall(KSPSetFromOptions(ksp)); 238*2b730f8bSJeremy L Thompson PetscCall(KSPSetOperators(ksp, mat_O, mat_O)); 239ed264d09SValeria Barra 240ed264d09SValeria Barra // First run, if benchmarking 241ed264d09SValeria Barra if (benchmark_mode) { 242*2b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1)); 243ed264d09SValeria Barra my_rt_start = MPI_Wtime(); 244*2b730f8bSJeremy L Thompson PetscCall(KSPSolve(ksp, rhs, X)); 245ed264d09SValeria Barra my_rt = MPI_Wtime() - my_rt_start; 246*2b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm)); 247ed264d09SValeria Barra // Set maxits based on first iteration timing 248ed264d09SValeria Barra if (my_rt > 0.02) { 249*2b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5)); 250ed264d09SValeria Barra } else { 251*2b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20)); 252ed264d09SValeria Barra } 253ed264d09SValeria Barra } 254ed264d09SValeria Barra 255ed264d09SValeria Barra // Timed solve 256*2b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X)); 257*2b730f8bSJeremy L Thompson PetscCall(PetscBarrier((PetscObject)ksp)); 25809a940d7Sjeremylt 25909a940d7Sjeremylt // -- Performance logging 260*2b730f8bSJeremy L Thompson PetscCall(PetscLogStageRegister("Solve Stage", &solve_stage)); 261*2b730f8bSJeremy L Thompson PetscCall(PetscLogStagePush(solve_stage)); 26209a940d7Sjeremylt 26309a940d7Sjeremylt // -- Solve 264ed264d09SValeria Barra my_rt_start = MPI_Wtime(); 265*2b730f8bSJeremy L Thompson PetscCall(KSPSolve(ksp, rhs, X)); 266ed264d09SValeria Barra my_rt = MPI_Wtime() - my_rt_start; 267ed264d09SValeria Barra 26809a940d7Sjeremylt // -- Performance logging 269*2b730f8bSJeremy L Thompson PetscCall(PetscLogStagePop()); 27009a940d7Sjeremylt 271ed264d09SValeria Barra // Output results 272ed264d09SValeria Barra { 2739b072555Sjeremylt KSPType ksp_type; 274ed264d09SValeria Barra KSPConvergedReason reason; 275ed264d09SValeria Barra PetscReal rnorm; 276ed264d09SValeria Barra PetscInt its; 277*2b730f8bSJeremy L Thompson PetscCall(KSPGetType(ksp, &ksp_type)); 278*2b730f8bSJeremy L Thompson PetscCall(KSPGetConvergedReason(ksp, &reason)); 279*2b730f8bSJeremy L Thompson PetscCall(KSPGetIterationNumber(ksp, &its)); 280*2b730f8bSJeremy L Thompson PetscCall(KSPGetResidualNorm(ksp, &rnorm)); 281ed264d09SValeria Barra if (!test_mode || reason < 0 || rnorm > 1e-8) { 282*2b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 283ed264d09SValeria Barra " KSP:\n" 284ed264d09SValeria Barra " KSP Type : %s\n" 285ed264d09SValeria Barra " KSP Convergence : %s\n" 286a9b2c5ddSrezgarshakeri " Total KSP Iterations : %" PetscInt_FMT "\n" 287ed264d09SValeria Barra " Final rnorm : %e\n", 288*2b730f8bSJeremy L Thompson ksp_type, KSPConvergedReasons[reason], its, (double)rnorm)); 289ed264d09SValeria Barra } 290ed264d09SValeria Barra if (!test_mode) { 291*2b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, " Performance:\n")); 292ed264d09SValeria Barra } 293ed264d09SValeria Barra { 2946c88e6a2Srezgarshakeri // Set up error operator context 295*2b730f8bSJeremy L Thompson PetscCall(SetupErrorOperatorCtx(comm, dm, ceed, ceed_data, X_loc, op_error, op_error_ctx)); 29638f32c05Srezgarshakeri PetscScalar l2_error; 297*2b730f8bSJeremy L Thompson PetscCall(ComputeL2Error(X, &l2_error, op_error_ctx)); 298587be3cdSvaleriabarra PetscReal tol = 5e-4; 29938f32c05Srezgarshakeri if (!test_mode || l2_error > tol) { 300*2b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm)); 301*2b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm)); 302*2b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 30338f32c05Srezgarshakeri " L2 Error : %e\n" 304ed264d09SValeria Barra " CG Solve Time : %g (%g) sec\n", 305*2b730f8bSJeremy L Thompson (double)l2_error, rt_max, rt_min)); 306ed264d09SValeria Barra } 307ed264d09SValeria Barra } 308ed264d09SValeria Barra if (benchmark_mode && (!test_mode)) { 309*2b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, " DoFs/Sec in CG : %g (%g) million\n", 1e-6 * g_size * its / rt_max, 310*2b730f8bSJeremy L Thompson 1e-6 * g_size * its / rt_min)); 311ed264d09SValeria Barra } 312ed264d09SValeria Barra } 313ed264d09SValeria Barra 314ed264d09SValeria Barra // Output solution 315ed264d09SValeria Barra if (write_solution) { 3169b072555Sjeremylt PetscViewer vtk_viewer_soln; 317ed264d09SValeria Barra 318*2b730f8bSJeremy L Thompson PetscCall(PetscViewerCreate(comm, &vtk_viewer_soln)); 319*2b730f8bSJeremy L Thompson PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK)); 320*2b730f8bSJeremy L Thompson PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu")); 321*2b730f8bSJeremy L Thompson PetscCall(VecView(X, vtk_viewer_soln)); 322*2b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&vtk_viewer_soln)); 323ed264d09SValeria Barra } 324ed264d09SValeria Barra 325ed264d09SValeria Barra // Cleanup 326*2b730f8bSJeremy L Thompson PetscCall(VecDestroy(&X)); 327*2b730f8bSJeremy L Thompson PetscCall(VecDestroy(&X_loc)); 328*2b730f8bSJeremy L Thompson PetscCall(VecDestroy(&op_apply_ctx->Y_loc)); 329*2b730f8bSJeremy L Thompson PetscCall(VecDestroy(&op_error_ctx->Y_loc)); 330*2b730f8bSJeremy L Thompson PetscCall(MatDestroy(&mat_O)); 331*2b730f8bSJeremy L Thompson PetscCall(PetscFree(op_apply_ctx)); 332*2b730f8bSJeremy L Thompson PetscCall(PetscFree(op_error_ctx)); 333*2b730f8bSJeremy L Thompson PetscCall(CeedDataDestroy(0, ceed_data)); 334*2b730f8bSJeremy L Thompson PetscCall(DMDestroy(&dm)); 335ed264d09SValeria Barra 336*2b730f8bSJeremy L Thompson PetscCall(VecDestroy(&rhs)); 337*2b730f8bSJeremy L Thompson PetscCall(VecDestroy(&rhs_loc)); 338*2b730f8bSJeremy L Thompson PetscCall(KSPDestroy(&ksp)); 339ed264d09SValeria Barra CeedVectorDestroy(&target); 3409b072555Sjeremylt CeedQFunctionDestroy(&qf_error); 3419b072555Sjeremylt CeedOperatorDestroy(&op_error); 342ed264d09SValeria Barra CeedDestroy(&ceed); 343ed264d09SValeria Barra return PetscFinalize(); 344ed264d09SValeria Barra } 345