Actual source code: ex56.c
petsc-3.8.4 2018-03-24
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>
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;
11: PetscMPIInt size,rank;
12: PetscScalar x[6][9],y[3][3],one=1.0;
13: PetscBool flg,testsbaij=PETSC_FALSE;
15: PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
16: MPI_Comm_size(PETSC_COMM_WORLD,&size);
17: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
19: PetscOptionsHasName(NULL,NULL,"-test_mat_sbaij",&testsbaij);
21: if (testsbaij) {
22: MatCreateSBAIJ(PETSC_COMM_WORLD,bs,m*bs,n*bs,PETSC_DECIDE,PETSC_DECIDE,1,NULL,1,NULL,&A);
23: } else {
24: MatCreateBAIJ(PETSC_COMM_WORLD,bs,m*bs,n*bs,PETSC_DECIDE,PETSC_DECIDE,1,NULL,1,NULL,&A);
25: }
26: MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);
27: eval = 9;
29: PetscOptionsHasName(NULL,NULL,"-ass_extern",&flg);
30: if (flg && (size != 1)) rstart = m*((rank+1)%size);
31: else rstart = m*(rank);
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: }
39: MatSetValuesBlocked(A,2,row,3,col,&x[0][0],INSERT_VALUES);
43: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
44: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
46: /*
47: This option does not work for rectangular matrices
48: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
49: */
51: MatSetValuesBlocked(A,2,row,3,col,&x[0][0],INSERT_VALUES);
53: /* Do another MatSetValues to test the case when only one local block is specified */
54: for (i=0; i<3; i++) {
55: for (j =0; j<3; j++) y[i][j] = (PetscScalar)(10 + i*eval + j);
56: }
57: MatSetValuesBlocked(A,1,row,1,col,&y[0][0],INSERT_VALUES);
58: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
59: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
62: PetscOptionsHasName(NULL,NULL,"-zero_rows",&flg);
63: if (flg) {
64: col[0] = rstart*bs+0;
65: col[1] = rstart*bs+1;
66: col[2] = rstart*bs+2;
67: MatZeroRows(A,3,col,one,0,0);
68: }
70: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
72: MatDestroy(&A);
73: PetscFinalize();
74: return ierr;
75: }