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