Actual source code: ex76.c
1: /*$Id: ex76.c,v 1.18 2001/08/07 21:30:08 bsmith Exp $*/
3: static char help[] = "Tests matrix permutation for factorization and solve on matrix with MatSBAIJ format. Modified from ex74.cn";
5: #include petscmat.h
7: int main(int argc,char **args)
8: {
9: Vec x,y,b;
10: Mat A; /* linear system matrix */
11: Mat sA,sC; /* symmetric part of the matrices */
12: int n,mbs=16,bs=1,nz=3,prob=1;
13: int ierr,i,j,col[3],size,block, row,I,J,n1,*ip_ptr;
14: int lf; /* level of fill for icc */
15: PetscReal norm1,norm2,tol=1.e-10,fill;
16: PetscScalar neg_one = -1.0,four=4.0,value[3];
17: IS perm;
18: PetscRandom rand;
19: PetscTruth reorder=PETSC_TRUE;
21: PetscInitialize(&argc,&args,(char *)0,help);
22: MPI_Comm_size(PETSC_COMM_WORLD,&size);
23: if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
24: PetscOptionsGetInt(PETSC_NULL,"-bs",&bs,PETSC_NULL);
25: PetscOptionsGetInt(PETSC_NULL,"-mbs",&mbs,PETSC_NULL);
27: n = mbs*bs;
28: ierr=MatCreateSeqBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,PETSC_NULL, &A);
29: ierr=MatCreateSeqSBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,PETSC_NULL, &sA);
31: /* Test MatGetOwnershipRange() */
32: MatGetOwnershipRange(A,&I,&J);
33: MatGetOwnershipRange(sA,&i,&j);
34: if (i-I || j-J){
35: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetOwnershipRange() in MatSBAIJ formatn");
36: }
38: /* Assemble matrix */
39: if (bs == 1){
40: PetscOptionsGetInt(PETSC_NULL,"-test_problem",&prob,PETSC_NULL);
41: if (prob == 1){ /* tridiagonal matrix */
42: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
43: for (i=1; i<n-1; i++) {
44: col[0] = i-1; col[1] = i; col[2] = i+1;
45: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
46: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
47: }
48: i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;
49: value[0]= 0.1; value[1]=-1; value[2]=2;
50: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
51: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
53: i = 0; col[0] = 0; col[1] = 1; col[2]=n-1;
54: value[0] = 2.0; value[1] = -1.0; value[2]=0.1;
55: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
56: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
57: }
58: else if (prob ==2){ /* matrix for the five point stencil */
59: n1 = (int) (sqrt((PetscReal)n) + 0.001);
60: if (n1*n1 - n) SETERRQ(PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive interger!");
61: for (i=0; i<n1; i++) {
62: for (j=0; j<n1; j++) {
63: I = j + n1*i;
64: if (i>0) {
65: J = I - n1;
66: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
67: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
68: }
69: if (i<n1-1) {
70: J = I + n1;
71: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
72: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
73: }
74: if (j>0) {
75: J = I - 1;
76: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
77: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
78: }
79: if (j<n1-1) {
80: J = I + 1;
81: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
82: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
83: }
84: MatSetValues(A,1,&I,1,&I,&four,INSERT_VALUES);
85: MatSetValues(sA,1,&I,1,&I,&four,INSERT_VALUES);
86: }
87: }
88: }
89: }
90: else { /* bs > 1 */
91: for (block=0; block<n/bs; block++){
92: /* diagonal blocks */
93: value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
94: for (i=1+block*bs; i<bs-1+block*bs; i++) {
95: col[0] = i-1; col[1] = i; col[2] = i+1;
96: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
97: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
98: }
99: i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
100: value[0]=-1.0; value[1]=4.0;
101: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
102: MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
104: i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
105: value[0]=4.0; value[1] = -1.0;
106: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
107: MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
108: }
109: /* off-diagonal blocks */
110: value[0]=-1.0;
111: for (i=0; i<(n/bs-1)*bs; i++){
112: col[0]=i+bs;
113: MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);
114: MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES);
115: col[0]=i; row=i+bs;
116: MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
117: MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);
118: }
119: if (bs == 2){
120: /* insert a value to off-diag blocks */
121: row = 2; col[0] = 5; value[0] = 0.01;
122: MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
123: MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);
124: row = 0; col[0] = 3; value[0] = 0.01;
125: MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
126: MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);
127: }
128: }
129: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
130: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
131: /* PetscPrintf(PETSC_COMM_SELF,"n The Matrix: n");
132: MatView(A, PETSC_VIEWER_DRAW_WORLD);
133: MatView(A, PETSC_VIEWER_STDOUT_WORLD); */
135: MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY);
136: MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY);
137: /* PetscPrintf(PETSC_COMM_SELF,"n Symmetric Part of Matrix: n"); */
138: /* MatView(sA, PETSC_VIEWER_DRAW_WORLD); */
139: /* MatView(sA, PETSC_VIEWER_STDOUT_WORLD); */
141: /* Vectors */
142: PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&rand);
143: VecCreateSeq(PETSC_COMM_SELF,n,&x);
144: VecDuplicate(x,&b);
145: VecDuplicate(x,&y);
146: VecSetRandom(rand,x);
148: /* Test MatReordering() */
149: PetscMalloc(mbs*sizeof(int),&ip_ptr);
150: for (i=0; i<mbs; i++) ip_ptr[i] = i;
151: if(reorder){
152: i = ip_ptr[1]; ip_ptr[1] = ip_ptr[mbs-2]; ip_ptr[mbs-2] = i;
153: /* i = ip_ptr[0]; ip_ptr[0] = ip_ptr[mbs-1]; ip_ptr[mbs-1] = i; */
154: /* i = ip_ptr[2]; ip_ptr[2] = ip_ptr[mbs-3]; ip_ptr[mbs-3] = i; */
155: }
156: ISCreateGeneral(PETSC_COMM_SELF,mbs,ip_ptr,&perm);
157: ISSetPermutation(perm);
158:
159: /* Test MatCholeskyFactor(), MatICCFactor() */
160: norm1 = tol;
161: for (lf=-1; lf<10*bs; lf += bs){
162: if (lf==-1) { /* Cholesky factor */
163: fill = 5.0;
164: MatCholeskyFactorSymbolic(sA,perm,fill,&sC);
165: } else { /* incomplete Cholesky factor */
166: fill = 5.0;
167: MatICCFactorSymbolic(sA,perm,fill,lf,&sC);
168: }
169: MatCholeskyFactorNumeric(sA,&sC);
170: /* MatView(sC, PETSC_VIEWER_DRAW_WORLD); */ /* view factored matrix */
171: /* MatView(sC, PETSC_VIEWER_STDOUT_WORLD); */
172:
173: MatMult(sA,x,b);
174: MatSolve(sC,b,y);
175: MatDestroy(sC);
177: /* Check the error */
178: VecAXPY(&neg_one,x,y);
179: VecNorm(y,NORM_2,&norm2);
180: /* printf("lf: %d, error: %gn", lf,norm2); */
181: if (10*norm1 < norm2 && lf-bs != -1){
182: PetscPrintf(PETSC_COMM_SELF,"lf=%d, %d, Norm of error=%g, %gn",lf-bs,lf,norm1,norm2);
183: }
184: norm1 = norm2;
185: if (norm2 < tol && lf != -1) break;
186: }
188: ISDestroy(perm);
189: PetscFree(ip_ptr);
190: MatDestroy(A);
191: MatDestroy(sA);
192: VecDestroy(x);
193: VecDestroy(y);
194: VecDestroy(b);
195: PetscRandomDestroy(rand);
197: PetscFinalize();
198: return 0;
199: }