xref: /petsc/src/mat/tests/ex86.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
1 static char help[] = "Testing MatCreateMPIMatConcatenateSeqMat().\n\n";
2 
3 #include <petscmat.h>
4 int main(int argc,char **argv)
5 {
6   Mat            seqmat,mpimat;
7   PetscMPIInt    rank;
8   PetscScalar    value[3],*vals;
9   PetscInt       i,col[3],n=5,bs=1;
10 
11   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
12   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
13   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
14 
15   /* Create seqaij matrices of size (n+rank) by n */
16   PetscCall(MatCreate(PETSC_COMM_SELF,&seqmat));
17   PetscCall(MatSetSizes(seqmat,(n+rank)*bs,PETSC_DECIDE,PETSC_DECIDE,n*bs));
18   PetscCall(MatSetFromOptions(seqmat));
19   PetscCall(MatSeqAIJSetPreallocation(seqmat,3*bs,NULL));
20   PetscCall(MatSeqBAIJSetPreallocation(seqmat,bs,3,NULL));
21 
22   if (bs == 1) {
23     value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
24     for (i=1; i<n-1; i++) {
25       col[0] = i-1; col[1] = i; col[2] = i+1;
26       PetscCall(MatSetValues(seqmat,1,&i,3,col,value,INSERT_VALUES));
27     }
28     i = n - 1; col[0] = n - 2; col[1] = n - 1;
29     PetscCall(MatSetValues(seqmat,1,&i,2,col,value,INSERT_VALUES));
30 
31     i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
32     PetscCall(MatSetValues(seqmat,1,&i,2,col,value,INSERT_VALUES));
33   } else {
34     PetscInt *rows,*cols,j;
35     PetscCall(PetscMalloc3(bs*bs,&vals,bs,&rows,bs,&cols));
36     /* diagonal blocks */
37     for (i=0; i<bs*bs; i++) vals[i] = 2.0;
38     for (i=0; i<n*bs; i+=bs) {
39       for (j=0; j<bs; j++) {rows[j] = i+j; cols[j] = i+j;}
40       PetscCall(MatSetValues(seqmat,bs,rows,bs,cols,vals,INSERT_VALUES));
41     }
42     /* off-diagonal blocks */
43     for (i=0; i<bs*bs; i++) vals[i] = -1.0;
44     for (i=0; i<(n-1)*bs; i+=bs) {
45       for (j=0; j<bs; j++) {rows[j] = i+j; cols[j] = i+bs+j;}
46       PetscCall(MatSetValues(seqmat,bs,rows,bs,cols,vals,INSERT_VALUES));
47     }
48 
49     PetscCall(PetscFree3(vals,rows,cols));
50   }
51   PetscCall(MatAssemblyBegin(seqmat,MAT_FINAL_ASSEMBLY));
52   PetscCall(MatAssemblyEnd(seqmat,MAT_FINAL_ASSEMBLY));
53   if (rank == 0) {
54     PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] seqmat:\n",rank));
55     PetscCall(MatView(seqmat,PETSC_VIEWER_STDOUT_SELF));
56   }
57 
58   PetscCall(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,seqmat,PETSC_DECIDE,MAT_INITIAL_MATRIX,&mpimat));
59   PetscCall(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,seqmat,PETSC_DECIDE,MAT_REUSE_MATRIX,&mpimat));
60   PetscCall(MatView(mpimat,PETSC_VIEWER_STDOUT_WORLD));
61 
62   PetscCall(MatDestroy(&seqmat));
63   PetscCall(MatDestroy(&mpimat));
64   PetscCall(PetscFinalize());
65   return 0;
66 }
67 
68 /*TEST
69 
70    test:
71       nsize: 3
72 
73    test:
74       suffix: 2
75       nsize: 3
76       args: -mat_type baij
77 
78    test:
79       suffix: 3
80       nsize: 3
81       args: -mat_type baij -bs 2
82 
83 TEST*/
84