xref: /petsc/src/mat/tests/ex52.c (revision 8fb5bd83c3955fefcf33a54e3bb66920a9fa884b)
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) PetscCall(MatSetOption(A,MAT_ROW_ORIENTED,PETSC_FALSE));
28 
29   /* inproc assembly */
30   for (row=start; row<end; row++) {
31     for (col=start; col<end; col++,data+=1) {
32       PetscCall(MatSetValues(A,1,&row,1,&col,&data,INSERT_VALUES));
33     }
34   }
35   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
36   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
37 
38   /* offproc assembly */
39   data = 5.0;
40   row  = (M+start-1)%M;
41   for (col=0; col<M; col++) {
42     PetscCall(MatSetValues(A,1,&row,1,&col,&data,ADD_VALUES));
43   }
44   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
45   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
46 
47   /* Test MatSetValuesBlocked() */
48   PetscCall(PetscOptionsHasName(NULL,NULL,"-test_setvaluesblocked",&flg));
49   if (flg) {
50     PetscScalar *bval;
51     row /= bs;
52     col  = start/bs;
53     PetscCall(PetscMalloc1(bs*bs,&bval));
54     k = 1;
55     /* row oriented - default */
56     for (i=0; i<bs; i++) {
57       for (j=0; j<bs; j++) {
58         bval[i*bs+j] = (PetscScalar)k; k++;
59       }
60     }
61     PetscCall(MatSetValuesBlocked(A,1,&row,1,&col,bval,INSERT_VALUES));
62     PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
63     PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
64     PetscCall(PetscFree(bval));
65   }
66 
67   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
68   PetscCall(MatDestroy(&A));
69   PetscCall(PetscFinalize());
70   return 0;
71 }
72 
73 /*TEST
74 
75    test:
76       suffix: 1
77       nsize: 3
78       args: -mat_block_size 2 -test_setvaluesblocked
79 
80    test:
81       suffix: 2
82       nsize: 3
83       args: -mat_block_size 2 -test_setvaluesblocked -column_oriented
84 
85    test:
86       suffix: 3
87       nsize: 3
88       args: -mat_block_size 1 -test_setvaluesblocked
89 
90    test:
91       suffix: 4
92       nsize: 3
93       args: -mat_block_size 1 -test_setvaluesblocked -column_oriented
94 
95 TEST*/
96