xref: /petsc/src/mat/tests/ex52.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
1 
2 static char help[] = "Tests various routines in MatMPIBAIJ format.\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc,char **args)
7 {
8   Mat            A;
9   PetscInt       m=2,bs=1,M,row,col,start,end,i,j,k;
10   PetscMPIInt    rank,size;
11   PetscScalar    data=100;
12   PetscBool      flg;
13 
14   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
15   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
16   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
17 
18   /* Test MatSetValues() and MatGetValues() */
19   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL));
20   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL));
21 
22   M    = m*bs*size;
23   PetscCall(MatCreateBAIJ(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,M,M,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&A));
24 
25   PetscCall(MatGetOwnershipRange(A,&start,&end));
26   PetscCall(PetscOptionsHasName(NULL,NULL,"-column_oriented",&flg));
27   if (flg) {
28     PetscCall(MatSetOption(A,MAT_ROW_ORIENTED,PETSC_FALSE));
29   }
30 
31   /* inproc assembly */
32   for (row=start; row<end; row++) {
33     for (col=start; col<end; col++,data+=1) {
34       PetscCall(MatSetValues(A,1,&row,1,&col,&data,INSERT_VALUES));
35     }
36   }
37   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
38   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
39 
40   /* offproc assembly */
41   data = 5.0;
42   row  = (M+start-1)%M;
43   for (col=0; col<M; col++) {
44     PetscCall(MatSetValues(A,1,&row,1,&col,&data,ADD_VALUES));
45   }
46   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
47   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
48 
49   /* Test MatSetValuesBlocked() */
50   PetscCall(PetscOptionsHasName(NULL,NULL,"-test_setvaluesblocked",&flg));
51   if (flg) {
52     PetscScalar *bval;
53     row /= bs;
54     col  = start/bs;
55     PetscCall(PetscMalloc1(bs*bs,&bval));
56     k = 1;
57     /* row oriented - default */
58     for (i=0; i<bs; i++) {
59       for (j=0; j<bs; j++) {
60         bval[i*bs+j] = (PetscScalar)k; k++;
61       }
62     }
63     PetscCall(MatSetValuesBlocked(A,1,&row,1,&col,bval,INSERT_VALUES));
64     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
65     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
66     PetscCall(PetscFree(bval));
67   }
68 
69   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
70   PetscCall(MatDestroy(&A));
71   PetscCall(PetscFinalize());
72   return 0;
73 }
74 
75 /*TEST
76 
77    test:
78       suffix: 1
79       nsize: 3
80       args: -mat_block_size 2 -test_setvaluesblocked
81 
82    test:
83       suffix: 2
84       nsize: 3
85       args: -mat_block_size 2 -test_setvaluesblocked -column_oriented
86 
87    test:
88       suffix: 3
89       nsize: 3
90       args: -mat_block_size 1 -test_setvaluesblocked
91 
92    test:
93       suffix: 4
94       nsize: 3
95       args: -mat_block_size 1 -test_setvaluesblocked -column_oriented
96 
97 TEST*/
98