xref: /libCEED/examples/petsc/bps.c (revision b6aaa75cb37e8d8b1128eb26ab4f9329fab4b3ba)
1 // Copyright (c) 2017-2022, Lawrence Livermore National Security, LLC and other CEED contributors.
2 // All Rights Reserved. See the top-level LICENSE and NOTICE files for details.
3 //
4 // SPDX-License-Identifier: BSD-2-Clause
5 //
6 // This file is part of CEED:  http://github.com/ceed
7 
8 //                        libCEED + PETSc Example: CEED BPs
9 //
10 // This example demonstrates a simple usage of libCEED with PETSc to solve the CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
11 //
12 // The code uses higher level communication protocols in DMPlex.
13 //
14 // Build with:
15 //
16 //     make bps [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
17 //
18 // Sample runs:
19 //
20 //     ./bps -problem bp1 -degree 3
21 //     ./bps -problem bp2 -degree 3
22 //     ./bps -problem bp3 -degree 3
23 //     ./bps -problem bp4 -degree 3
24 //     ./bps -problem bp5 -degree 3 -ceed /cpu/self
25 //     ./bps -problem bp6 -degree 3 -ceed /gpu/cuda
26 //
27 //TESTARGS(name="BP3, tet elements") -ceed {ceed_resource} -test -problem bp3 -degree 3 -ksp_max_it_clip 50,50 -simplex
28 //TESTARGS(name="BP5, hex elements") -ceed {ceed_resource} -test -problem bp5 -degree 3 -ksp_max_it_clip 15,15
29 
30 /// @file
31 /// CEED BPs example using PETSc with DMPlex
32 /// See bpsraw.c for a "raw" implementation using a structured grid.
33 const char help[] = "Solve CEED BPs using PETSc with DMPlex\n";
34 
35 #include "bps.h"
36 
37 #include <ceed.h>
38 #include <petscdmplex.h>
39 #include <petscksp.h>
40 #include <stdbool.h>
41 #include <string.h>
42 
43 #include "include/bpsproblemdata.h"
44 #include "include/libceedsetup.h"
45 #include "include/matops.h"
46 #include "include/petscutils.h"
47 #include "include/petscversion.h"
48 #include "include/structs.h"
49 
50 // -----------------------------------------------------------------------------
51 // Main body of program, called in a loop for performance benchmarking purposes
52 // -----------------------------------------------------------------------------
53 static PetscErrorCode RunWithDM(RunParams rp, DM dm, const char *ceed_resource) {
54   double               my_rt_start, my_rt, rt_min, rt_max;
55   PetscInt             xl_size, l_size, g_size;
56   Vec                  X, X_loc, rhs, rhs_loc;
57   Mat                  mat_O;
58   KSP                  ksp;
59   OperatorApplyContext op_apply_ctx, op_error_ctx;
60   Ceed                 ceed;
61   CeedData             ceed_data;
62   CeedQFunction        qf_error;
63   CeedOperator         op_error;
64   CeedVector           rhs_ceed, target;
65   VecType              vec_type;
66   PetscMemType         mem_type;
67 
68   PetscFunctionBeginUser;
69   // Set up libCEED
70   CeedInit(ceed_resource, &ceed);
71   CeedMemType mem_type_backend;
72   CeedGetPreferredMemType(ceed, &mem_type_backend);
73 
74   PetscCall(DMGetVecType(dm, &vec_type));
75   if (!vec_type) {  // Not yet set by user -dm_vec_type
76     switch (mem_type_backend) {
77       case CEED_MEM_HOST:
78         vec_type = VECSTANDARD;
79         break;
80       case CEED_MEM_DEVICE: {
81         const char *resolved;
82         CeedGetResource(ceed, &resolved);
83         if (strstr(resolved, "/gpu/cuda")) vec_type = VECCUDA;
84         else if (strstr(resolved, "/gpu/hip/occa")) vec_type = VECSTANDARD;  // https://github.com/CEED/libCEED/issues/678
85         else if (strstr(resolved, "/gpu/hip")) vec_type = VECHIP;
86         else vec_type = VECSTANDARD;
87       }
88     }
89     PetscCall(DMSetVecType(dm, vec_type));
90   }
91 
92   // Create global and local solution vectors
93   PetscCall(DMCreateGlobalVector(dm, &X));
94   PetscCall(VecGetLocalSize(X, &l_size));
95   PetscCall(VecGetSize(X, &g_size));
96   PetscCall(DMCreateLocalVector(dm, &X_loc));
97   PetscCall(VecGetSize(X_loc, &xl_size));
98   PetscCall(VecDuplicate(X, &rhs));
99 
100   // Operator
101   PetscCall(PetscMalloc1(1, &op_apply_ctx));
102   PetscCall(PetscMalloc1(1, &op_error_ctx));
103   PetscCall(MatCreateShell(rp->comm, l_size, l_size, g_size, g_size, op_apply_ctx, &mat_O));
104   PetscCall(MatShellSetOperation(mat_O, MATOP_MULT, (void (*)(void))MatMult_Ceed));
105   PetscCall(MatShellSetOperation(mat_O, MATOP_GET_DIAGONAL, (void (*)(void))MatGetDiag));
106   PetscCall(MatShellSetVecType(mat_O, vec_type));
107 
108   // Print summary
109   if (!rp->test_mode) {
110     PetscInt P = rp->degree + 1, Q = P + rp->q_extra;
111 
112     const char *used_resource;
113     CeedGetResource(ceed, &used_resource);
114 
115     VecType vec_type;
116     PetscCall(VecGetType(X, &vec_type));
117 
118     PetscInt c_start, c_end;
119     PetscCall(DMPlexGetHeightStratum(dm, 0, &c_start, &c_end));
120     DMPolytopeType cell_type;
121     PetscCall(DMPlexGetCellType(dm, c_start, &cell_type));
122     CeedElemTopology elem_topo = ElemTopologyP2C(cell_type);
123     PetscMPIInt      comm_size;
124     PetscCall(MPI_Comm_size(rp->comm, &comm_size));
125     PetscCall(PetscPrintf(rp->comm,
126                           "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc --\n"
127                           "  MPI:\n"
128                           "    Hostname                                : %s\n"
129                           "    Total ranks                             : %d\n"
130                           "    Ranks per compute node                  : %d\n"
131                           "  PETSc:\n"
132                           "    PETSc Vec Type                          : %s\n"
133                           "  libCEED:\n"
134                           "    libCEED Backend                         : %s\n"
135                           "    libCEED Backend MemType                 : %s\n"
136                           "  Mesh:\n"
137                           "    Solution Order (P)                      : %" CeedInt_FMT "\n"
138                           "    Quadrature  Order (Q)                   : %" CeedInt_FMT "\n"
139                           "    Additional quadrature points (q_extra)  : %" CeedInt_FMT "\n"
140                           "    Global nodes                            : %" PetscInt_FMT "\n"
141                           "    Local Elements                          : %" PetscInt_FMT "\n"
142                           "    Element topology                        : %s\n"
143                           "    Owned nodes                             : %" PetscInt_FMT "\n"
144                           "    DoF per node                            : %" PetscInt_FMT "\n",
145                           rp->bp_choice + 1, rp->hostname, comm_size, rp->ranks_per_node, vec_type, used_resource, CeedMemTypes[mem_type_backend], P,
146                           Q, rp->q_extra, g_size / rp->num_comp_u, c_end - c_start, CeedElemTopologies[elem_topo], l_size / rp->num_comp_u,
147                           rp->num_comp_u));
148   }
149 
150   // Create RHS vector
151   PetscCall(VecDuplicate(X_loc, &rhs_loc));
152   PetscCall(VecZeroEntries(rhs_loc));
153   CeedVectorCreate(ceed, xl_size, &rhs_ceed);
154   PetscCall(VecP2C(rhs_loc, &mem_type, rhs_ceed));
155 
156   PetscCall(PetscMalloc1(1, &ceed_data));
157   PetscCall(SetupLibceedByDegree(dm, ceed, rp->degree, rp->dim, rp->q_extra, rp->dim, rp->num_comp_u, g_size, xl_size, bp_options[rp->bp_choice],
158                                  ceed_data, true, rhs_ceed, &target));
159 
160   // Gather RHS
161   PetscCall(VecC2P(rhs_ceed, mem_type, rhs_loc));
162   PetscCall(VecZeroEntries(rhs));
163   PetscCall(DMLocalToGlobal(dm, rhs_loc, ADD_VALUES, rhs));
164   CeedVectorDestroy(&rhs_ceed);
165 
166   // Create the error QFunction
167   CeedQFunctionCreateInterior(ceed, 1, bp_options[rp->bp_choice].error, bp_options[rp->bp_choice].error_loc, &qf_error);
168   CeedQFunctionAddInput(qf_error, "u", rp->num_comp_u, CEED_EVAL_INTERP);
169   CeedQFunctionAddInput(qf_error, "true_soln", rp->num_comp_u, CEED_EVAL_NONE);
170   CeedQFunctionAddInput(qf_error, "qdata", ceed_data->q_data_size, CEED_EVAL_NONE);
171   CeedQFunctionAddOutput(qf_error, "error", rp->num_comp_u, CEED_EVAL_INTERP);
172 
173   // Create the error operator
174   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, &op_error);
175   CeedOperatorSetField(op_error, "u", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE);
176   CeedOperatorSetField(op_error, "true_soln", ceed_data->elem_restr_u_i, CEED_BASIS_NONE, target);
177   CeedOperatorSetField(op_error, "qdata", ceed_data->elem_restr_qd_i, CEED_BASIS_NONE, ceed_data->q_data);
178   CeedOperatorSetField(op_error, "error", ceed_data->elem_restr_u, ceed_data->basis_u, CEED_VECTOR_ACTIVE);
179 
180   // Set up apply operator context
181   PetscCall(SetupApplyOperatorCtx(rp->comm, dm, ceed, ceed_data, X_loc, op_apply_ctx));
182   PetscCall(KSPCreate(rp->comm, &ksp));
183   {
184     PC pc;
185     PetscCall(KSPGetPC(ksp, &pc));
186     if (rp->bp_choice == CEED_BP1 || rp->bp_choice == CEED_BP2) {
187       PetscCall(PCSetType(pc, PCJACOBI));
188       if (rp->simplex) {
189         PetscCall(PCJacobiSetType(pc, PC_JACOBI_DIAGONAL));
190       } else {
191         PetscCall(PCJacobiSetType(pc, PC_JACOBI_ROWSUM));
192       }
193     } else {
194       PetscCall(PCSetType(pc, PCNONE));
195     }
196     PetscCall(KSPSetType(ksp, KSPCG));
197     PetscCall(KSPSetNormType(ksp, KSP_NORM_NATURAL));
198     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT));
199   }
200   PetscCall(KSPSetOperators(ksp, mat_O, mat_O));
201 
202   // First run's performance log is not considered for benchmarking purposes
203   PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1));
204   my_rt_start = MPI_Wtime();
205   PetscCall(KSPSolve(ksp, rhs, X));
206   my_rt = MPI_Wtime() - my_rt_start;
207   PetscCall(MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, rp->comm));
208   // Set maxits based on first iteration timing
209   if (my_rt > 0.02) {
210     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, rp->ksp_max_it_clip[0]));
211   } else {
212     PetscCall(KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, rp->ksp_max_it_clip[1]));
213   }
214   PetscCall(KSPSetFromOptions(ksp));
215 
216   // Timed solve
217   PetscCall(VecZeroEntries(X));
218   PetscCall(PetscBarrier((PetscObject)ksp));
219 
220   // -- Performance logging
221   PetscCall(PetscLogStagePush(rp->solve_stage));
222 
223   // -- Solve
224   my_rt_start = MPI_Wtime();
225   PetscCall(KSPSolve(ksp, rhs, X));
226   my_rt = MPI_Wtime() - my_rt_start;
227 
228   // -- Performance logging
229   PetscCall(PetscLogStagePop());
230 
231   // Output results
232   {
233     KSPType            ksp_type;
234     KSPConvergedReason reason;
235     PetscReal          rnorm;
236     PetscInt           its;
237     PetscCall(KSPGetType(ksp, &ksp_type));
238     PetscCall(KSPGetConvergedReason(ksp, &reason));
239     PetscCall(KSPGetIterationNumber(ksp, &its));
240     PetscCall(KSPGetResidualNorm(ksp, &rnorm));
241     if (!rp->test_mode || reason < 0 || rnorm > 1e-8) {
242       PetscCall(PetscPrintf(rp->comm,
243                             "  KSP:\n"
244                             "    KSP Type                                : %s\n"
245                             "    KSP Convergence                         : %s\n"
246                             "    Total KSP Iterations                    : %" PetscInt_FMT "\n"
247                             "    Final rnorm                             : %e\n",
248                             ksp_type, KSPConvergedReasons[reason], its, (double)rnorm));
249     }
250     if (!rp->test_mode) {
251       PetscCall(PetscPrintf(rp->comm, "  Performance:\n"));
252     }
253     {
254       // Set up error operator context
255       PetscCall(SetupErrorOperatorCtx(rp->comm, dm, ceed, ceed_data, X_loc, op_error, op_error_ctx));
256       PetscScalar l2_error;
257       PetscCall(ComputeL2Error(X, &l2_error, op_error_ctx));
258       PetscReal tol = 5e-2;
259       if (!rp->test_mode || l2_error > tol) {
260         PetscCall(MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, rp->comm));
261         PetscCall(MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, rp->comm));
262         PetscCall(PetscPrintf(rp->comm,
263                               "    L2 Error                                : %e\n"
264                               "    CG Solve Time                           : %g (%g) sec\n",
265                               (double)l2_error, rt_max, rt_min));
266       }
267     }
268     if (!rp->test_mode) {
269       PetscCall(PetscPrintf(rp->comm, "    DoFs/Sec in CG                          : %g (%g) million\n", 1e-6 * g_size * its / rt_max,
270                             1e-6 * g_size * its / rt_min));
271     }
272   }
273 
274   if (rp->write_solution) {
275     PetscViewer vtk_viewer_soln;
276 
277     PetscCall(PetscViewerCreate(rp->comm, &vtk_viewer_soln));
278     PetscCall(PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK));
279     PetscCall(PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"));
280     PetscCall(VecView(X, vtk_viewer_soln));
281     PetscCall(PetscViewerDestroy(&vtk_viewer_soln));
282   }
283 
284   // Cleanup
285   PetscCall(VecDestroy(&X));
286   PetscCall(VecDestroy(&X_loc));
287   PetscCall(VecDestroy(&op_apply_ctx->Y_loc));
288   PetscCall(VecDestroy(&op_error_ctx->Y_loc));
289   PetscCall(MatDestroy(&mat_O));
290   PetscCall(PetscFree(op_apply_ctx));
291   PetscCall(PetscFree(op_error_ctx));
292   PetscCall(CeedDataDestroy(0, ceed_data));
293 
294   PetscCall(VecDestroy(&rhs));
295   PetscCall(VecDestroy(&rhs_loc));
296   PetscCall(KSPDestroy(&ksp));
297   CeedVectorDestroy(&target);
298   CeedQFunctionDestroy(&qf_error);
299   CeedOperatorDestroy(&op_error);
300   CeedDestroy(&ceed);
301   PetscFunctionReturn(PETSC_SUCCESS);
302 }
303 
304 static PetscErrorCode Run(RunParams rp, PetscInt num_resources, char *const *ceed_resources, PetscInt num_bp_choices, const BPType *bp_choices) {
305   DM dm;
306 
307   PetscFunctionBeginUser;
308   // Setup DM
309   PetscCall(CreateDistributedDM(rp, &dm));
310 
311   for (PetscInt b = 0; b < num_bp_choices; b++) {
312     DM       dm_deg;
313     VecType  vec_type;
314     PetscInt q_extra = rp->q_extra;
315     rp->bp_choice    = bp_choices[b];
316     rp->num_comp_u   = bp_options[rp->bp_choice].num_comp_u;
317     rp->q_extra      = q_extra < 0 ? bp_options[rp->bp_choice].q_extra : q_extra;
318     PetscCall(DMClone(dm, &dm_deg));
319     PetscCall(DMGetVecType(dm, &vec_type));
320     PetscCall(DMSetVecType(dm_deg, vec_type));
321     // Create DM
322     PetscInt dim;
323     PetscCall(DMGetDimension(dm_deg, &dim));
324     PetscCall(SetupDMByDegree(dm_deg, rp->degree, rp->q_extra, rp->num_comp_u, dim, bp_options[rp->bp_choice].enforce_bc));
325     for (PetscInt r = 0; r < num_resources; r++) {
326       PetscCall(RunWithDM(rp, dm_deg, ceed_resources[r]));
327     }
328     PetscCall(DMDestroy(&dm_deg));
329     rp->q_extra = q_extra;
330   }
331 
332   PetscCall(DMDestroy(&dm));
333   PetscFunctionReturn(PETSC_SUCCESS);
334 }
335 
336 int main(int argc, char **argv) {
337   PetscInt  comm_size;
338   RunParams rp;
339   MPI_Comm  comm;
340   char      filename[PETSC_MAX_PATH_LEN];
341   char     *ceed_resources[30];
342   PetscInt  num_ceed_resources = 30;
343   char      hostname[PETSC_MAX_PATH_LEN];
344 
345   PetscInt    dim = 3, mesh_elem[3] = {3, 3, 3};
346   PetscInt    num_degrees = 30, degree[30] = {0}, num_local_nodes = 2, local_nodes[2] = {0};
347   PetscMPIInt ranks_per_node;
348   PetscBool   degree_set;
349   BPType      bp_choices[10];
350   PetscInt    num_bp_choices = 10;
351 
352   // Initialize PETSc
353   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
354   comm = PETSC_COMM_WORLD;
355   PetscCall(MPI_Comm_size(comm, &comm_size));
356 #if defined(PETSC_HAVE_MPI_PROCESS_SHARED_MEMORY)
357   {
358     MPI_Comm splitcomm;
359     PetscCall(MPI_Comm_split_type(comm, MPI_COMM_TYPE_SHARED, 0, MPI_INFO_NULL, &splitcomm));
360     PetscCall(MPI_Comm_size(splitcomm, &ranks_per_node));
361     PetscCall(MPI_Comm_free(&splitcomm));
362   }
363 #else
364   ranks_per_node = -1;  // Unknown
365 #endif
366 
367   // Setup all parameters needed in Run()
368   PetscCall(PetscMalloc1(1, &rp));
369   rp->comm = comm;
370 
371   // Read command line options
372   PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
373   {
374     PetscBool set;
375     PetscCall(PetscOptionsEnumArray("-problem", "CEED benchmark problem to solve", NULL, bp_types, (PetscEnum *)bp_choices, &num_bp_choices, &set));
376     if (!set) {
377       bp_choices[0]  = CEED_BP1;
378       num_bp_choices = 1;
379     }
380   }
381   rp->test_mode = PETSC_FALSE;
382   PetscCall(PetscOptionsBool("-test", "Testing mode (do not print unless error is large)", NULL, rp->test_mode, &rp->test_mode, NULL));
383   rp->write_solution = PETSC_FALSE;
384   PetscCall(PetscOptionsBool("-write_solution", "Write solution for visualization", NULL, rp->write_solution, &rp->write_solution, NULL));
385   rp->simplex = PETSC_FALSE;
386   PetscCall(PetscOptionsBool("-simplex", "Element topology (default:hex)", NULL, rp->simplex, &rp->simplex, NULL));
387   if ((bp_choices[0] == CEED_BP5 || bp_choices[0] == CEED_BP6) && (rp->simplex)) {
388     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "BP5/6 is not supported with simplex");
389   }
390   degree[0] = rp->test_mode ? 3 : 2;
391   PetscCall(PetscOptionsIntArray("-degree", "Polynomial degree of tensor product basis", NULL, degree, &num_degrees, &degree_set));
392   if (!degree_set) num_degrees = 1;
393   rp->q_extra = PETSC_DECIDE;
394   PetscCall(PetscOptionsInt("-q_extra", "Number of extra quadrature points (-1 for auto)", NULL, rp->q_extra, &rp->q_extra, NULL));
395   {
396     PetscBool set;
397     PetscCall(PetscOptionsStringArray("-ceed", "CEED resource specifier (comma-separated list)", NULL, ceed_resources, &num_ceed_resources, &set));
398     if (!set) {
399       PetscCall(PetscStrallocpy("/cpu/self", &ceed_resources[0]));
400       num_ceed_resources = 1;
401     }
402   }
403   PetscCall(PetscGetHostName(hostname, sizeof hostname));
404   PetscCall(PetscOptionsString("-hostname", "Hostname for output", NULL, hostname, hostname, sizeof(hostname), NULL));
405   rp->read_mesh = PETSC_FALSE;
406   PetscCall(PetscOptionsString("-mesh", "Read mesh from file", NULL, filename, filename, sizeof(filename), &rp->read_mesh));
407   rp->filename = filename;
408   if (!rp->read_mesh) {
409     PetscInt tmp = dim;
410     PetscCall(PetscOptionsIntArray("-cells", "Number of cells per dimension", NULL, mesh_elem, &tmp, NULL));
411   }
412   local_nodes[0] = 1000;
413   PetscCall(PetscOptionsIntArray("-local_nodes",
414                                  "Target number of locally owned nodes per "
415                                  "process (single value or min,max)",
416                                  NULL, local_nodes, &num_local_nodes, &rp->user_l_nodes));
417   if (num_local_nodes < 2) local_nodes[1] = 2 * local_nodes[0];
418   {
419     PetscInt two           = 2;
420     rp->ksp_max_it_clip[0] = 5;
421     rp->ksp_max_it_clip[1] = 20;
422     PetscCall(PetscOptionsIntArray("-ksp_max_it_clip", "Min and max number of iterations to use during benchmarking", NULL, rp->ksp_max_it_clip, &two,
423                                    NULL));
424   }
425   if (!degree_set) {
426     PetscInt max_degree = 8;
427     PetscCall(PetscOptionsInt("-max_degree", "Range of degrees [1, max_degree] to run with", NULL, max_degree, &max_degree, NULL));
428     for (PetscInt i = 0; i < max_degree; i++) degree[i] = i + 1;
429     num_degrees = max_degree;
430   }
431   {
432     PetscBool flg;
433     PetscInt  p = ranks_per_node;
434     PetscCall(PetscOptionsInt("-p", "Number of MPI ranks per node", NULL, p, &p, &flg));
435     if (flg) ranks_per_node = p;
436   }
437 
438   PetscOptionsEnd();
439 
440   // Register PETSc logging stage
441   PetscCall(PetscLogStageRegister("Solve Stage", &rp->solve_stage));
442 
443   rp->hostname       = hostname;
444   rp->dim            = dim;
445   rp->mesh_elem      = mesh_elem;
446   rp->ranks_per_node = ranks_per_node;
447 
448   for (PetscInt d = 0; d < num_degrees; d++) {
449     PetscInt deg = degree[d];
450     for (PetscInt n = local_nodes[0]; n < local_nodes[1]; n *= 2) {
451       rp->degree      = deg;
452       rp->local_nodes = n;
453       PetscCall(Run(rp, num_ceed_resources, ceed_resources, num_bp_choices, bp_choices));
454     }
455   }
456   // Clear memory
457   PetscCall(PetscFree(rp));
458   for (PetscInt i = 0; i < num_ceed_resources; i++) {
459     PetscCall(PetscFree(ceed_resources[i]));
460   }
461   return PetscFinalize();
462 }
463