Actual source code: ex56.c
2: static char help[] = "Test the use of MatSetValuesBlocked(), MatZeroRows() for rectangular MatBAIJ matrix, test MatSetValuesBlocked() for MatSBAIJ matrix (-test_mat_sbaij).";
4: #include petscmat.h
8: int main(int argc,char **args)
9: {
10: Mat A;
11: PetscInt bs=3,m=4,n=6,i,j,val = 10,row[2],col[3],eval,rstart;
13: PetscMPIInt size,rank;
14: PetscScalar x[6][9],y[3][3],one=1.0;
15: PetscTruth flg,testsbaij=PETSC_FALSE;
17: PetscInitialize(&argc,&args,(char *)0,help);
19: MPI_Comm_size(PETSC_COMM_WORLD,&size);
20: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
22: PetscOptionsHasName(PETSC_NULL,"-test_mat_sbaij",&testsbaij);
23:
24: if (testsbaij){
25: if (size == 1) {
26: MatCreateSeqSBAIJ(PETSC_COMM_SELF,bs,m*bs,n*bs,1,PETSC_NULL,&A);
27: } else {
28: MatCreateMPISBAIJ(PETSC_COMM_WORLD,bs,m*bs,n*bs,PETSC_DECIDE,PETSC_DECIDE,1,
29: PETSC_NULL,1,PETSC_NULL,&A);
30: }
31: } else {
32: if (size == 1) {
33: MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,m*bs,n*bs,1,PETSC_NULL,&A);
34: } else {
35: MatCreateMPIBAIJ(PETSC_COMM_WORLD,bs,m*bs,n*bs,PETSC_DECIDE,PETSC_DECIDE,1,
36: PETSC_NULL,1,PETSC_NULL,&A);
37: }
38: }
39: PetscOptionsHasName(PETSC_NULL,"-column_oriented",&flg);
40: if (flg) {
41: MatSetOption(A,MAT_COLUMN_ORIENTED);
42: eval = 6;
43: } else {
44: eval = 9;
45: }
47: PetscOptionsHasName(PETSC_NULL,"-ass_extern",&flg);
48: if (flg && (size != 1)) rstart = m*((rank+1)%size);
49: else rstart = m*(rank);
51: row[0] =rstart+0; row[1] =rstart+2;
52: col[0] =rstart+0; col[1] =rstart+1; col[2] =rstart+3;
53: for (i=0; i<6; i++) {
54: for (j =0; j< 9; j++) x[i][j] = (PetscScalar)val++;
55: }
57: MatSetValuesBlocked(A,2,row,3,col,&x[0][0],INSERT_VALUES);
61: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
62: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
64: /*
65: This option does not work for rectangular matrices
66: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR);
67: */
68:
69: MatSetValuesBlocked(A,2,row,3,col,&x[0][0],INSERT_VALUES);
71: /* Do another MatSetValues to test the case when only one local block is specified */
72: for (i=0; i<3; i++) {
73: for (j =0; j<3 ; j++) y[i][j] = (PetscScalar)(10 + i*eval + j);
74: }
75: MatSetValuesBlocked(A,1,row,1,col,&y[0][0],INSERT_VALUES);
76: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
77: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
79:
80: PetscOptionsHasName(PETSC_NULL,"-zero_rows",&flg);
81: if (flg) {
82: col[0] = rstart*bs+0;
83: col[1] = rstart*bs+1;
84: col[2] = rstart*bs+2;
85: MatZeroRows(A,3,col,one);
86: }
88: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
90: MatDestroy(A);
91: PetscFinalize();
92: return 0;
93: }