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