xref: /petsc/src/mat/tests/ex86.c (revision 1e1ea65d8de51fde77ce8a787efbef25e407badc)
1 static char help[] = "Testing MatCreateMPIMatConcatenateSeqMat().\n\n";
2 
3 #include <petscmat.h>
4 int main(int argc,char **argv)
5 {
6   PetscErrorCode ierr;
7   Mat            seqmat,mpimat;
8   PetscMPIInt    rank;
9   PetscScalar    value[3],*vals;
10   PetscInt       i,col[3],n=5,bs=1;
11 
12   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
13   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
14   ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr);
15 
16   /* Create seqaij matrices of size (n+rank) by n */
17   ierr = MatCreate(PETSC_COMM_SELF,&seqmat);CHKERRQ(ierr);
18   ierr = MatSetSizes(seqmat,(n+rank)*bs,PETSC_DECIDE,PETSC_DECIDE,n*bs);CHKERRQ(ierr);
19   ierr = MatSetFromOptions(seqmat);CHKERRQ(ierr);
20   ierr = MatSeqAIJSetPreallocation(seqmat,3*bs,NULL);CHKERRQ(ierr);
21   ierr = MatSeqBAIJSetPreallocation(seqmat,bs,3,NULL);CHKERRQ(ierr);
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       ierr   = MatSetValues(seqmat,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr);
28     }
29     i = n - 1; col[0] = n - 2; col[1] = n - 1;
30     ierr = MatSetValues(seqmat,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
31 
32     i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
33     ierr = MatSetValues(seqmat,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
34   } else {
35     PetscInt *rows,*cols,j;
36     ierr = PetscMalloc3(bs*bs,&vals,bs,&rows,bs,&cols);CHKERRQ(ierr);
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       ierr = MatSetValues(seqmat,bs,rows,bs,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
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       ierr = MatSetValues(seqmat,bs,rows,bs,cols,vals,INSERT_VALUES);CHKERRQ(ierr);
48     }
49 
50     ierr = PetscFree3(vals,rows,cols);CHKERRQ(ierr);
51   }
52   ierr = MatAssemblyBegin(seqmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
53   ierr = MatAssemblyEnd(seqmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
54   if (!rank) {
55     ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] seqmat:\n",rank);CHKERRQ(ierr);
56     ierr = MatView(seqmat,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
57   }
58 
59   ierr = MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,seqmat,PETSC_DECIDE,MAT_INITIAL_MATRIX,&mpimat);CHKERRQ(ierr);
60   ierr = MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,seqmat,PETSC_DECIDE,MAT_REUSE_MATRIX,&mpimat);CHKERRQ(ierr);
61   ierr = MatView(mpimat,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
62 
63   ierr = MatDestroy(&seqmat);CHKERRQ(ierr);
64   ierr = MatDestroy(&mpimat);CHKERRQ(ierr);
65   ierr = PetscFinalize();
66   return ierr;
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