1c4762a1bSJed Brown static char help[] = "Benchmark for MatMatMult() of AIJ matrices using different 2d finite-difference stencils.\n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown /* Converts 3d grid coordinates (i,j,k) for a grid of size m \times n to global indexing. Pass k = 0 for a 2d grid. */ 6c4762a1bSJed Brown int global_index(PetscInt i,PetscInt j,PetscInt k, PetscInt m, PetscInt n) { return i + j * m + k * m * n; } 7c4762a1bSJed Brown 8c4762a1bSJed Brown int main(int argc,char **argv) 9c4762a1bSJed Brown { 10c4762a1bSJed Brown Mat A,B,C,PtAP,PtAP_copy,PtAP_squared; 11c4762a1bSJed Brown PetscInt i,M,N,Istart,Iend,n=7,j,J,Ii,m=8,k,o=1; 12c4762a1bSJed Brown PetscScalar v; 13c4762a1bSJed Brown PetscBool equal=PETSC_FALSE,mat_view=PETSC_FALSE; 14c4762a1bSJed Brown char stencil[PETSC_MAX_PATH_LEN]; 15c4762a1bSJed Brown #if defined(PETSC_USE_LOG) 16c4762a1bSJed Brown PetscLogStage fullMatMatMultStage; 17c4762a1bSJed Brown #endif 18c4762a1bSJed Brown 19*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 20*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 21*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 22*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-o",&o,NULL)); 23*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-result_view",&mat_view)); 24*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-stencil",stencil,sizeof(stencil),NULL)); 25c4762a1bSJed Brown 26c4762a1bSJed Brown /* Create a aij matrix A */ 27c4762a1bSJed Brown M = N = m*n*o; 28*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 29*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N)); 30*9566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATAIJ)); 31*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Consistency checks */ 342c71b3e2SJacob Faibussowitsch PetscCheckFalse(o < 1 || m < 1 || n < 1,PETSC_COMM_WORLD,PETSC_ERR_USER,"Dimensions need to be larger than zero!"); 35c4762a1bSJed Brown 36c4762a1bSJed Brown /************ 2D stencils ***************/ 37*9566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(stencil, "2d5point", &equal)); 38c4762a1bSJed Brown if (equal) { /* 5-point stencil, 2D */ 39*9566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL)); 40*9566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,5,NULL)); 41*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); 42c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 43c4762a1bSJed Brown v = -1.0; k = Ii / (m*n); j = (Ii - k * m * n) / m; i = (Ii - k * m * n - j * m); 44*9566063dSJacob Faibussowitsch if (i>0) {J = global_index(i-1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 45*9566063dSJacob Faibussowitsch if (i<m-1) {J = global_index(i+1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 46*9566063dSJacob Faibussowitsch if (j>0) {J = global_index(i,j-1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 47*9566063dSJacob Faibussowitsch if (j<n-1) {J = global_index(i,j+1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 48*9566063dSJacob Faibussowitsch v = 4.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown } 51*9566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(stencil, "2d9point", &equal)); 52c4762a1bSJed Brown if (equal) { /* 9-point stencil, 2D */ 53*9566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,9,NULL,9,NULL)); 54*9566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,9,NULL)); 55*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); 56c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 57c4762a1bSJed Brown v = -1.0; k = Ii / (m*n); j = (Ii - k * m * n) / m; i = (Ii - k * m * n - j * m); 58*9566063dSJacob Faibussowitsch if (i>0) {J = global_index(i-1,j, k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 59*9566063dSJacob Faibussowitsch if (i>0 && j>0) {J = global_index(i-1,j-1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 60*9566063dSJacob Faibussowitsch if (j>0) {J = global_index(i, j-1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 61*9566063dSJacob Faibussowitsch if (i<m-1 && j>0) {J = global_index(i+1,j-1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 62*9566063dSJacob Faibussowitsch if (i<m-1) {J = global_index(i+1,j, k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 63*9566063dSJacob Faibussowitsch if (i<m-1 && j<n-1) {J = global_index(i+1,j+1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 64*9566063dSJacob Faibussowitsch if (j<n-1) {J = global_index(i, j+1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 65*9566063dSJacob Faibussowitsch if (i>0 && j<n-1) {J = global_index(i-1,j+1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 66*9566063dSJacob Faibussowitsch v = 8.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 67c4762a1bSJed Brown } 68c4762a1bSJed Brown } 69*9566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(stencil, "2d9point2", &equal)); 70c4762a1bSJed Brown if (equal) { /* 9-point Cartesian stencil (width 2 per coordinate), 2D */ 71*9566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,9,NULL,9,NULL)); 72*9566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,9,NULL)); 73*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); 74c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 75c4762a1bSJed Brown v = -1.0; k = Ii / (m*n); j = (Ii - k * m * n) / m; i = (Ii - k * m * n - j * m); 76*9566063dSJacob Faibussowitsch if (i>0) {J = global_index(i-1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 77*9566063dSJacob Faibussowitsch if (i>1) {J = global_index(i-2,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 78*9566063dSJacob Faibussowitsch if (i<m-1) {J = global_index(i+1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 79*9566063dSJacob Faibussowitsch if (i<m-2) {J = global_index(i+2,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 80*9566063dSJacob Faibussowitsch if (j>0) {J = global_index(i,j-1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 81*9566063dSJacob Faibussowitsch if (j>1) {J = global_index(i,j-2,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 82*9566063dSJacob Faibussowitsch if (j<n-1) {J = global_index(i,j+1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 83*9566063dSJacob Faibussowitsch if (j<n-2) {J = global_index(i,j+2,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 84*9566063dSJacob Faibussowitsch v = 8.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 85c4762a1bSJed Brown } 86c4762a1bSJed Brown } 87*9566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(stencil, "2d13point", &equal)); 88c4762a1bSJed Brown if (equal) { /* 13-point Cartesian stencil (width 3 per coordinate), 2D */ 89*9566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,13,NULL,13,NULL)); 90*9566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,13,NULL)); 91*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); 92c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 93c4762a1bSJed Brown v = -1.0; k = Ii / (m*n); j = (Ii - k * m * n) / m; i = (Ii - k * m * n - j * m); 94*9566063dSJacob Faibussowitsch if (i>0) {J = global_index(i-1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 95*9566063dSJacob Faibussowitsch if (i>1) {J = global_index(i-2,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 96*9566063dSJacob Faibussowitsch if (i>2) {J = global_index(i-3,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 97*9566063dSJacob Faibussowitsch if (i<m-1) {J = global_index(i+1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 98*9566063dSJacob Faibussowitsch if (i<m-2) {J = global_index(i+2,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 99*9566063dSJacob Faibussowitsch if (i<m-3) {J = global_index(i+3,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 100*9566063dSJacob Faibussowitsch if (j>0) {J = global_index(i,j-1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 101*9566063dSJacob Faibussowitsch if (j>1) {J = global_index(i,j-2,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 102*9566063dSJacob Faibussowitsch if (j>2) {J = global_index(i,j-3,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 103*9566063dSJacob Faibussowitsch if (j<n-1) {J = global_index(i,j+1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 104*9566063dSJacob Faibussowitsch if (j<n-2) {J = global_index(i,j+2,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 105*9566063dSJacob Faibussowitsch if (j<n-3) {J = global_index(i,j+3,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 106*9566063dSJacob Faibussowitsch v = 12.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 107c4762a1bSJed Brown } 108c4762a1bSJed Brown } 109c4762a1bSJed Brown /************ 3D stencils ***************/ 110*9566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(stencil, "3d7point", &equal)); 111c4762a1bSJed Brown if (equal) { /* 7-point stencil, 3D */ 112*9566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,7,NULL,7,NULL)); 113*9566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,7,NULL)); 114*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); 115c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 116c4762a1bSJed Brown v = -1.0; k = Ii / (m*n); j = (Ii - k * m * n) / m; i = (Ii - k * m * n - j * m); 117*9566063dSJacob Faibussowitsch if (i>0) {J = global_index(i-1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 118*9566063dSJacob Faibussowitsch if (i<m-1) {J = global_index(i+1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 119*9566063dSJacob Faibussowitsch if (j>0) {J = global_index(i,j-1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 120*9566063dSJacob Faibussowitsch if (j<n-1) {J = global_index(i,j+1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 121*9566063dSJacob Faibussowitsch if (k>0) {J = global_index(i,j,k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 122*9566063dSJacob Faibussowitsch if (k<o-1) {J = global_index(i,j,k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 123*9566063dSJacob Faibussowitsch v = 6.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 124c4762a1bSJed Brown } 125c4762a1bSJed Brown } 126*9566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(stencil, "3d13point", &equal)); 127c4762a1bSJed Brown if (equal) { /* 13-point stencil, 3D */ 128*9566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,13,NULL,13,NULL)); 129*9566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,13,NULL)); 130*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); 131c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 132c4762a1bSJed Brown v = -1.0; k = Ii / (m*n); j = (Ii - k * m * n) / m; i = (Ii - k * m * n - j * m); 133*9566063dSJacob Faibussowitsch if (i>0) {J = global_index(i-1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 134*9566063dSJacob Faibussowitsch if (i>1) {J = global_index(i-2,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 135*9566063dSJacob Faibussowitsch if (i<m-1) {J = global_index(i+1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 136*9566063dSJacob Faibussowitsch if (i<m-2) {J = global_index(i+2,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 137*9566063dSJacob Faibussowitsch if (j>0) {J = global_index(i,j-1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 138*9566063dSJacob Faibussowitsch if (j>1) {J = global_index(i,j-2,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 139*9566063dSJacob Faibussowitsch if (j<n-1) {J = global_index(i,j+1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 140*9566063dSJacob Faibussowitsch if (j<n-2) {J = global_index(i,j+2,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 141*9566063dSJacob Faibussowitsch if (k>0) {J = global_index(i,j,k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 142*9566063dSJacob Faibussowitsch if (k>1) {J = global_index(i,j,k-2,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 143*9566063dSJacob Faibussowitsch if (k<o-1) {J = global_index(i,j,k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 144*9566063dSJacob Faibussowitsch if (k<o-2) {J = global_index(i,j,k+2,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 145*9566063dSJacob Faibussowitsch v = 12.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 146c4762a1bSJed Brown } 147c4762a1bSJed Brown } 148*9566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(stencil, "3d19point", &equal)); 149c4762a1bSJed Brown if (equal) { /* 19-point stencil, 3D */ 150*9566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,19,NULL,19,NULL)); 151*9566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,19,NULL)); 152*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); 153c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 154c4762a1bSJed Brown v = -1.0; k = Ii / (m*n); j = (Ii - k * m * n) / m; i = (Ii - k * m * n - j * m); 155c4762a1bSJed Brown /* one hop */ 156*9566063dSJacob Faibussowitsch if (i>0) {J = global_index(i-1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 157*9566063dSJacob Faibussowitsch if (i<m-1) {J = global_index(i+1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 158*9566063dSJacob Faibussowitsch if (j>0) {J = global_index(i,j-1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 159*9566063dSJacob Faibussowitsch if (j<n-1) {J = global_index(i,j+1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 160*9566063dSJacob Faibussowitsch if (k>0) {J = global_index(i,j,k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 161*9566063dSJacob Faibussowitsch if (k<o-1) {J = global_index(i,j,k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 162c4762a1bSJed Brown /* two hops */ 163*9566063dSJacob Faibussowitsch if (i>0 && j>0) {J = global_index(i-1,j-1,k ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 164*9566063dSJacob Faibussowitsch if (i>0 && k>0) {J = global_index(i-1,j, k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 165*9566063dSJacob Faibussowitsch if (i>0 && j<n-1) {J = global_index(i-1,j+1,k ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 166*9566063dSJacob Faibussowitsch if (i>0 && k<o-1) {J = global_index(i-1,j, k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 167*9566063dSJacob Faibussowitsch if (i<m-1 && j>0) {J = global_index(i+1,j-1,k ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 168*9566063dSJacob Faibussowitsch if (i<m-1 && k>0) {J = global_index(i+1,j, k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 169*9566063dSJacob Faibussowitsch if (i<m-1 && j<n-1) {J = global_index(i+1,j+1,k ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 170*9566063dSJacob Faibussowitsch if (i<m-1 && k<o-1) {J = global_index(i+1,j, k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 171*9566063dSJacob Faibussowitsch if (j>0 && k>0) {J = global_index(i, j-1,k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 172*9566063dSJacob Faibussowitsch if (j>0 && k<o-1) {J = global_index(i, j-1,k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 173*9566063dSJacob Faibussowitsch if (j<n-1 && k>0) {J = global_index(i, j+1,k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 174*9566063dSJacob Faibussowitsch if (j<n-1 && k<o-1) {J = global_index(i, j+1,k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 175*9566063dSJacob Faibussowitsch v = 18.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 176c4762a1bSJed Brown } 177c4762a1bSJed Brown } 178*9566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(stencil, "3d27point", &equal)); 179c4762a1bSJed Brown if (equal) { /* 27-point stencil, 3D */ 180*9566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(A,27,NULL,27,NULL)); 181*9566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(A,27,NULL)); 182*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&Istart,&Iend)); 183c4762a1bSJed Brown for (Ii=Istart; Ii<Iend; Ii++) { 184c4762a1bSJed Brown v = -1.0; k = Ii / (m*n); j = (Ii - k * m * n) / m; i = (Ii - k * m * n - j * m); 185c4762a1bSJed Brown if (k>0) { 186c4762a1bSJed Brown if (j>0) { 187*9566063dSJacob Faibussowitsch if (i>0) {J = global_index(i-1,j-1,k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 188*9566063dSJacob Faibussowitsch J = global_index(i, j-1,k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES)); 189*9566063dSJacob Faibussowitsch if (i<m-1) {J = global_index(i+1,j-1,k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 190c4762a1bSJed Brown } 191c4762a1bSJed Brown { 192*9566063dSJacob Faibussowitsch if (i>0) {J = global_index(i-1,j, k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 193*9566063dSJacob Faibussowitsch J = global_index(i, j, k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES)); 194*9566063dSJacob Faibussowitsch if (i<m-1) {J = global_index(i+1,j, k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 195c4762a1bSJed Brown } 196c4762a1bSJed Brown if (j<n-1) { 197*9566063dSJacob Faibussowitsch if (i>0) {J = global_index(i-1,j+1,k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 198*9566063dSJacob Faibussowitsch J = global_index(i, j+1,k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES)); 199*9566063dSJacob Faibussowitsch if (i<m-1) {J = global_index(i+1,j+1,k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 200c4762a1bSJed Brown } 201c4762a1bSJed Brown } 202c4762a1bSJed Brown { 203c4762a1bSJed Brown if (j>0) { 204*9566063dSJacob Faibussowitsch if (i>0) {J = global_index(i-1,j-1,k ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 205*9566063dSJacob Faibussowitsch J = global_index(i, j-1,k ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES)); 206*9566063dSJacob Faibussowitsch if (i<m-1) {J = global_index(i+1,j-1,k ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 207c4762a1bSJed Brown } 208c4762a1bSJed Brown { 209*9566063dSJacob Faibussowitsch if (i>0) {J = global_index(i-1,j, k ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 210*9566063dSJacob Faibussowitsch J = global_index(i, j, k ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES)); 211*9566063dSJacob Faibussowitsch if (i<m-1) {J = global_index(i+1,j, k ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 212c4762a1bSJed Brown } 213c4762a1bSJed Brown if (j<n-1) { 214*9566063dSJacob Faibussowitsch if (i>0) {J = global_index(i-1,j+1,k ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 215*9566063dSJacob Faibussowitsch J = global_index(i, j+1,k ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES)); 216*9566063dSJacob Faibussowitsch if (i<m-1) {J = global_index(i+1,j+1,k ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 217c4762a1bSJed Brown } 218c4762a1bSJed Brown } 219c4762a1bSJed Brown if (k<o-1) { 220c4762a1bSJed Brown if (j>0) { 221*9566063dSJacob Faibussowitsch if (i>0) {J = global_index(i-1,j-1,k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 222*9566063dSJacob Faibussowitsch J = global_index(i, j-1,k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES)); 223*9566063dSJacob Faibussowitsch if (i<m-1) {J = global_index(i+1,j-1,k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 224c4762a1bSJed Brown } 225c4762a1bSJed Brown { 226*9566063dSJacob Faibussowitsch if (i>0) {J = global_index(i-1,j, k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 227*9566063dSJacob Faibussowitsch J = global_index(i, j, k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES)); 228*9566063dSJacob Faibussowitsch if (i<m-1) {J = global_index(i+1,j, k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 229c4762a1bSJed Brown } 230c4762a1bSJed Brown if (j<n-1) { 231*9566063dSJacob Faibussowitsch if (i>0) {J = global_index(i-1,j+1,k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 232*9566063dSJacob Faibussowitsch J = global_index(i, j+1,k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES)); 233*9566063dSJacob Faibussowitsch if (i<m-1) {J = global_index(i+1,j+1,k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 234c4762a1bSJed Brown } 235c4762a1bSJed Brown } 236*9566063dSJacob Faibussowitsch v = 26.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 237c4762a1bSJed Brown } 238c4762a1bSJed Brown } 239*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 240*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* Copy A into B in order to have a more representative benchmark (A*A has more cache hits than A*B) */ 243*9566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&B)); 244c4762a1bSJed Brown 245*9566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Full MatMatMult",&fullMatMatMultStage)); 246c4762a1bSJed Brown 247c4762a1bSJed Brown /* Test C = A*B */ 248*9566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(fullMatMatMultStage)); 249*9566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); 250c4762a1bSJed Brown 251c4762a1bSJed Brown /* Test PtAP_squared = PtAP(C,C)*PtAP(C,C) */ 252*9566063dSJacob Faibussowitsch PetscCall(MatPtAP(C,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&PtAP)); 253*9566063dSJacob Faibussowitsch PetscCall(MatDuplicate(PtAP,MAT_COPY_VALUES,&PtAP_copy)); 254*9566063dSJacob Faibussowitsch PetscCall(MatMatMult(PtAP,PtAP_copy,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&PtAP_squared)); 255c4762a1bSJed Brown 256*9566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 257*9566063dSJacob Faibussowitsch PetscCall(MatView(PtAP_squared,PETSC_VIEWER_STDOUT_WORLD)); 258c4762a1bSJed Brown 259*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&PtAP_squared)); 260*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&PtAP_copy)); 261*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&PtAP)); 262*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 263*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 264*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 265*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 266b122ec5aSJacob Faibussowitsch return 0; 267c4762a1bSJed Brown } 268c4762a1bSJed Brown 269c4762a1bSJed Brown /*TEST 270c4762a1bSJed Brown 271c4762a1bSJed Brown test: 272c4762a1bSJed Brown suffix: 1 273c4762a1bSJed Brown nsize: 1 274c20d7725SJed Brown args: -m 8 -n 8 -stencil 2d5point -matmatmult_via sorted 275c4762a1bSJed Brown 276c4762a1bSJed Brown test: 277c4762a1bSJed Brown suffix: 2 278c4762a1bSJed Brown nsize: 1 279c4762a1bSJed Brown args: -m 5 -n 5 -o 5 -stencil 3d27point -matmatmult_via rowmerge 280c4762a1bSJed Brown 281c4762a1bSJed Brown test: 282c4762a1bSJed Brown suffix: 3 283c4762a1bSJed Brown nsize: 4 284c4762a1bSJed Brown args: -m 6 -n 6 -stencil 2d5point -matmatmult_via seqmpi 285c4762a1bSJed Brown 286c4762a1bSJed Brown TEST*/ 287