Actual source code: ex56.c
petsc-3.6.1 2015-08-06
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: PetscBool 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(NULL,"-test_mat_sbaij",&testsbaij);
24: if (testsbaij) {
25: MatCreateSBAIJ(PETSC_COMM_WORLD,bs,m*bs,n*bs,PETSC_DECIDE,PETSC_DECIDE,1,NULL,1,NULL,&A);
26: } else {
27: MatCreateBAIJ(PETSC_COMM_WORLD,bs,m*bs,n*bs,PETSC_DECIDE,PETSC_DECIDE,1,NULL,1,NULL,&A);
28: }
29: MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
30: eval = 9;
32: PetscOptionsHasName(NULL,"-ass_extern",&flg);
33: if (flg && (size != 1)) rstart = m*((rank+1)%size);
34: else rstart = m*(rank);
36: row[0] =rstart+0; row[1] =rstart+2;
37: col[0] =rstart+0; col[1] =rstart+1; col[2] =rstart+3;
38: for (i=0; i<6; i++) {
39: for (j =0; j< 9; j++) x[i][j] = (PetscScalar)val++;
40: }
42: MatSetValuesBlocked(A,2,row,3,col,&x[0][0],INSERT_VALUES);
46: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
47: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
49: /*
50: This option does not work for rectangular matrices
51: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
52: */
54: MatSetValuesBlocked(A,2,row,3,col,&x[0][0],INSERT_VALUES);
56: /* Do another MatSetValues to test the case when only one local block is specified */
57: for (i=0; i<3; i++) {
58: for (j =0; j<3; j++) y[i][j] = (PetscScalar)(10 + i*eval + j);
59: }
60: MatSetValuesBlocked(A,1,row,1,col,&y[0][0],INSERT_VALUES);
61: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
62: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
65: PetscOptionsHasName(NULL,"-zero_rows",&flg);
66: if (flg) {
67: col[0] = rstart*bs+0;
68: col[1] = rstart*bs+1;
69: col[2] = rstart*bs+2;
70: MatZeroRows(A,3,col,one,0,0);
71: }
73: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
75: MatDestroy(&A);
76: PetscFinalize();
77: return 0;
78: }