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. 36c5df90dSjeremylt // 43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause 56c5df90dSjeremylt // 63d8e8822SJeremy L Thompson // This file is part of CEED: http://github.com/ceed 76c5df90dSjeremylt 86c5df90dSjeremylt // libCEED + PETSc Example: CEED BPs 3-6 with Multigrid 96c5df90dSjeremylt // 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. 116c5df90dSjeremylt // 126c5df90dSjeremylt // The code uses higher level communication protocols in DMPlex. 136c5df90dSjeremylt // 146c5df90dSjeremylt // Build with: 156c5df90dSjeremylt // 166c5df90dSjeremylt // make multigrid [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>] 176c5df90dSjeremylt // 186c5df90dSjeremylt // Sample runs: 196c5df90dSjeremylt // 206c5df90dSjeremylt // multigrid -problem bp3 2128688798Sjeremylt // multigrid -problem bp4 2228688798Sjeremylt // multigrid -problem bp5 -ceed /cpu/self 236c5df90dSjeremylt // multigrid -problem bp6 -ceed /gpu/cuda 246c5df90dSjeremylt // 25fa5adaf5SJeremy L Thompson //TESTARGS(name="BP3, hex elements") -ceed {ceed_resource} -test -problem bp3 -degree 3 26fa5adaf5SJeremy L Thompson //TESTARGS(name="BP3, tet elements") -ceed {ceed_resource} -test -problem bp3 -degree 3 -simplex 276c5df90dSjeremylt 286c5df90dSjeremylt /// @file 296c5df90dSjeremylt /// CEED BPs 1-6 multigrid example using PETSc 306c5df90dSjeremylt const char help[] = "Solve CEED BPs using p-multigrid with PETSc and DMPlex\n"; 316c5df90dSjeremylt 32636cccdbSjeremylt #include <ceed.h> 33636cccdbSjeremylt #include <petsc.h> 34636cccdbSjeremylt #include <petscdmplex.h> 35636cccdbSjeremylt #include <petscksp.h> 36636cccdbSjeremylt #include <petscsys.h> 372b730f8bSJeremy L Thompson #include <stdbool.h> 382b730f8bSJeremy L Thompson #include <string.h> 39636cccdbSjeremylt 40e83e87a5Sjeremylt #include "bps.h" 41636cccdbSjeremylt #include "include/bpsproblemdata.h" 422b730f8bSJeremy L Thompson #include "include/libceedsetup.h" 432b730f8bSJeremy L Thompson #include "include/matops.h" 44636cccdbSjeremylt #include "include/petscutils.h" 45b8962995SJeremy L Thompson #include "include/petscversion.h" 46636cccdbSjeremylt #include "include/structs.h" 47636cccdbSjeremylt 486c5df90dSjeremylt int main(int argc, char **argv) { 496c5df90dSjeremylt MPI_Comm comm; 502b730f8bSJeremy L Thompson char filename[PETSC_MAX_PATH_LEN], ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self"; 516c5df90dSjeremylt double my_rt_start, my_rt, rt_min, rt_max; 522b730f8bSJeremy L Thompson PetscInt degree = 3, q_extra, *l_size, *xl_size, *g_size, dim = 3, fine_level, mesh_elem[3] = {3, 3, 3}, num_comp_u = 1, num_levels = degree, 532b730f8bSJeremy L Thompson *level_degrees; 546c5df90dSjeremylt PetscScalar *r; 55cfa59c5bSRey PetscScalar eps = 1.0; 56de1229c5Srezgarshakeri PetscBool test_mode, benchmark_mode, read_mesh, write_solution, simplex; 579b072555Sjeremylt PetscLogStage solve_stage; 5805b9c820SJed Brown PetscLogEvent assemble_event; 599b072555Sjeremylt DM *dm, dm_orig; 606c5df90dSjeremylt KSP ksp; 616c5df90dSjeremylt PC pc; 629b072555Sjeremylt Mat *mat_O, *mat_pr, mat_coarse; 639b072555Sjeremylt Vec *X, *X_loc, *mult, rhs, rhs_loc; 649b072555Sjeremylt PetscMemType mem_type; 656c88e6a2Srezgarshakeri OperatorApplyContext *op_apply_ctx, op_error_ctx; 66d4d45553Srezgarshakeri ProlongRestrContext *pr_restr_ctx; 676c5df90dSjeremylt Ceed ceed; 689b072555Sjeremylt CeedData *ceed_data; 699b072555Sjeremylt CeedVector rhs_ceed, target; 701c9a79dbSrezgarshakeri CeedQFunction qf_error; 719b072555Sjeremylt CeedOperator op_error; 729b072555Sjeremylt BPType bp_choice; 739b072555Sjeremylt CoarsenType coarsen; 746c5df90dSjeremylt 752b730f8bSJeremy L Thompson PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 766c5df90dSjeremylt comm = PETSC_COMM_WORLD; 776c5df90dSjeremylt 786c5df90dSjeremylt // Parse command line options 7967490bc6SJeremy L Thompson PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); 809b072555Sjeremylt bp_choice = CEED_BP3; 812b730f8bSJeremy L Thompson PetscCall(PetscOptionsEnum("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice, NULL)); 829b072555Sjeremylt num_comp_u = bp_options[bp_choice].num_comp_u; 836c5df90dSjeremylt test_mode = PETSC_FALSE; 842b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL)); 856c5df90dSjeremylt benchmark_mode = PETSC_FALSE; 862b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-benchmark", "Benchmarking mode (prints benchmark statistics)", NULL, benchmark_mode, &benchmark_mode, NULL)); 876c5df90dSjeremylt write_solution = PETSC_FALSE; 882b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, write_solution, &write_solution, NULL)); 89de1229c5Srezgarshakeri simplex = PETSC_FALSE; 902b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-simplex", "Element topology (default:hex)", NULL, simplex, &simplex, NULL)); 91eb6a3b92Srezgarshakeri if ((bp_choice == CEED_BP5 || bp_choice == CEED_BP6) && (simplex)) { 922b730f8bSJeremy L Thompson SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "BP5/6 is not supported with simplex"); 93eb6a3b92Srezgarshakeri } 942b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-eps", "Epsilon parameter for Kershaw mesh transformation", NULL, eps, &eps, NULL)); 952b730f8bSJeremy L Thompson if (eps > 1 || eps <= 0) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "-eps %g must be (0,1]", (double)PetscRealPart(eps)); 966c5df90dSjeremylt degree = test_mode ? 3 : 2; 972b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, degree, °ree, NULL)); 982b730f8bSJeremy L Thompson if (degree < 1) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "-degree %" PetscInt_FMT " must be at least 1", degree); 999b072555Sjeremylt q_extra = bp_options[bp_choice].q_extra; 1002b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra, NULL)); 1012b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resource, sizeof(ceed_resource), NULL)); 1026c5df90dSjeremylt coarsen = COARSEN_UNIFORM; 1032b730f8bSJeremy L Thompson PetscCall(PetscOptionsEnum("-coarsen", "Coarsening strategy to use", NULL, coarsen_types, (PetscEnum)coarsen, (PetscEnum *)&coarsen, NULL)); 104cb32e2e7SValeria Barra read_mesh = PETSC_FALSE; 1052b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &read_mesh)); 1066c5df90dSjeremylt if (!read_mesh) { 1076c5df90dSjeremylt PetscInt tmp = dim; 1082b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL, mesh_elem, &tmp, NULL)); 1096c5df90dSjeremylt } 11067490bc6SJeremy L Thompson PetscOptionsEnd(); 1116c5df90dSjeremylt 1129396343dSjeremylt // Set up libCEED 1139b072555Sjeremylt CeedInit(ceed_resource, &ceed); 1149b072555Sjeremylt CeedMemType mem_type_backend; 1159b072555Sjeremylt CeedGetPreferredMemType(ceed, &mem_type_backend); 1169396343dSjeremylt 1176c5df90dSjeremylt // Setup DM 1186c5df90dSjeremylt if (read_mesh) { 1192b730f8bSJeremy L Thompson PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, NULL, PETSC_TRUE, &dm_orig)); 1206c5df90dSjeremylt } else { 1212b730f8bSJeremy L Thompson PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, simplex, mesh_elem, NULL, NULL, NULL, PETSC_TRUE, &dm_orig)); 1226c5df90dSjeremylt } 1236c5df90dSjeremylt 1249b072555Sjeremylt VecType vec_type; 1259b072555Sjeremylt switch (mem_type_backend) { 1262b730f8bSJeremy L Thompson case CEED_MEM_HOST: 1272b730f8bSJeremy L Thompson vec_type = VECSTANDARD; 1282b730f8bSJeremy L Thompson break; 129b68a8d79SJed Brown case CEED_MEM_DEVICE: { 130b68a8d79SJed Brown const char *resolved; 131b68a8d79SJed Brown CeedGetResource(ceed, &resolved); 1329b072555Sjeremylt if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA; 1332b730f8bSJeremy L Thompson else if (strstr(resolved, "/gpu/hip/occa")) vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678 1349b072555Sjeremylt else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP; 1359b072555Sjeremylt else vec_type = VECSTANDARD; 136b68a8d79SJed Brown } 137b68a8d79SJed Brown } 1382b730f8bSJeremy L Thompson PetscCall(DMSetVecType(dm_orig, vec_type)); 1392b730f8bSJeremy L Thompson PetscCall(DMSetFromOptions(dm_orig)); 1402b730f8bSJeremy L Thompson PetscCall(DMViewFromOptions(dm_orig, NULL, "-dm_view")); 1413fc8a154SJed Brown 1423fc8a154SJed Brown // Apply Kershaw mesh transformation 1432b730f8bSJeremy L Thompson PetscCall(Kershaw(dm_orig, eps)); 144b68a8d79SJed Brown 1456c5df90dSjeremylt // Allocate arrays for PETSc objects for each level 1466c5df90dSjeremylt switch (coarsen) { 1476c5df90dSjeremylt case COARSEN_UNIFORM: 1489b072555Sjeremylt num_levels = degree; 1496c5df90dSjeremylt break; 150dc7d240cSValeria Barra case COARSEN_LOGARITHMIC: 1519b072555Sjeremylt num_levels = ceil(log(degree) / log(2)) + 1; 1526c5df90dSjeremylt break; 1536c5df90dSjeremylt } 1542b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &level_degrees)); 1559b072555Sjeremylt fine_level = num_levels - 1; 15661608365Sjeremylt 1576c5df90dSjeremylt switch (coarsen) { 1586c5df90dSjeremylt case COARSEN_UNIFORM: 1599b072555Sjeremylt for (int i = 0; i < num_levels; i++) level_degrees[i] = i + 1; 1606c5df90dSjeremylt break; 161dc7d240cSValeria Barra case COARSEN_LOGARITHMIC: 1629b072555Sjeremylt for (int i = 0; i < num_levels - 1; i++) level_degrees[i] = pow(2, i); 1639b072555Sjeremylt level_degrees[fine_level] = degree; 1646c5df90dSjeremylt break; 1656c5df90dSjeremylt } 1662b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &dm)); 1672b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &X)); 1682b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &X_loc)); 1692b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &mult)); 1702b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &op_apply_ctx)); 1712b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &pr_restr_ctx)); 1722b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &mat_O)); 1732b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &mat_pr)); 1742b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &l_size)); 1752b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &xl_size)); 1762b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &g_size)); 1776c5df90dSjeremylt 1781c9a79dbSrezgarshakeri PetscInt c_start, c_end; 1792b730f8bSJeremy L Thompson PetscCall(DMPlexGetHeightStratum(dm_orig, 0, &c_start, &c_end)); 1801c9a79dbSrezgarshakeri DMPolytopeType cell_type; 1812b730f8bSJeremy L Thompson PetscCall(DMPlexGetCellType(dm_orig, c_start, &cell_type)); 1821c9a79dbSrezgarshakeri CeedElemTopology elem_topo = ElemTopologyP2C(cell_type); 1831c9a79dbSrezgarshakeri 1846c5df90dSjeremylt // Setup DM and Operator Mat Shells for each level 1859b072555Sjeremylt for (CeedInt i = 0; i < num_levels; i++) { 1866c5df90dSjeremylt // Create DM 1872b730f8bSJeremy L Thompson PetscCall(DMClone(dm_orig, &dm[i])); 1882b730f8bSJeremy L Thompson PetscCall(DMGetVecType(dm_orig, &vec_type)); 1892b730f8bSJeremy L Thompson PetscCall(DMSetVecType(dm[i], vec_type)); 190e83e87a5Sjeremylt PetscInt dim; 1912b730f8bSJeremy L Thompson PetscCall(DMGetDimension(dm[i], &dim)); 1922b730f8bSJeremy L Thompson PetscCall(SetupDMByDegree(dm[i], level_degrees[fine_level], q_extra, num_comp_u, dim, bp_options[bp_choice].enforce_bc)); 1936c5df90dSjeremylt 1946c5df90dSjeremylt // Create vectors 1952b730f8bSJeremy L Thompson PetscCall(DMCreateGlobalVector(dm[i], &X[i])); 1962b730f8bSJeremy L Thompson PetscCall(VecGetLocalSize(X[i], &l_size[i])); 1972b730f8bSJeremy L Thompson PetscCall(VecGetSize(X[i], &g_size[i])); 1982b730f8bSJeremy L Thompson PetscCall(DMCreateLocalVector(dm[i], &X_loc[i])); 1992b730f8bSJeremy L Thompson PetscCall(VecGetSize(X_loc[i], &xl_size[i])); 2006c5df90dSjeremylt 2016c5df90dSjeremylt // Operator 2022b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &op_apply_ctx[i])); 2032b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &op_error_ctx)); 2042b730f8bSJeremy L Thompson PetscCall(MatCreateShell(comm, l_size[i], l_size[i], g_size[i], g_size[i], op_apply_ctx[i], &mat_O[i])); 2052b730f8bSJeremy L Thompson PetscCall(MatShellSetOperation(mat_O[i], MATOP_MULT, (void (*)(void))MatMult_Ceed)); 2062b730f8bSJeremy L Thompson PetscCall(MatShellSetOperation(mat_O[i], MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag)); 2072b730f8bSJeremy L Thompson PetscCall(MatShellSetVecType(mat_O[i], vec_type)); 2086c5df90dSjeremylt 2096c5df90dSjeremylt // Level transfers 2106c5df90dSjeremylt if (i > 0) { 2116c5df90dSjeremylt // Interp 2122b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &pr_restr_ctx[i])); 2132b730f8bSJeremy L Thompson PetscCall(MatCreateShell(comm, l_size[i], l_size[i - 1], g_size[i], g_size[i - 1], pr_restr_ctx[i], &mat_pr[i])); 2142b730f8bSJeremy L Thompson PetscCall(MatShellSetOperation(mat_pr[i], MATOP_MULT, (void (*)(void))MatMult_Prolong)); 2152b730f8bSJeremy L Thompson PetscCall(MatShellSetOperation(mat_pr[i], MATOP_MULT_TRANSPOSE, (void (*)(void))MatMult_Restrict)); 2162b730f8bSJeremy L Thompson PetscCall(MatShellSetVecType(mat_pr[i], vec_type)); 2176c5df90dSjeremylt } 2186c5df90dSjeremylt } 2192b730f8bSJeremy L Thompson PetscCall(VecDuplicate(X[fine_level], &rhs)); 2206c5df90dSjeremylt 2216c5df90dSjeremylt // Print global grid information 2226c5df90dSjeremylt if (!test_mode) { 2239b072555Sjeremylt PetscInt P = degree + 1, Q = P + q_extra; 2249396343dSjeremylt 2259b072555Sjeremylt const char *used_resource; 2269b072555Sjeremylt CeedGetResource(ceed, &used_resource); 2279396343dSjeremylt 2282b730f8bSJeremy L Thompson PetscCall(VecGetType(X[0], &vec_type)); 2299396343dSjeremylt 2302b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 231990fdeb6SJeremy L Thompson "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc + PCMG --\n" 2329396343dSjeremylt " PETSc:\n" 2339396343dSjeremylt " PETSc Vec Type : %s\n" 2346c5df90dSjeremylt " libCEED:\n" 2356c5df90dSjeremylt " libCEED Backend : %s\n" 2369396343dSjeremylt " libCEED Backend MemType : %s\n" 2376c5df90dSjeremylt " Mesh:\n" 238751eb813Srezgarshakeri " Solution Order (P) : %" CeedInt_FMT "\n" 239751eb813Srezgarshakeri " Quadrature Order (Q) : %" CeedInt_FMT "\n" 240de1229c5Srezgarshakeri " Additional quadrature points (q_extra) : %" CeedInt_FMT "\n" 24108140895SJed Brown " Global Nodes : %" PetscInt_FMT "\n" 24208140895SJed Brown " Owned Nodes : %" PetscInt_FMT "\n" 24308140895SJed Brown " DoF per node : %" PetscInt_FMT "\n" 24451ad7d5bSrezgarshakeri " Element topology : %s\n" 2456c5df90dSjeremylt " Multigrid:\n" 246990fdeb6SJeremy L Thompson " Number of Levels : %" CeedInt_FMT "\n", 2472b730f8bSJeremy L Thompson bp_choice + 1, vec_type, used_resource, CeedMemTypes[mem_type_backend], P, Q, q_extra, g_size[fine_level] / num_comp_u, 2482b730f8bSJeremy L Thompson l_size[fine_level] / num_comp_u, num_comp_u, CeedElemTopologies[elem_topo], num_levels)); 2496c5df90dSjeremylt } 2506c5df90dSjeremylt 2516c5df90dSjeremylt // Create RHS vector 2522b730f8bSJeremy L Thompson PetscCall(VecDuplicate(X_loc[fine_level], &rhs_loc)); 2532b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(rhs_loc)); 2542b730f8bSJeremy L Thompson PetscCall(VecGetArrayAndMemType(rhs_loc, &r, &mem_type)); 2559b072555Sjeremylt CeedVectorCreate(ceed, xl_size[fine_level], &rhs_ceed); 2569b072555Sjeremylt CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r); 2576c5df90dSjeremylt 2586c5df90dSjeremylt // Set up libCEED operators on each level 2592b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &ceed_data)); 260990fdeb6SJeremy L Thompson for (PetscInt i = 0; i < num_levels; i++) { 2616c5df90dSjeremylt // Print level information 2629b072555Sjeremylt if (!test_mode && (i == 0 || i == fine_level)) { 2632b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 2642b730f8bSJeremy L Thompson " Level %" PetscInt_FMT " (%s):\n" 265751eb813Srezgarshakeri " Solution Order (P) : %" CeedInt_FMT "\n" 26608140895SJed Brown " Global Nodes : %" PetscInt_FMT "\n" 26708140895SJed Brown " Owned Nodes : %" PetscInt_FMT "\n", 2682b730f8bSJeremy L Thompson i, (i ? "fine" : "coarse"), level_degrees[i] + 1, g_size[i] / num_comp_u, l_size[i] / num_comp_u)); 2696c5df90dSjeremylt } 2702b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &ceed_data[i])); 2712b730f8bSJeremy L Thompson PetscCall(SetupLibceedByDegree(dm[i], ceed, level_degrees[i], dim, q_extra, dim, num_comp_u, g_size[i], xl_size[i], bp_options[bp_choice], 2722b730f8bSJeremy L Thompson ceed_data[i], i == (fine_level), rhs_ceed, &target)); 2736c5df90dSjeremylt } 2746c5df90dSjeremylt 2756c5df90dSjeremylt // Gather RHS 2769b072555Sjeremylt CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL); 2772b730f8bSJeremy L Thompson PetscCall(VecRestoreArrayAndMemType(rhs_loc, &r)); 2782b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(rhs)); 2792b730f8bSJeremy L Thompson PetscCall(DMLocalToGlobal(dm[fine_level], rhs_loc, ADD_VALUES, rhs)); 2809b072555Sjeremylt CeedVectorDestroy(&rhs_ceed); 2816c5df90dSjeremylt 28243eb8658SJeremy L Thompson // Create the error QFunction 2832b730f8bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error, bp_options[bp_choice].error_loc, &qf_error); 2849b072555Sjeremylt CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP); 2859b072555Sjeremylt CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE); 2862b730f8bSJeremy L Thompson CeedQFunctionAddInput(qf_error, "qdata", ceed_data[fine_level]->q_data_size, CEED_EVAL_NONE); 28738f32c05Srezgarshakeri CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_INTERP); 2886c5df90dSjeremylt 2896c5df90dSjeremylt // Create the error operator 2902b730f8bSJeremy L Thompson CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_error); 2912b730f8bSJeremy L Thompson CeedOperatorSetField(op_error, "u", ceed_data[fine_level]->elem_restr_u, ceed_data[fine_level]->basis_u, CEED_VECTOR_ACTIVE); 292*356036faSJeremy L Thompson CeedOperatorSetField(op_error, "true_soln", ceed_data[fine_level]->elem_restr_u_i, CEED_BASIS_NONE, target); 293*356036faSJeremy L Thompson CeedOperatorSetField(op_error, "qdata", ceed_data[fine_level]->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data[fine_level]->q_data); 2942b730f8bSJeremy L Thompson CeedOperatorSetField(op_error, "error", ceed_data[fine_level]->elem_restr_u, ceed_data[fine_level]->basis_u, CEED_VECTOR_ACTIVE); 2956c5df90dSjeremylt 2966c5df90dSjeremylt // Calculate multiplicity 2979b072555Sjeremylt for (int i = 0; i < num_levels; i++) { 2986c5df90dSjeremylt PetscScalar *x; 2996c5df90dSjeremylt 3006c5df90dSjeremylt // CEED vector 3012b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X_loc[i])); 3022b730f8bSJeremy L Thompson PetscCall(VecGetArray(X_loc[i], &x)); 3039b072555Sjeremylt CeedVectorSetArray(ceed_data[i]->x_ceed, CEED_MEM_HOST, CEED_USE_POINTER, x); 3046c5df90dSjeremylt 3056c5df90dSjeremylt // Multiplicity 3062b730f8bSJeremy L Thompson CeedElemRestrictionGetMultiplicity(ceed_data[i]->elem_restr_u, ceed_data[i]->x_ceed); 3079b072555Sjeremylt CeedVectorSyncArray(ceed_data[i]->x_ceed, CEED_MEM_HOST); 3086c5df90dSjeremylt 3096c5df90dSjeremylt // Restore vector 3102b730f8bSJeremy L Thompson PetscCall(VecRestoreArray(X_loc[i], &x)); 3116c5df90dSjeremylt 3126c5df90dSjeremylt // Creat mult vector 3132b730f8bSJeremy L Thompson PetscCall(VecDuplicate(X_loc[i], &mult[i])); 3146c5df90dSjeremylt 3156c5df90dSjeremylt // Local-to-global 3162b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X[i])); 3172b730f8bSJeremy L Thompson PetscCall(DMLocalToGlobal(dm[i], X_loc[i], ADD_VALUES, X[i])); 3182b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X_loc[i])); 3196c5df90dSjeremylt 3206c5df90dSjeremylt // Global-to-local 3212b730f8bSJeremy L Thompson PetscCall(DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i])); 3222b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X[i])); 3236c5df90dSjeremylt 3246c5df90dSjeremylt // Multiplicity scaling 3252b730f8bSJeremy L Thompson PetscCall(VecReciprocal(mult[i])); 3266c5df90dSjeremylt } 3276c5df90dSjeremylt 3286c5df90dSjeremylt // Set up Mat 3299b072555Sjeremylt for (int i = 0; i < num_levels; i++) { 3306c88e6a2Srezgarshakeri // Set up apply operator context 3312b730f8bSJeremy L Thompson PetscCall(SetupApplyOperatorCtx(comm, dm[i], ceed, ceed_data[i], X_loc[i], op_apply_ctx[i])); 3326c5df90dSjeremylt 3336c5df90dSjeremylt if (i > 0) { 334a97643b0Sjeremylt // Prolongation/Restriction Operator 3352b730f8bSJeremy L Thompson PetscCall(CeedLevelTransferSetup(dm[i - 1], ceed, i, num_comp_u, ceed_data, bp_options[bp_choice], mult[i])); 336d4d45553Srezgarshakeri pr_restr_ctx[i]->comm = comm; 337d4d45553Srezgarshakeri pr_restr_ctx[i]->dmf = dm[i]; 338d4d45553Srezgarshakeri pr_restr_ctx[i]->dmc = dm[i - 1]; 339d4d45553Srezgarshakeri pr_restr_ctx[i]->loc_vec_c = X_loc[i - 1]; 340d4d45553Srezgarshakeri pr_restr_ctx[i]->loc_vec_f = op_apply_ctx[i]->Y_loc; 341d4d45553Srezgarshakeri pr_restr_ctx[i]->mult_vec = mult[i]; 342d4d45553Srezgarshakeri pr_restr_ctx[i]->ceed_vec_c = op_apply_ctx[i - 1]->x_ceed; 343d4d45553Srezgarshakeri pr_restr_ctx[i]->ceed_vec_f = op_apply_ctx[i]->y_ceed; 344d4d45553Srezgarshakeri pr_restr_ctx[i]->op_prolong = ceed_data[i]->op_prolong; 345d4d45553Srezgarshakeri pr_restr_ctx[i]->op_restrict = ceed_data[i]->op_restrict; 346d4d45553Srezgarshakeri pr_restr_ctx[i]->ceed = ceed; 3476c5df90dSjeremylt } 3486c5df90dSjeremylt } 3496c5df90dSjeremylt 35053b04fa6SJed Brown // Assemble coarse grid Jacobian for AMG (or other sparse matrix) solve 3512b730f8bSJeremy L Thompson PetscCall(DMCreateMatrix(dm[0], &mat_coarse)); 35253b04fa6SJed Brown 3532b730f8bSJeremy L Thompson PetscCall(PetscLogEventRegister("AssembleMatrix", MAT_CLASSID, &assemble_event)); 354cffe6a52SJeremy L Thompson { 355cffe6a52SJeremy L Thompson // Assemble matrix analytically 3563047f789SJeremy L Thompson PetscCount num_entries; 3573047f789SJeremy L Thompson CeedInt *rows, *cols; 35853b04fa6SJed Brown CeedVector coo_values; 3592b730f8bSJeremy L Thompson CeedOperatorLinearAssembleSymbolic(op_apply_ctx[0]->op, &num_entries, &rows, &cols); 36053b04fa6SJed Brown ISLocalToGlobalMapping ltog_row, ltog_col; 3612b730f8bSJeremy L Thompson PetscCall(MatGetLocalToGlobalMapping(mat_coarse, <og_row, <og_col)); 3622b730f8bSJeremy L Thompson PetscCall(ISLocalToGlobalMappingApply(ltog_row, num_entries, rows, rows)); 3632b730f8bSJeremy L Thompson PetscCall(ISLocalToGlobalMappingApply(ltog_col, num_entries, cols, cols)); 3642b730f8bSJeremy L Thompson PetscCall(MatSetPreallocationCOO(mat_coarse, num_entries, rows, cols)); 36553b04fa6SJed Brown free(rows); 36653b04fa6SJed Brown free(cols); 36753b04fa6SJed Brown CeedVectorCreate(ceed, num_entries, &coo_values); 3682b730f8bSJeremy L Thompson PetscCall(PetscLogEventBegin(assemble_event, mat_coarse, 0, 0, 0)); 369d4d45553Srezgarshakeri CeedOperatorLinearAssemble(op_apply_ctx[0]->op, coo_values); 37053b04fa6SJed Brown const CeedScalar *values; 37153b04fa6SJed Brown CeedVectorGetArrayRead(coo_values, CEED_MEM_HOST, &values); 3722b730f8bSJeremy L Thompson PetscCall(MatSetValuesCOO(mat_coarse, values, ADD_VALUES)); 37353b04fa6SJed Brown CeedVectorRestoreArrayRead(coo_values, &values); 3742b730f8bSJeremy L Thompson PetscCall(PetscLogEventEnd(assemble_event, mat_coarse, 0, 0, 0)); 37553b04fa6SJed Brown CeedVectorDestroy(&coo_values); 37653b04fa6SJed Brown } 37715ce0ef0Sjeremylt 3786c5df90dSjeremylt // Set up KSP 3792b730f8bSJeremy L Thompson PetscCall(KSPCreate(comm, &ksp)); 3806c5df90dSjeremylt { 3812b730f8bSJeremy L Thompson PetscCall(KSPSetType(ksp, KSPCG)); 3822b730f8bSJeremy L Thompson PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL)); 3832b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT)); 3846c5df90dSjeremylt } 3852b730f8bSJeremy L Thompson PetscCall(KSPSetFromOptions(ksp)); 3862b730f8bSJeremy L Thompson PetscCall(KSPSetOperators(ksp, mat_O[fine_level], mat_O[fine_level])); 3876c5df90dSjeremylt 3886c5df90dSjeremylt // Set up PCMG 3892b730f8bSJeremy L Thompson PetscCall(KSPGetPC(ksp, &pc)); 3909b072555Sjeremylt PCMGCycleType pcmg_cycle_type = PC_MG_CYCLE_V; 3916c5df90dSjeremylt { 3922b730f8bSJeremy L Thompson PetscCall(PCSetType(pc, PCMG)); 3936c5df90dSjeremylt 3946c5df90dSjeremylt // PCMG levels 3952b730f8bSJeremy L Thompson PetscCall(PCMGSetLevels(pc, num_levels, NULL)); 3969b072555Sjeremylt for (int i = 0; i < num_levels; i++) { 3976c5df90dSjeremylt // Smoother 3986c5df90dSjeremylt KSP smoother; 3996c5df90dSjeremylt PC smoother_pc; 4002b730f8bSJeremy L Thompson PetscCall(PCMGGetSmoother(pc, i, &smoother)); 4012b730f8bSJeremy L Thompson PetscCall(KSPSetType(smoother, KSPCHEBYSHEV)); 4022b730f8bSJeremy L Thompson PetscCall(KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1)); 4032b730f8bSJeremy L Thompson PetscCall(KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE)); 4042b730f8bSJeremy L Thompson PetscCall(KSPSetOperators(smoother, mat_O[i], mat_O[i])); 4052b730f8bSJeremy L Thompson PetscCall(KSPGetPC(smoother, &smoother_pc)); 4062b730f8bSJeremy L Thompson PetscCall(PCSetType(smoother_pc, PCJACOBI)); 4072b730f8bSJeremy L Thompson PetscCall(PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL)); 4086c5df90dSjeremylt 4096c5df90dSjeremylt // Work vector 4109b072555Sjeremylt if (i < num_levels - 1) { 4112b730f8bSJeremy L Thompson PetscCall(PCMGSetX(pc, i, X[i])); 4126c5df90dSjeremylt } 4136c5df90dSjeremylt 4146c5df90dSjeremylt // Level transfers 4156c5df90dSjeremylt if (i > 0) { 4166c5df90dSjeremylt // Interpolation 4172b730f8bSJeremy L Thompson PetscCall(PCMGSetInterpolation(pc, i, mat_pr[i])); 4186c5df90dSjeremylt } 4196c5df90dSjeremylt 4206c5df90dSjeremylt // Coarse solve 4216c5df90dSjeremylt KSP coarse; 4226c5df90dSjeremylt PC coarse_pc; 4232b730f8bSJeremy L Thompson PetscCall(PCMGGetCoarseSolve(pc, &coarse)); 4242b730f8bSJeremy L Thompson PetscCall(KSPSetType(coarse, KSPPREONLY)); 4252b730f8bSJeremy L Thompson PetscCall(KSPSetOperators(coarse, mat_coarse, mat_coarse)); 42615ce0ef0Sjeremylt 4272b730f8bSJeremy L Thompson PetscCall(KSPGetPC(coarse, &coarse_pc)); 4282b730f8bSJeremy L Thompson PetscCall(PCSetType(coarse_pc, PCGAMG)); 42915ce0ef0Sjeremylt 4302b730f8bSJeremy L Thompson PetscCall(KSPSetOptionsPrefix(coarse, "coarse_")); 4312b730f8bSJeremy L Thompson PetscCall(PCSetOptionsPrefix(coarse_pc, "coarse_")); 4322b730f8bSJeremy L Thompson PetscCall(KSPSetFromOptions(coarse)); 4332b730f8bSJeremy L Thompson PetscCall(PCSetFromOptions(coarse_pc)); 4346c5df90dSjeremylt } 4356c5df90dSjeremylt 4366c5df90dSjeremylt // PCMG options 4372b730f8bSJeremy L Thompson PetscCall(PCMGSetType(pc, PC_MG_MULTIPLICATIVE)); 4382b730f8bSJeremy L Thompson PetscCall(PCMGSetNumberSmooth(pc, 3)); 4392b730f8bSJeremy L Thompson PetscCall(PCMGSetCycleType(pc, pcmg_cycle_type)); 4406c5df90dSjeremylt } 4416c5df90dSjeremylt 4426c5df90dSjeremylt // First run, if benchmarking 4436c5df90dSjeremylt if (benchmark_mode) { 4442b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1)); 4452b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X[fine_level])); 4466c5df90dSjeremylt my_rt_start = MPI_Wtime(); 4472b730f8bSJeremy L Thompson PetscCall(KSPSolve(ksp, rhs, X[fine_level])); 4486c5df90dSjeremylt my_rt = MPI_Wtime() - my_rt_start; 4492b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm)); 4506c5df90dSjeremylt // Set maxits based on first iteration timing 4516c5df90dSjeremylt if (my_rt > 0.02) { 4522b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5)); 4536c5df90dSjeremylt } else { 4542b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20)); 4556c5df90dSjeremylt } 4566c5df90dSjeremylt } 4576c5df90dSjeremylt 4586c5df90dSjeremylt // Timed solve 4592b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X[fine_level])); 4602b730f8bSJeremy L Thompson PetscCall(PetscBarrier((PetscObject)ksp)); 46109a940d7Sjeremylt 46209a940d7Sjeremylt // -- Performance logging 4632b730f8bSJeremy L Thompson PetscCall(PetscLogStageRegister("Solve Stage", &solve_stage)); 4642b730f8bSJeremy L Thompson PetscCall(PetscLogStagePush(solve_stage)); 46509a940d7Sjeremylt 46609a940d7Sjeremylt // -- Solve 4676c5df90dSjeremylt my_rt_start = MPI_Wtime(); 4682b730f8bSJeremy L Thompson PetscCall(KSPSolve(ksp, rhs, X[fine_level])); 4696c5df90dSjeremylt my_rt = MPI_Wtime() - my_rt_start; 4706c5df90dSjeremylt 47109a940d7Sjeremylt // -- Performance logging 4722b730f8bSJeremy L Thompson PetscCall(PetscLogStagePop()); 47309a940d7Sjeremylt 4746c5df90dSjeremylt // Output results 4756c5df90dSjeremylt { 4769b072555Sjeremylt KSPType ksp_type; 4779b072555Sjeremylt PCMGType pcmg_type; 4786c5df90dSjeremylt KSPConvergedReason reason; 4796c5df90dSjeremylt PetscReal rnorm; 4806c5df90dSjeremylt PetscInt its; 4812b730f8bSJeremy L Thompson PetscCall(KSPGetType(ksp, &ksp_type)); 4822b730f8bSJeremy L Thompson PetscCall(KSPGetConvergedReason(ksp, &reason)); 4832b730f8bSJeremy L Thompson PetscCall(KSPGetIterationNumber(ksp, &its)); 4842b730f8bSJeremy L Thompson PetscCall(KSPGetResidualNorm(ksp, &rnorm)); 4852b730f8bSJeremy L Thompson PetscCall(PCMGGetType(pc, &pcmg_type)); 4866c5df90dSjeremylt if (!test_mode || reason < 0 || rnorm > 1e-8) { 4872b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 4886c5df90dSjeremylt " KSP:\n" 4896c5df90dSjeremylt " KSP Type : %s\n" 4906c5df90dSjeremylt " KSP Convergence : %s\n" 491a9b2c5ddSrezgarshakeri " Total KSP Iterations : %" PetscInt_FMT "\n" 4926c5df90dSjeremylt " Final rnorm : %e\n", 4932b730f8bSJeremy L Thompson ksp_type, KSPConvergedReasons[reason], its, (double)rnorm)); 4942b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 4956c5df90dSjeremylt " PCMG:\n" 4966c5df90dSjeremylt " PCMG Type : %s\n" 4976c5df90dSjeremylt " PCMG Cycle Type : %s\n", 4982b730f8bSJeremy L Thompson PCMGTypes[pcmg_type], PCMGCycleTypes[pcmg_cycle_type])); 4996c5df90dSjeremylt } 5006c5df90dSjeremylt if (!test_mode) { 5012b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, " Performance:\n")); 5026c5df90dSjeremylt } 5036c5df90dSjeremylt { 5046c88e6a2Srezgarshakeri // Set up error operator context 5052b730f8bSJeremy L Thompson PetscCall(SetupErrorOperatorCtx(comm, dm[fine_level], ceed, ceed_data[fine_level], X_loc[fine_level], op_error, op_error_ctx)); 50638f32c05Srezgarshakeri PetscScalar l2_error; 5072b730f8bSJeremy L Thompson PetscCall(ComputeL2Error(X[fine_level], &l2_error, op_error_ctx)); 5080a8fc04aSrezgarshakeri PetscReal tol = 5e-2; 50938f32c05Srezgarshakeri if (!test_mode || l2_error > tol) { 5102b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm)); 5112b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm)); 5122b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, 51338f32c05Srezgarshakeri " L2 Error : %e\n" 5146c5df90dSjeremylt " CG Solve Time : %g (%g) sec\n", 5152b730f8bSJeremy L Thompson (double)l2_error, rt_max, rt_min)); 5166c5df90dSjeremylt } 5176c5df90dSjeremylt } 5186c5df90dSjeremylt if (benchmark_mode && (!test_mode)) { 5192b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, " DoFs/Sec in CG : %g (%g) million\n", 1e-6 * g_size[fine_level] * its / rt_max, 5202b730f8bSJeremy L Thompson 1e-6 * g_size[fine_level] * its / rt_min)); 5216c5df90dSjeremylt } 5226c5df90dSjeremylt } 5236c5df90dSjeremylt 5246c5df90dSjeremylt if (write_solution) { 5259b072555Sjeremylt PetscViewer vtk_viewer_soln; 5266c5df90dSjeremylt 5272b730f8bSJeremy L Thompson PetscCall(PetscViewerCreate(comm, &vtk_viewer_soln)); 5282b730f8bSJeremy L Thompson PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK)); 5292b730f8bSJeremy L Thompson PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu")); 5302b730f8bSJeremy L Thompson PetscCall(VecView(X[fine_level], vtk_viewer_soln)); 5312b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&vtk_viewer_soln)); 5326c5df90dSjeremylt } 5336c5df90dSjeremylt 5346c5df90dSjeremylt // Cleanup 5359b072555Sjeremylt for (int i = 0; i < num_levels; i++) { 5362b730f8bSJeremy L Thompson PetscCall(VecDestroy(&X[i])); 5372b730f8bSJeremy L Thompson PetscCall(VecDestroy(&X_loc[i])); 5382b730f8bSJeremy L Thompson PetscCall(VecDestroy(&mult[i])); 5392b730f8bSJeremy L Thompson PetscCall(VecDestroy(&op_apply_ctx[i]->Y_loc)); 5402b730f8bSJeremy L Thompson PetscCall(MatDestroy(&mat_O[i])); 5412b730f8bSJeremy L Thompson PetscCall(PetscFree(op_apply_ctx[i])); 5426c5df90dSjeremylt if (i > 0) { 5432b730f8bSJeremy L Thompson PetscCall(MatDestroy(&mat_pr[i])); 5442b730f8bSJeremy L Thompson PetscCall(PetscFree(pr_restr_ctx[i])); 5456c5df90dSjeremylt } 5462b730f8bSJeremy L Thompson PetscCall(CeedDataDestroy(i, ceed_data[i])); 5472b730f8bSJeremy L Thompson PetscCall(DMDestroy(&dm[i])); 5486c5df90dSjeremylt } 5492b730f8bSJeremy L Thompson PetscCall(PetscFree(level_degrees)); 5502b730f8bSJeremy L Thompson PetscCall(PetscFree(dm)); 5512b730f8bSJeremy L Thompson PetscCall(PetscFree(X)); 5522b730f8bSJeremy L Thompson PetscCall(PetscFree(X_loc)); 5532b730f8bSJeremy L Thompson PetscCall(VecDestroy(&op_error_ctx->Y_loc)); 5542b730f8bSJeremy L Thompson PetscCall(PetscFree(mult)); 5552b730f8bSJeremy L Thompson PetscCall(PetscFree(mat_O)); 5562b730f8bSJeremy L Thompson PetscCall(PetscFree(mat_pr)); 5572b730f8bSJeremy L Thompson PetscCall(PetscFree(ceed_data)); 5582b730f8bSJeremy L Thompson PetscCall(PetscFree(op_apply_ctx)); 5592b730f8bSJeremy L Thompson PetscCall(PetscFree(op_error_ctx)); 5602b730f8bSJeremy L Thompson PetscCall(PetscFree(pr_restr_ctx)); 5612b730f8bSJeremy L Thompson PetscCall(PetscFree(l_size)); 5622b730f8bSJeremy L Thompson PetscCall(PetscFree(xl_size)); 5632b730f8bSJeremy L Thompson PetscCall(PetscFree(g_size)); 5642b730f8bSJeremy L Thompson PetscCall(VecDestroy(&rhs)); 5652b730f8bSJeremy L Thompson PetscCall(VecDestroy(&rhs_loc)); 5662b730f8bSJeremy L Thompson PetscCall(MatDestroy(&mat_coarse)); 5672b730f8bSJeremy L Thompson PetscCall(KSPDestroy(&ksp)); 5682b730f8bSJeremy L Thompson PetscCall(DMDestroy(&dm_orig)); 5696c5df90dSjeremylt CeedVectorDestroy(&target); 5709b072555Sjeremylt CeedQFunctionDestroy(&qf_error); 5719b072555Sjeremylt CeedOperatorDestroy(&op_error); 5726c5df90dSjeremylt CeedDestroy(&ceed); 5736c5df90dSjeremylt return PetscFinalize(); 5746c5df90dSjeremylt } 575