xref: /petsc/src/mat/tests/ex95.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1 static char help[] = "Testing MatCreateMPIAIJSumSeqAIJ().\n\n";
2 
3 #include <petscmat.h>
4 
5 int main(int argc,char **argv)
6 {
7   Mat            A,B;
8   MatScalar      a[1],alpha;
9   PetscMPIInt    size,rank;
10   PetscInt       m,n,i,col, prid;
11 
12   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
13   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
14   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
15   prid = size;
16   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-prid",&prid,NULL));
17 
18   m    = n = 10*size;
19   CHKERRQ(MatCreate(PETSC_COMM_SELF,&A));
20   CHKERRQ(MatSetSizes(A,PETSC_DETERMINE,PETSC_DETERMINE,m,n));
21   CHKERRQ(MatSetType(A,MATSEQAIJ));
22   CHKERRQ(MatSetUp(A));
23 
24   a[0] = rank+1;
25   for (i=0; i<m-rank; i++) {
26     col  = i+rank;
27     CHKERRQ(MatSetValues(A,1,&i,1,&col,a,INSERT_VALUES));
28   }
29   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
30   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
31 
32   if (rank == prid) {
33     CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d] A: \n",rank));
34     CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_SELF));
35   }
36 
37   /* Test MatCreateMPIAIJSumSeqAIJ */
38   CHKERRQ(MatCreateMPIAIJSumSeqAIJ(PETSC_COMM_WORLD,A,PETSC_DECIDE,PETSC_DECIDE,MAT_INITIAL_MATRIX,&B));
39 
40   /* Test MAT_REUSE_MATRIX */
41   alpha = 0.1;
42   for (i=0; i<3; i++) {
43     CHKERRQ(MatScale(A,alpha));
44     CHKERRQ(MatCreateMPIAIJSumSeqAIJ(PETSC_COMM_WORLD,A,PETSC_DECIDE,PETSC_DECIDE,MAT_REUSE_MATRIX,&B));
45   }
46   CHKERRQ(MatView(B, PETSC_VIEWER_STDOUT_WORLD));
47   CHKERRQ(MatDestroy(&B));
48   CHKERRQ(MatDestroy(&A));
49   CHKERRQ(PetscFinalize());
50   return 0;
51 }
52 
53 /*TEST
54 
55    test:
56       nsize: 3
57       filter: grep -v "MPI processes"
58 
59    test:
60       suffix: 2
61       filter: grep -v "MPI processes"
62 
63 TEST*/
64