Actual source code: ex56.c
1: /*$Id: ex56.c,v 1.32 2001/08/07 03:03:07 balay Exp $*/
2: static char help[] = "Test the use of MatSetValuesBlocked(), MatZeroRows() for rectangular MatBAIJ matrix.";
4: #include petscmat.h
6: int main(int argc,char **args)
7: {
8: Mat A;
9: int bs=3,m=4,n=6,i,j,val = 10,ierr,size,rank,rstart;
10: PetscScalar x[6][9],y[3][3],one=1.0;
11: int row[2],col[3],eval;
12: IS is;
13: PetscTruth flg;
15: PetscInitialize(&argc,&args,(char *)0,help);
17: MPI_Comm_size(PETSC_COMM_WORLD,&size);
18: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
19:
20: if (size == 1) {
21: MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,m*bs,n*bs,1,PETSC_NULL,&A);
22: } else {
23: MatCreateMPIBAIJ(PETSC_COMM_WORLD,bs,m*bs,n*bs,PETSC_DECIDE,PETSC_DECIDE,1,
24: PETSC_NULL,1,PETSC_NULL,&A);
25: }
27: PetscOptionsHasName(PETSC_NULL,"-column_oriented",&flg);
28: if (flg) {
29: MatSetOption(A,MAT_COLUMN_ORIENTED);
30: eval = 6;
31: } else {
32: eval = 9;
33: }
35: PetscOptionsHasName(PETSC_NULL,"-ass_extern",&flg);
36: if (flg && (size != 1)) rstart = m*((rank+1)%size);
37: else rstart = m*(rank);
39: row[0] =rstart+0; row[1] =rstart+2;
40: col[0] =rstart+0; col[1] =rstart+1; col[2] =rstart+3;
41: for (i=0; i<6; i++) {
42: for (j =0; j< 9; j++) x[i][j] = (PetscScalar)val++;
43: }
45: MatSetValuesBlocked(A,2,row,3,col,&x[0][0],INSERT_VALUES);
49: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
50: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
52: /*
53: This option does not work for rectangular matrices
54: MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR);
55: */
56:
57: MatSetValuesBlocked(A,2,row,3,col,&x[0][0],INSERT_VALUES);
59: /* Do another MatSetValues to test the case when only one local block is specified */
60: for (i=0; i<3; i++) {
61: for (j =0; j<3 ; j++) y[i][j] = (PetscScalar)(10 + i*eval + j);
62: }
63: MatSetValuesBlocked(A,1,row,1,col,&y[0][0],INSERT_VALUES);
64: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
65: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
67:
68: PetscOptionsHasName(PETSC_NULL,"-zero_rows",&flg);
69: if (flg) {
70: col[0] = rstart*bs+0;
71: col[1] = rstart*bs+1;
72: col[2] = rstart*bs+2;
73: ISCreateGeneral(MPI_COMM_SELF,3,col,&is);
74: MatZeroRows(A,is,&one);
75: ISDestroy(is);
76: }
78: MatView(A,PETSC_VIEWER_STDOUT_WORLD);
80: MatDestroy(A);
81: PetscFinalize();
82: return 0;
83: }