xref: /petsc/src/mat/tests/ex56.c (revision 6ffe77eaecce1557e50d00ca5347a7f48e598865)
1 
2 static char help[] = "Test the use of MatSetValuesBlocked(), MatZeroRows() for rectangular MatBAIJ matrix, test MatSetValuesBlocked() for MatSBAIJ matrix (-test_mat_sbaij).";
3 
4 #include <petscmat.h>
5 
6 int main(int argc,char **args)
7 {
8   Mat            A;
9   PetscInt       bs=3,m=4,n=6,i,j,val = 10,row[2],col[3],eval,rstart;
10   PetscMPIInt    size,rank;
11   PetscScalar    x[6][9],y[3][3],one=1.0;
12   PetscBool      flg,testsbaij=PETSC_FALSE;
13 
14   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
15   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
16   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
17 
18   PetscCall(PetscOptionsHasName(NULL,NULL,"-test_mat_sbaij",&testsbaij));
19 
20   if (testsbaij) {
21     PetscCall(MatCreateSBAIJ(PETSC_COMM_WORLD,bs,m*bs,n*bs,PETSC_DECIDE,PETSC_DECIDE,1,NULL,1,NULL,&A));
22   } else {
23     PetscCall(MatCreateBAIJ(PETSC_COMM_WORLD,bs,m*bs,n*bs,PETSC_DECIDE,PETSC_DECIDE,1,NULL,1,NULL,&A));
24   }
25   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
26   eval = 9;
27 
28   PetscCall(PetscOptionsHasName(NULL,NULL,"-ass_extern",&flg));
29   if (flg && (size != 1)) rstart = m*((rank+1)%size);
30   else rstart = m*(rank);
31 
32   row[0] =rstart+0;  row[1] =rstart+2;
33   col[0] =rstart+0;  col[1] =rstart+1;  col[2] =rstart+3;
34   for (i=0; i<6; i++) {
35     for (j =0; j< 9; j++) x[i][j] = (PetscScalar)val++;
36   }
37 
38   PetscCall(MatSetValuesBlocked(A,2,row,3,col,&x[0][0],INSERT_VALUES));
39 
40   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
41   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
42 
43   /*
44   This option does not work for rectangular matrices
45   PetscCall(MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
46   */
47 
48   PetscCall(MatSetValuesBlocked(A,2,row,3,col,&x[0][0],INSERT_VALUES));
49 
50   /* Do another MatSetValues to test the case when only one local block is specified */
51   for (i=0; i<3; i++) {
52     for (j =0; j<3; j++) y[i][j] = (PetscScalar)(10 + i*eval + j);
53   }
54   PetscCall(MatSetValuesBlocked(A,1,row,1,col,&y[0][0],INSERT_VALUES));
55   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
56   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
57 
58   PetscCall(PetscOptionsHasName(NULL,NULL,"-zero_rows",&flg));
59   if (flg) {
60     col[0] = rstart*bs+0;
61     col[1] = rstart*bs+1;
62     col[2] = rstart*bs+2;
63     PetscCall(MatZeroRows(A,3,col,one,0,0));
64   }
65 
66   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD));
67 
68   PetscCall(MatDestroy(&A));
69   PetscCall(PetscFinalize());
70   return 0;
71 }
72 
73 /*TEST
74 
75    test:
76       filter: grep -v " MPI process"
77 
78    test:
79       suffix: 4
80       nsize: 3
81       args: -ass_extern
82       filter: grep -v " MPI process"
83 
84    test:
85       suffix: 5
86       nsize: 3
87       args: -ass_extern -zero_rows
88       filter: grep -v " MPI process"
89 
90 TEST*/
91