xref: /petsc/src/mat/tests/ex226.c (revision bc30f867be46b9f027b3139a042f37469b892957)
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 
199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-o",&o,NULL));
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-result_view",&mat_view));
249566063dSJacob 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;
289566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
299566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,M,N));
309566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATAIJ));
319566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   /* Consistency checks */
34*bc30f867SBarry Smith   PetscCheck(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 ***************/
379566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(stencil, "2d5point", &equal));
38c4762a1bSJed Brown   if (equal) {   /* 5-point stencil, 2D */
399566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(A,5,NULL,5,NULL));
409566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(A,5,NULL));
419566063dSJacob 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);
449566063dSJacob Faibussowitsch       if (i>0)   {J = global_index(i-1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
459566063dSJacob Faibussowitsch       if (i<m-1) {J = global_index(i+1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
469566063dSJacob Faibussowitsch       if (j>0)   {J = global_index(i,j-1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
479566063dSJacob Faibussowitsch       if (j<n-1) {J = global_index(i,j+1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
489566063dSJacob Faibussowitsch       v = 4.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
49c4762a1bSJed Brown     }
50c4762a1bSJed Brown   }
519566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(stencil, "2d9point", &equal));
52c4762a1bSJed Brown   if (equal) {      /* 9-point stencil, 2D */
539566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(A,9,NULL,9,NULL));
549566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(A,9,NULL));
559566063dSJacob 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);
589566063dSJacob Faibussowitsch       if (i>0)            {J = global_index(i-1,j,  k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
599566063dSJacob 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));}
609566063dSJacob Faibussowitsch       if (j>0)            {J = global_index(i,  j-1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
619566063dSJacob 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));}
629566063dSJacob Faibussowitsch       if (i<m-1)          {J = global_index(i+1,j,  k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
639566063dSJacob 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));}
649566063dSJacob Faibussowitsch       if (j<n-1)          {J = global_index(i,  j+1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
659566063dSJacob 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));}
669566063dSJacob Faibussowitsch       v = 8.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
67c4762a1bSJed Brown     }
68c4762a1bSJed Brown   }
699566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(stencil, "2d9point2", &equal));
70c4762a1bSJed Brown   if (equal) {      /* 9-point Cartesian stencil (width 2 per coordinate), 2D */
719566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(A,9,NULL,9,NULL));
729566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(A,9,NULL));
739566063dSJacob 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);
769566063dSJacob Faibussowitsch       if (i>0)   {J = global_index(i-1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
779566063dSJacob Faibussowitsch       if (i>1)   {J = global_index(i-2,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
789566063dSJacob Faibussowitsch       if (i<m-1) {J = global_index(i+1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
799566063dSJacob Faibussowitsch       if (i<m-2) {J = global_index(i+2,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
809566063dSJacob Faibussowitsch       if (j>0)   {J = global_index(i,j-1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
819566063dSJacob Faibussowitsch       if (j>1)   {J = global_index(i,j-2,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
829566063dSJacob Faibussowitsch       if (j<n-1) {J = global_index(i,j+1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
839566063dSJacob Faibussowitsch       if (j<n-2) {J = global_index(i,j+2,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
849566063dSJacob Faibussowitsch       v = 8.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
85c4762a1bSJed Brown     }
86c4762a1bSJed Brown   }
879566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(stencil, "2d13point", &equal));
88c4762a1bSJed Brown   if (equal) {      /* 13-point Cartesian stencil (width 3 per coordinate), 2D */
899566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(A,13,NULL,13,NULL));
909566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(A,13,NULL));
919566063dSJacob 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);
949566063dSJacob Faibussowitsch       if (i>0)   {J = global_index(i-1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
959566063dSJacob Faibussowitsch       if (i>1)   {J = global_index(i-2,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
969566063dSJacob Faibussowitsch       if (i>2)   {J = global_index(i-3,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
979566063dSJacob Faibussowitsch       if (i<m-1) {J = global_index(i+1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
989566063dSJacob Faibussowitsch       if (i<m-2) {J = global_index(i+2,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
999566063dSJacob Faibussowitsch       if (i<m-3) {J = global_index(i+3,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1009566063dSJacob Faibussowitsch       if (j>0)   {J = global_index(i,j-1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1019566063dSJacob Faibussowitsch       if (j>1)   {J = global_index(i,j-2,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1029566063dSJacob Faibussowitsch       if (j>2)   {J = global_index(i,j-3,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1039566063dSJacob Faibussowitsch       if (j<n-1) {J = global_index(i,j+1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1049566063dSJacob Faibussowitsch       if (j<n-2) {J = global_index(i,j+2,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1059566063dSJacob Faibussowitsch       if (j<n-3) {J = global_index(i,j+3,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1069566063dSJacob Faibussowitsch       v = 12.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
107c4762a1bSJed Brown     }
108c4762a1bSJed Brown   }
109c4762a1bSJed Brown   /************ 3D stencils ***************/
1109566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(stencil, "3d7point", &equal));
111c4762a1bSJed Brown   if (equal) {      /* 7-point stencil, 3D */
1129566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(A,7,NULL,7,NULL));
1139566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(A,7,NULL));
1149566063dSJacob 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);
1179566063dSJacob Faibussowitsch       if (i>0)   {J = global_index(i-1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1189566063dSJacob Faibussowitsch       if (i<m-1) {J = global_index(i+1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1199566063dSJacob Faibussowitsch       if (j>0)   {J = global_index(i,j-1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1209566063dSJacob Faibussowitsch       if (j<n-1) {J = global_index(i,j+1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1219566063dSJacob Faibussowitsch       if (k>0)   {J = global_index(i,j,k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1229566063dSJacob Faibussowitsch       if (k<o-1) {J = global_index(i,j,k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1239566063dSJacob Faibussowitsch       v = 6.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
124c4762a1bSJed Brown     }
125c4762a1bSJed Brown   }
1269566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(stencil, "3d13point", &equal));
127c4762a1bSJed Brown   if (equal) {      /* 13-point stencil, 3D */
1289566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(A,13,NULL,13,NULL));
1299566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(A,13,NULL));
1309566063dSJacob 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);
1339566063dSJacob Faibussowitsch       if (i>0)   {J = global_index(i-1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1349566063dSJacob Faibussowitsch       if (i>1)   {J = global_index(i-2,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1359566063dSJacob Faibussowitsch       if (i<m-1) {J = global_index(i+1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1369566063dSJacob Faibussowitsch       if (i<m-2) {J = global_index(i+2,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1379566063dSJacob Faibussowitsch       if (j>0)   {J = global_index(i,j-1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1389566063dSJacob Faibussowitsch       if (j>1)   {J = global_index(i,j-2,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1399566063dSJacob Faibussowitsch       if (j<n-1) {J = global_index(i,j+1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1409566063dSJacob Faibussowitsch       if (j<n-2) {J = global_index(i,j+2,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1419566063dSJacob Faibussowitsch       if (k>0)   {J = global_index(i,j,k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1429566063dSJacob Faibussowitsch       if (k>1)   {J = global_index(i,j,k-2,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1439566063dSJacob Faibussowitsch       if (k<o-1) {J = global_index(i,j,k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1449566063dSJacob Faibussowitsch       if (k<o-2) {J = global_index(i,j,k+2,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1459566063dSJacob Faibussowitsch       v = 12.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
146c4762a1bSJed Brown     }
147c4762a1bSJed Brown   }
1489566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(stencil, "3d19point", &equal));
149c4762a1bSJed Brown   if (equal) {      /* 19-point stencil, 3D */
1509566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(A,19,NULL,19,NULL));
1519566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(A,19,NULL));
1529566063dSJacob 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 */
1569566063dSJacob Faibussowitsch       if (i>0)   {J = global_index(i-1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1579566063dSJacob Faibussowitsch       if (i<m-1) {J = global_index(i+1,j,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1589566063dSJacob Faibussowitsch       if (j>0)   {J = global_index(i,j-1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1599566063dSJacob Faibussowitsch       if (j<n-1) {J = global_index(i,j+1,k,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1609566063dSJacob Faibussowitsch       if (k>0)   {J = global_index(i,j,k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1619566063dSJacob 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 */
1639566063dSJacob 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));}
1649566063dSJacob 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));}
1659566063dSJacob 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));}
1669566063dSJacob 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));}
1679566063dSJacob 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));}
1689566063dSJacob 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));}
1699566063dSJacob 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));}
1709566063dSJacob 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));}
1719566063dSJacob 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));}
1729566063dSJacob 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));}
1739566063dSJacob 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));}
1749566063dSJacob 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));}
1759566063dSJacob Faibussowitsch       v = 18.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
176c4762a1bSJed Brown     }
177c4762a1bSJed Brown   }
1789566063dSJacob Faibussowitsch   PetscCall(PetscStrcmp(stencil, "3d27point", &equal));
179c4762a1bSJed Brown   if (equal) {      /* 27-point stencil, 3D */
1809566063dSJacob Faibussowitsch     PetscCall(MatMPIAIJSetPreallocation(A,27,NULL,27,NULL));
1819566063dSJacob Faibussowitsch     PetscCall(MatSeqAIJSetPreallocation(A,27,NULL));
1829566063dSJacob 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) {
1879566063dSJacob 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));}
1889566063dSJacob Faibussowitsch                       J = global_index(i,  j-1,k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
1899566063dSJacob 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         {
1929566063dSJacob Faibussowitsch           if (i>0)   {J = global_index(i-1,j,  k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
1939566063dSJacob Faibussowitsch                       J = global_index(i,  j,  k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
1949566063dSJacob 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) {
1979566063dSJacob 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));}
1989566063dSJacob Faibussowitsch                       J = global_index(i,  j+1,k-1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
1999566063dSJacob 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) {
2049566063dSJacob Faibussowitsch           if (i>0)   {J = global_index(i-1,j-1,k  ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
2059566063dSJacob Faibussowitsch                       J = global_index(i,  j-1,k  ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
2069566063dSJacob 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         {
2099566063dSJacob Faibussowitsch           if (i>0)   {J = global_index(i-1,j,  k  ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
2109566063dSJacob Faibussowitsch                       J = global_index(i,  j,  k  ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
2119566063dSJacob 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) {
2149566063dSJacob Faibussowitsch           if (i>0)   {J = global_index(i-1,j+1,k  ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
2159566063dSJacob Faibussowitsch                       J = global_index(i,  j+1,k  ,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
2169566063dSJacob 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) {
2219566063dSJacob 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));}
2229566063dSJacob Faibussowitsch                       J = global_index(i,  j-1,k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
2239566063dSJacob 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         {
2269566063dSJacob Faibussowitsch           if (i>0)   {J = global_index(i-1,j,  k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
2279566063dSJacob Faibussowitsch                       J = global_index(i,  j,  k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
2289566063dSJacob 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) {
2319566063dSJacob 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));}
2329566063dSJacob Faibussowitsch                       J = global_index(i,  j+1,k+1,m,n); PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));
2339566063dSJacob 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       }
2369566063dSJacob Faibussowitsch       v = 26.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
237c4762a1bSJed Brown     }
238c4762a1bSJed Brown   }
2399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
2409566063dSJacob 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) */
2439566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&B));
244c4762a1bSJed Brown 
2459566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("Full MatMatMult",&fullMatMatMultStage));
246c4762a1bSJed Brown 
247c4762a1bSJed Brown   /* Test C = A*B */
2489566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(fullMatMatMultStage));
2499566063dSJacob Faibussowitsch   PetscCall(MatMatMult(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C));
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   /* Test PtAP_squared = PtAP(C,C)*PtAP(C,C)  */
2529566063dSJacob Faibussowitsch   PetscCall(MatPtAP(C,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&PtAP));
2539566063dSJacob Faibussowitsch   PetscCall(MatDuplicate(PtAP,MAT_COPY_VALUES,&PtAP_copy));
2549566063dSJacob Faibussowitsch   PetscCall(MatMatMult(PtAP,PtAP_copy,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&PtAP_squared));
255c4762a1bSJed Brown 
2569566063dSJacob Faibussowitsch   PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD));
2579566063dSJacob Faibussowitsch   PetscCall(MatView(PtAP_squared,PETSC_VIEWER_STDOUT_WORLD));
258c4762a1bSJed Brown 
2599566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&PtAP_squared));
2609566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&PtAP_copy));
2619566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&PtAP));
2629566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
2639566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&B));
2649566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2659566063dSJacob 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