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 } 969b072555Sjeremylt static int CreateRestriction(Ceed ceed, const CeedInt mesh_elem[3], CeedInt P, 979b072555Sjeremylt CeedInt num_comp, CeedElemRestriction *elem_restr) { 989b072555Sjeremylt const PetscInt num_elem = mesh_elem[0]*mesh_elem[1]*mesh_elem[2]; 999b072555Sjeremylt PetscInt m_nodes[3], *idx, *idx_p; 100cb32e2e7SValeria Barra 101cb32e2e7SValeria Barra // Get indicies 1029b072555Sjeremylt for (int d=0; d<3; d++) m_nodes[d] = mesh_elem[d]*(P-1) + 1; 1039b072555Sjeremylt idx_p = idx = malloc(num_elem*P*P*P*sizeof idx[0]); 1049b072555Sjeremylt for (CeedInt i=0; i<mesh_elem[0]; i++) 1059b072555Sjeremylt for (CeedInt j=0; j<mesh_elem[1]; j++) 1069b072555Sjeremylt for (CeedInt k=0; k<mesh_elem[2]; k++,idx_p += P*P*P) 107cb32e2e7SValeria Barra for (CeedInt ii=0; ii<P; ii++) 108cb32e2e7SValeria Barra for (CeedInt jj=0; jj<P; jj++) 109cb32e2e7SValeria Barra for (CeedInt kk=0; kk<P; kk++) { 110cb32e2e7SValeria Barra if (0) { // This is the C-style (i,j,k) ordering that I prefer 1119b072555Sjeremylt idx_p[(ii*P+jj)*P+kk] = num_comp*(((i*(P-1)+ii)*m_nodes[1] 1129b072555Sjeremylt + (j*(P-1)+jj))*m_nodes[2] 113cb32e2e7SValeria Barra + (k*(P-1)+kk)); 114cb32e2e7SValeria Barra } else { // (k,j,i) ordering for consistency with MFEM example 1159b072555Sjeremylt idx_p[ii+P*(jj+P*kk)] = num_comp*(((i*(P-1)+ii)*m_nodes[1] 1169b072555Sjeremylt + (j*(P-1)+jj))*m_nodes[2] 117cb32e2e7SValeria Barra + (k*(P-1)+kk)); 118cb32e2e7SValeria Barra } 119cb32e2e7SValeria Barra } 120cb32e2e7SValeria Barra 121cb32e2e7SValeria Barra // Setup CEED restriction 1229b072555Sjeremylt CeedElemRestrictionCreate(ceed, num_elem, P*P*P, num_comp, 1, 1239b072555Sjeremylt m_nodes[0]*m_nodes[1]*m_nodes[2]*num_comp, 1249b072555Sjeremylt CEED_MEM_HOST, CEED_OWN_POINTER, idx, elem_restr); 125cb32e2e7SValeria Barra 126cb32e2e7SValeria Barra PetscFunctionReturn(0); 127cb32e2e7SValeria Barra } 128cb32e2e7SValeria Barra 129cb32e2e7SValeria Barra // Data for PETSc 130cb32e2e7SValeria Barra typedef struct User_ *User; 131cb32e2e7SValeria Barra struct User_ { 132cb32e2e7SValeria Barra MPI_Comm comm; 1339b072555Sjeremylt VecScatter l_to_g; // Scatter for all entries 1349b072555Sjeremylt VecScatter l_to_g_0; // Skip Dirichlet values 1359b072555Sjeremylt VecScatter g_to_g_D; // global-to-global; only Dirichlet values 1369b072555Sjeremylt Vec X_loc, Y_loc; 1379b072555Sjeremylt CeedVector x_ceed, y_ceed; 138cb32e2e7SValeria Barra CeedOperator op; 1399b072555Sjeremylt CeedVector q_data; 140cb32e2e7SValeria Barra Ceed ceed; 141cb32e2e7SValeria Barra }; 142cb32e2e7SValeria Barra 143cb32e2e7SValeria Barra // BP Options 144cb32e2e7SValeria Barra typedef enum { 145cb32e2e7SValeria Barra CEED_BP1 = 0, CEED_BP2 = 1, CEED_BP3 = 2, 146cb32e2e7SValeria Barra CEED_BP4 = 3, CEED_BP5 = 4, CEED_BP6 = 5 1479b072555Sjeremylt } BPType; 1489b072555Sjeremylt static const char *const bp_types[] = {"bp1","bp2","bp3","bp4","bp5","bp6", 1499b072555Sjeremylt "BPType","CEED_BP",0 150cb32e2e7SValeria Barra }; 151cb32e2e7SValeria Barra 152cb32e2e7SValeria Barra // BP specific data 153cb32e2e7SValeria Barra typedef struct { 1549b072555Sjeremylt CeedInt num_comp_u, q_data_size, q_extra; 1559b072555Sjeremylt CeedQFunctionUser setup_geo, setup_rhs, apply, error; 1569b072555Sjeremylt const char *setup_geo_loc, *setup_rhs_loc, *apply_loc, *error_loc; 1579b072555Sjeremylt CeedEvalMode in_mode, out_mode; 1589b072555Sjeremylt CeedQuadMode q_mode; 1599b072555Sjeremylt } BPData; 160cb32e2e7SValeria Barra 1619b072555Sjeremylt BPData bp_options[6] = { 162cb32e2e7SValeria Barra [CEED_BP1] = { 1639b072555Sjeremylt .num_comp_u = 1, 1649b072555Sjeremylt .q_data_size = 1, 1659b072555Sjeremylt .q_extra = 1, 1669b072555Sjeremylt .setup_geo = SetupMassGeo, 1679b072555Sjeremylt .setup_rhs = SetupMassRhs, 168cb32e2e7SValeria Barra .apply = Mass, 169cb32e2e7SValeria Barra .error = Error, 1709b072555Sjeremylt .setup_geo_loc = SetupMassGeo_loc, 1719b072555Sjeremylt .setup_rhs_loc = SetupMassRhs_loc, 1729b072555Sjeremylt .apply_loc = Mass_loc, 1739b072555Sjeremylt .error_loc = Error_loc, 1749b072555Sjeremylt .in_mode = CEED_EVAL_INTERP, 1759b072555Sjeremylt .out_mode = CEED_EVAL_INTERP, 1769b072555Sjeremylt .q_mode = CEED_GAUSS 177cb32e2e7SValeria Barra }, 178cb32e2e7SValeria Barra [CEED_BP2] = { 1799b072555Sjeremylt .num_comp_u = 3, 1809b072555Sjeremylt .q_data_size = 1, 1819b072555Sjeremylt .q_extra = 1, 1829b072555Sjeremylt .setup_geo = SetupMassGeo, 1839b072555Sjeremylt .setup_rhs = SetupMassRhs3, 184cb32e2e7SValeria Barra .apply = Mass3, 185cb32e2e7SValeria Barra .error = Error3, 1869b072555Sjeremylt .setup_geo_loc = SetupMassGeo_loc, 1879b072555Sjeremylt .setup_rhs_loc = SetupMassRhs3_loc, 1889b072555Sjeremylt .apply_loc = Mass3_loc, 1899b072555Sjeremylt .error_loc = Error3_loc, 1909b072555Sjeremylt .in_mode = CEED_EVAL_INTERP, 1919b072555Sjeremylt .out_mode = CEED_EVAL_INTERP, 1929b072555Sjeremylt .q_mode = CEED_GAUSS 193cb32e2e7SValeria Barra }, 194cb32e2e7SValeria Barra [CEED_BP3] = { 1959b072555Sjeremylt .num_comp_u = 1, 1969b072555Sjeremylt .q_data_size = 7, 1979b072555Sjeremylt .q_extra = 1, 1989b072555Sjeremylt .setup_geo = SetupDiffGeo, 1999b072555Sjeremylt .setup_rhs = SetupDiffRhs, 200cb32e2e7SValeria Barra .apply = Diff, 201cb32e2e7SValeria Barra .error = Error, 2029b072555Sjeremylt .setup_geo_loc = SetupDiffGeo_loc, 2039b072555Sjeremylt .setup_rhs_loc = SetupDiffRhs_loc, 2049b072555Sjeremylt .apply_loc = Diff_loc, 2059b072555Sjeremylt .error_loc = Error_loc, 2069b072555Sjeremylt .in_mode = CEED_EVAL_GRAD, 2079b072555Sjeremylt .out_mode = CEED_EVAL_GRAD, 2089b072555Sjeremylt .q_mode = CEED_GAUSS 209cb32e2e7SValeria Barra }, 210cb32e2e7SValeria Barra [CEED_BP4] = { 2119b072555Sjeremylt .num_comp_u = 3, 2129b072555Sjeremylt .q_data_size = 7, 2139b072555Sjeremylt .q_extra = 1, 2149b072555Sjeremylt .setup_geo = SetupDiffGeo, 2159b072555Sjeremylt .setup_rhs = SetupDiffRhs3, 216cb32e2e7SValeria Barra .apply = Diff3, 217cb32e2e7SValeria Barra .error = Error3, 2189b072555Sjeremylt .setup_geo_loc = SetupDiffGeo_loc, 2199b072555Sjeremylt .setup_rhs_loc = SetupDiffRhs3_loc, 2209b072555Sjeremylt .apply_loc = Diff3_loc, 2219b072555Sjeremylt .error_loc = Error3_loc, 2229b072555Sjeremylt .in_mode = CEED_EVAL_GRAD, 2239b072555Sjeremylt .out_mode = CEED_EVAL_GRAD, 2249b072555Sjeremylt .q_mode = CEED_GAUSS 225cb32e2e7SValeria Barra }, 226cb32e2e7SValeria Barra [CEED_BP5] = { 2279b072555Sjeremylt .num_comp_u = 1, 2289b072555Sjeremylt .q_data_size = 7, 2299b072555Sjeremylt .q_extra = 0, 2309b072555Sjeremylt .setup_geo = SetupDiffGeo, 2319b072555Sjeremylt .setup_rhs = SetupDiffRhs, 232cb32e2e7SValeria Barra .apply = Diff, 233cb32e2e7SValeria Barra .error = Error, 2349b072555Sjeremylt .setup_geo_loc = SetupDiffGeo_loc, 2359b072555Sjeremylt .setup_rhs_loc = SetupDiffRhs_loc, 2369b072555Sjeremylt .apply_loc = Diff_loc, 2379b072555Sjeremylt .error_loc = Error_loc, 2389b072555Sjeremylt .in_mode = CEED_EVAL_GRAD, 2399b072555Sjeremylt .out_mode = CEED_EVAL_GRAD, 2409b072555Sjeremylt .q_mode = CEED_GAUSS_LOBATTO 241cb32e2e7SValeria Barra }, 242cb32e2e7SValeria Barra [CEED_BP6] = { 2439b072555Sjeremylt .num_comp_u = 3, 2449b072555Sjeremylt .q_data_size = 7, 2459b072555Sjeremylt .q_extra = 0, 2469b072555Sjeremylt .setup_geo = SetupDiffGeo, 2479b072555Sjeremylt .setup_rhs = SetupDiffRhs3, 248cb32e2e7SValeria Barra .apply = Diff3, 249cb32e2e7SValeria Barra .error = Error3, 2509b072555Sjeremylt .setup_geo_loc = SetupDiffGeo_loc, 2519b072555Sjeremylt .setup_rhs_loc = SetupDiffRhs3_loc, 2529b072555Sjeremylt .apply_loc = Diff3_loc, 2539b072555Sjeremylt .error_loc = Error3_loc, 2549b072555Sjeremylt .in_mode = CEED_EVAL_GRAD, 2559b072555Sjeremylt .out_mode = CEED_EVAL_GRAD, 2569b072555Sjeremylt .q_mode = CEED_GAUSS_LOBATTO 257cb32e2e7SValeria Barra } 258cb32e2e7SValeria Barra }; 259cb32e2e7SValeria Barra 260cb32e2e7SValeria Barra // This function uses libCEED to compute the action of the mass matrix 261cb32e2e7SValeria Barra static PetscErrorCode MatMult_Mass(Mat A, Vec X, Vec Y) { 262cb32e2e7SValeria Barra PetscErrorCode ierr; 263cb32e2e7SValeria Barra User user; 264cb32e2e7SValeria Barra PetscScalar *x, *y; 2659b072555Sjeremylt PetscMemType x_mem_type, y_mem_type; 266cb32e2e7SValeria Barra 267cb32e2e7SValeria Barra PetscFunctionBeginUser; 2689396343dSjeremylt 269cb32e2e7SValeria Barra ierr = MatShellGetContext(A, &user); CHKERRQ(ierr); 2709396343dSjeremylt 2719396343dSjeremylt // Global-to-local 2729b072555Sjeremylt ierr = VecScatterBegin(user->l_to_g, X, user->X_loc, INSERT_VALUES, 273cb32e2e7SValeria Barra SCATTER_REVERSE); CHKERRQ(ierr); 2749b072555Sjeremylt ierr = VecScatterEnd(user->l_to_g, X, user->X_loc, INSERT_VALUES, 2759396343dSjeremylt SCATTER_REVERSE); CHKERRQ(ierr); 276cb32e2e7SValeria Barra 2779396343dSjeremylt // Setup libCEED vectors 2789b072555Sjeremylt ierr = VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x, 2799b072555Sjeremylt &x_mem_type); CHKERRQ(ierr); 2809b072555Sjeremylt ierr = VecGetArrayAndMemType(user->Y_loc, &y, &y_mem_type); CHKERRQ(ierr); 2819b072555Sjeremylt CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x); 2829b072555Sjeremylt CeedVectorSetArray(user->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y); 283cb32e2e7SValeria Barra 2849396343dSjeremylt // Apply libCEED operator 2859b072555Sjeremylt CeedOperatorApply(user->op, user->x_ceed, user->y_ceed, 286cb32e2e7SValeria Barra CEED_REQUEST_IMMEDIATE); 287cb32e2e7SValeria Barra 2889396343dSjeremylt // Restore PETSc vectors 2899b072555Sjeremylt CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL); 2909b072555Sjeremylt CeedVectorTakeArray(user->y_ceed, MemTypeP2C(y_mem_type), NULL); 2919b072555Sjeremylt ierr = VecRestoreArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x); 2929396343dSjeremylt CHKERRQ(ierr); 2939b072555Sjeremylt ierr = VecRestoreArrayAndMemType(user->Y_loc, &y); CHKERRQ(ierr); 294cb32e2e7SValeria Barra 2959396343dSjeremylt // Local-to-global 296cb32e2e7SValeria Barra if (Y) { 297cb32e2e7SValeria Barra ierr = VecZeroEntries(Y); CHKERRQ(ierr); 2989b072555Sjeremylt ierr = VecScatterBegin(user->l_to_g, user->Y_loc, Y, ADD_VALUES, 2999396343dSjeremylt SCATTER_FORWARD); CHKERRQ(ierr); 3009b072555Sjeremylt ierr = VecScatterEnd(user->l_to_g, user->Y_loc, Y, ADD_VALUES, 3019396343dSjeremylt SCATTER_FORWARD); CHKERRQ(ierr); 302cb32e2e7SValeria Barra } 303cb32e2e7SValeria Barra PetscFunctionReturn(0); 304cb32e2e7SValeria Barra } 305cb32e2e7SValeria Barra 306cb32e2e7SValeria Barra // This function uses libCEED to compute the action of the Laplacian with 307cb32e2e7SValeria Barra // Dirichlet boundary conditions 308cb32e2e7SValeria Barra static PetscErrorCode MatMult_Diff(Mat A, Vec X, Vec Y) { 309cb32e2e7SValeria Barra PetscErrorCode ierr; 310cb32e2e7SValeria Barra User user; 311cb32e2e7SValeria Barra PetscScalar *x, *y; 3129b072555Sjeremylt PetscMemType x_mem_type, y_mem_type; 313cb32e2e7SValeria Barra 314cb32e2e7SValeria Barra PetscFunctionBeginUser; 3159396343dSjeremylt 316cb32e2e7SValeria Barra ierr = MatShellGetContext(A, &user); CHKERRQ(ierr); 317cb32e2e7SValeria Barra 318cb32e2e7SValeria Barra // Global-to-local 3199b072555Sjeremylt ierr = VecScatterBegin(user->l_to_g_0, X, user->X_loc, INSERT_VALUES, 320cb32e2e7SValeria Barra SCATTER_REVERSE); CHKERRQ(ierr); 3219b072555Sjeremylt ierr = VecScatterEnd(user->l_to_g_0, X, user->X_loc, INSERT_VALUES, 322cb32e2e7SValeria Barra SCATTER_REVERSE); 323cb32e2e7SValeria Barra CHKERRQ(ierr); 324cb32e2e7SValeria Barra 3259396343dSjeremylt // Setup libCEED vectors 3269b072555Sjeremylt ierr = VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x, 3279b072555Sjeremylt &x_mem_type); CHKERRQ(ierr); 3289b072555Sjeremylt ierr = VecGetArrayAndMemType(user->Y_loc, &y, &y_mem_type); CHKERRQ(ierr); 3299b072555Sjeremylt CeedVectorSetArray(user->x_ceed, MemTypeP2C(x_mem_type), CEED_USE_POINTER, x); 3309b072555Sjeremylt CeedVectorSetArray(user->y_ceed, MemTypeP2C(y_mem_type), CEED_USE_POINTER, y); 331cb32e2e7SValeria Barra 3329396343dSjeremylt // Apply libCEED operator 3339b072555Sjeremylt CeedOperatorApply(user->op, user->x_ceed, user->y_ceed, 334cb32e2e7SValeria Barra CEED_REQUEST_IMMEDIATE); 335cb32e2e7SValeria Barra 336cb32e2e7SValeria Barra // Restore PETSc vectors 3379b072555Sjeremylt CeedVectorTakeArray(user->x_ceed, MemTypeP2C(x_mem_type), NULL); 3389b072555Sjeremylt CeedVectorTakeArray(user->y_ceed, MemTypeP2C(y_mem_type), NULL); 3399b072555Sjeremylt ierr = VecRestoreArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x); 3409396343dSjeremylt CHKERRQ(ierr); 3419b072555Sjeremylt ierr = VecRestoreArrayAndMemType(user->Y_loc, &y); CHKERRQ(ierr); 342cb32e2e7SValeria Barra 343cb32e2e7SValeria Barra // Local-to-global 344cb32e2e7SValeria Barra ierr = VecZeroEntries(Y); CHKERRQ(ierr); 3459b072555Sjeremylt ierr = VecScatterBegin(user->g_to_g_D, X, Y, INSERT_VALUES, SCATTER_FORWARD); 346cb32e2e7SValeria Barra CHKERRQ(ierr); 3479b072555Sjeremylt ierr = VecScatterEnd(user->g_to_g_D, X, Y, INSERT_VALUES, SCATTER_FORWARD); 348cb32e2e7SValeria Barra CHKERRQ(ierr); 3499b072555Sjeremylt ierr = VecScatterBegin(user->l_to_g_0, user->Y_loc, Y, ADD_VALUES, 3509396343dSjeremylt SCATTER_FORWARD); CHKERRQ(ierr); 3519b072555Sjeremylt ierr = VecScatterEnd(user->l_to_g_0, user->Y_loc, Y, ADD_VALUES, 3529b072555Sjeremylt SCATTER_FORWARD); 353cb32e2e7SValeria Barra CHKERRQ(ierr); 354cb32e2e7SValeria Barra 355cb32e2e7SValeria Barra PetscFunctionReturn(0); 356cb32e2e7SValeria Barra } 357cb32e2e7SValeria Barra 358cb32e2e7SValeria Barra // This function calculates the error in the final solution 3599b072555Sjeremylt static PetscErrorCode ComputeErrorMax(User user, CeedOperator op_error, Vec X, 3609b072555Sjeremylt CeedVector target, PetscReal *max_error) { 361cb32e2e7SValeria Barra PetscErrorCode ierr; 362cb32e2e7SValeria Barra PetscScalar *x; 3639b072555Sjeremylt PetscMemType mem_type; 364cb32e2e7SValeria Barra CeedVector collocated_error; 3651f9221feSJeremy L Thompson CeedSize length; 366cb32e2e7SValeria Barra 367cb32e2e7SValeria Barra PetscFunctionBeginUser; 3689396343dSjeremylt 369cb32e2e7SValeria Barra CeedVectorGetLength(target, &length); 370cb32e2e7SValeria Barra CeedVectorCreate(user->ceed, length, &collocated_error); 371cb32e2e7SValeria Barra 372cb32e2e7SValeria Barra // Global-to-local 3739b072555Sjeremylt ierr = VecScatterBegin(user->l_to_g, X, user->X_loc, INSERT_VALUES, 374cb32e2e7SValeria Barra SCATTER_REVERSE); CHKERRQ(ierr); 3759b072555Sjeremylt ierr = VecScatterEnd(user->l_to_g, X, user->X_loc, INSERT_VALUES, 3769396343dSjeremylt SCATTER_REVERSE); CHKERRQ(ierr); 3779396343dSjeremylt 3789396343dSjeremylt // Setup libCEED vector 3799b072555Sjeremylt ierr = VecGetArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x, 3809b072555Sjeremylt &mem_type); CHKERRQ(ierr); 3819b072555Sjeremylt CeedVectorSetArray(user->x_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, x); 382cb32e2e7SValeria Barra 3839396343dSjeremylt // Apply libCEED operator 3849b072555Sjeremylt CeedOperatorApply(op_error, user->x_ceed, collocated_error, 385cb32e2e7SValeria Barra CEED_REQUEST_IMMEDIATE); 386cb32e2e7SValeria Barra 387cb32e2e7SValeria Barra // Restore PETSc vector 3889b072555Sjeremylt CeedVectorTakeArray(user->x_ceed, MemTypeP2C(mem_type), NULL); 3899b072555Sjeremylt ierr = VecRestoreArrayReadAndMemType(user->X_loc, (const PetscScalar **)&x); 3909396343dSjeremylt CHKERRQ(ierr); 391cb32e2e7SValeria Barra 392cb32e2e7SValeria Barra // Reduce max error 3939b072555Sjeremylt *max_error = 0; 394cb32e2e7SValeria Barra const CeedScalar *e; 395cb32e2e7SValeria Barra CeedVectorGetArrayRead(collocated_error, CEED_MEM_HOST, &e); 396cb32e2e7SValeria Barra for (CeedInt i=0; i<length; i++) { 3979b072555Sjeremylt *max_error = PetscMax(*max_error, PetscAbsScalar(e[i])); 398cb32e2e7SValeria Barra } 399cb32e2e7SValeria Barra CeedVectorRestoreArrayRead(collocated_error, &e); 4009b072555Sjeremylt ierr = MPI_Allreduce(MPI_IN_PLACE, max_error, 1, MPIU_REAL, MPIU_MAX, 4019396343dSjeremylt user->comm); CHKERRQ(ierr); 402cb32e2e7SValeria Barra 403cb32e2e7SValeria Barra // Cleanup 404cb32e2e7SValeria Barra CeedVectorDestroy(&collocated_error); 405cb32e2e7SValeria Barra 406cb32e2e7SValeria Barra PetscFunctionReturn(0); 407cb32e2e7SValeria Barra } 408cb32e2e7SValeria Barra 409cb32e2e7SValeria Barra int main(int argc, char **argv) { 410cb32e2e7SValeria Barra PetscInt ierr; 411cb32e2e7SValeria Barra MPI_Comm comm; 4129b072555Sjeremylt char ceed_resource[PETSC_MAX_PATH_LEN] = "/cpu/self"; 413cb32e2e7SValeria Barra double my_rt_start, my_rt, rt_min, rt_max; 4149b072555Sjeremylt PetscInt degree, q_extra, local_nodes, local_elem, mesh_elem[3], m_nodes[3], 4159b072555Sjeremylt p[3], 4169b072555Sjeremylt i_rank[3], l_nodes[3], l_size, num_comp_u = 1, ksp_max_it_clip[2]; 417cb32e2e7SValeria Barra PetscScalar *r; 418cb32e2e7SValeria Barra PetscBool test_mode, benchmark_mode, write_solution; 419cb32e2e7SValeria Barra PetscMPIInt size, rank; 4209b072555Sjeremylt PetscLogStage solve_stage; 4219b072555Sjeremylt Vec X, X_loc, rhs, rhs_loc; 422cb32e2e7SValeria Barra Mat mat; 423cb32e2e7SValeria Barra KSP ksp; 4249b072555Sjeremylt VecScatter l_to_g, l_to_g_0, g_to_g_D; 4259b072555Sjeremylt PetscMemType mem_type; 426cb32e2e7SValeria Barra User user; 427cb32e2e7SValeria Barra Ceed ceed; 4289b072555Sjeremylt CeedBasis basis_x, basis_u; 4299b072555Sjeremylt CeedElemRestriction elem_restr_x, elem_restr_u, elem_restr_u_i, elem_restr_qd_i; 4309b072555Sjeremylt CeedQFunction qf_setup_geo, qf_setup_rhs, qf_apply, qf_error; 4319b072555Sjeremylt CeedOperator op_setup_geo, op_setup_rhs, op_apply, op_error; 4329b072555Sjeremylt CeedVector x_coord, q_data, rhs_ceed, target; 433cb32e2e7SValeria Barra CeedInt P, Q; 4349b072555Sjeremylt const CeedInt dim = 3, num_comp_x = 3; 4359b072555Sjeremylt BPType bp_choice; 436cb32e2e7SValeria Barra 437cb32e2e7SValeria Barra ierr = PetscInitialize(&argc, &argv, NULL, help); 438cb32e2e7SValeria Barra if (ierr) return ierr; 439cb32e2e7SValeria Barra comm = PETSC_COMM_WORLD; 44032d2ee49SValeria Barra 44132d2ee49SValeria Barra // Read command line options 442*67490bc6SJeremy L Thompson PetscOptionsBegin(comm, NULL, "CEED BPs in PETSc", NULL); 4439b072555Sjeremylt bp_choice = CEED_BP1; 444cb32e2e7SValeria Barra ierr = PetscOptionsEnum("-problem", 445cb32e2e7SValeria Barra "CEED benchmark problem to solve", NULL, 4469b072555Sjeremylt bp_types, (PetscEnum)bp_choice, (PetscEnum *)&bp_choice, 447cb32e2e7SValeria Barra NULL); CHKERRQ(ierr); 4489b072555Sjeremylt num_comp_u = bp_options[bp_choice].num_comp_u; 449cb32e2e7SValeria Barra test_mode = PETSC_FALSE; 450cb32e2e7SValeria Barra ierr = PetscOptionsBool("-test", 451cb32e2e7SValeria Barra "Testing mode (do not print unless error is large)", 452cb32e2e7SValeria Barra NULL, test_mode, &test_mode, NULL); CHKERRQ(ierr); 453cb32e2e7SValeria Barra benchmark_mode = PETSC_FALSE; 454cb32e2e7SValeria Barra ierr = PetscOptionsBool("-benchmark", 455cb32e2e7SValeria Barra "Benchmarking mode (prints benchmark statistics)", 456cb32e2e7SValeria Barra NULL, benchmark_mode, &benchmark_mode, NULL); 457cb32e2e7SValeria Barra CHKERRQ(ierr); 458cb32e2e7SValeria Barra write_solution = PETSC_FALSE; 459cb32e2e7SValeria Barra ierr = PetscOptionsBool("-write_solution", 460cb32e2e7SValeria Barra "Write solution for visualization", 461cb32e2e7SValeria Barra NULL, write_solution, &write_solution, NULL); 462cb32e2e7SValeria Barra CHKERRQ(ierr); 463cb32e2e7SValeria Barra degree = test_mode ? 3 : 1; 464cb32e2e7SValeria Barra ierr = PetscOptionsInt("-degree", "Polynomial degree of tensor product basis", 465cb32e2e7SValeria Barra NULL, degree, °ree, NULL); CHKERRQ(ierr); 4669b072555Sjeremylt q_extra = bp_options[bp_choice].q_extra; 4679b072555Sjeremylt ierr = PetscOptionsInt("-q_extra", "Number of extra quadrature points", 4689b072555Sjeremylt NULL, q_extra, &q_extra, NULL); CHKERRQ(ierr); 469cb32e2e7SValeria Barra ierr = PetscOptionsString("-ceed", "CEED resource specifier", 4709b072555Sjeremylt NULL, ceed_resource, ceed_resource, 4719b072555Sjeremylt sizeof(ceed_resource), NULL); CHKERRQ(ierr); 4729b072555Sjeremylt local_nodes = 1000; 473cb32e2e7SValeria Barra ierr = PetscOptionsInt("-local", 474cb32e2e7SValeria Barra "Target number of locally owned nodes per process", 4759b072555Sjeremylt NULL, local_nodes, &local_nodes, NULL); CHKERRQ(ierr); 4762fbc6e41SJeremy L Thompson PetscInt two = 2; 4772fbc6e41SJeremy L Thompson ksp_max_it_clip[0] = 5; 4782fbc6e41SJeremy L Thompson ksp_max_it_clip[1] = 20; 4792fbc6e41SJeremy L Thompson ierr = PetscOptionsIntArray("-ksp_max_it_clip", 4802fbc6e41SJeremy L Thompson "Min and max number of iterations to use during benchmarking", 4812fbc6e41SJeremy L Thompson NULL, ksp_max_it_clip, &two, NULL); 4822fbc6e41SJeremy L Thompson CHKERRQ(ierr); 483*67490bc6SJeremy L Thompson PetscOptionsEnd(); 484cb32e2e7SValeria Barra P = degree + 1; 4859b072555Sjeremylt Q = P + q_extra; 486cb32e2e7SValeria Barra 4879396343dSjeremylt // Set up libCEED 4889b072555Sjeremylt CeedInit(ceed_resource, &ceed); 4899b072555Sjeremylt CeedMemType mem_type_backend; 4909b072555Sjeremylt CeedGetPreferredMemType(ceed, &mem_type_backend); 4919396343dSjeremylt 4929b072555Sjeremylt VecType default_vec_type = NULL, vec_type; 4939b072555Sjeremylt switch (mem_type_backend) { 4949b072555Sjeremylt case CEED_MEM_HOST: default_vec_type = VECSTANDARD; break; 495b68a8d79SJed Brown case CEED_MEM_DEVICE: { 496b68a8d79SJed Brown const char *resolved; 497b68a8d79SJed Brown CeedGetResource(ceed, &resolved); 4989b072555Sjeremylt if (strstr(resolved, "/gpu/cuda")) default_vec_type = VECCUDA; 499b68a8d79SJed Brown else if (strstr(resolved, "/gpu/hip/occa")) 5009b072555Sjeremylt default_vec_type = VECSTANDARD; // https://github.com/CEED/libCEED/issues/678 5019b072555Sjeremylt else if (strstr(resolved, "/gpu/hip")) default_vec_type = VECHIP; 5029b072555Sjeremylt else default_vec_type = VECSTANDARD; 503b68a8d79SJed Brown } 504b68a8d79SJed Brown } 5059396343dSjeremylt 506cb32e2e7SValeria Barra // Determine size of process grid 507cb32e2e7SValeria Barra ierr = MPI_Comm_size(comm, &size); CHKERRQ(ierr); 508cb32e2e7SValeria Barra Split3(size, p, false); 509cb32e2e7SValeria Barra 5109b072555Sjeremylt // Find a nicely composite number of elements no less than local_nodes 5119b072555Sjeremylt for (local_elem = PetscMax(1, local_nodes / (degree*degree*degree)); ; 5129b072555Sjeremylt local_elem++) { 5139b072555Sjeremylt Split3(local_elem, mesh_elem, true); 5149b072555Sjeremylt if (Max3(mesh_elem) / Min3(mesh_elem) <= 2) break; 515cb32e2e7SValeria Barra } 516cb32e2e7SValeria Barra 517cb32e2e7SValeria Barra // Find my location in the process grid 518cb32e2e7SValeria Barra ierr = MPI_Comm_rank(comm, &rank); CHKERRQ(ierr); 5199b072555Sjeremylt for (int d=0, rank_left=rank; d<dim; d++) { 520cb32e2e7SValeria Barra const int pstride[3] = {p[1] *p[2], p[2], 1}; 5219b072555Sjeremylt i_rank[d] = rank_left / pstride[d]; 5229b072555Sjeremylt rank_left -= i_rank[d] * pstride[d]; 523cb32e2e7SValeria Barra } 524cb32e2e7SValeria Barra 5259b072555Sjeremylt GlobalNodes(p, i_rank, degree, mesh_elem, m_nodes); 526cb32e2e7SValeria Barra 527cb32e2e7SValeria Barra // Setup global vector 5284a5da23aSJed Brown ierr = VecCreate(comm, &X); CHKERRQ(ierr); 5299b072555Sjeremylt ierr = VecSetType(X, default_vec_type); CHKERRQ(ierr); 5309b072555Sjeremylt ierr = VecSetSizes(X, m_nodes[0]*m_nodes[1]*m_nodes[2]*num_comp_u, 5319b072555Sjeremylt PETSC_DECIDE); 532cb32e2e7SValeria Barra CHKERRQ(ierr); 533b68a8d79SJed Brown ierr = VecSetFromOptions(X); CHKERRQ(ierr); 534cb32e2e7SValeria Barra ierr = VecSetUp(X); CHKERRQ(ierr); 535cb32e2e7SValeria Barra 536cb32e2e7SValeria Barra // Set up libCEED 5379b072555Sjeremylt CeedInit(ceed_resource, &ceed); 538cb32e2e7SValeria Barra 539cb32e2e7SValeria Barra // Print summary 540cb32e2e7SValeria Barra CeedInt gsize; 541cb32e2e7SValeria Barra ierr = VecGetSize(X, &gsize); CHKERRQ(ierr); 542dc7d240cSValeria Barra if (!test_mode) { 5439b072555Sjeremylt const char *used_resource; 5449b072555Sjeremylt CeedGetResource(ceed, &used_resource); 5459396343dSjeremylt 5469b072555Sjeremylt ierr = VecGetType(X, &vec_type); CHKERRQ(ierr); 5479396343dSjeremylt 548cb32e2e7SValeria Barra ierr = PetscPrintf(comm, 549cb32e2e7SValeria Barra "\n-- CEED Benchmark Problem %d -- libCEED + PETSc --\n" 5509396343dSjeremylt " PETSc:\n" 5519396343dSjeremylt " PETSc Vec Type : %s\n" 552cb32e2e7SValeria Barra " libCEED:\n" 553cb32e2e7SValeria Barra " libCEED Backend : %s\n" 5549396343dSjeremylt " libCEED Backend MemType : %s\n" 555cb32e2e7SValeria Barra " Mesh:\n" 5568d0bb2bbSvaleriabarra " Number of 1D Basis Nodes (P) : %d\n" 5578d0bb2bbSvaleriabarra " Number of 1D Quadrature Points (Q) : %d\n" 558cb32e2e7SValeria Barra " Global nodes : %D\n" 559cb32e2e7SValeria Barra " Process Decomposition : %D %D %D\n" 560cb32e2e7SValeria Barra " Local Elements : %D = %D %D %D\n" 561db419314Sjeremylt " Owned nodes : %D = %D %D %D\n" 562db419314Sjeremylt " DoF per node : %D\n", 5639b072555Sjeremylt bp_choice+1, vec_type, used_resource, 5649b072555Sjeremylt CeedMemTypes[mem_type_backend], 5659b072555Sjeremylt P, Q, gsize/num_comp_u, p[0], p[1], p[2], local_elem, 5669b072555Sjeremylt mesh_elem[0], mesh_elem[1], mesh_elem[2], 5679b072555Sjeremylt m_nodes[0]*m_nodes[1]*m_nodes[2], m_nodes[0], m_nodes[1], 5689b072555Sjeremylt m_nodes[2], num_comp_u); CHKERRQ(ierr); 569cb32e2e7SValeria Barra } 570cb32e2e7SValeria Barra 571cb32e2e7SValeria Barra { 5729b072555Sjeremylt l_size = 1; 573cb32e2e7SValeria Barra for (int d=0; d<dim; d++) { 5749b072555Sjeremylt l_nodes[d] = mesh_elem[d]*degree + 1; 5759b072555Sjeremylt l_size *= l_nodes[d]; 576cb32e2e7SValeria Barra } 5779b072555Sjeremylt ierr = VecCreate(PETSC_COMM_SELF, &X_loc); CHKERRQ(ierr); 5789b072555Sjeremylt ierr = VecSetType(X_loc, default_vec_type); CHKERRQ(ierr); 5799b072555Sjeremylt ierr = VecSetSizes(X_loc, l_size*num_comp_u, PETSC_DECIDE); CHKERRQ(ierr); 5809b072555Sjeremylt ierr = VecSetFromOptions(X_loc); CHKERRQ(ierr); 5819b072555Sjeremylt ierr = VecSetUp(X_loc); CHKERRQ(ierr); 582cb32e2e7SValeria Barra 583cb32e2e7SValeria Barra // Create local-to-global scatter 5849b072555Sjeremylt PetscInt *l_to_g_ind, *l_to_g_ind_0, *loc_ind, l_0_count; 5859b072555Sjeremylt IS l_to_g_is, l_to_g_is_0, loc_is; 5869b072555Sjeremylt PetscInt g_start[2][2][2], g_m_nodes[2][2][2][dim]; 587cb32e2e7SValeria Barra 588cb32e2e7SValeria Barra for (int i=0; i<2; i++) { 589cb32e2e7SValeria Barra for (int j=0; j<2; j++) { 590cb32e2e7SValeria Barra for (int k=0; k<2; k++) { 5919b072555Sjeremylt PetscInt ijk_rank[3] = {i_rank[0]+i, i_rank[1]+j, i_rank[2]+k}; 5929b072555Sjeremylt g_start[i][j][k] = GlobalStart(p, ijk_rank, degree, mesh_elem); 5939b072555Sjeremylt GlobalNodes(p, ijk_rank, degree, mesh_elem, g_m_nodes[i][j][k]); 594cb32e2e7SValeria Barra } 595cb32e2e7SValeria Barra } 596cb32e2e7SValeria Barra } 597cb32e2e7SValeria Barra 5989b072555Sjeremylt ierr = PetscMalloc1(l_size, &l_to_g_ind); CHKERRQ(ierr); 5999b072555Sjeremylt ierr = PetscMalloc1(l_size, &l_to_g_ind_0); CHKERRQ(ierr); 6009b072555Sjeremylt ierr = PetscMalloc1(l_size, &loc_ind); CHKERRQ(ierr); 6019b072555Sjeremylt l_0_count = 0; 6029b072555Sjeremylt for (PetscInt i=0,ir,ii; ir=i>=m_nodes[0], ii=i-ir*m_nodes[0], i<l_nodes[0]; 6039b072555Sjeremylt i++) 6049b072555Sjeremylt for (PetscInt j=0,jr,jj; jr=j>=m_nodes[1], jj=j-jr*m_nodes[1], j<l_nodes[1]; 6059b072555Sjeremylt j++) 6069b072555Sjeremylt for (PetscInt k=0,kr,kk; kr=k>=m_nodes[2], kk=k-kr*m_nodes[2], k<l_nodes[2]; 6079b072555Sjeremylt k++) { 6089b072555Sjeremylt PetscInt here = (i*l_nodes[1]+j)*l_nodes[2]+k; 6099b072555Sjeremylt l_to_g_ind[here] = 6109b072555Sjeremylt g_start[ir][jr][kr] + (ii*g_m_nodes[ir][jr][kr][1]+jj)*g_m_nodes[ir][jr][kr][2] 6119b072555Sjeremylt +kk; 6129b072555Sjeremylt if ((i_rank[0] == 0 && i == 0) 6139b072555Sjeremylt || (i_rank[1] == 0 && j == 0) 6149b072555Sjeremylt || (i_rank[2] == 0 && k == 0) 6159b072555Sjeremylt || (i_rank[0]+1 == p[0] && i+1 == l_nodes[0]) 6169b072555Sjeremylt || (i_rank[1]+1 == p[1] && j+1 == l_nodes[1]) 6179b072555Sjeremylt || (i_rank[2]+1 == p[2] && k+1 == l_nodes[2])) 618cb32e2e7SValeria Barra continue; 6199b072555Sjeremylt l_to_g_ind_0[l_0_count] = l_to_g_ind[here]; 6209b072555Sjeremylt loc_ind[l_0_count++] = here; 621cb32e2e7SValeria Barra } 6229b072555Sjeremylt ierr = ISCreateBlock(comm, num_comp_u, l_size, l_to_g_ind, PETSC_OWN_POINTER, 6239b072555Sjeremylt &l_to_g_is); CHKERRQ(ierr); 6249b072555Sjeremylt ierr = VecScatterCreate(X_loc, NULL, X, l_to_g_is, &l_to_g); CHKERRQ(ierr); 625cb32e2e7SValeria Barra CHKERRQ(ierr); 6269b072555Sjeremylt ierr = ISCreateBlock(comm, num_comp_u, l_0_count, l_to_g_ind_0, 6279b072555Sjeremylt PETSC_OWN_POINTER, 6289b072555Sjeremylt &l_to_g_is_0); CHKERRQ(ierr); 6299b072555Sjeremylt ierr = ISCreateBlock(comm, num_comp_u, l_0_count, loc_ind, PETSC_OWN_POINTER, 6309b072555Sjeremylt &loc_is); CHKERRQ(ierr); 6319b072555Sjeremylt ierr = VecScatterCreate(X_loc, loc_is, X, l_to_g_is_0, &l_to_g_0); 6329b072555Sjeremylt CHKERRQ(ierr); 633cb32e2e7SValeria Barra { 634cb32e2e7SValeria Barra // Create global-to-global scatter for Dirichlet values (everything not in 6359b072555Sjeremylt // l_to_g_is_0, which is the range of l_to_g_0) 6369b072555Sjeremylt PetscInt x_start, x_end, *ind_D, count_D = 0; 6379b072555Sjeremylt IS is_D; 638cb32e2e7SValeria Barra const PetscScalar *x; 6399b072555Sjeremylt ierr = VecZeroEntries(X_loc); CHKERRQ(ierr); 640cb32e2e7SValeria Barra ierr = VecSet(X, 1.0); CHKERRQ(ierr); 6419b072555Sjeremylt ierr = VecScatterBegin(l_to_g_0, X_loc, X, INSERT_VALUES, SCATTER_FORWARD); 642cb32e2e7SValeria Barra CHKERRQ(ierr); 6439b072555Sjeremylt ierr = VecScatterEnd(l_to_g_0, X_loc, X, INSERT_VALUES, SCATTER_FORWARD); 644cb32e2e7SValeria Barra CHKERRQ(ierr); 6459b072555Sjeremylt ierr = VecGetOwnershipRange(X, &x_start, &x_end); CHKERRQ(ierr); 6469b072555Sjeremylt ierr = PetscMalloc1(x_end-x_start, &ind_D); CHKERRQ(ierr); 647cb32e2e7SValeria Barra ierr = VecGetArrayRead(X, &x); CHKERRQ(ierr); 6489b072555Sjeremylt for (PetscInt i=0; i<x_end-x_start; i++) { 6499b072555Sjeremylt if (x[i] == 1.) ind_D[count_D++] = x_start + i; 650cb32e2e7SValeria Barra } 651cb32e2e7SValeria Barra ierr = VecRestoreArrayRead(X, &x); CHKERRQ(ierr); 6529b072555Sjeremylt ierr = ISCreateGeneral(comm, count_D, ind_D, PETSC_COPY_VALUES, &is_D); 653cb32e2e7SValeria Barra CHKERRQ(ierr); 6549b072555Sjeremylt ierr = PetscFree(ind_D); CHKERRQ(ierr); 6559b072555Sjeremylt ierr = VecScatterCreate(X, is_D, X, is_D, &g_to_g_D); CHKERRQ(ierr); 6569b072555Sjeremylt ierr = ISDestroy(&is_D); CHKERRQ(ierr); 657cb32e2e7SValeria Barra } 6589b072555Sjeremylt ierr = ISDestroy(&l_to_g_is); CHKERRQ(ierr); 6599b072555Sjeremylt ierr = ISDestroy(&l_to_g_is_0); CHKERRQ(ierr); 6609b072555Sjeremylt ierr = ISDestroy(&loc_is); CHKERRQ(ierr); 661cb32e2e7SValeria Barra } 662cb32e2e7SValeria Barra 663cb32e2e7SValeria Barra // CEED bases 6649b072555Sjeremylt CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_u, P, Q, 6659b072555Sjeremylt bp_options[bp_choice].q_mode, &basis_u); 6669b072555Sjeremylt CeedBasisCreateTensorH1Lagrange(ceed, dim, num_comp_x, 2, Q, 6679b072555Sjeremylt bp_options[bp_choice].q_mode, &basis_x); 668cb32e2e7SValeria Barra 669cb32e2e7SValeria Barra // CEED restrictions 6709b072555Sjeremylt CreateRestriction(ceed, mesh_elem, P, num_comp_u, &elem_restr_u); 6719b072555Sjeremylt CreateRestriction(ceed, mesh_elem, 2, dim, &elem_restr_x); 6729b072555Sjeremylt CeedInt num_elem = mesh_elem[0]*mesh_elem[1]*mesh_elem[2]; 6739b072555Sjeremylt CeedElemRestrictionCreateStrided(ceed, num_elem, Q*Q*Q, num_comp_u, 6749b072555Sjeremylt num_comp_u*num_elem*Q*Q*Q, 6759b072555Sjeremylt CEED_STRIDES_BACKEND, &elem_restr_u_i); 6769b072555Sjeremylt CeedElemRestrictionCreateStrided(ceed, num_elem, Q*Q*Q, 6779b072555Sjeremylt bp_options[bp_choice].q_data_size, 6789b072555Sjeremylt bp_options[bp_choice].q_data_size*num_elem*Q*Q*Q, 6799b072555Sjeremylt CEED_STRIDES_BACKEND, &elem_restr_qd_i); 680cb32e2e7SValeria Barra { 6819b072555Sjeremylt CeedScalar *x_loc; 6829b072555Sjeremylt CeedInt shape[3] = {mesh_elem[0]+1, mesh_elem[1]+1, mesh_elem[2]+1}, len = 683cb32e2e7SValeria Barra shape[0]*shape[1]*shape[2]; 6849b072555Sjeremylt x_loc = malloc(len*num_comp_x*sizeof x_loc[0]); 685cb32e2e7SValeria Barra for (CeedInt i=0; i<shape[0]; i++) { 686cb32e2e7SValeria Barra for (CeedInt j=0; j<shape[1]; j++) { 687cb32e2e7SValeria Barra for (CeedInt k=0; k<shape[2]; k++) { 6889b072555Sjeremylt x_loc[dim*((i*shape[1]+j)*shape[2]+k) + 0] = 1.*(i_rank[0]*mesh_elem[0]+i) / 6899b072555Sjeremylt (p[0]*mesh_elem[0]); 6909b072555Sjeremylt x_loc[dim*((i*shape[1]+j)*shape[2]+k) + 1] = 1.*(i_rank[1]*mesh_elem[1]+j) / 6919b072555Sjeremylt (p[1]*mesh_elem[1]); 6929b072555Sjeremylt x_loc[dim*((i*shape[1]+j)*shape[2]+k) + 2] = 1.*(i_rank[2]*mesh_elem[2]+k) / 6939b072555Sjeremylt (p[2]*mesh_elem[2]); 694cb32e2e7SValeria Barra } 695cb32e2e7SValeria Barra } 696cb32e2e7SValeria Barra } 6979b072555Sjeremylt CeedVectorCreate(ceed, len*num_comp_x, &x_coord); 6989b072555Sjeremylt CeedVectorSetArray(x_coord, CEED_MEM_HOST, CEED_OWN_POINTER, x_loc); 699cb32e2e7SValeria Barra } 700cb32e2e7SValeria Barra 7019b072555Sjeremylt // Create the QFunction that builds the operator quadrature data 7029b072555Sjeremylt CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].setup_geo, 7039b072555Sjeremylt bp_options[bp_choice].setup_geo_loc, &qf_setup_geo); 7049b072555Sjeremylt CeedQFunctionAddInput(qf_setup_geo, "x", num_comp_x, CEED_EVAL_INTERP); 7059b072555Sjeremylt CeedQFunctionAddInput(qf_setup_geo, "dx", num_comp_x*dim, CEED_EVAL_GRAD); 7069b072555Sjeremylt CeedQFunctionAddInput(qf_setup_geo, "weight", 1, CEED_EVAL_WEIGHT); 7079b072555Sjeremylt CeedQFunctionAddOutput(qf_setup_geo, "q_data", 7089b072555Sjeremylt bp_options[bp_choice].q_data_size, 709cb32e2e7SValeria Barra CEED_EVAL_NONE); 710cb32e2e7SValeria Barra 7119b072555Sjeremylt // Create the QFunction that sets up the RHS and true solution 7129b072555Sjeremylt CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].setup_rhs, 7139b072555Sjeremylt bp_options[bp_choice].setup_rhs_loc, &qf_setup_rhs); 7149b072555Sjeremylt CeedQFunctionAddInput(qf_setup_rhs, "x", num_comp_x, CEED_EVAL_INTERP); 7159b072555Sjeremylt CeedQFunctionAddInput(qf_setup_rhs, "q_data", bp_options[bp_choice].q_data_size, 716e83e87a5Sjeremylt CEED_EVAL_NONE); 7179b072555Sjeremylt CeedQFunctionAddOutput(qf_setup_rhs, "true_soln", num_comp_u, CEED_EVAL_NONE); 7189b072555Sjeremylt CeedQFunctionAddOutput(qf_setup_rhs, "rhs", num_comp_u, CEED_EVAL_INTERP); 719cb32e2e7SValeria Barra 720cb32e2e7SValeria Barra // Set up PDE operator 7219b072555Sjeremylt CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].apply, 7229b072555Sjeremylt bp_options[bp_choice].apply_loc, &qf_apply); 723cb32e2e7SValeria Barra // Add inputs and outputs 7249b072555Sjeremylt CeedInt in_scale = bp_options[bp_choice].in_mode==CEED_EVAL_GRAD ? 3 : 1; 7259b072555Sjeremylt CeedInt out_scale = bp_options[bp_choice].out_mode==CEED_EVAL_GRAD ? 3 : 1; 7269b072555Sjeremylt CeedQFunctionAddInput(qf_apply, "u", num_comp_u*in_scale, 7279b072555Sjeremylt bp_options[bp_choice].in_mode); 7289b072555Sjeremylt CeedQFunctionAddInput(qf_apply, "q_data", bp_options[bp_choice].q_data_size, 729cb32e2e7SValeria Barra CEED_EVAL_NONE); 7309b072555Sjeremylt CeedQFunctionAddOutput(qf_apply, "v", num_comp_u*out_scale, 7319b072555Sjeremylt bp_options[bp_choice].out_mode); 732cb32e2e7SValeria Barra 733cb32e2e7SValeria Barra // Create the error qfunction 7349b072555Sjeremylt CeedQFunctionCreateInterior(ceed, 1, bp_options[bp_choice].error, 7359b072555Sjeremylt bp_options[bp_choice].error_loc, &qf_error); 7369b072555Sjeremylt CeedQFunctionAddInput(qf_error, "u", num_comp_u, CEED_EVAL_INTERP); 7379b072555Sjeremylt CeedQFunctionAddInput(qf_error, "true_soln", num_comp_u, CEED_EVAL_NONE); 7389b072555Sjeremylt CeedQFunctionAddOutput(qf_error, "error", num_comp_u, CEED_EVAL_NONE); 739cb32e2e7SValeria Barra 740cb32e2e7SValeria Barra // Create the persistent vectors that will be needed in setup 7419b072555Sjeremylt CeedInt num_qpts; 7429b072555Sjeremylt CeedBasisGetNumQuadraturePoints(basis_u, &num_qpts); 7439b072555Sjeremylt CeedVectorCreate(ceed, bp_options[bp_choice].q_data_size*num_elem*num_qpts, 7449b072555Sjeremylt &q_data); 7459b072555Sjeremylt CeedVectorCreate(ceed, num_elem*num_qpts*num_comp_u, &target); 7469b072555Sjeremylt CeedVectorCreate(ceed, l_size*num_comp_u, &rhs_ceed); 747cb32e2e7SValeria Barra 748cb32e2e7SValeria Barra // Create the operator that builds the quadrature data for the ceed operator 7499b072555Sjeremylt CeedOperatorCreate(ceed, qf_setup_geo, CEED_QFUNCTION_NONE, 7509b072555Sjeremylt CEED_QFUNCTION_NONE, &op_setup_geo); 7519b072555Sjeremylt CeedOperatorSetField(op_setup_geo, "x", elem_restr_x, basis_x, 752e83e87a5Sjeremylt CEED_VECTOR_ACTIVE); 7539b072555Sjeremylt CeedOperatorSetField(op_setup_geo, "dx", elem_restr_x, basis_x, 754a8d32208Sjeremylt CEED_VECTOR_ACTIVE); 7559b072555Sjeremylt CeedOperatorSetField(op_setup_geo, "weight", CEED_ELEMRESTRICTION_NONE, basis_x, 756a8d32208Sjeremylt CEED_VECTOR_NONE); 7579b072555Sjeremylt CeedOperatorSetField(op_setup_geo, "q_data", elem_restr_qd_i, 758cb32e2e7SValeria Barra CEED_BASIS_COLLOCATED, CEED_VECTOR_ACTIVE); 759cb32e2e7SValeria Barra 760cb32e2e7SValeria Barra // Create the operator that builds the RHS and true solution 7619b072555Sjeremylt CeedOperatorCreate(ceed, qf_setup_rhs, CEED_QFUNCTION_NONE, 7629b072555Sjeremylt CEED_QFUNCTION_NONE, &op_setup_rhs); 7639b072555Sjeremylt CeedOperatorSetField(op_setup_rhs, "x", elem_restr_x, basis_x, 764a8d32208Sjeremylt CEED_VECTOR_ACTIVE); 7659b072555Sjeremylt CeedOperatorSetField(op_setup_rhs, "q_data", elem_restr_qd_i, 7669b072555Sjeremylt CEED_BASIS_COLLOCATED, 7679b072555Sjeremylt q_data); 7689b072555Sjeremylt CeedOperatorSetField(op_setup_rhs, "true_soln", elem_restr_u_i, 769cb32e2e7SValeria Barra CEED_BASIS_COLLOCATED, target); 7709b072555Sjeremylt CeedOperatorSetField(op_setup_rhs, "rhs", elem_restr_u, basis_u, 771a8d32208Sjeremylt CEED_VECTOR_ACTIVE); 772cb32e2e7SValeria Barra 773cb32e2e7SValeria Barra // Create the mass or diff operator 7749b072555Sjeremylt CeedOperatorCreate(ceed, qf_apply, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 7759b072555Sjeremylt &op_apply); 7769b072555Sjeremylt CeedOperatorSetField(op_apply, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE); 7779b072555Sjeremylt CeedOperatorSetField(op_apply, "q_data", elem_restr_qd_i, CEED_BASIS_COLLOCATED, 7789b072555Sjeremylt q_data); 7799b072555Sjeremylt CeedOperatorSetField(op_apply, "v", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE); 780cb32e2e7SValeria Barra 781cb32e2e7SValeria Barra // Create the error operator 7829b072555Sjeremylt CeedOperatorCreate(ceed, qf_error, CEED_QFUNCTION_NONE, CEED_QFUNCTION_NONE, 7839b072555Sjeremylt &op_error); 7849b072555Sjeremylt CeedOperatorSetField(op_error, "u", elem_restr_u, basis_u, CEED_VECTOR_ACTIVE); 7859b072555Sjeremylt CeedOperatorSetField(op_error, "true_soln", elem_restr_u_i, 786cb32e2e7SValeria Barra CEED_BASIS_COLLOCATED, target); 7879b072555Sjeremylt CeedOperatorSetField(op_error, "error", elem_restr_u_i, CEED_BASIS_COLLOCATED, 788a8d32208Sjeremylt CEED_VECTOR_ACTIVE); 789cb32e2e7SValeria Barra 790cb32e2e7SValeria Barra // Set up Mat 791cb32e2e7SValeria Barra ierr = PetscMalloc1(1, &user); CHKERRQ(ierr); 792cb32e2e7SValeria Barra user->comm = comm; 7939b072555Sjeremylt user->l_to_g = l_to_g; 7949b072555Sjeremylt if (bp_choice != CEED_BP1 && bp_choice != CEED_BP2) { 7959b072555Sjeremylt user->l_to_g_0 = l_to_g_0; 7969b072555Sjeremylt user->g_to_g_D = g_to_g_D; 797cb32e2e7SValeria Barra } 7989b072555Sjeremylt user->X_loc = X_loc; 7999b072555Sjeremylt ierr = VecDuplicate(X_loc, &user->Y_loc); CHKERRQ(ierr); 8009b072555Sjeremylt CeedVectorCreate(ceed, l_size*num_comp_u, &user->x_ceed); 8019b072555Sjeremylt CeedVectorCreate(ceed, l_size*num_comp_u, &user->y_ceed); 8029b072555Sjeremylt user->op = op_apply; 8039b072555Sjeremylt user->q_data = q_data; 804cb32e2e7SValeria Barra user->ceed = ceed; 805cb32e2e7SValeria Barra 8069b072555Sjeremylt ierr = MatCreateShell(comm, m_nodes[0]*m_nodes[1]*m_nodes[2]*num_comp_u, 8079b072555Sjeremylt m_nodes[0]*m_nodes[1]*m_nodes[2]*num_comp_u, 808cb32e2e7SValeria Barra PETSC_DECIDE, PETSC_DECIDE, user, &mat); CHKERRQ(ierr); 8099b072555Sjeremylt if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) { 810cb32e2e7SValeria Barra ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Mass); 811cb32e2e7SValeria Barra CHKERRQ(ierr); 812cb32e2e7SValeria Barra } else { 813cb32e2e7SValeria Barra ierr = MatShellSetOperation(mat, MATOP_MULT, (void(*)(void))MatMult_Diff); 814cb32e2e7SValeria Barra CHKERRQ(ierr); 815cb32e2e7SValeria Barra } 8169b072555Sjeremylt ierr = VecGetType(X, &vec_type); CHKERRQ(ierr); 8179b072555Sjeremylt ierr = MatShellSetVecType(mat, vec_type); CHKERRQ(ierr); 818cb32e2e7SValeria Barra 819cb32e2e7SValeria Barra // Get RHS vector 8209396343dSjeremylt ierr = VecDuplicate(X, &rhs); CHKERRQ(ierr); 8219b072555Sjeremylt ierr = VecDuplicate(X_loc, &rhs_loc); CHKERRQ(ierr); 8229b072555Sjeremylt ierr = VecZeroEntries(rhs_loc); CHKERRQ(ierr); 8239b072555Sjeremylt ierr = VecGetArrayAndMemType(rhs_loc, &r, &mem_type); CHKERRQ(ierr); 8249b072555Sjeremylt CeedVectorSetArray(rhs_ceed, MemTypeP2C(mem_type), CEED_USE_POINTER, r); 825cb32e2e7SValeria Barra 8269b072555Sjeremylt // Setup q_data, rhs, and target 8279b072555Sjeremylt CeedOperatorApply(op_setup_geo, x_coord, q_data, CEED_REQUEST_IMMEDIATE); 8289b072555Sjeremylt CeedOperatorApply(op_setup_rhs, x_coord, rhs_ceed, CEED_REQUEST_IMMEDIATE); 8299b072555Sjeremylt CeedVectorDestroy(&x_coord); 830cb32e2e7SValeria Barra 831cb32e2e7SValeria Barra // Gather RHS 8329b072555Sjeremylt ierr = CeedVectorTakeArray(rhs_ceed, MemTypeP2C(mem_type), NULL); CHKERRQ(ierr); 8339b072555Sjeremylt ierr = VecRestoreArrayAndMemType(rhs_loc, &r); CHKERRQ(ierr); 834cb32e2e7SValeria Barra ierr = VecZeroEntries(rhs); CHKERRQ(ierr); 8359b072555Sjeremylt ierr = VecScatterBegin(l_to_g, rhs_loc, rhs, ADD_VALUES, SCATTER_FORWARD); 836cb32e2e7SValeria Barra CHKERRQ(ierr); 8379b072555Sjeremylt ierr = VecScatterEnd(l_to_g, rhs_loc, rhs, ADD_VALUES, SCATTER_FORWARD); 838cb32e2e7SValeria Barra CHKERRQ(ierr); 8399b072555Sjeremylt CeedVectorDestroy(&rhs_ceed); 840cb32e2e7SValeria Barra 841cb32e2e7SValeria Barra ierr = KSPCreate(comm, &ksp); CHKERRQ(ierr); 842cb32e2e7SValeria Barra { 843cb32e2e7SValeria Barra PC pc; 844cb32e2e7SValeria Barra ierr = KSPGetPC(ksp, &pc); CHKERRQ(ierr); 8459b072555Sjeremylt if (bp_choice == CEED_BP1 || bp_choice == CEED_BP2) { 846cb32e2e7SValeria Barra ierr = PCSetType(pc, PCJACOBI); CHKERRQ(ierr); 847cb32e2e7SValeria Barra ierr = PCJacobiSetType(pc, PC_JACOBI_ROWSUM); CHKERRQ(ierr); 848cb32e2e7SValeria Barra } else { 849cb32e2e7SValeria Barra ierr = PCSetType(pc, PCNONE); CHKERRQ(ierr); 850cb32e2e7SValeria Barra } 851cb32e2e7SValeria Barra ierr = KSPSetType(ksp, KSPCG); CHKERRQ(ierr); 852cb32e2e7SValeria Barra ierr = KSPSetNormType(ksp, KSP_NORM_NATURAL); CHKERRQ(ierr); 853cb32e2e7SValeria Barra ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 854cb32e2e7SValeria Barra PETSC_DEFAULT); CHKERRQ(ierr); 855cb32e2e7SValeria Barra } 856cb32e2e7SValeria Barra ierr = KSPSetOperators(ksp, mat, mat); CHKERRQ(ierr); 857da9108adSvaleriabarra // First run's performance log is not considered for benchmarking purposes 858cb32e2e7SValeria Barra ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 1); 859cb32e2e7SValeria Barra CHKERRQ(ierr); 860cb32e2e7SValeria Barra my_rt_start = MPI_Wtime(); 861cb32e2e7SValeria Barra ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 862cb32e2e7SValeria Barra my_rt = MPI_Wtime() - my_rt_start; 863cb32e2e7SValeria Barra ierr = MPI_Allreduce(MPI_IN_PLACE, &my_rt, 1, MPI_DOUBLE, MPI_MIN, comm); 864cb32e2e7SValeria Barra CHKERRQ(ierr); 865cb32e2e7SValeria Barra // Set maxits based on first iteration timing 866cb32e2e7SValeria Barra if (my_rt > 0.02) { 8672fbc6e41SJeremy L Thompson ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 8682fbc6e41SJeremy L Thompson ksp_max_it_clip[0]); CHKERRQ(ierr); 869cb32e2e7SValeria Barra } else { 8702fbc6e41SJeremy L Thompson ierr = KSPSetTolerances(ksp, 1e-10, PETSC_DEFAULT, PETSC_DEFAULT, 8712fbc6e41SJeremy L Thompson ksp_max_it_clip[1]); CHKERRQ(ierr); 872cb32e2e7SValeria Barra } 873debcf919SJed Brown ierr = KSPSetFromOptions(ksp); CHKERRQ(ierr); 87409a940d7Sjeremylt 875cb32e2e7SValeria Barra // Timed solve 87609a940d7Sjeremylt ierr = VecZeroEntries(X); CHKERRQ(ierr); 877cb32e2e7SValeria Barra ierr = PetscBarrier((PetscObject)ksp); CHKERRQ(ierr); 87809a940d7Sjeremylt 87909a940d7Sjeremylt // -- Performance logging 8809b072555Sjeremylt ierr = PetscLogStageRegister("Solve Stage", &solve_stage); CHKERRQ(ierr); 8819b072555Sjeremylt ierr = PetscLogStagePush(solve_stage); CHKERRQ(ierr); 88209a940d7Sjeremylt 88309a940d7Sjeremylt // -- Solve 884cb32e2e7SValeria Barra my_rt_start = MPI_Wtime(); 885cb32e2e7SValeria Barra ierr = KSPSolve(ksp, rhs, X); CHKERRQ(ierr); 886cb32e2e7SValeria Barra my_rt = MPI_Wtime() - my_rt_start; 88709a940d7Sjeremylt 88809a940d7Sjeremylt // -- Performance logging 88909a940d7Sjeremylt ierr = PetscLogStagePop(); 89009a940d7Sjeremylt 8918e87e98bSJed Brown // Output results 892cb32e2e7SValeria Barra { 8939b072555Sjeremylt KSPType ksp_type; 894cb32e2e7SValeria Barra KSPConvergedReason reason; 895cb32e2e7SValeria Barra PetscReal rnorm; 896cb32e2e7SValeria Barra PetscInt its; 8979b072555Sjeremylt ierr = KSPGetType(ksp, &ksp_type); CHKERRQ(ierr); 898cb32e2e7SValeria Barra ierr = KSPGetConvergedReason(ksp, &reason); CHKERRQ(ierr); 899cb32e2e7SValeria Barra ierr = KSPGetIterationNumber(ksp, &its); CHKERRQ(ierr); 900cb32e2e7SValeria Barra ierr = KSPGetResidualNorm(ksp, &rnorm); CHKERRQ(ierr); 901a22c4fb5SJed Brown if (!test_mode || reason < 0 || rnorm > 1e-8) { 902cb32e2e7SValeria Barra ierr = PetscPrintf(comm, 903cb32e2e7SValeria Barra " KSP:\n" 904cb32e2e7SValeria Barra " KSP Type : %s\n" 905cb32e2e7SValeria Barra " KSP Convergence : %s\n" 906cb32e2e7SValeria Barra " Total KSP Iterations : %D\n" 907cb32e2e7SValeria Barra " Final rnorm : %e\n", 9089b072555Sjeremylt ksp_type, KSPConvergedReasons[reason], its, 909cb32e2e7SValeria Barra (double)rnorm); CHKERRQ(ierr); 910cb32e2e7SValeria Barra } 9118e87e98bSJed Brown if (!test_mode) { 9128e87e98bSJed Brown ierr = PetscPrintf(comm," Performance:\n"); CHKERRQ(ierr); 9138e87e98bSJed Brown } 9148e87e98bSJed Brown { 9159b072555Sjeremylt PetscReal max_error; 9169b072555Sjeremylt ierr = ComputeErrorMax(user, op_error, X, target, &max_error); 9178e87e98bSJed Brown CHKERRQ(ierr); 9188e87e98bSJed Brown PetscReal tol = 5e-2; 9199b072555Sjeremylt if (!test_mode || max_error > tol) { 920cb32e2e7SValeria Barra ierr = MPI_Allreduce(&my_rt, &rt_min, 1, MPI_DOUBLE, MPI_MIN, comm); 921cb32e2e7SValeria Barra CHKERRQ(ierr); 922cb32e2e7SValeria Barra ierr = MPI_Allreduce(&my_rt, &rt_max, 1, MPI_DOUBLE, MPI_MAX, comm); 923cb32e2e7SValeria Barra CHKERRQ(ierr); 924cb32e2e7SValeria Barra ierr = PetscPrintf(comm, 9258e87e98bSJed Brown " Pointwise Error (max) : %e\n" 9268e87e98bSJed Brown " CG Solve Time : %g (%g) sec\n", 9279b072555Sjeremylt (double)max_error, rt_max, rt_min); CHKERRQ(ierr); 928cb32e2e7SValeria Barra } 929cb32e2e7SValeria Barra } 9308d0bb2bbSvaleriabarra if (!test_mode) { 931cb32e2e7SValeria Barra ierr = PetscPrintf(comm, 9328e87e98bSJed Brown " DoFs/Sec in CG : %g (%g) million\n", 9338e87e98bSJed Brown 1e-6*gsize*its/rt_max, 9348e87e98bSJed Brown 1e-6*gsize*its/rt_min); CHKERRQ(ierr); 935cb32e2e7SValeria Barra } 936cb32e2e7SValeria Barra } 937cb32e2e7SValeria Barra 938cb32e2e7SValeria Barra if (write_solution) { 9399b072555Sjeremylt PetscViewer vtk_viewer_soln; 940cb32e2e7SValeria Barra 9419b072555Sjeremylt ierr = PetscViewerCreate(comm, &vtk_viewer_soln); CHKERRQ(ierr); 9429b072555Sjeremylt ierr = PetscViewerSetType(vtk_viewer_soln, PETSCVIEWERVTK); CHKERRQ(ierr); 9439b072555Sjeremylt ierr = PetscViewerFileSetName(vtk_viewer_soln, "solution.vtu"); CHKERRQ(ierr); 9449b072555Sjeremylt ierr = VecView(X, vtk_viewer_soln); CHKERRQ(ierr); 9459b072555Sjeremylt ierr = PetscViewerDestroy(&vtk_viewer_soln); CHKERRQ(ierr); 946cb32e2e7SValeria Barra } 947cb32e2e7SValeria Barra 948cb32e2e7SValeria Barra ierr = VecDestroy(&rhs); CHKERRQ(ierr); 9499b072555Sjeremylt ierr = VecDestroy(&rhs_loc); CHKERRQ(ierr); 950cb32e2e7SValeria Barra ierr = VecDestroy(&X); CHKERRQ(ierr); 9519b072555Sjeremylt ierr = VecDestroy(&user->X_loc); CHKERRQ(ierr); 9529b072555Sjeremylt ierr = VecDestroy(&user->Y_loc); CHKERRQ(ierr); 9539b072555Sjeremylt ierr = VecScatterDestroy(&l_to_g); CHKERRQ(ierr); 9549b072555Sjeremylt ierr = VecScatterDestroy(&l_to_g_0); CHKERRQ(ierr); 9559b072555Sjeremylt ierr = VecScatterDestroy(&g_to_g_D); CHKERRQ(ierr); 956cb32e2e7SValeria Barra ierr = MatDestroy(&mat); CHKERRQ(ierr); 957cb32e2e7SValeria Barra ierr = KSPDestroy(&ksp); CHKERRQ(ierr); 958cb32e2e7SValeria Barra 9599b072555Sjeremylt CeedVectorDestroy(&user->x_ceed); 9609b072555Sjeremylt CeedVectorDestroy(&user->y_ceed); 9619b072555Sjeremylt CeedVectorDestroy(&user->q_data); 962cb32e2e7SValeria Barra CeedVectorDestroy(&target); 9639b072555Sjeremylt CeedOperatorDestroy(&op_setup_geo); 9649b072555Sjeremylt CeedOperatorDestroy(&op_setup_rhs); 9659b072555Sjeremylt CeedOperatorDestroy(&op_apply); 9669b072555Sjeremylt CeedOperatorDestroy(&op_error); 9679b072555Sjeremylt CeedElemRestrictionDestroy(&elem_restr_u); 9689b072555Sjeremylt CeedElemRestrictionDestroy(&elem_restr_x); 9699b072555Sjeremylt CeedElemRestrictionDestroy(&elem_restr_u_i); 9709b072555Sjeremylt CeedElemRestrictionDestroy(&elem_restr_qd_i); 9719b072555Sjeremylt CeedQFunctionDestroy(&qf_setup_geo); 9729b072555Sjeremylt CeedQFunctionDestroy(&qf_setup_rhs); 9739b072555Sjeremylt CeedQFunctionDestroy(&qf_apply); 9749b072555Sjeremylt CeedQFunctionDestroy(&qf_error); 9759b072555Sjeremylt CeedBasisDestroy(&basis_u); 9769b072555Sjeremylt CeedBasisDestroy(&basis_x); 977cb32e2e7SValeria Barra CeedDestroy(&ceed); 978cb32e2e7SValeria Barra ierr = PetscFree(user); CHKERRQ(ierr); 979cb32e2e7SValeria Barra return PetscFinalize(); 980cb32e2e7SValeria Barra } 981