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