xref: /libCEED/examples/petsc/bpsraw.c (revision 751eb813a262a2147c23ac2e55685e7f0d6b8af0)
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.
3cb32e2e7SValeria Barra //
43d8e8822SJeremy L Thompson // SPDX-License-Identifier: BSD-2-Clause
5cb32e2e7SValeria Barra //
63d8e8822SJeremy L Thompson // This file is part of CEED:  http://github.com/ceed
7cb32e2e7SValeria Barra 
8cb32e2e7SValeria Barra //                        libCEED + PETSc Example: CEED BPs
9cb32e2e7SValeria Barra //
10cb32e2e7SValeria Barra // This example demonstrates a simple usage of libCEED with PETSc to solve the
11cb32e2e7SValeria Barra // CEED BP benchmark problems, see http://ceed.exascaleproject.org/bps.
12cb32e2e7SValeria Barra //
13cb32e2e7SValeria Barra // The code is intentionally "raw", using only low-level communication
14cb32e2e7SValeria Barra // primitives.
15cb32e2e7SValeria Barra //
16cb32e2e7SValeria Barra // Build with:
17cb32e2e7SValeria Barra //
18cb32e2e7SValeria Barra //     make bpsraw [PETSC_DIR=</path/to/petsc>] [CEED_DIR=</path/to/libceed>]
19cb32e2e7SValeria Barra //
20cb32e2e7SValeria Barra // Sample runs:
21cb32e2e7SValeria Barra //
2232d2ee49SValeria Barra //     ./bpsraw -problem bp1
2328688798Sjeremylt //     ./bpsraw -problem bp2
2428688798Sjeremylt //     ./bpsraw -problem bp3
2528688798Sjeremylt //     ./bpsraw -problem bp4
2628688798Sjeremylt //     ./bpsraw -problem bp5 -ceed /cpu/self
2728688798Sjeremylt //     ./bpsraw -problem bp6 -ceed /gpu/cuda
28cb32e2e7SValeria Barra //
299b072555Sjeremylt //TESTARGS -ceed {ceed_resource} -test -problem bp2 -degree 5 -q_extra 1 -ksp_max_it_clip 15,15
30cb32e2e7SValeria Barra 
31cb32e2e7SValeria Barra /// @file
32cb32e2e7SValeria Barra /// CEED BPs example using PETSc
33cb32e2e7SValeria Barra /// See bps.c for an implementation using DMPlex unstructured grids.
34cb32e2e7SValeria Barra const char help[] = "Solve CEED BPs using PETSc\n";
35cb32e2e7SValeria Barra 
363d576824SJeremy L Thompson #include <ceed.h>
373d576824SJeremy L Thompson #include <petscksp.h>
383d576824SJeremy L Thompson #include <petscsys.h>
39cb32e2e7SValeria Barra #include <stdbool.h>
40cb32e2e7SValeria Barra #include <string.h>
41cb32e2e7SValeria Barra #include "qfunctions/bps/bp1.h"
42cb32e2e7SValeria Barra #include "qfunctions/bps/bp2.h"
43cb32e2e7SValeria Barra #include "qfunctions/bps/bp3.h"
44cb32e2e7SValeria Barra #include "qfunctions/bps/bp4.h"
453d576824SJeremy L Thompson #include "qfunctions/bps/common.h"
46cb32e2e7SValeria Barra 
479396343dSjeremylt #if PETSC_VERSION_LT(3,12,0)
489396343dSjeremylt #ifdef PETSC_HAVE_CUDA
499396343dSjeremylt #include <petsccuda.h>
509396343dSjeremylt // Note: With PETSc prior to version 3.12.0, providing the source path to
519396343dSjeremylt //       include 'cublas_v2.h' will be needed to use 'petsccuda.h'.
529396343dSjeremylt #endif
539396343dSjeremylt #endif
549396343dSjeremylt 
559b072555Sjeremylt static CeedMemType MemTypeP2C(PetscMemType mem_type) {
569b072555Sjeremylt   return PetscMemTypeDevice(mem_type) ? CEED_MEM_DEVICE : CEED_MEM_HOST;
57b68a8d79SJed Brown }
58b68a8d79SJed Brown 
59cb32e2e7SValeria Barra static void Split3(PetscInt size, PetscInt m[3], bool reverse) {
609b072555Sjeremylt   for (PetscInt d=0, size_left=size; d<3; d++) {
619b072555Sjeremylt     PetscInt try = (PetscInt)PetscCeilReal(PetscPowReal(size_left, 1./(3 - d)));
629b072555Sjeremylt     while (try * (size_left / try) != size_left) try++;
63cb32e2e7SValeria Barra     m[reverse ? 2-d : d] = try;
649b072555Sjeremylt     size_left /= try;
65cb32e2e7SValeria Barra   }
66cb32e2e7SValeria Barra }
67cb32e2e7SValeria Barra 
68cb32e2e7SValeria Barra static PetscInt Max3(const PetscInt a[3]) {
69cb32e2e7SValeria Barra   return PetscMax(a[0], PetscMax(a[1], a[2]));
70cb32e2e7SValeria Barra }
71cb32e2e7SValeria Barra static PetscInt Min3(const PetscInt a[3]) {
72cb32e2e7SValeria Barra   return PetscMin(a[0], PetscMin(a[1], a[2]));
73cb32e2e7SValeria Barra }
749b072555Sjeremylt static void GlobalNodes(const PetscInt p[3], const PetscInt i_rank[3],
759b072555Sjeremylt                         PetscInt degree, const PetscInt mesh_elem[3],
769b072555Sjeremylt                         PetscInt m_nodes[3]) {
77cb32e2e7SValeria Barra   for (int d=0; d<3; d++)
789b072555Sjeremylt     m_nodes[d] = degree*mesh_elem[d] + (i_rank[d] == p[d]-1);
79cb32e2e7SValeria Barra }
809b072555Sjeremylt static PetscInt GlobalStart(const PetscInt p[3], const PetscInt i_rank[3],
819b072555Sjeremylt                             PetscInt degree, const PetscInt mesh_elem[3]) {
82cb32e2e7SValeria Barra   PetscInt start = 0;
83cb32e2e7SValeria Barra   // Dumb brute-force is easier to read
84cb32e2e7SValeria Barra   for (PetscInt i=0; i<p[0]; i++) {
85cb32e2e7SValeria Barra     for (PetscInt j=0; j<p[1]; j++) {
86cb32e2e7SValeria Barra       for (PetscInt k=0; k<p[2]; k++) {
879b072555Sjeremylt         PetscInt m_nodes[3], ijk_rank[] = {i,j,k};
889b072555Sjeremylt         if (i == i_rank[0] && j == i_rank[1] && k == i_rank[2]) return start;
899b072555Sjeremylt         GlobalNodes(p, ijk_rank, degree, mesh_elem, m_nodes);
909b072555Sjeremylt         start += m_nodes[0] * m_nodes[1] * m_nodes[2];
91cb32e2e7SValeria Barra       }
92cb32e2e7SValeria Barra     }
93cb32e2e7SValeria Barra   }
94cb32e2e7SValeria Barra   return -1;
95cb32e2e7SValeria Barra }
965dfaedb8SJed Brown static PetscErrorCode CreateRestriction(Ceed ceed, const CeedInt mesh_elem[3],
975dfaedb8SJed Brown                                         CeedInt P,
989b072555Sjeremylt                                         CeedInt num_comp, CeedElemRestriction *elem_restr) {
999b072555Sjeremylt   const PetscInt num_elem = mesh_elem[0]*mesh_elem[1]*mesh_elem[2];
1009b072555Sjeremylt   PetscInt m_nodes[3], *idx, *idx_p;
101cb32e2e7SValeria Barra 
1025dfaedb8SJed Brown   PetscFunctionBeginUser;
103cb32e2e7SValeria Barra   // Get indicies
1049b072555Sjeremylt   for (int d=0; d<3; d++) m_nodes[d] = mesh_elem[d]*(P-1) + 1;
1059b072555Sjeremylt   idx_p = idx = malloc(num_elem*P*P*P*sizeof idx[0]);
1069b072555Sjeremylt   for (CeedInt i=0; i<mesh_elem[0]; i++)
1079b072555Sjeremylt     for (CeedInt j=0; j<mesh_elem[1]; j++)
1089b072555Sjeremylt       for (CeedInt k=0; k<mesh_elem[2]; k++,idx_p += P*P*P)
109cb32e2e7SValeria Barra         for (CeedInt ii=0; ii<P; ii++)
110cb32e2e7SValeria Barra           for (CeedInt jj=0; jj<P; jj++)
111cb32e2e7SValeria Barra             for (CeedInt kk=0; kk<P; kk++) {
112cb32e2e7SValeria Barra               if (0) { // This is the C-style (i,j,k) ordering that I prefer
1139b072555Sjeremylt                 idx_p[(ii*P+jj)*P+kk] = num_comp*(((i*(P-1)+ii)*m_nodes[1]
1149b072555Sjeremylt                                                    + (j*(P-1)+jj))*m_nodes[2]
115cb32e2e7SValeria Barra                                                   + (k*(P-1)+kk));
116cb32e2e7SValeria Barra               } else { // (k,j,i) ordering for consistency with MFEM example
1179b072555Sjeremylt                 idx_p[ii+P*(jj+P*kk)] = num_comp*(((i*(P-1)+ii)*m_nodes[1]
1189b072555Sjeremylt                                                    + (j*(P-1)+jj))*m_nodes[2]
119cb32e2e7SValeria Barra                                                   + (k*(P-1)+kk));
120cb32e2e7SValeria Barra               }
121cb32e2e7SValeria Barra             }
122cb32e2e7SValeria Barra 
123cb32e2e7SValeria Barra   // Setup CEED restriction
1249b072555Sjeremylt   CeedElemRestrictionCreate(ceed, num_elem, P*P*P, num_comp, 1,
1259b072555Sjeremylt                             m_nodes[0]*m_nodes[1]*m_nodes[2]*num_comp,
1269b072555Sjeremylt                             CEED_MEM_HOST, CEED_OWN_POINTER, idx, elem_restr);
127cb32e2e7SValeria Barra 
128cb32e2e7SValeria Barra   PetscFunctionReturn(0);
129cb32e2e7SValeria Barra }
130cb32e2e7SValeria Barra 
131cb32e2e7SValeria Barra // Data for PETSc
132cb32e2e7SValeria Barra typedef struct User_ *User;
133cb32e2e7SValeria Barra struct User_ {
134cb32e2e7SValeria Barra   MPI_Comm comm;
1359b072555Sjeremylt   VecScatter l_to_g;              // Scatter for all entries
1369b072555Sjeremylt   VecScatter l_to_g_0;            // Skip Dirichlet values
1379b072555Sjeremylt   VecScatter g_to_g_D;            // global-to-global; only Dirichlet values
1389b072555Sjeremylt   Vec X_loc, Y_loc;
1399b072555Sjeremylt   CeedVector x_ceed, y_ceed;
140cb32e2e7SValeria Barra   CeedOperator op;
1419b072555Sjeremylt   CeedVector q_data;
142cb32e2e7SValeria Barra   Ceed ceed;
143cb32e2e7SValeria Barra };
144cb32e2e7SValeria Barra 
145cb32e2e7SValeria Barra // BP Options
146cb32e2e7SValeria Barra typedef enum {
147cb32e2e7SValeria Barra   CEED_BP1 = 0, CEED_BP2 = 1, CEED_BP3 = 2,
148cb32e2e7SValeria Barra   CEED_BP4 = 3, CEED_BP5 = 4, CEED_BP6 = 5
1499b072555Sjeremylt } BPType;
1509b072555Sjeremylt static const char *const bp_types[] = {"bp1","bp2","bp3","bp4","bp5","bp6",
1519b072555Sjeremylt                                        "BPType","CEED_BP",0
152cb32e2e7SValeria Barra                                       };
153cb32e2e7SValeria Barra 
154cb32e2e7SValeria Barra // BP specific data
155cb32e2e7SValeria Barra typedef struct {
1569b072555Sjeremylt   CeedInt num_comp_u, q_data_size, q_extra;
1579b072555Sjeremylt   CeedQFunctionUser setup_geo, setup_rhs, apply, error;
1589b072555Sjeremylt   const char *setup_geo_loc, *setup_rhs_loc, *apply_loc, *error_loc;
1599b072555Sjeremylt   CeedEvalMode in_mode, out_mode;
1609b072555Sjeremylt   CeedQuadMode q_mode;
1619b072555Sjeremylt } BPData;
162cb32e2e7SValeria Barra 
1639b072555Sjeremylt BPData bp_options[6] = {
164cb32e2e7SValeria Barra   [CEED_BP1] = {
1659b072555Sjeremylt     .num_comp_u = 1,
1669b072555Sjeremylt     .q_data_size = 1,
1679b072555Sjeremylt     .q_extra = 1,
1689b072555Sjeremylt     .setup_geo = SetupMassGeo,
1699b072555Sjeremylt     .setup_rhs = SetupMassRhs,
170cb32e2e7SValeria Barra     .apply = Mass,
171cb32e2e7SValeria Barra     .error = Error,
1729b072555Sjeremylt     .setup_geo_loc = SetupMassGeo_loc,
1739b072555Sjeremylt     .setup_rhs_loc = SetupMassRhs_loc,
1749b072555Sjeremylt     .apply_loc = Mass_loc,
1759b072555Sjeremylt     .error_loc = Error_loc,
1769b072555Sjeremylt     .in_mode = CEED_EVAL_INTERP,
1779b072555Sjeremylt     .out_mode = CEED_EVAL_INTERP,
1789b072555Sjeremylt     .q_mode = CEED_GAUSS
179cb32e2e7SValeria Barra   },
180cb32e2e7SValeria Barra   [CEED_BP2] = {
1819b072555Sjeremylt     .num_comp_u = 3,
1829b072555Sjeremylt     .q_data_size = 1,
1839b072555Sjeremylt     .q_extra = 1,
1849b072555Sjeremylt     .setup_geo = SetupMassGeo,
1859b072555Sjeremylt     .setup_rhs = SetupMassRhs3,
186cb32e2e7SValeria Barra     .apply = Mass3,
187cb32e2e7SValeria Barra     .error = Error3,
1889b072555Sjeremylt     .setup_geo_loc = SetupMassGeo_loc,
1899b072555Sjeremylt     .setup_rhs_loc = SetupMassRhs3_loc,
1909b072555Sjeremylt     .apply_loc = Mass3_loc,
1919b072555Sjeremylt     .error_loc = Error3_loc,
1929b072555Sjeremylt     .in_mode = CEED_EVAL_INTERP,
1939b072555Sjeremylt     .out_mode = CEED_EVAL_INTERP,
1949b072555Sjeremylt     .q_mode = CEED_GAUSS
195cb32e2e7SValeria Barra   },
196cb32e2e7SValeria Barra   [CEED_BP3] = {
1979b072555Sjeremylt     .num_comp_u = 1,
1989b072555Sjeremylt     .q_data_size = 7,
1999b072555Sjeremylt     .q_extra = 1,
2009b072555Sjeremylt     .setup_geo = SetupDiffGeo,
2019b072555Sjeremylt     .setup_rhs = SetupDiffRhs,
202cb32e2e7SValeria Barra     .apply = Diff,
203cb32e2e7SValeria Barra     .error = Error,
2049b072555Sjeremylt     .setup_geo_loc = SetupDiffGeo_loc,
2059b072555Sjeremylt     .setup_rhs_loc = SetupDiffRhs_loc,
2069b072555Sjeremylt     .apply_loc = Diff_loc,
2079b072555Sjeremylt     .error_loc = Error_loc,
2089b072555Sjeremylt     .in_mode = CEED_EVAL_GRAD,
2099b072555Sjeremylt     .out_mode = CEED_EVAL_GRAD,
2109b072555Sjeremylt     .q_mode = CEED_GAUSS
211cb32e2e7SValeria Barra   },
212cb32e2e7SValeria Barra   [CEED_BP4] = {
2139b072555Sjeremylt     .num_comp_u = 3,
2149b072555Sjeremylt     .q_data_size = 7,
2159b072555Sjeremylt     .q_extra = 1,
2169b072555Sjeremylt     .setup_geo = SetupDiffGeo,
2179b072555Sjeremylt     .setup_rhs = SetupDiffRhs3,
218cb32e2e7SValeria Barra     .apply = Diff3,
219cb32e2e7SValeria Barra     .error = Error3,
2209b072555Sjeremylt     .setup_geo_loc = SetupDiffGeo_loc,
2219b072555Sjeremylt     .setup_rhs_loc = SetupDiffRhs3_loc,
2229b072555Sjeremylt     .apply_loc = Diff3_loc,
2239b072555Sjeremylt     .error_loc = Error3_loc,
2249b072555Sjeremylt     .in_mode = CEED_EVAL_GRAD,
2259b072555Sjeremylt     .out_mode = CEED_EVAL_GRAD,
2269b072555Sjeremylt     .q_mode = CEED_GAUSS
227cb32e2e7SValeria Barra   },
228cb32e2e7SValeria Barra   [CEED_BP5] = {
2299b072555Sjeremylt     .num_comp_u = 1,
2309b072555Sjeremylt     .q_data_size = 7,
2319b072555Sjeremylt     .q_extra = 0,
2329b072555Sjeremylt     .setup_geo = SetupDiffGeo,
2339b072555Sjeremylt     .setup_rhs = SetupDiffRhs,
234cb32e2e7SValeria Barra     .apply = Diff,
235cb32e2e7SValeria Barra     .error = Error,
2369b072555Sjeremylt     .setup_geo_loc = SetupDiffGeo_loc,
2379b072555Sjeremylt     .setup_rhs_loc = SetupDiffRhs_loc,
2389b072555Sjeremylt     .apply_loc = Diff_loc,
2399b072555Sjeremylt     .error_loc = Error_loc,
2409b072555Sjeremylt     .in_mode = CEED_EVAL_GRAD,
2419b072555Sjeremylt     .out_mode = CEED_EVAL_GRAD,
2429b072555Sjeremylt     .q_mode = CEED_GAUSS_LOBATTO
243cb32e2e7SValeria Barra   },
244cb32e2e7SValeria Barra   [CEED_BP6] = {
2459b072555Sjeremylt     .num_comp_u = 3,
2469b072555Sjeremylt     .q_data_size = 7,
2479b072555Sjeremylt     .q_extra = 0,
2489b072555Sjeremylt     .setup_geo = SetupDiffGeo,
2499b072555Sjeremylt     .setup_rhs = SetupDiffRhs3,
250cb32e2e7SValeria Barra     .apply = Diff3,
251cb32e2e7SValeria Barra     .error = Error3,
2529b072555Sjeremylt     .setup_geo_loc = SetupDiffGeo_loc,
2539b072555Sjeremylt     .setup_rhs_loc = SetupDiffRhs3_loc,
2549b072555Sjeremylt     .apply_loc = Diff3_loc,
2559b072555Sjeremylt     .error_loc = Error3_loc,
2569b072555Sjeremylt     .in_mode = CEED_EVAL_GRAD,
2579b072555Sjeremylt     .out_mode = CEED_EVAL_GRAD,
2589b072555Sjeremylt     .q_mode = CEED_GAUSS_LOBATTO
259cb32e2e7SValeria Barra   }
260cb32e2e7SValeria Barra };
261cb32e2e7SValeria Barra 
262cb32e2e7SValeria Barra // This function uses libCEED to compute the action of the mass matrix
263cb32e2e7SValeria Barra static PetscErrorCode MatMult_Mass(Mat A, Vec X, Vec Y) {
264cb32e2e7SValeria Barra   PetscErrorCode ierr;
265cb32e2e7SValeria Barra   User user;
266cb32e2e7SValeria Barra   PetscScalar *x, *y;
2679b072555Sjeremylt   PetscMemType x_mem_type, y_mem_type;
268cb32e2e7SValeria Barra 
269cb32e2e7SValeria Barra   PetscFunctionBeginUser;
2709396343dSjeremylt 
271cb32e2e7SValeria Barra   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
2729396343dSjeremylt 
2739396343dSjeremylt   // Global-to-local
2749b072555Sjeremylt   ierr = VecScatterBegin(user->l_to_g, X, user->X_loc, INSERT_VALUES,
275cb32e2e7SValeria Barra                          SCATTER_REVERSE); CHKERRQ(ierr);
2769b072555Sjeremylt   ierr = VecScatterEnd(user->l_to_g, X, user->X_loc, INSERT_VALUES,
2779396343dSjeremylt                        SCATTER_REVERSE); CHKERRQ(ierr);
278cb32e2e7SValeria Barra 
2799396343dSjeremylt   // Setup libCEED vectors
2809b072555Sjeremylt   ierr = VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x,
2819b072555Sjeremylt                                    &x_mem_type); CHKERRQ(ierr);
2829b072555Sjeremylt   ierr = VecGetArrayAndMemType(user->Y_loc, &y, &y_mem_type); CHKERRQ(ierr);
2839b072555Sjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
2849b072555Sjeremylt   CeedVectorSetArray(user->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y);
285cb32e2e7SValeria Barra 
2869396343dSjeremylt   // Apply libCEED operator
2879b072555Sjeremylt   CeedOperatorApply(user->op, user->x_ceed, user->y_ceed,
288cb32e2e7SValeria Barra                     CEED_REQUEST_IMMEDIATE);
289cb32e2e7SValeria Barra 
2909396343dSjeremylt   // Restore PETSc vectors
2919b072555Sjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL);
2929b072555Sjeremylt   CeedVectorTakeArray(user->y_ceed, MemTypeP2C(y_mem_type), NULL);
2939b072555Sjeremylt   ierr = VecRestoreArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x);
2949396343dSjeremylt   CHKERRQ(ierr);
2959b072555Sjeremylt   ierr = VecRestoreArrayAndMemType(user->Y_loc, &y); CHKERRQ(ierr);
296cb32e2e7SValeria Barra 
2979396343dSjeremylt   // Local-to-global
298cb32e2e7SValeria Barra   if (Y) {
299cb32e2e7SValeria Barra     ierr = VecZeroEntries(Y); CHKERRQ(ierr);
3009b072555Sjeremylt     ierr = VecScatterBegin(user->l_to_g, user->Y_loc, Y, ADD_VALUES,
3019396343dSjeremylt                            SCATTER_FORWARD); CHKERRQ(ierr);
3029b072555Sjeremylt     ierr = VecScatterEnd(user->l_to_g, user->Y_loc, Y, ADD_VALUES,
3039396343dSjeremylt                          SCATTER_FORWARD); CHKERRQ(ierr);
304cb32e2e7SValeria Barra   }
305cb32e2e7SValeria Barra   PetscFunctionReturn(0);
306cb32e2e7SValeria Barra }
307cb32e2e7SValeria Barra 
308cb32e2e7SValeria Barra // This function uses libCEED to compute the action of the Laplacian with
309cb32e2e7SValeria Barra // Dirichlet boundary conditions
310cb32e2e7SValeria Barra static PetscErrorCode MatMult_Diff(Mat A, Vec X, Vec Y) {
311cb32e2e7SValeria Barra   PetscErrorCode ierr;
312cb32e2e7SValeria Barra   User user;
313cb32e2e7SValeria Barra   PetscScalar *x, *y;
3149b072555Sjeremylt   PetscMemType x_mem_type, y_mem_type;
315cb32e2e7SValeria Barra 
316cb32e2e7SValeria Barra   PetscFunctionBeginUser;
3179396343dSjeremylt 
318cb32e2e7SValeria Barra   ierr = MatShellGetContext(A, &user); CHKERRQ(ierr);
319cb32e2e7SValeria Barra 
320cb32e2e7SValeria Barra   // Global-to-local
3219b072555Sjeremylt   ierr = VecScatterBegin(user->l_to_g_0, X, user->X_loc, INSERT_VALUES,
322cb32e2e7SValeria Barra                          SCATTER_REVERSE); CHKERRQ(ierr);
3239b072555Sjeremylt   ierr = VecScatterEnd(user->l_to_g_0, X, user->X_loc, INSERT_VALUES,
324cb32e2e7SValeria Barra                        SCATTER_REVERSE);
325cb32e2e7SValeria Barra   CHKERRQ(ierr);
326cb32e2e7SValeria Barra 
3279396343dSjeremylt   // Setup libCEED vectors
3289b072555Sjeremylt   ierr = VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x,
3299b072555Sjeremylt                                    &x_mem_type); CHKERRQ(ierr);
3309b072555Sjeremylt   ierr = VecGetArrayAndMemType(user->Y_loc, &y, &y_mem_type); CHKERRQ(ierr);
3319b072555Sjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x);
3329b072555Sjeremylt   CeedVectorSetArray(user->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y);
333cb32e2e7SValeria Barra 
3349396343dSjeremylt   // Apply libCEED operator
3359b072555Sjeremylt   CeedOperatorApply(user->op, user->x_ceed, user->y_ceed,
336cb32e2e7SValeria Barra                     CEED_REQUEST_IMMEDIATE);
337cb32e2e7SValeria Barra 
338cb32e2e7SValeria Barra   // Restore PETSc vectors
3399b072555Sjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL);
3409b072555Sjeremylt   CeedVectorTakeArray(user->y_ceed, MemTypeP2C(y_mem_type), NULL);
3419b072555Sjeremylt   ierr = VecRestoreArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x);
3429396343dSjeremylt   CHKERRQ(ierr);
3439b072555Sjeremylt   ierr = VecRestoreArrayAndMemType(user->Y_loc, &y); CHKERRQ(ierr);
344cb32e2e7SValeria Barra 
345cb32e2e7SValeria Barra   // Local-to-global
346cb32e2e7SValeria Barra   ierr = VecZeroEntries(Y); CHKERRQ(ierr);
3479b072555Sjeremylt   ierr = VecScatterBegin(user->g_to_g_D, X, Y, INSERT_VALUES, SCATTER_FORWARD);
348cb32e2e7SValeria Barra   CHKERRQ(ierr);
3499b072555Sjeremylt   ierr = VecScatterEnd(user->g_to_g_D, X, Y, INSERT_VALUES, SCATTER_FORWARD);
350cb32e2e7SValeria Barra   CHKERRQ(ierr);
3519b072555Sjeremylt   ierr = VecScatterBegin(user->l_to_g_0, user->Y_loc, Y, ADD_VALUES,
3529396343dSjeremylt                          SCATTER_FORWARD); CHKERRQ(ierr);
3539b072555Sjeremylt   ierr = VecScatterEnd(user->l_to_g_0, user->Y_loc, Y, ADD_VALUES,
3549b072555Sjeremylt                        SCATTER_FORWARD);
355cb32e2e7SValeria Barra   CHKERRQ(ierr);
356cb32e2e7SValeria Barra 
357cb32e2e7SValeria Barra   PetscFunctionReturn(0);
358cb32e2e7SValeria Barra }
359cb32e2e7SValeria Barra 
360cb32e2e7SValeria Barra // This function calculates the error in the final solution
3619b072555Sjeremylt static PetscErrorCode ComputeErrorMax(User user, CeedOperator op_error, Vec X,
3629b072555Sjeremylt                                       CeedVector target, PetscReal *max_error) {
363cb32e2e7SValeria Barra   PetscErrorCode ierr;
364cb32e2e7SValeria Barra   PetscScalar *x;
3659b072555Sjeremylt   PetscMemType mem_type;
366cb32e2e7SValeria Barra   CeedVector collocated_error;
3671f9221feSJeremy L Thompson   CeedSize length;
368cb32e2e7SValeria Barra 
369cb32e2e7SValeria Barra   PetscFunctionBeginUser;
3709396343dSjeremylt 
371cb32e2e7SValeria Barra   CeedVectorGetLength(target, &length);
372cb32e2e7SValeria Barra   CeedVectorCreate(user->ceed, length, &collocated_error);
373cb32e2e7SValeria Barra 
374cb32e2e7SValeria Barra   // Global-to-local
3759b072555Sjeremylt   ierr = VecScatterBegin(user->l_to_g, X, user->X_loc, INSERT_VALUES,
376cb32e2e7SValeria Barra                          SCATTER_REVERSE); CHKERRQ(ierr);
3779b072555Sjeremylt   ierr = VecScatterEnd(user->l_to_g, X, user->X_loc, INSERT_VALUES,
3789396343dSjeremylt                        SCATTER_REVERSE); CHKERRQ(ierr);
3799396343dSjeremylt 
3809396343dSjeremylt   // Setup libCEED vector
3819b072555Sjeremylt   ierr = VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x,
3829b072555Sjeremylt                                    &mem_type); CHKERRQ(ierr);
3839b072555Sjeremylt   CeedVectorSetArray(user->x_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, x);
384cb32e2e7SValeria Barra 
3859396343dSjeremylt   // Apply libCEED operator
3869b072555Sjeremylt   CeedOperatorApply(op_error, user->x_ceed, collocated_error,
387cb32e2e7SValeria Barra                     CEED_REQUEST_IMMEDIATE);
388cb32e2e7SValeria Barra 
389cb32e2e7SValeria Barra   // Restore PETSc vector
3909b072555Sjeremylt   CeedVectorTakeArray(user->x_ceed, MemTypeP2C(mem_type), NULL);
3919b072555Sjeremylt   ierr = VecRestoreArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x);
3929396343dSjeremylt   CHKERRQ(ierr);
393cb32e2e7SValeria Barra 
394cb32e2e7SValeria Barra   // Reduce max error
3959b072555Sjeremylt   *max_error = 0;
396cb32e2e7SValeria Barra   const CeedScalar *e;
397cb32e2e7SValeria Barra   CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e);
398cb32e2e7SValeria Barra   for (CeedInt i=0; i<length; i++) {
3999b072555Sjeremylt     *max_error = PetscMax(*max_error, PetscAbsScalar(e[i]));
400cb32e2e7SValeria Barra   }
401cb32e2e7SValeria Barra   CeedVectorRestoreArrayRead(collocated_error, &e);
4029b072555Sjeremylt   ierr = MPI_Allreduce(MPI_IN_PLACE, max_error, 1, MPIU_REAL, MPIU_MAX,
4039396343dSjeremylt                        user->comm); CHKERRQ(ierr);
404cb32e2e7SValeria Barra 
405cb32e2e7SValeria Barra   // Cleanup
406cb32e2e7SValeria Barra   CeedVectorDestroy(&collocated_error);
407cb32e2e7SValeria Barra 
408cb32e2e7SValeria Barra   PetscFunctionReturn(0);
409cb32e2e7SValeria Barra }
410cb32e2e7SValeria Barra 
411cb32e2e7SValeria Barra int main(int argc, char **argv) {
412cb32e2e7SValeria Barra   PetscInt ierr;
413cb32e2e7SValeria Barra   MPI_Comm comm;
4149b072555Sjeremylt   char ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self";
415cb32e2e7SValeria Barra   double my_rt_start, my_rt, rt_min, rt_max;
4169b072555Sjeremylt   PetscInt degree, q_extra, local_nodes, local_elem, mesh_elem[3], m_nodes[3],
4179b072555Sjeremylt            p[3],
4189b072555Sjeremylt            i_rank[3], l_nodes[3], l_size, num_comp_u = 1, ksp_max_it_clip[2];
419cb32e2e7SValeria Barra   PetscScalar *r;
420cb32e2e7SValeria Barra   PetscBool test_mode, benchmark_mode, write_solution;
421cb32e2e7SValeria Barra   PetscMPIInt size, rank;
4229b072555Sjeremylt   PetscLogStage solve_stage;
4239b072555Sjeremylt   Vec X, X_loc, rhs, rhs_loc;
424cb32e2e7SValeria Barra   Mat mat;
425cb32e2e7SValeria Barra   KSP ksp;
4269b072555Sjeremylt   VecScatter l_to_g, l_to_g_0, g_to_g_D;
4279b072555Sjeremylt   PetscMemType mem_type;
428cb32e2e7SValeria Barra   User user;
429cb32e2e7SValeria Barra   Ceed ceed;
4309b072555Sjeremylt   CeedBasis basis_x, basis_u;
4319b072555Sjeremylt   CeedElemRestriction elem_restr_x, elem_restr_u, elem_restr_u_i, elem_restr_qd_i;
4329b072555Sjeremylt   CeedQFunction qf_setup_geo, qf_setup_rhs, qf_apply, qf_error;
4339b072555Sjeremylt   CeedOperator op_setup_geo, op_setup_rhs, op_apply, op_error;
4349b072555Sjeremylt   CeedVector x_coord, q_data, rhs_ceed, target;
435cb32e2e7SValeria Barra   CeedInt P, Q;
4369b072555Sjeremylt   const CeedInt dim = 3, num_comp_x = 3;
4379b072555Sjeremylt   BPType bp_choice;
438cb32e2e7SValeria Barra 
439cb32e2e7SValeria Barra   ierr = PetscInitialize(&argc, &argv, NULL, help);
440cb32e2e7SValeria Barra   if (ierr) return ierr;
441cb32e2e7SValeria Barra   comm = PETSC_COMM_WORLD;
44232d2ee49SValeria Barra 
44332d2ee49SValeria Barra   // Read command line options
44467490bc6SJeremy L Thompson   PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL);
4459b072555Sjeremylt   bp_choice = CEED_BP1;
446cb32e2e7SValeria Barra   ierr = PetscOptionsEnum("-problem",
447cb32e2e7SValeria Barra                           "CEED benchmark problem to solve", NULL,
4489b072555Sjeremylt                           bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice,
449cb32e2e7SValeria Barra                           NULL); CHKERRQ(ierr);
4509b072555Sjeremylt   num_comp_u = bp_options[bp_choice].num_comp_u;
451cb32e2e7SValeria Barra   test_mode = PETSC_FALSE;
452cb32e2e7SValeria Barra   ierr = PetscOptionsBool("-test",
453cb32e2e7SValeria Barra                           "Testing mode (do not print unless error is large)",
454cb32e2e7SValeria Barra                           NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr);
455cb32e2e7SValeria Barra   benchmark_mode = PETSC_FALSE;
456cb32e2e7SValeria Barra   ierr = PetscOptionsBool("-benchmark",
457cb32e2e7SValeria Barra                           "Benchmarking mode (prints benchmark statistics)",
458cb32e2e7SValeria Barra                           NULL, benchmark_mode, &benchmark_mode, NULL);
459cb32e2e7SValeria Barra   CHKERRQ(ierr);
460cb32e2e7SValeria Barra   write_solution = PETSC_FALSE;
461cb32e2e7SValeria Barra   ierr = PetscOptionsBool("-write_solution",
462cb32e2e7SValeria Barra                           "Write solution for visualization",
463cb32e2e7SValeria Barra                           NULL, write_solution, &write_solution, NULL);
464cb32e2e7SValeria Barra   CHKERRQ(ierr);
465cb32e2e7SValeria Barra   degree = test_mode ? 3 : 1;
466cb32e2e7SValeria Barra   ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis",
467cb32e2e7SValeria Barra                          NULL, degree, &degree, NULL); CHKERRQ(ierr);
4689b072555Sjeremylt   q_extra = bp_options[bp_choice].q_extra;
4699b072555Sjeremylt   ierr = PetscOptionsInt("-q_extra", "Number of extra quadrature points",
4709b072555Sjeremylt                          NULL, q_extra, &q_extra, NULL); CHKERRQ(ierr);
471cb32e2e7SValeria Barra   ierr = PetscOptionsString("-ceed", "CEED resource specifier",
4729b072555Sjeremylt                             NULL, ceed_resource, ceed_resource,
4739b072555Sjeremylt                             sizeof(ceed_resource), NULL); CHKERRQ(ierr);
4749b072555Sjeremylt   local_nodes = 1000;
475cb32e2e7SValeria Barra   ierr = PetscOptionsInt("-local",
476cb32e2e7SValeria Barra                          "Target number of locally owned nodes per process",
4779b072555Sjeremylt                          NULL, local_nodes, &local_nodes, NULL); CHKERRQ(ierr);
4782fbc6e41SJeremy L Thompson   PetscInt two = 2;
4792fbc6e41SJeremy L Thompson   ksp_max_it_clip[0] = 5;
4802fbc6e41SJeremy L Thompson   ksp_max_it_clip[1] = 20;
4812fbc6e41SJeremy L Thompson   ierr = PetscOptionsIntArray("-ksp_max_it_clip",
4822fbc6e41SJeremy L Thompson                               "Min and max number of iterations to use during benchmarking",
4832fbc6e41SJeremy L Thompson                               NULL, ksp_max_it_clip, &two, NULL);
4842fbc6e41SJeremy L Thompson   CHKERRQ(ierr);
48567490bc6SJeremy L Thompson   PetscOptionsEnd();
486cb32e2e7SValeria Barra   P = degree + 1;
4879b072555Sjeremylt   Q = P + q_extra;
488cb32e2e7SValeria Barra 
4899396343dSjeremylt   // Set up libCEED
4909b072555Sjeremylt   CeedInit(ceed_resource, &ceed);
4919b072555Sjeremylt   CeedMemType mem_type_backend;
4929b072555Sjeremylt   CeedGetPreferredMemType(ceed, &mem_type_backend);
4939396343dSjeremylt 
4949b072555Sjeremylt   VecType default_vec_type = NULL, vec_type;
4959b072555Sjeremylt   switch (mem_type_backend) {
4969b072555Sjeremylt   case CEED_MEM_HOST: default_vec_type = VECSTANDARD; break;
497b68a8d79SJed Brown   case CEED_MEM_DEVICE: {
498b68a8d79SJed Brown     const char *resolved;
499b68a8d79SJed Brown     CeedGetResource(ceed, &resolved);
5009b072555Sjeremylt     if (strstr(resolved, "/gpu/cuda")) default_vec_type = VECCUDA;
501b68a8d79SJed Brown     else if (strstr(resolved, "/gpu/hip/occa"))
5029b072555Sjeremylt       default_vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678
5039b072555Sjeremylt     else if (strstr(resolved, "/gpu/hip")) default_vec_type = VECHIP;
5049b072555Sjeremylt     else default_vec_type = VECSTANDARD;
505b68a8d79SJed Brown   }
506b68a8d79SJed Brown   }
5079396343dSjeremylt 
508cb32e2e7SValeria Barra   // Determine size of process grid
509cb32e2e7SValeria Barra   ierr = MPI_Comm_size(comm, &size); CHKERRQ(ierr);
510cb32e2e7SValeria Barra   Split3(size, p, false);
511cb32e2e7SValeria Barra 
5129b072555Sjeremylt   // Find a nicely composite number of elements no less than local_nodes
5139b072555Sjeremylt   for (local_elem = PetscMax(1, local_nodes / (degree*degree*degree)); ;
5149b072555Sjeremylt        local_elem++) {
5159b072555Sjeremylt     Split3(local_elem, mesh_elem, true);
5169b072555Sjeremylt     if (Max3(mesh_elem) / Min3(mesh_elem) <= 2) break;
517cb32e2e7SValeria Barra   }
518cb32e2e7SValeria Barra 
519cb32e2e7SValeria Barra   // Find my location in the process grid
520cb32e2e7SValeria Barra   ierr = MPI_Comm_rank(comm, &rank); CHKERRQ(ierr);
5219b072555Sjeremylt   for (int d=0, rank_left=rank; d<dim; d++) {
522cb32e2e7SValeria Barra     const int pstride[3] = {p[1] *p[2], p[2], 1};
5239b072555Sjeremylt     i_rank[d] = rank_left / pstride[d];
5249b072555Sjeremylt     rank_left -= i_rank[d] * pstride[d];
525cb32e2e7SValeria Barra   }
526cb32e2e7SValeria Barra 
5279b072555Sjeremylt   GlobalNodes(p, i_rank, degree, mesh_elem, m_nodes);
528cb32e2e7SValeria Barra 
529cb32e2e7SValeria Barra   // Setup global vector
5304a5da23aSJed Brown   ierr = VecCreate(comm, &X); CHKERRQ(ierr);
5319b072555Sjeremylt   ierr = VecSetType(X, default_vec_type); CHKERRQ(ierr);
5329b072555Sjeremylt   ierr = VecSetSizes(X, m_nodes[0]*m_nodes[1]*m_nodes[2]*num_comp_u,
5339b072555Sjeremylt                      PETSC_DECIDE);
534cb32e2e7SValeria Barra   CHKERRQ(ierr);
535b68a8d79SJed Brown   ierr = VecSetFromOptions(X); CHKERRQ(ierr);
536cb32e2e7SValeria Barra   ierr = VecSetUp(X); CHKERRQ(ierr);
537cb32e2e7SValeria Barra 
538cb32e2e7SValeria Barra   // Set up libCEED
5399b072555Sjeremylt   CeedInit(ceed_resource, &ceed);
540cb32e2e7SValeria Barra 
541cb32e2e7SValeria Barra   // Print summary
542cb32e2e7SValeria Barra   CeedInt gsize;
543cb32e2e7SValeria Barra   ierr = VecGetSize(X, &gsize); CHKERRQ(ierr);
544dc7d240cSValeria Barra   if (!test_mode) {
5459b072555Sjeremylt     const char *used_resource;
5469b072555Sjeremylt     CeedGetResource(ceed, &used_resource);
5479396343dSjeremylt 
5489b072555Sjeremylt     ierr = VecGetType(X, &vec_type); CHKERRQ(ierr);
5499396343dSjeremylt 
550cb32e2e7SValeria Barra     ierr = PetscPrintf(comm,
551990fdeb6SJeremy L Thompson                        "\n-- CEED Benchmark Problem %" CeedInt_FMT " -- libCEED + PETSc --\n"
5529396343dSjeremylt                        "  PETSc:\n"
5539396343dSjeremylt                        "    PETSc Vec Type                     : %s\n"
554cb32e2e7SValeria Barra                        "  libCEED:\n"
555cb32e2e7SValeria Barra                        "    libCEED Backend                    : %s\n"
5569396343dSjeremylt                        "    libCEED Backend MemType            : %s\n"
557cb32e2e7SValeria Barra                        "  Mesh:\n"
558*751eb813Srezgarshakeri                        "    Solution Order (P)                 : %" CeedInt_FMT "\n"
559*751eb813Srezgarshakeri                        "    Quadrature  Order (Q)              : %" CeedInt_FMT "\n"
56008140895SJed Brown                        "    Global nodes                       : %" PetscInt_FMT "\n"
56151ad7d5bSrezgarshakeri                        "    Process Decomposition              : %" PetscInt_FMT
56251ad7d5bSrezgarshakeri                        " %" PetscInt_FMT " %" PetscInt_FMT "\n"
56351ad7d5bSrezgarshakeri                        "    Local Elements                     : %" PetscInt_FMT
56451ad7d5bSrezgarshakeri                        " = %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n"
56551ad7d5bSrezgarshakeri                        "    Owned nodes                        : %" PetscInt_FMT
56651ad7d5bSrezgarshakeri                        " = %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n"
56708140895SJed Brown                        "    DoF per node                       : %" PetscInt_FMT "\n",
5689b072555Sjeremylt                        bp_choice+1, vec_type, used_resource,
5699b072555Sjeremylt                        CeedMemTypes[mem_type_backend],
5709b072555Sjeremylt                        P, Q,  gsize/num_comp_u, p[0], p[1], p[2], local_elem,
5719b072555Sjeremylt                        mesh_elem[0], mesh_elem[1], mesh_elem[2],
5729b072555Sjeremylt                        m_nodes[0]*m_nodes[1]*m_nodes[2], m_nodes[0], m_nodes[1],
5739b072555Sjeremylt                        m_nodes[2], num_comp_u); CHKERRQ(ierr);
574cb32e2e7SValeria Barra   }
575cb32e2e7SValeria Barra 
576cb32e2e7SValeria Barra   {
5779b072555Sjeremylt     l_size = 1;
578cb32e2e7SValeria Barra     for (int d=0; d<dim; d++) {
5799b072555Sjeremylt       l_nodes[d] = mesh_elem[d]*degree + 1;
5809b072555Sjeremylt       l_size *= l_nodes[d];
581cb32e2e7SValeria Barra     }
5829b072555Sjeremylt     ierr = VecCreate(PETSC_COMM_SELF, &X_loc); CHKERRQ(ierr);
5839b072555Sjeremylt     ierr = VecSetType(X_loc, default_vec_type); CHKERRQ(ierr);
5849b072555Sjeremylt     ierr = VecSetSizes(X_loc, l_size*num_comp_u, PETSC_DECIDE); CHKERRQ(ierr);
5859b072555Sjeremylt     ierr = VecSetFromOptions(X_loc); CHKERRQ(ierr);
5869b072555Sjeremylt     ierr = VecSetUp(X_loc); CHKERRQ(ierr);
587cb32e2e7SValeria Barra 
588cb32e2e7SValeria Barra     // Create local-to-global scatter
5899b072555Sjeremylt     PetscInt *l_to_g_ind, *l_to_g_ind_0, *loc_ind, l_0_count;
5909b072555Sjeremylt     IS l_to_g_is, l_to_g_is_0, loc_is;
5919b072555Sjeremylt     PetscInt g_start[2][2][2], g_m_nodes[2][2][2][dim];
592cb32e2e7SValeria Barra 
593cb32e2e7SValeria Barra     for (int i=0; i<2; i++) {
594cb32e2e7SValeria Barra       for (int j=0; j<2; j++) {
595cb32e2e7SValeria Barra         for (int k=0; k<2; k++) {
5969b072555Sjeremylt           PetscInt ijk_rank[3] = {i_rank[0]+i, i_rank[1]+j, i_rank[2]+k};
5979b072555Sjeremylt           g_start[i][j][k] = GlobalStart(p, ijk_rank, degree, mesh_elem);
5989b072555Sjeremylt           GlobalNodes(p, ijk_rank, degree, mesh_elem, g_m_nodes[i][j][k]);
599cb32e2e7SValeria Barra         }
600cb32e2e7SValeria Barra       }
601cb32e2e7SValeria Barra     }
602cb32e2e7SValeria Barra 
6039b072555Sjeremylt     ierr = PetscMalloc1(l_size, &l_to_g_ind); CHKERRQ(ierr);
6049b072555Sjeremylt     ierr = PetscMalloc1(l_size, &l_to_g_ind_0); CHKERRQ(ierr);
6059b072555Sjeremylt     ierr = PetscMalloc1(l_size, &loc_ind); CHKERRQ(ierr);
6069b072555Sjeremylt     l_0_count = 0;
6079b072555Sjeremylt     for (PetscInt i=0,ir,ii; ir=i>=m_nodes[0], ii=i-ir*m_nodes[0], i<l_nodes[0];
6089b072555Sjeremylt          i++)
6099b072555Sjeremylt       for (PetscInt j=0,jr,jj; jr=j>=m_nodes[1], jj=j-jr*m_nodes[1], j<l_nodes[1];
6109b072555Sjeremylt            j++)
6119b072555Sjeremylt         for (PetscInt k=0,kr,kk; kr=k>=m_nodes[2], kk=k-kr*m_nodes[2], k<l_nodes[2];
6129b072555Sjeremylt              k++) {
6139b072555Sjeremylt           PetscInt here = (i*l_nodes[1]+j)*l_nodes[2]+k;
6149b072555Sjeremylt           l_to_g_ind[here] =
6159b072555Sjeremylt             g_start[ir][jr][kr] + (ii*g_m_nodes[ir][jr][kr][1]+jj)*g_m_nodes[ir][jr][kr][2]
6169b072555Sjeremylt             +kk;
6179b072555Sjeremylt           if ((i_rank[0] == 0 && i == 0)
6189b072555Sjeremylt               || (i_rank[1] == 0 && j == 0)
6199b072555Sjeremylt               || (i_rank[2] == 0 && k == 0)
6209b072555Sjeremylt               || (i_rank[0]+1 == p[0] && i+1 == l_nodes[0])
6219b072555Sjeremylt               || (i_rank[1]+1 == p[1] && j+1 == l_nodes[1])
6229b072555Sjeremylt               || (i_rank[2]+1 == p[2] && k+1 == l_nodes[2]))
623cb32e2e7SValeria Barra             continue;
6249b072555Sjeremylt           l_to_g_ind_0[l_0_count] = l_to_g_ind[here];
6259b072555Sjeremylt           loc_ind[l_0_count++] = here;
626cb32e2e7SValeria Barra         }
6279b072555Sjeremylt     ierr = ISCreateBlock(comm, num_comp_u, l_size, l_to_g_ind, PETSC_OWN_POINTER,
6289b072555Sjeremylt                          &l_to_g_is); CHKERRQ(ierr);
6299b072555Sjeremylt     ierr = VecScatterCreate(X_loc, NULL, X, l_to_g_is, &l_to_g); CHKERRQ(ierr);
630cb32e2e7SValeria Barra     CHKERRQ(ierr);
6319b072555Sjeremylt     ierr = ISCreateBlock(comm, num_comp_u, l_0_count, l_to_g_ind_0,
6329b072555Sjeremylt                          PETSC_OWN_POINTER,
6339b072555Sjeremylt                          &l_to_g_is_0); CHKERRQ(ierr);
6349b072555Sjeremylt     ierr = ISCreateBlock(comm, num_comp_u, l_0_count, loc_ind, PETSC_OWN_POINTER,
6359b072555Sjeremylt                          &loc_is); CHKERRQ(ierr);
6369b072555Sjeremylt     ierr = VecScatterCreate(X_loc, loc_is, X, l_to_g_is_0, &l_to_g_0);
6379b072555Sjeremylt     CHKERRQ(ierr);
638cb32e2e7SValeria Barra     {
639cb32e2e7SValeria Barra       // Create global-to-global scatter for Dirichlet values (everything not in
6409b072555Sjeremylt       // l_to_g_is_0, which is the range of l_to_g_0)
6419b072555Sjeremylt       PetscInt x_start, x_end, *ind_D, count_D = 0;
6429b072555Sjeremylt       IS is_D;
643cb32e2e7SValeria Barra       const PetscScalar *x;
6449b072555Sjeremylt       ierr = VecZeroEntries(X_loc); CHKERRQ(ierr);
645cb32e2e7SValeria Barra       ierr = VecSet(X, 1.0); CHKERRQ(ierr);
6469b072555Sjeremylt       ierr = VecScatterBegin(l_to_g_0, X_loc, X, INSERT_VALUES, SCATTER_FORWARD);
647cb32e2e7SValeria Barra       CHKERRQ(ierr);
6489b072555Sjeremylt       ierr = VecScatterEnd(l_to_g_0, X_loc, X, INSERT_VALUES, SCATTER_FORWARD);
649cb32e2e7SValeria Barra       CHKERRQ(ierr);
6509b072555Sjeremylt       ierr = VecGetOwnershipRange(X, &x_start, &x_end); CHKERRQ(ierr);
6519b072555Sjeremylt       ierr = PetscMalloc1(x_end-x_start, &ind_D); CHKERRQ(ierr);
652cb32e2e7SValeria Barra       ierr = VecGetArrayRead(X, &x); CHKERRQ(ierr);
6539b072555Sjeremylt       for (PetscInt i=0; i<x_end-x_start; i++) {
6549b072555Sjeremylt         if (x[i] == 1.) ind_D[count_D++] = x_start + i;
655cb32e2e7SValeria Barra       }
656cb32e2e7SValeria Barra       ierr = VecRestoreArrayRead(X, &x); CHKERRQ(ierr);
6579b072555Sjeremylt       ierr = ISCreateGeneral(comm, count_D, ind_D, PETSC_COPY_VALUES, &is_D);
658cb32e2e7SValeria Barra       CHKERRQ(ierr);
6599b072555Sjeremylt       ierr = PetscFree(ind_D); CHKERRQ(ierr);
6609b072555Sjeremylt       ierr = VecScatterCreate(X, is_D, X, is_D, &g_to_g_D); CHKERRQ(ierr);
6619b072555Sjeremylt       ierr = ISDestroy(&is_D); CHKERRQ(ierr);
662cb32e2e7SValeria Barra     }
6639b072555Sjeremylt     ierr = ISDestroy(&l_to_g_is); CHKERRQ(ierr);
6649b072555Sjeremylt     ierr = ISDestroy(&l_to_g_is_0); CHKERRQ(ierr);
6659b072555Sjeremylt     ierr = ISDestroy(&loc_is); CHKERRQ(ierr);
666cb32e2e7SValeria Barra   }
667cb32e2e7SValeria Barra 
668cb32e2e7SValeria Barra   // CEED bases
6699b072555Sjeremylt   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_u, P, Q,
6709b072555Sjeremylt                                   bp_options[bp_choice].q_mode, &basis_u);
6719b072555Sjeremylt   CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q,
6729b072555Sjeremylt                                   bp_options[bp_choice].q_mode, &basis_x);
673cb32e2e7SValeria Barra 
674cb32e2e7SValeria Barra   // CEED restrictions
6759b072555Sjeremylt   CreateRestriction(ceed, mesh_elem, P, num_comp_u, &elem_restr_u);
6769b072555Sjeremylt   CreateRestriction(ceed, mesh_elem, 2, dim, &elem_restr_x);
6779b072555Sjeremylt   CeedInt num_elem = mesh_elem[0]*mesh_elem[1]*mesh_elem[2];
6789b072555Sjeremylt   CeedElemRestrictionCreateStrided(ceed, num_elem, Q*Q*Q, num_comp_u,
6799b072555Sjeremylt                                    num_comp_u*num_elem*Q*Q*Q,
6809b072555Sjeremylt                                    CEED_STRIDES_BACKEND, &elem_restr_u_i);
6819b072555Sjeremylt   CeedElemRestrictionCreateStrided(ceed, num_elem, Q*Q*Q,
6829b072555Sjeremylt                                    bp_options[bp_choice].q_data_size,
6839b072555Sjeremylt                                    bp_options[bp_choice].q_data_size*num_elem*Q*Q*Q,
6849b072555Sjeremylt                                    CEED_STRIDES_BACKEND, &elem_restr_qd_i);
685cb32e2e7SValeria Barra   {
6869b072555Sjeremylt     CeedScalar *x_loc;
6879b072555Sjeremylt     CeedInt shape[3] = {mesh_elem[0]+1, mesh_elem[1]+1, mesh_elem[2]+1}, len =
688cb32e2e7SValeria Barra                          shape[0]*shape[1]*shape[2];
6899b072555Sjeremylt     x_loc = malloc(len*num_comp_x*sizeof x_loc[0]);
690cb32e2e7SValeria Barra     for (CeedInt i=0; i<shape[0]; i++) {
691cb32e2e7SValeria Barra       for (CeedInt j=0; j<shape[1]; j++) {
692cb32e2e7SValeria Barra         for (CeedInt k=0; k<shape[2]; k++) {
6939b072555Sjeremylt           x_loc[dim*((i*shape[1]+j)*shape[2]+k) + 0] = 1.*(i_rank[0]*mesh_elem[0]+i) /
6949b072555Sjeremylt               (p[0]*mesh_elem[0]);
6959b072555Sjeremylt           x_loc[dim*((i*shape[1]+j)*shape[2]+k) + 1] = 1.*(i_rank[1]*mesh_elem[1]+j) /
6969b072555Sjeremylt               (p[1]*mesh_elem[1]);
6979b072555Sjeremylt           x_loc[dim*((i*shape[1]+j)*shape[2]+k) + 2] = 1.*(i_rank[2]*mesh_elem[2]+k) /
6989b072555Sjeremylt               (p[2]*mesh_elem[2]);
699cb32e2e7SValeria Barra         }
700cb32e2e7SValeria Barra       }
701cb32e2e7SValeria Barra     }
7029b072555Sjeremylt     CeedVectorCreate(ceed, len*num_comp_x, &x_coord);
7039b072555Sjeremylt     CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_OWN_POINTER, x_loc);
704cb32e2e7SValeria Barra   }
705cb32e2e7SValeria Barra 
7069b072555Sjeremylt   // Create the QFunction that builds the operator quadrature data
7079b072555Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].setup_geo,
7089b072555Sjeremylt                               bp_options[bp_choice].setup_geo_loc, &qf_setup_geo);
7099b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP);
7109b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x*dim, CEED_EVAL_GRAD);
7119b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT);
7129b072555Sjeremylt   CeedQFunctionAddOutput(qf_setup_geo, "q_data",
7139b072555Sjeremylt                          bp_options[bp_choice].q_data_size,
714cb32e2e7SValeria Barra                          CEED_EVAL_NONE);
715cb32e2e7SValeria Barra 
7169b072555Sjeremylt   // Create the QFunction that sets up the RHS and true solution
7179b072555Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].setup_rhs,
7189b072555Sjeremylt                               bp_options[bp_choice].setup_rhs_loc, &qf_setup_rhs);
7199b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP);
7209b072555Sjeremylt   CeedQFunctionAddInput(qf_setup_rhs, "q_data", bp_options[bp_choice].q_data_size,
721e83e87a5Sjeremylt                         CEED_EVAL_NONE);
7229b072555Sjeremylt   CeedQFunctionAddOutput(qf_setup_rhs, "true_soln", num_comp_u, CEED_EVAL_NONE);
7239b072555Sjeremylt   CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP);
724cb32e2e7SValeria Barra 
725cb32e2e7SValeria Barra   // Set up PDE operator
7269b072555Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].apply,
7279b072555Sjeremylt                               bp_options[bp_choice].apply_loc, &qf_apply);
728cb32e2e7SValeria Barra   // Add inputs and outputs
7299b072555Sjeremylt   CeedInt in_scale = bp_options[bp_choice].in_mode==CEED_EVAL_GRAD ? 3 : 1;
7309b072555Sjeremylt   CeedInt out_scale = bp_options[bp_choice].out_mode==CEED_EVAL_GRAD ? 3 : 1;
7319b072555Sjeremylt   CeedQFunctionAddInput(qf_apply, "u", num_comp_u*in_scale,
7329b072555Sjeremylt                         bp_options[bp_choice].in_mode);
7339b072555Sjeremylt   CeedQFunctionAddInput(qf_apply, "q_data", bp_options[bp_choice].q_data_size,
734cb32e2e7SValeria Barra                         CEED_EVAL_NONE);
7359b072555Sjeremylt   CeedQFunctionAddOutput(qf_apply, "v", num_comp_u*out_scale,
7369b072555Sjeremylt                          bp_options[bp_choice].out_mode);
737cb32e2e7SValeria Barra 
738cb32e2e7SValeria Barra   // Create the error qfunction
7399b072555Sjeremylt   CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error,
7409b072555Sjeremylt                               bp_options[bp_choice].error_loc, &qf_error);
7419b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP);
7429b072555Sjeremylt   CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE);
7439b072555Sjeremylt   CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_NONE);
744cb32e2e7SValeria Barra 
745cb32e2e7SValeria Barra   // Create the persistent vectors that will be needed in setup
7469b072555Sjeremylt   CeedInt num_qpts;
7479b072555Sjeremylt   CeedBasisGetNumQuadraturePoints(basis_u, &num_qpts);
7489b072555Sjeremylt   CeedVectorCreate(ceed, bp_options[bp_choice].q_data_size*num_elem*num_qpts,
7499b072555Sjeremylt                    &q_data);
7509b072555Sjeremylt   CeedVectorCreate(ceed, num_elem*num_qpts*num_comp_u, &target);
7519b072555Sjeremylt   CeedVectorCreate(ceed, l_size*num_comp_u, &rhs_ceed);
752cb32e2e7SValeria Barra 
753cb32e2e7SValeria Barra   // Create the operator that builds the quadrature data for the ceed operator
7549b072555Sjeremylt   CeedOperatorCreate(ceed, qf_setup_geo, CEED_QFUNCTION_NONE,
7559b072555Sjeremylt                      CEED_QFUNCTION_NONE, &op_setup_geo);
7569b072555Sjeremylt   CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x,
757e83e87a5Sjeremylt                        CEED_VECTOR_ACTIVE);
7589b072555Sjeremylt   CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x,
759a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
7609b072555Sjeremylt   CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x,
761a8d32208Sjeremylt                        CEED_VECTOR_NONE);
7629b072555Sjeremylt   CeedOperatorSetField(op_setup_geo, "q_data", elem_restr_qd_i,
763cb32e2e7SValeria Barra                        CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE);
764cb32e2e7SValeria Barra 
765cb32e2e7SValeria Barra   // Create the operator that builds the RHS and true solution
7669b072555Sjeremylt   CeedOperatorCreate(ceed, qf_setup_rhs, CEED_QFUNCTION_NONE,
7679b072555Sjeremylt                      CEED_QFUNCTION_NONE, &op_setup_rhs);
7689b072555Sjeremylt   CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x,
769a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
7709b072555Sjeremylt   CeedOperatorSetField(op_setup_rhs, "q_data", elem_restr_qd_i,
7719b072555Sjeremylt                        CEED_BASIS_COLLOCATED,
7729b072555Sjeremylt                        q_data);
7739b072555Sjeremylt   CeedOperatorSetField(op_setup_rhs, "true_soln", elem_restr_u_i,
774cb32e2e7SValeria Barra                        CEED_BASIS_COLLOCATED, target);
7759b072555Sjeremylt   CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u,
776a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
777cb32e2e7SValeria Barra 
778cb32e2e7SValeria Barra   // Create the mass or diff operator
7799b072555Sjeremylt   CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
7809b072555Sjeremylt                      &op_apply);
7819b072555Sjeremylt   CeedOperatorSetField(op_apply, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
7829b072555Sjeremylt   CeedOperatorSetField(op_apply, "q_data", elem_restr_qd_i, CEED_BASIS_COLLOCATED,
7839b072555Sjeremylt                        q_data);
7849b072555Sjeremylt   CeedOperatorSetField(op_apply, "v", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
785cb32e2e7SValeria Barra 
786cb32e2e7SValeria Barra   // Create the error operator
7879b072555Sjeremylt   CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE,
7889b072555Sjeremylt                      &op_error);
7899b072555Sjeremylt   CeedOperatorSetField(op_error, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE);
7909b072555Sjeremylt   CeedOperatorSetField(op_error, "true_soln", elem_restr_u_i,
791cb32e2e7SValeria Barra                        CEED_BASIS_COLLOCATED, target);
7929b072555Sjeremylt   CeedOperatorSetField(op_error, "error", elem_restr_u_i, CEED_BASIS_COLLOCATED,
793a8d32208Sjeremylt                        CEED_VECTOR_ACTIVE);
794cb32e2e7SValeria Barra 
795cb32e2e7SValeria Barra   // Set up Mat
796cb32e2e7SValeria Barra   ierr = PetscMalloc1(1, &user); CHKERRQ(ierr);
797cb32e2e7SValeria Barra   user->comm = comm;
7989b072555Sjeremylt   user->l_to_g = l_to_g;
7999b072555Sjeremylt   if (bp_choice != CEED_BP1 && bp_choice != CEED_BP2) {
8009b072555Sjeremylt     user->l_to_g_0 = l_to_g_0;
8019b072555Sjeremylt     user->g_to_g_D = g_to_g_D;
802cb32e2e7SValeria Barra   }
8039b072555Sjeremylt   user->X_loc = X_loc;
8049b072555Sjeremylt   ierr = VecDuplicate(X_loc, &user->Y_loc); CHKERRQ(ierr);
8059b072555Sjeremylt   CeedVectorCreate(ceed, l_size*num_comp_u, &user->x_ceed);
8069b072555Sjeremylt   CeedVectorCreate(ceed, l_size*num_comp_u, &user->y_ceed);
8079b072555Sjeremylt   user->op = op_apply;
8089b072555Sjeremylt   user->q_data = q_data;
809cb32e2e7SValeria Barra   user->ceed = ceed;
810cb32e2e7SValeria Barra 
8119b072555Sjeremylt   ierr = MatCreateShell(comm, m_nodes[0]*m_nodes[1]*m_nodes[2]*num_comp_u,
8129b072555Sjeremylt                         m_nodes[0]*m_nodes[1]*m_nodes[2]*num_comp_u,
813cb32e2e7SValeria Barra                         PETSC_DECIDE, PETSC_DECIDE, user, &mat); CHKERRQ(ierr);
8149b072555Sjeremylt   if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) {
815cb32e2e7SValeria Barra     ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Mass);
816cb32e2e7SValeria Barra     CHKERRQ(ierr);
817cb32e2e7SValeria Barra   } else {
818cb32e2e7SValeria Barra     ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Diff);
819cb32e2e7SValeria Barra     CHKERRQ(ierr);
820cb32e2e7SValeria Barra   }
8219b072555Sjeremylt   ierr = VecGetType(X, &vec_type); CHKERRQ(ierr);
8229b072555Sjeremylt   ierr = MatShellSetVecType(mat, vec_type); CHKERRQ(ierr);
823cb32e2e7SValeria Barra 
824cb32e2e7SValeria Barra   // Get RHS vector
8259396343dSjeremylt   ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr);
8269b072555Sjeremylt   ierr = VecDuplicate(X_loc, &rhs_loc); CHKERRQ(ierr);
8279b072555Sjeremylt   ierr = VecZeroEntries(rhs_loc); CHKERRQ(ierr);
8289b072555Sjeremylt   ierr = VecGetArrayAndMemType(rhs_loc, &r, &mem_type); CHKERRQ(ierr);
8299b072555Sjeremylt   CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r);
830cb32e2e7SValeria Barra 
8319b072555Sjeremylt   // Setup q_data, rhs, and target
8329b072555Sjeremylt   CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE);
8339b072555Sjeremylt   CeedOperatorApply(op_setup_rhs, x_coord, rhs_ceed, CEED_REQUEST_IMMEDIATE);
8349b072555Sjeremylt   CeedVectorDestroy(&x_coord);
835cb32e2e7SValeria Barra 
836cb32e2e7SValeria Barra   // Gather RHS
8379b072555Sjeremylt   ierr = CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL); CHKERRQ(ierr);
8389b072555Sjeremylt   ierr = VecRestoreArrayAndMemType(rhs_loc, &r); CHKERRQ(ierr);
839cb32e2e7SValeria Barra   ierr = VecZeroEntries(rhs); CHKERRQ(ierr);
8409b072555Sjeremylt   ierr = VecScatterBegin(l_to_g, rhs_loc, rhs, ADD_VALUES, SCATTER_FORWARD);
841cb32e2e7SValeria Barra   CHKERRQ(ierr);
8429b072555Sjeremylt   ierr = VecScatterEnd(l_to_g, rhs_loc, rhs, ADD_VALUES, SCATTER_FORWARD);
843cb32e2e7SValeria Barra   CHKERRQ(ierr);
8449b072555Sjeremylt   CeedVectorDestroy(&rhs_ceed);
845cb32e2e7SValeria Barra 
846cb32e2e7SValeria Barra   ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr);
847cb32e2e7SValeria Barra   {
848cb32e2e7SValeria Barra     PC pc;
849cb32e2e7SValeria Barra     ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr);
8509b072555Sjeremylt     if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) {
851cb32e2e7SValeria Barra       ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr);
852cb32e2e7SValeria Barra       ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr);
853cb32e2e7SValeria Barra     } else {
854cb32e2e7SValeria Barra       ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr);
855cb32e2e7SValeria Barra     }
856cb32e2e7SValeria Barra     ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr);
857cb32e2e7SValeria Barra     ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr);
858cb32e2e7SValeria Barra     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
859cb32e2e7SValeria Barra                             PETSC_DEFAULT); CHKERRQ(ierr);
860cb32e2e7SValeria Barra   }
861cb32e2e7SValeria Barra   ierr = KSPSetOperators(ksp, mat, mat); CHKERRQ(ierr);
862da9108adSvaleriabarra   // First run's performance log is not considered for benchmarking purposes
863cb32e2e7SValeria Barra   ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1);
864cb32e2e7SValeria Barra   CHKERRQ(ierr);
865cb32e2e7SValeria Barra   my_rt_start = MPI_Wtime();
866cb32e2e7SValeria Barra   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
867cb32e2e7SValeria Barra   my_rt = MPI_Wtime() - my_rt_start;
868cb32e2e7SValeria Barra   ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm);
869cb32e2e7SValeria Barra   CHKERRQ(ierr);
870cb32e2e7SValeria Barra   // Set maxits based on first iteration timing
871cb32e2e7SValeria Barra   if (my_rt > 0.02) {
8722fbc6e41SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
8732fbc6e41SJeremy L Thompson                             ksp_max_it_clip[0]); CHKERRQ(ierr);
874cb32e2e7SValeria Barra   } else {
8752fbc6e41SJeremy L Thompson     ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT,
8762fbc6e41SJeremy L Thompson                             ksp_max_it_clip[1]); CHKERRQ(ierr);
877cb32e2e7SValeria Barra   }
878debcf919SJed Brown   ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr);
87909a940d7Sjeremylt 
880cb32e2e7SValeria Barra   // Timed solve
88109a940d7Sjeremylt   ierr = VecZeroEntries(X); CHKERRQ(ierr);
882cb32e2e7SValeria Barra   ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr);
88309a940d7Sjeremylt 
88409a940d7Sjeremylt   // -- Performance logging
8859b072555Sjeremylt   ierr = PetscLogStageRegister("Solve Stage", &solve_stage); CHKERRQ(ierr);
8869b072555Sjeremylt   ierr = PetscLogStagePush(solve_stage); CHKERRQ(ierr);
88709a940d7Sjeremylt 
88809a940d7Sjeremylt   // -- Solve
889cb32e2e7SValeria Barra   my_rt_start = MPI_Wtime();
890cb32e2e7SValeria Barra   ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr);
891cb32e2e7SValeria Barra   my_rt = MPI_Wtime() - my_rt_start;
89209a940d7Sjeremylt 
89309a940d7Sjeremylt   // -- Performance logging
89409a940d7Sjeremylt   ierr = PetscLogStagePop();
89509a940d7Sjeremylt 
8968e87e98bSJed Brown   // Output results
897cb32e2e7SValeria Barra   {
8989b072555Sjeremylt     KSPType ksp_type;
899cb32e2e7SValeria Barra     KSPConvergedReason reason;
900cb32e2e7SValeria Barra     PetscReal rnorm;
901cb32e2e7SValeria Barra     PetscInt its;
9029b072555Sjeremylt     ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr);
903cb32e2e7SValeria Barra     ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr);
904cb32e2e7SValeria Barra     ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr);
905cb32e2e7SValeria Barra     ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr);
906a22c4fb5SJed Brown     if (!test_mode || reason < 0 || rnorm > 1e-8) {
907cb32e2e7SValeria Barra       ierr = PetscPrintf(comm,
908cb32e2e7SValeria Barra                          "  KSP:\n"
909cb32e2e7SValeria Barra                          "    KSP Type                           : %s\n"
910cb32e2e7SValeria Barra                          "    KSP Convergence                    : %s\n"
91108140895SJed Brown                          "    Total KSP Iterations               : %" PetscInt_FMT "\n"
912cb32e2e7SValeria Barra                          "    Final rnorm                        : %e\n",
9139b072555Sjeremylt                          ksp_type, KSPConvergedReasons[reason], its,
914cb32e2e7SValeria Barra                          (double)rnorm); CHKERRQ(ierr);
915cb32e2e7SValeria Barra     }
9168e87e98bSJed Brown     if (!test_mode) {
9178e87e98bSJed Brown       ierr = PetscPrintf(comm,"  Performance:\n"); CHKERRQ(ierr);
9188e87e98bSJed Brown     }
9198e87e98bSJed Brown     {
9209b072555Sjeremylt       PetscReal max_error;
9219b072555Sjeremylt       ierr = ComputeErrorMax(user, op_error, X, target, &max_error);
9228e87e98bSJed Brown       CHKERRQ(ierr);
9238e87e98bSJed Brown       PetscReal tol = 5e-2;
9249b072555Sjeremylt       if (!test_mode || max_error > tol) {
925cb32e2e7SValeria Barra         ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm);
926cb32e2e7SValeria Barra         CHKERRQ(ierr);
927cb32e2e7SValeria Barra         ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm);
928cb32e2e7SValeria Barra         CHKERRQ(ierr);
929cb32e2e7SValeria Barra         ierr = PetscPrintf(comm,
9308e87e98bSJed Brown                            "    Pointwise Error (max)              : %e\n"
9318e87e98bSJed Brown                            "    CG Solve Time                      : %g (%g) sec\n",
9329b072555Sjeremylt                            (double)max_error, rt_max, rt_min); CHKERRQ(ierr);
933cb32e2e7SValeria Barra       }
934cb32e2e7SValeria Barra     }
9358d0bb2bbSvaleriabarra     if (!test_mode) {
936cb32e2e7SValeria Barra       ierr = PetscPrintf(comm,
9378e87e98bSJed Brown                          "    DoFs/Sec in CG                     : %g (%g) million\n",
9388e87e98bSJed Brown                          1e-6*gsize*its/rt_max,
9398e87e98bSJed Brown                          1e-6*gsize*its/rt_min); CHKERRQ(ierr);
940cb32e2e7SValeria Barra     }
941cb32e2e7SValeria Barra   }
942cb32e2e7SValeria Barra 
943cb32e2e7SValeria Barra   if (write_solution) {
9449b072555Sjeremylt     PetscViewer vtk_viewer_soln;
945cb32e2e7SValeria Barra 
9469b072555Sjeremylt     ierr = PetscViewerCreate(comm, &vtk_viewer_soln); CHKERRQ(ierr);
9479b072555Sjeremylt     ierr = PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK); CHKERRQ(ierr);
9489b072555Sjeremylt     ierr = PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"); CHKERRQ(ierr);
9499b072555Sjeremylt     ierr = VecView(X, vtk_viewer_soln); CHKERRQ(ierr);
9509b072555Sjeremylt     ierr = PetscViewerDestroy(&vtk_viewer_soln); CHKERRQ(ierr);
951cb32e2e7SValeria Barra   }
952cb32e2e7SValeria Barra 
953cb32e2e7SValeria Barra   ierr = VecDestroy(&rhs); CHKERRQ(ierr);
9549b072555Sjeremylt   ierr = VecDestroy(&rhs_loc); CHKERRQ(ierr);
955cb32e2e7SValeria Barra   ierr = VecDestroy(&X); CHKERRQ(ierr);
9569b072555Sjeremylt   ierr = VecDestroy(&user->X_loc); CHKERRQ(ierr);
9579b072555Sjeremylt   ierr = VecDestroy(&user->Y_loc); CHKERRQ(ierr);
9589b072555Sjeremylt   ierr = VecScatterDestroy(&l_to_g); CHKERRQ(ierr);
9599b072555Sjeremylt   ierr = VecScatterDestroy(&l_to_g_0); CHKERRQ(ierr);
9609b072555Sjeremylt   ierr = VecScatterDestroy(&g_to_g_D); CHKERRQ(ierr);
961cb32e2e7SValeria Barra   ierr = MatDestroy(&mat); CHKERRQ(ierr);
962cb32e2e7SValeria Barra   ierr = KSPDestroy(&ksp); CHKERRQ(ierr);
963cb32e2e7SValeria Barra 
9649b072555Sjeremylt   CeedVectorDestroy(&user->x_ceed);
9659b072555Sjeremylt   CeedVectorDestroy(&user->y_ceed);
9669b072555Sjeremylt   CeedVectorDestroy(&user->q_data);
967cb32e2e7SValeria Barra   CeedVectorDestroy(&target);
9689b072555Sjeremylt   CeedOperatorDestroy(&op_setup_geo);
9699b072555Sjeremylt   CeedOperatorDestroy(&op_setup_rhs);
9709b072555Sjeremylt   CeedOperatorDestroy(&op_apply);
9719b072555Sjeremylt   CeedOperatorDestroy(&op_error);
9729b072555Sjeremylt   CeedElemRestrictionDestroy(&elem_restr_u);
9739b072555Sjeremylt   CeedElemRestrictionDestroy(&elem_restr_x);
9749b072555Sjeremylt   CeedElemRestrictionDestroy(&elem_restr_u_i);
9759b072555Sjeremylt   CeedElemRestrictionDestroy(&elem_restr_qd_i);
9769b072555Sjeremylt   CeedQFunctionDestroy(&qf_setup_geo);
9779b072555Sjeremylt   CeedQFunctionDestroy(&qf_setup_rhs);
9789b072555Sjeremylt   CeedQFunctionDestroy(&qf_apply);
9799b072555Sjeremylt   CeedQFunctionDestroy(&qf_error);
9809b072555Sjeremylt   CeedBasisDestroy(&basis_u);
9819b072555Sjeremylt   CeedBasisDestroy(&basis_x);
982cb32e2e7SValeria Barra   CeedDestroy(&ceed);
983cb32e2e7SValeria Barra   ierr = PetscFree(user); CHKERRQ(ierr);
984cb32e2e7SValeria Barra   return PetscFinalize();
985cb32e2e7SValeria Barra }
986