1*9ba83ac0SJeremy L Thompson // Copyright (c) 2017-2026, 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
main(int argc,char ** argv)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;
54cfa59c5bSRey PetscScalar eps = 1.0;
55de1229c5Srezgarshakeri PetscBool test_mode, benchmark_mode, read_mesh, write_solution, simplex;
569b072555Sjeremylt PetscLogStage solve_stage;
5705b9c820SJed Brown PetscLogEvent assemble_event;
589b072555Sjeremylt DM *dm, dm_orig;
596c5df90dSjeremylt KSP ksp;
606c5df90dSjeremylt PC pc;
619b072555Sjeremylt Mat *mat_O, *mat_pr, mat_coarse;
629b072555Sjeremylt Vec *X, *X_loc, *mult, rhs, rhs_loc;
639b072555Sjeremylt PetscMemType mem_type;
646c88e6a2Srezgarshakeri OperatorApplyContext *op_apply_ctx, op_error_ctx;
65d4d45553Srezgarshakeri ProlongRestrContext *pr_restr_ctx;
666c5df90dSjeremylt Ceed ceed;
679b072555Sjeremylt CeedData *ceed_data;
689b072555Sjeremylt CeedVector rhs_ceed, target;
691c9a79dbSrezgarshakeri CeedQFunction qf_error;
709b072555Sjeremylt CeedOperator op_error;
719b072555Sjeremylt BPType bp_choice;
729b072555Sjeremylt CoarsenType coarsen;
736c5df90dSjeremylt
742b730f8bSJeremy L Thompson PetscCall(PetscInitialize(&argc, &argv, NULL, help));
756c5df90dSjeremylt comm = PETSC_COMM_WORLD;
766c5df90dSjeremylt
776c5df90dSjeremylt // Parse command line options
7867490bc6SJeremy L Thompson PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
799b072555Sjeremylt bp_choice = CEED_BP3;
802b730f8bSJeremy L Thompson PetscCall(PetscOptionsEnum("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice, NULL));
819b072555Sjeremylt num_comp_u = bp_options[bp_choice].num_comp_u;
826c5df90dSjeremylt test_mode = PETSC_FALSE;
832b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, test_mode, &test_mode, NULL));
846c5df90dSjeremylt benchmark_mode = PETSC_FALSE;
852b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-benchmark", "Benchmarking mode (prints benchmark statistics)", NULL, benchmark_mode, &benchmark_mode, NULL));
866c5df90dSjeremylt write_solution = PETSC_FALSE;
872b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, write_solution, &write_solution, NULL));
88de1229c5Srezgarshakeri simplex = PETSC_FALSE;
892b730f8bSJeremy L Thompson PetscCall(PetscOptionsBool("-simplex", "Element topology (default:hex)", NULL, simplex, &simplex, NULL));
90eb6a3b92Srezgarshakeri if ((bp_choice == CEED_BP5 || bp_choice == CEED_BP6) && (simplex)) {
912b730f8bSJeremy L Thompson SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "BP5/6 is not supported with simplex");
92eb6a3b92Srezgarshakeri }
932b730f8bSJeremy L Thompson PetscCall(PetscOptionsScalar("-eps", "Epsilon parameter for Kershaw mesh transformation", NULL, eps, &eps, NULL));
942b730f8bSJeremy L Thompson if (eps > 1 || eps <= 0) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "-eps %g must be (0,1]", (double)PetscRealPart(eps));
956c5df90dSjeremylt degree = test_mode ? 3 : 2;
962b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", NULL, degree, °ree, NULL));
972b730f8bSJeremy L Thompson if (degree < 1) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "-degree %" PetscInt_FMT " must be at least 1", degree);
989b072555Sjeremylt q_extra = bp_options[bp_choice].q_extra;
992b730f8bSJeremy L Thompson PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points", NULL, q_extra, &q_extra, NULL));
1002b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-ceed", "CEED resource specifier", NULL, ceed_resource, ceed_resource, sizeof(ceed_resource), NULL));
1016c5df90dSjeremylt coarsen = COARSEN_UNIFORM;
1022b730f8bSJeremy L Thompson PetscCall(PetscOptionsEnum("-coarsen", "Coarsening strategy to use", NULL, coarsen_types, (PetscEnum)coarsen, (PetscEnum *)&coarsen, NULL));
103cb32e2e7SValeria Barra read_mesh = PETSC_FALSE;
1042b730f8bSJeremy L Thompson PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &read_mesh));
1056c5df90dSjeremylt if (!read_mesh) {
1066c5df90dSjeremylt PetscInt tmp = dim;
1072b730f8bSJeremy L Thompson PetscCall(PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL, mesh_elem, &tmp, NULL));
1086c5df90dSjeremylt }
10967490bc6SJeremy L Thompson PetscOptionsEnd();
1106c5df90dSjeremylt
1119396343dSjeremylt // Set up libCEED
1129b072555Sjeremylt CeedInit(ceed_resource, &ceed);
1139b072555Sjeremylt CeedMemType mem_type_backend;
1149b072555Sjeremylt CeedGetPreferredMemType(ceed, &mem_type_backend);
1159396343dSjeremylt
1166c5df90dSjeremylt // Setup DM
1176c5df90dSjeremylt if (read_mesh) {
1182b730f8bSJeremy L Thompson PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, filename, NULL, PETSC_TRUE, &dm_orig));
1196c5df90dSjeremylt } else {
1207cf95199SJeremy L Thompson PetscCall(DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, simplex, mesh_elem, NULL, NULL, NULL, PETSC_TRUE, 0, PETSC_FALSE, &dm_orig));
1216c5df90dSjeremylt }
1226c5df90dSjeremylt
1238c03e814SJeremy L Thompson // Set mesh vec_type
1248c03e814SJeremy L Thompson VecType vec_type = VECSTANDARD;
1258c03e814SJeremy L Thompson
1269b072555Sjeremylt switch (mem_type_backend) {
1272b730f8bSJeremy L Thompson case CEED_MEM_HOST:
1282b730f8bSJeremy L Thompson vec_type = VECSTANDARD;
1292b730f8bSJeremy L Thompson break;
130b68a8d79SJed Brown case CEED_MEM_DEVICE: {
131b68a8d79SJed Brown const char *resolved;
1328c03e814SJeremy L Thompson
133b68a8d79SJed Brown CeedGetResource(ceed, &resolved);
1349b072555Sjeremylt if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
1359b072555Sjeremylt else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
1369b072555Sjeremylt else vec_type = VECSTANDARD;
137b68a8d79SJed Brown }
138b68a8d79SJed Brown }
1392b730f8bSJeremy L Thompson PetscCall(DMSetVecType(dm_orig, vec_type));
1402b730f8bSJeremy L Thompson PetscCall(DMSetFromOptions(dm_orig));
1412b730f8bSJeremy L Thompson PetscCall(DMViewFromOptions(dm_orig, NULL, "-dm_view"));
1423fc8a154SJed Brown
1433fc8a154SJed Brown // Apply Kershaw mesh transformation
1442b730f8bSJeremy L Thompson PetscCall(Kershaw(dm_orig, eps));
145b68a8d79SJed Brown
1466c5df90dSjeremylt // Allocate arrays for PETSc objects for each level
1476c5df90dSjeremylt switch (coarsen) {
1486c5df90dSjeremylt case COARSEN_UNIFORM:
1499b072555Sjeremylt num_levels = degree;
1506c5df90dSjeremylt break;
151dc7d240cSValeria Barra case COARSEN_LOGARITHMIC:
1529b072555Sjeremylt num_levels = ceil(log(degree) / log(2)) + 1;
1536c5df90dSjeremylt break;
1546c5df90dSjeremylt }
1552b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &level_degrees));
1569b072555Sjeremylt fine_level = num_levels - 1;
15761608365Sjeremylt
1586c5df90dSjeremylt switch (coarsen) {
1596c5df90dSjeremylt case COARSEN_UNIFORM:
1604dbe2ad5SJeremy L Thompson for (PetscInt i = 0; i < num_levels; i++) level_degrees[i] = i + 1;
1616c5df90dSjeremylt break;
162dc7d240cSValeria Barra case COARSEN_LOGARITHMIC:
1634dbe2ad5SJeremy L Thompson for (PetscInt i = 0; i < num_levels - 1; i++) level_degrees[i] = pow(2, i);
1649b072555Sjeremylt level_degrees[fine_level] = degree;
1656c5df90dSjeremylt break;
1666c5df90dSjeremylt }
1672b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &dm));
1682b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &X));
1692b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &X_loc));
1702b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &mult));
1712b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &op_apply_ctx));
1722b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &pr_restr_ctx));
1732b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &mat_O));
1742b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &mat_pr));
1752b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &l_size));
1762b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &xl_size));
1772b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &g_size));
1786c5df90dSjeremylt
1791c9a79dbSrezgarshakeri PetscInt c_start, c_end;
1802b730f8bSJeremy L Thompson PetscCall(DMPlexGetHeightStratum(dm_orig, 0, &c_start, &c_end));
1811c9a79dbSrezgarshakeri DMPolytopeType cell_type;
1822b730f8bSJeremy L Thompson PetscCall(DMPlexGetCellType(dm_orig, c_start, &cell_type));
1831c9a79dbSrezgarshakeri CeedElemTopology elem_topo = ElemTopologyP2C(cell_type);
1841c9a79dbSrezgarshakeri
1856c5df90dSjeremylt // Setup DM and Operator Mat Shells for each level
186d2fd4e27SJeremy L Thompson for (PetscInt i = 0; i < num_levels; i++) {
1876c5df90dSjeremylt // Create DM
1882b730f8bSJeremy L Thompson PetscCall(DMClone(dm_orig, &dm[i]));
1892b730f8bSJeremy L Thompson PetscCall(DMGetVecType(dm_orig, &vec_type));
1902b730f8bSJeremy L Thompson PetscCall(DMSetVecType(dm[i], vec_type));
191e83e87a5Sjeremylt PetscInt dim;
1922b730f8bSJeremy L Thompson PetscCall(DMGetDimension(dm[i], &dim));
1932b730f8bSJeremy L Thompson PetscCall(SetupDMByDegree(dm[i], level_degrees[fine_level], q_extra, num_comp_u, dim, bp_options[bp_choice].enforce_bc));
1946c5df90dSjeremylt
1956c5df90dSjeremylt // Create vectors
1962b730f8bSJeremy L Thompson PetscCall(DMCreateGlobalVector(dm[i], &X[i]));
1972b730f8bSJeremy L Thompson PetscCall(VecGetLocalSize(X[i], &l_size[i]));
1982b730f8bSJeremy L Thompson PetscCall(VecGetSize(X[i], &g_size[i]));
1992b730f8bSJeremy L Thompson PetscCall(DMCreateLocalVector(dm[i], &X_loc[i]));
2002b730f8bSJeremy L Thompson PetscCall(VecGetSize(X_loc[i], &xl_size[i]));
2016c5df90dSjeremylt
2026c5df90dSjeremylt // Operator
2032b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &op_apply_ctx[i]));
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"
2384d00b080SJeremy L Thompson " Solution Order (P) : %" PetscInt_FMT "\n"
2394d00b080SJeremy L Thompson " Quadrature Order (Q) : %" PetscInt_FMT "\n"
2404d00b080SJeremy L Thompson " Additional quadrature points (q_extra) : %" PetscInt_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"
2464d00b080SJeremy L Thompson " Number of Levels : %" PetscInt_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));
2549b072555Sjeremylt CeedVectorCreate(ceed, xl_size[fine_level], &rhs_ceed);
255179e5961SZach Atkins PetscCall(VecP2C(rhs_loc, &mem_type, rhs_ceed));
2566c5df90dSjeremylt
2576c5df90dSjeremylt // Set up libCEED operators on each level
2582b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(num_levels, &ceed_data));
259990fdeb6SJeremy L Thompson for (PetscInt i = 0; i < num_levels; i++) {
2606c5df90dSjeremylt // Print level information
2619b072555Sjeremylt if (!test_mode && (i == 0 || i == fine_level)) {
2622b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm,
2632b730f8bSJeremy L Thompson " Level %" PetscInt_FMT " (%s):\n"
2644d00b080SJeremy L Thompson " Solution Order (P) : %" PetscInt_FMT "\n"
26508140895SJed Brown " Global Nodes : %" PetscInt_FMT "\n"
26608140895SJed Brown " Owned Nodes : %" PetscInt_FMT "\n",
2672b730f8bSJeremy L Thompson i, (i ? "fine" : "coarse"), level_degrees[i] + 1, g_size[i] / num_comp_u, l_size[i] / num_comp_u));
2686c5df90dSjeremylt }
2692b730f8bSJeremy L Thompson PetscCall(PetscMalloc1(1, &ceed_data[i]));
2702b730f8bSJeremy 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],
2714dbe2ad5SJeremy L Thompson ceed_data[i], i == fine_level, i == fine_level, rhs_ceed, &target));
2726c5df90dSjeremylt }
2736c5df90dSjeremylt
2746c5df90dSjeremylt // Gather RHS
275179e5961SZach Atkins PetscCall(VecC2P(rhs_ceed, mem_type, rhs_loc));
2762b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(rhs));
2772b730f8bSJeremy L Thompson PetscCall(DMLocalToGlobal(dm[fine_level], rhs_loc, ADD_VALUES, rhs));
2789b072555Sjeremylt CeedVectorDestroy(&rhs_ceed);
2796c5df90dSjeremylt
28043eb8658SJeremy L Thompson // Create the error QFunction
2812b730f8bSJeremy L Thompson CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error, bp_options[bp_choice].error_loc, &qf_error);
2829b072555Sjeremylt CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
2839b072555Sjeremylt CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
2842b730f8bSJeremy L Thompson CeedQFunctionAddInput(qf_error, "qdata", ceed_data[fine_level]->q_data_size, CEED_EVAL_NONE);
28538f32c05Srezgarshakeri CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_INTERP);
2866c5df90dSjeremylt
2876c5df90dSjeremylt // Create the error operator
2882b730f8bSJeremy L Thompson CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_error);
2892b730f8bSJeremy L Thompson CeedOperatorSetField(op_error, "u", ceed_data[fine_level]->elem_restr_u, ceed_data[fine_level]->basis_u, CEED_VECTOR_ACTIVE);
290356036faSJeremy L Thompson CeedOperatorSetField(op_error, "true_soln", ceed_data[fine_level]->elem_restr_u_i, CEED_BASIS_NONE, target);
291356036faSJeremy L Thompson CeedOperatorSetField(op_error, "qdata", ceed_data[fine_level]->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data[fine_level]->q_data);
2922b730f8bSJeremy L Thompson CeedOperatorSetField(op_error, "error", ceed_data[fine_level]->elem_restr_u, ceed_data[fine_level]->basis_u, CEED_VECTOR_ACTIVE);
2936c5df90dSjeremylt
2946c5df90dSjeremylt // Calculate multiplicity
2954dbe2ad5SJeremy L Thompson for (PetscInt i = 0; i < num_levels; i++) {
296179e5961SZach Atkins PetscMemType mem_type;
2976c5df90dSjeremylt
2986c5df90dSjeremylt // CEED vector
2992b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X_loc[i]));
300179e5961SZach Atkins PetscCall(VecP2C(X_loc[i], &mem_type, ceed_data[i]->x_ceed));
3016c5df90dSjeremylt
3026c5df90dSjeremylt // Multiplicity
3032b730f8bSJeremy L Thompson CeedElemRestrictionGetMultiplicity(ceed_data[i]->elem_restr_u, ceed_data[i]->x_ceed);
3049b072555Sjeremylt CeedVectorSyncArray(ceed_data[i]->x_ceed, CEED_MEM_HOST);
3056c5df90dSjeremylt
3066c5df90dSjeremylt // Restore vector
307179e5961SZach Atkins PetscCall(VecC2P(ceed_data[i]->x_ceed, mem_type, X_loc[i]));
3086c5df90dSjeremylt
3096c5df90dSjeremylt // Creat mult vector
3102b730f8bSJeremy L Thompson PetscCall(VecDuplicate(X_loc[i], &mult[i]));
3116c5df90dSjeremylt
3126c5df90dSjeremylt // Local-to-global
3132b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X[i]));
3142b730f8bSJeremy L Thompson PetscCall(DMLocalToGlobal(dm[i], X_loc[i], ADD_VALUES, X[i]));
3152b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X_loc[i]));
3166c5df90dSjeremylt
3176c5df90dSjeremylt // Global-to-local
3182b730f8bSJeremy L Thompson PetscCall(DMGlobalToLocal(dm[i], X[i], INSERT_VALUES, mult[i]));
3192b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X[i]));
3206c5df90dSjeremylt
3216c5df90dSjeremylt // Multiplicity scaling
3222b730f8bSJeremy L Thompson PetscCall(VecReciprocal(mult[i]));
3236c5df90dSjeremylt }
3246c5df90dSjeremylt
3256c5df90dSjeremylt // Set up Mat
3264dbe2ad5SJeremy L Thompson for (PetscInt i = fine_level; i >= 0; i--) {
3276c88e6a2Srezgarshakeri // Set up apply operator context
3282b730f8bSJeremy L Thompson PetscCall(SetupApplyOperatorCtx(comm, dm[i], ceed, ceed_data[i], X_loc[i], op_apply_ctx[i]));
3296c5df90dSjeremylt
3306c5df90dSjeremylt if (i > 0) {
331a97643b0Sjeremylt // Prolongation/Restriction Operator
3322b730f8bSJeremy L Thompson PetscCall(CeedLevelTransferSetup(dm[i - 1], ceed, i, num_comp_u, ceed_data, bp_options[bp_choice], mult[i]));
333d4d45553Srezgarshakeri pr_restr_ctx[i]->comm = comm;
334d4d45553Srezgarshakeri pr_restr_ctx[i]->dmf = dm[i];
335d4d45553Srezgarshakeri pr_restr_ctx[i]->dmc = dm[i - 1];
336d4d45553Srezgarshakeri pr_restr_ctx[i]->loc_vec_c = X_loc[i - 1];
337d4d45553Srezgarshakeri pr_restr_ctx[i]->loc_vec_f = op_apply_ctx[i]->Y_loc;
338d4d45553Srezgarshakeri pr_restr_ctx[i]->mult_vec = mult[i];
3394dbe2ad5SJeremy L Thompson pr_restr_ctx[i]->ceed_vec_c = ceed_data[i - 1]->x_ceed;
3404dbe2ad5SJeremy L Thompson pr_restr_ctx[i]->ceed_vec_f = ceed_data[i]->y_ceed;
341d4d45553Srezgarshakeri pr_restr_ctx[i]->op_prolong = ceed_data[i]->op_prolong;
342d4d45553Srezgarshakeri pr_restr_ctx[i]->op_restrict = ceed_data[i]->op_restrict;
343d4d45553Srezgarshakeri pr_restr_ctx[i]->ceed = ceed;
3446c5df90dSjeremylt }
3456c5df90dSjeremylt }
3466c5df90dSjeremylt
34753b04fa6SJed Brown // Assemble coarse grid Jacobian for AMG (or other sparse matrix) solve
3482b730f8bSJeremy L Thompson PetscCall(DMCreateMatrix(dm[0], &mat_coarse));
34953b04fa6SJed Brown
3502b730f8bSJeremy L Thompson PetscCall(PetscLogEventRegister("AssembleMatrix", MAT_CLASSID, &assemble_event));
351cffe6a52SJeremy L Thompson {
352cffe6a52SJeremy L Thompson // Assemble matrix analytically
3533047f789SJeremy L Thompson PetscCount num_entries;
3544d00b080SJeremy L Thompson CeedInt *rows_ceed, *cols_ceed;
3554d00b080SJeremy L Thompson PetscInt *rows_petsc, *cols_petsc;
35653b04fa6SJed Brown ISLocalToGlobalMapping ltog_row, ltog_col;
3574d00b080SJeremy L Thompson CeedVector coo_values;
3584d00b080SJeremy L Thompson
3594d00b080SJeremy L Thompson CeedOperatorLinearAssembleSymbolic(op_apply_ctx[0]->op, &num_entries, &rows_ceed, &cols_ceed);
3604d00b080SJeremy L Thompson PetscCall(IntArrayCeedToPetsc(num_entries, &rows_ceed, &rows_petsc));
3614d00b080SJeremy L Thompson PetscCall(IntArrayCeedToPetsc(num_entries, &cols_ceed, &cols_petsc));
3622b730f8bSJeremy L Thompson PetscCall(MatGetLocalToGlobalMapping(mat_coarse, <og_row, <og_col));
3634d00b080SJeremy L Thompson PetscCall(ISLocalToGlobalMappingApply(ltog_row, num_entries, rows_petsc, rows_petsc));
3644d00b080SJeremy L Thompson PetscCall(ISLocalToGlobalMappingApply(ltog_col, num_entries, cols_petsc, cols_petsc));
3654d00b080SJeremy L Thompson PetscCall(MatSetPreallocationCOO(mat_coarse, num_entries, rows_petsc, cols_petsc));
3664d00b080SJeremy L Thompson free(rows_petsc);
3674d00b080SJeremy L Thompson free(cols_petsc);
36853b04fa6SJed Brown CeedVectorCreate(ceed, num_entries, &coo_values);
3692b730f8bSJeremy L Thompson PetscCall(PetscLogEventBegin(assemble_event, mat_coarse, 0, 0, 0));
370d4d45553Srezgarshakeri CeedOperatorLinearAssemble(op_apply_ctx[0]->op, coo_values);
37153b04fa6SJed Brown const CeedScalar *values;
37253b04fa6SJed Brown CeedVectorGetArrayRead(coo_values, CEED_MEM_HOST, &values);
3732b730f8bSJeremy L Thompson PetscCall(MatSetValuesCOO(mat_coarse, values, ADD_VALUES));
37453b04fa6SJed Brown CeedVectorRestoreArrayRead(coo_values, &values);
3752b730f8bSJeremy L Thompson PetscCall(PetscLogEventEnd(assemble_event, mat_coarse, 0, 0, 0));
37653b04fa6SJed Brown CeedVectorDestroy(&coo_values);
37753b04fa6SJed Brown }
37815ce0ef0Sjeremylt
3796c5df90dSjeremylt // Set up KSP
3802b730f8bSJeremy L Thompson PetscCall(KSPCreate(comm, &ksp));
3816c5df90dSjeremylt {
3822b730f8bSJeremy L Thompson PetscCall(KSPSetType(ksp, KSPCG));
3832b730f8bSJeremy L Thompson PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
3842b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
3856c5df90dSjeremylt }
3862b730f8bSJeremy L Thompson PetscCall(KSPSetFromOptions(ksp));
3872b730f8bSJeremy L Thompson PetscCall(KSPSetOperators(ksp, mat_O[fine_level], mat_O[fine_level]));
3886c5df90dSjeremylt
3896c5df90dSjeremylt // Set up PCMG
3902b730f8bSJeremy L Thompson PetscCall(KSPGetPC(ksp, &pc));
3919b072555Sjeremylt PCMGCycleType pcmg_cycle_type = PC_MG_CYCLE_V;
3926c5df90dSjeremylt {
3932b730f8bSJeremy L Thompson PetscCall(PCSetType(pc, PCMG));
3946c5df90dSjeremylt
3956c5df90dSjeremylt // PCMG levels
3962b730f8bSJeremy L Thompson PetscCall(PCMGSetLevels(pc, num_levels, NULL));
3974dbe2ad5SJeremy L Thompson for (PetscInt i = 0; i < num_levels; i++) {
3986c5df90dSjeremylt // Smoother
3996c5df90dSjeremylt KSP smoother;
4006c5df90dSjeremylt PC smoother_pc;
4012b730f8bSJeremy L Thompson PetscCall(PCMGGetSmoother(pc, i, &smoother));
4022b730f8bSJeremy L Thompson PetscCall(KSPSetType(smoother, KSPCHEBYSHEV));
4032b730f8bSJeremy L Thompson PetscCall(KSPChebyshevEstEigSet(smoother, 0, 0.1, 0, 1.1));
4042b730f8bSJeremy L Thompson PetscCall(KSPChebyshevEstEigSetUseNoisy(smoother, PETSC_TRUE));
4052b730f8bSJeremy L Thompson PetscCall(KSPSetOperators(smoother, mat_O[i], mat_O[i]));
4062b730f8bSJeremy L Thompson PetscCall(KSPGetPC(smoother, &smoother_pc));
4072b730f8bSJeremy L Thompson PetscCall(PCSetType(smoother_pc, PCJACOBI));
4082b730f8bSJeremy L Thompson PetscCall(PCJacobiSetType(smoother_pc, PC_JACOBI_DIAGONAL));
4096c5df90dSjeremylt
4106c5df90dSjeremylt // Work vector
4119b072555Sjeremylt if (i < num_levels - 1) {
4122b730f8bSJeremy L Thompson PetscCall(PCMGSetX(pc, i, X[i]));
4136c5df90dSjeremylt }
4146c5df90dSjeremylt
4156c5df90dSjeremylt // Level transfers
4166c5df90dSjeremylt if (i > 0) {
4176c5df90dSjeremylt // Interpolation
4182b730f8bSJeremy L Thompson PetscCall(PCMGSetInterpolation(pc, i, mat_pr[i]));
4196c5df90dSjeremylt }
4206c5df90dSjeremylt
4216c5df90dSjeremylt // Coarse solve
4226c5df90dSjeremylt KSP coarse;
4236c5df90dSjeremylt PC coarse_pc;
4242b730f8bSJeremy L Thompson PetscCall(PCMGGetCoarseSolve(pc, &coarse));
4252b730f8bSJeremy L Thompson PetscCall(KSPSetType(coarse, KSPPREONLY));
4262b730f8bSJeremy L Thompson PetscCall(KSPSetOperators(coarse, mat_coarse, mat_coarse));
42715ce0ef0Sjeremylt
4282b730f8bSJeremy L Thompson PetscCall(KSPGetPC(coarse, &coarse_pc));
4292b730f8bSJeremy L Thompson PetscCall(PCSetType(coarse_pc, PCGAMG));
43015ce0ef0Sjeremylt
4312b730f8bSJeremy L Thompson PetscCall(KSPSetOptionsPrefix(coarse, "coarse_"));
4322b730f8bSJeremy L Thompson PetscCall(PCSetOptionsPrefix(coarse_pc, "coarse_"));
4332b730f8bSJeremy L Thompson PetscCall(KSPSetFromOptions(coarse));
4342b730f8bSJeremy L Thompson PetscCall(PCSetFromOptions(coarse_pc));
4356c5df90dSjeremylt }
4366c5df90dSjeremylt
4376c5df90dSjeremylt // PCMG options
4382b730f8bSJeremy L Thompson PetscCall(PCMGSetType(pc, PC_MG_MULTIPLICATIVE));
4392b730f8bSJeremy L Thompson PetscCall(PCMGSetNumberSmooth(pc, 3));
4402b730f8bSJeremy L Thompson PetscCall(PCMGSetCycleType(pc, pcmg_cycle_type));
4416c5df90dSjeremylt }
4426c5df90dSjeremylt
4436c5df90dSjeremylt // First run, if benchmarking
4446c5df90dSjeremylt if (benchmark_mode) {
4452b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1));
4462b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X[fine_level]));
4476c5df90dSjeremylt my_rt_start = MPI_Wtime();
4482b730f8bSJeremy L Thompson PetscCall(KSPSolve(ksp, rhs, X[fine_level]));
4496c5df90dSjeremylt my_rt = MPI_Wtime() - my_rt_start;
4502b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm));
4516c5df90dSjeremylt // Set maxits based on first iteration timing
4526c5df90dSjeremylt if (my_rt > 0.02) {
4532b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 5));
4546c5df90dSjeremylt } else {
4552b730f8bSJeremy L Thompson PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 20));
4566c5df90dSjeremylt }
4576c5df90dSjeremylt }
4586c5df90dSjeremylt
4596c5df90dSjeremylt // Timed solve
4602b730f8bSJeremy L Thompson PetscCall(VecZeroEntries(X[fine_level]));
4612b730f8bSJeremy L Thompson PetscCall(PetscBarrier((PetscObject)ksp));
46209a940d7Sjeremylt
46309a940d7Sjeremylt // -- Performance logging
4642b730f8bSJeremy L Thompson PetscCall(PetscLogStageRegister("Solve Stage", &solve_stage));
4652b730f8bSJeremy L Thompson PetscCall(PetscLogStagePush(solve_stage));
46609a940d7Sjeremylt
46709a940d7Sjeremylt // -- Solve
4686c5df90dSjeremylt my_rt_start = MPI_Wtime();
4692b730f8bSJeremy L Thompson PetscCall(KSPSolve(ksp, rhs, X[fine_level]));
4706c5df90dSjeremylt my_rt = MPI_Wtime() - my_rt_start;
4716c5df90dSjeremylt
47209a940d7Sjeremylt // -- Performance logging
4732b730f8bSJeremy L Thompson PetscCall(PetscLogStagePop());
47409a940d7Sjeremylt
4756c5df90dSjeremylt // Output results
4766c5df90dSjeremylt {
4779b072555Sjeremylt KSPType ksp_type;
4789b072555Sjeremylt PCMGType pcmg_type;
4796c5df90dSjeremylt KSPConvergedReason reason;
4806c5df90dSjeremylt PetscReal rnorm;
4816c5df90dSjeremylt PetscInt its;
4822b730f8bSJeremy L Thompson PetscCall(KSPGetType(ksp, &ksp_type));
4832b730f8bSJeremy L Thompson PetscCall(KSPGetConvergedReason(ksp, &reason));
4842b730f8bSJeremy L Thompson PetscCall(KSPGetIterationNumber(ksp, &its));
4852b730f8bSJeremy L Thompson PetscCall(KSPGetResidualNorm(ksp, &rnorm));
4862b730f8bSJeremy L Thompson PetscCall(PCMGGetType(pc, &pcmg_type));
4876c5df90dSjeremylt if (!test_mode || reason < 0 || rnorm > 1e-8) {
4882b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm,
4896c5df90dSjeremylt " KSP:\n"
4906c5df90dSjeremylt " KSP Type : %s\n"
4916c5df90dSjeremylt " KSP Convergence : %s\n"
492a9b2c5ddSrezgarshakeri " Total KSP Iterations : %" PetscInt_FMT "\n"
4936c5df90dSjeremylt " Final rnorm : %e\n",
4942b730f8bSJeremy L Thompson ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
4952b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm,
4966c5df90dSjeremylt " PCMG:\n"
4976c5df90dSjeremylt " PCMG Type : %s\n"
4986c5df90dSjeremylt " PCMG Cycle Type : %s\n",
4992b730f8bSJeremy L Thompson PCMGTypes[pcmg_type], PCMGCycleTypes[pcmg_cycle_type]));
5006c5df90dSjeremylt }
5016c5df90dSjeremylt if (!test_mode) {
5022b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, " Performance:\n"));
5036c5df90dSjeremylt }
5046c5df90dSjeremylt {
5056c88e6a2Srezgarshakeri // Set up error operator context
5064dbe2ad5SJeremy L Thompson PetscCall(PetscMalloc1(1, &op_error_ctx));
5072b730f8bSJeremy L Thompson PetscCall(SetupErrorOperatorCtx(comm, dm[fine_level], ceed, ceed_data[fine_level], X_loc[fine_level], op_error, op_error_ctx));
50838f32c05Srezgarshakeri PetscScalar l2_error;
5092b730f8bSJeremy L Thompson PetscCall(ComputeL2Error(X[fine_level], &l2_error, op_error_ctx));
5100a8fc04aSrezgarshakeri PetscReal tol = 5e-2;
51138f32c05Srezgarshakeri if (!test_mode || l2_error > tol) {
5122b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm));
5132b730f8bSJeremy L Thompson PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm));
5142b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm,
51538f32c05Srezgarshakeri " L2 Error : %e\n"
5166c5df90dSjeremylt " CG Solve Time : %g (%g) sec\n",
5172b730f8bSJeremy L Thompson (double)l2_error, rt_max, rt_min));
5186c5df90dSjeremylt }
5196c5df90dSjeremylt }
5206c5df90dSjeremylt if (benchmark_mode && (!test_mode)) {
5212b730f8bSJeremy L Thompson PetscCall(PetscPrintf(comm, " DoFs/Sec in CG : %g (%g) million\n", 1e-6 * g_size[fine_level] * its / rt_max,
5222b730f8bSJeremy L Thompson 1e-6 * g_size[fine_level] * its / rt_min));
5236c5df90dSjeremylt }
5246c5df90dSjeremylt }
5256c5df90dSjeremylt
5266c5df90dSjeremylt if (write_solution) {
5279b072555Sjeremylt PetscViewer vtk_viewer_soln;
5286c5df90dSjeremylt
5292b730f8bSJeremy L Thompson PetscCall(PetscViewerCreate(comm, &vtk_viewer_soln));
5302b730f8bSJeremy L Thompson PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK));
5312b730f8bSJeremy L Thompson PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"));
5322b730f8bSJeremy L Thompson PetscCall(VecView(X[fine_level], vtk_viewer_soln));
5332b730f8bSJeremy L Thompson PetscCall(PetscViewerDestroy(&vtk_viewer_soln));
5346c5df90dSjeremylt }
5356c5df90dSjeremylt
5366c5df90dSjeremylt // Cleanup
5374dbe2ad5SJeremy L Thompson for (PetscInt i = 0; i < num_levels; i++) {
5382b730f8bSJeremy L Thompson PetscCall(VecDestroy(&X[i]));
5392b730f8bSJeremy L Thompson PetscCall(VecDestroy(&X_loc[i]));
5402b730f8bSJeremy L Thompson PetscCall(VecDestroy(&mult[i]));
5412b730f8bSJeremy L Thompson PetscCall(VecDestroy(&op_apply_ctx[i]->Y_loc));
5422b730f8bSJeremy L Thompson PetscCall(MatDestroy(&mat_O[i]));
5432b730f8bSJeremy L Thompson PetscCall(PetscFree(op_apply_ctx[i]));
5446c5df90dSjeremylt if (i > 0) {
5452b730f8bSJeremy L Thompson PetscCall(MatDestroy(&mat_pr[i]));
5462b730f8bSJeremy L Thompson PetscCall(PetscFree(pr_restr_ctx[i]));
5476c5df90dSjeremylt }
5482b730f8bSJeremy L Thompson PetscCall(CeedDataDestroy(i, ceed_data[i]));
5492b730f8bSJeremy L Thompson PetscCall(DMDestroy(&dm[i]));
5506c5df90dSjeremylt }
5512b730f8bSJeremy L Thompson PetscCall(PetscFree(level_degrees));
5522b730f8bSJeremy L Thompson PetscCall(PetscFree(dm));
5532b730f8bSJeremy L Thompson PetscCall(PetscFree(X));
5542b730f8bSJeremy L Thompson PetscCall(PetscFree(X_loc));
5552b730f8bSJeremy L Thompson PetscCall(VecDestroy(&op_error_ctx->Y_loc));
5562b730f8bSJeremy L Thompson PetscCall(PetscFree(mult));
5572b730f8bSJeremy L Thompson PetscCall(PetscFree(mat_O));
5582b730f8bSJeremy L Thompson PetscCall(PetscFree(mat_pr));
5592b730f8bSJeremy L Thompson PetscCall(PetscFree(ceed_data));
5602b730f8bSJeremy L Thompson PetscCall(PetscFree(op_apply_ctx));
5612b730f8bSJeremy L Thompson PetscCall(PetscFree(op_error_ctx));
5622b730f8bSJeremy L Thompson PetscCall(PetscFree(pr_restr_ctx));
5632b730f8bSJeremy L Thompson PetscCall(PetscFree(l_size));
5642b730f8bSJeremy L Thompson PetscCall(PetscFree(xl_size));
5652b730f8bSJeremy L Thompson PetscCall(PetscFree(g_size));
5662b730f8bSJeremy L Thompson PetscCall(VecDestroy(&rhs));
5672b730f8bSJeremy L Thompson PetscCall(VecDestroy(&rhs_loc));
5682b730f8bSJeremy L Thompson PetscCall(MatDestroy(&mat_coarse));
5692b730f8bSJeremy L Thompson PetscCall(KSPDestroy(&ksp));
5702b730f8bSJeremy L Thompson PetscCall(DMDestroy(&dm_orig));
5716c5df90dSjeremylt CeedVectorDestroy(&target);
5729b072555Sjeremylt CeedQFunctionDestroy(&qf_error);
5739b072555Sjeremylt CeedOperatorDestroy(&op_error);
5746c5df90dSjeremylt CeedDestroy(&ceed);
5756c5df90dSjeremylt return PetscFinalize();
5766c5df90dSjeremylt }
577