xref: /petsc/src/mat/tests/ex86.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown static char help[] = "Testing MatCreateMPIMatConcatenateSeqMat().\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown int main(int argc,char **argv)
5c4762a1bSJed Brown {
6c4762a1bSJed Brown   Mat            seqmat,mpimat;
7c4762a1bSJed Brown   PetscMPIInt    rank;
8c4762a1bSJed Brown   PetscScalar    value[3],*vals;
9c4762a1bSJed Brown   PetscInt       i,col[3],n=5,bs=1;
10c4762a1bSJed Brown 
11*327415f7SBarry Smith   PetscFunctionBeginUser;
129566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
139566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL));
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   /* Create seqaij matrices of size (n+rank) by n */
179566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF,&seqmat));
189566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(seqmat,(n+rank)*bs,PETSC_DECIDE,PETSC_DECIDE,n*bs));
199566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(seqmat));
209566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(seqmat,3*bs,NULL));
219566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(seqmat,bs,3,NULL));
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   if (bs == 1) {
24c4762a1bSJed Brown     value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
25c4762a1bSJed Brown     for (i=1; i<n-1; i++) {
26c4762a1bSJed Brown       col[0] = i-1; col[1] = i; col[2] = i+1;
279566063dSJacob Faibussowitsch       PetscCall(MatSetValues(seqmat,1,&i,3,col,value,INSERT_VALUES));
28c4762a1bSJed Brown     }
29c4762a1bSJed Brown     i = n - 1; col[0] = n - 2; col[1] = n - 1;
309566063dSJacob Faibussowitsch     PetscCall(MatSetValues(seqmat,1,&i,2,col,value,INSERT_VALUES));
31c4762a1bSJed Brown 
32c4762a1bSJed Brown     i = 0; col[0] = 0; col[1] = 1; value[0] = 2.0; value[1] = -1.0;
339566063dSJacob Faibussowitsch     PetscCall(MatSetValues(seqmat,1,&i,2,col,value,INSERT_VALUES));
34c4762a1bSJed Brown   } else {
35c4762a1bSJed Brown     PetscInt *rows,*cols,j;
369566063dSJacob Faibussowitsch     PetscCall(PetscMalloc3(bs*bs,&vals,bs,&rows,bs,&cols));
37c4762a1bSJed Brown     /* diagonal blocks */
38c4762a1bSJed Brown     for (i=0; i<bs*bs; i++) vals[i] = 2.0;
39c4762a1bSJed Brown     for (i=0; i<n*bs; i+=bs) {
40c4762a1bSJed Brown       for (j=0; j<bs; j++) {rows[j] = i+j; cols[j] = i+j;}
419566063dSJacob Faibussowitsch       PetscCall(MatSetValues(seqmat,bs,rows,bs,cols,vals,INSERT_VALUES));
42c4762a1bSJed Brown     }
43c4762a1bSJed Brown     /* off-diagonal blocks */
44c4762a1bSJed Brown     for (i=0; i<bs*bs; i++) vals[i] = -1.0;
45c4762a1bSJed Brown     for (i=0; i<(n-1)*bs; i+=bs) {
46c4762a1bSJed Brown       for (j=0; j<bs; j++) {rows[j] = i+j; cols[j] = i+bs+j;}
479566063dSJacob Faibussowitsch       PetscCall(MatSetValues(seqmat,bs,rows,bs,cols,vals,INSERT_VALUES));
48c4762a1bSJed Brown     }
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch     PetscCall(PetscFree3(vals,rows,cols));
51c4762a1bSJed Brown   }
529566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(seqmat,MAT_FINAL_ASSEMBLY));
539566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(seqmat,MAT_FINAL_ASSEMBLY));
54dd400576SPatrick Sanan   if (rank == 0) {
559566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] seqmat:\n",rank));
569566063dSJacob Faibussowitsch     PetscCall(MatView(seqmat,PETSC_VIEWER_STDOUT_SELF));
57c4762a1bSJed Brown   }
58c4762a1bSJed Brown 
599566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,seqmat,PETSC_DECIDE,MAT_INITIAL_MATRIX,&mpimat));
609566063dSJacob Faibussowitsch   PetscCall(MatCreateMPIMatConcatenateSeqMat(PETSC_COMM_WORLD,seqmat,PETSC_DECIDE,MAT_REUSE_MATRIX,&mpimat));
619566063dSJacob Faibussowitsch   PetscCall(MatView(mpimat,PETSC_VIEWER_STDOUT_WORLD));
62c4762a1bSJed Brown 
639566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&seqmat));
649566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&mpimat));
659566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
66b122ec5aSJacob Faibussowitsch   return 0;
67c4762a1bSJed Brown }
68c4762a1bSJed Brown 
69c4762a1bSJed Brown /*TEST
70c4762a1bSJed Brown 
71c4762a1bSJed Brown    test:
72c4762a1bSJed Brown       nsize: 3
73c4762a1bSJed Brown 
74c4762a1bSJed Brown    test:
75c4762a1bSJed Brown       suffix: 2
76c4762a1bSJed Brown       nsize: 3
77c4762a1bSJed Brown       args: -mat_type baij
78c4762a1bSJed Brown 
79c4762a1bSJed Brown    test:
80c4762a1bSJed Brown       suffix: 3
81c4762a1bSJed Brown       nsize: 3
82c4762a1bSJed Brown       args: -mat_type baij -bs 2
83c4762a1bSJed Brown 
84c4762a1bSJed Brown TEST*/
85