Actual source code: ex76.c
2: static char help[] = "Tests cholesky/icc factorization and solve on sequential aij, baij and sbaij matrices. \n";
4: #include petscmat.h
8: int main(int argc,char **args)
9: {
10: Vec x,y,b;
11: Mat A; /* linear system matrix */
12: Mat sA,sC; /* symmetric part of the matrices */
13: PetscInt n,mbs=16,bs=1,nz=3,prob=1,i,j,col[3],block, row,I,J,n1,*ip_ptr,lvl;
15: PetscMPIInt size;
16: PetscReal norm2,tol=1.e-10,err[10];
17: PetscScalar neg_one = -1.0,four=4.0,value[3];
18: IS perm;
19: PetscRandom rdm;
20: PetscInt reorder=0,displ=0;
21: MatFactorInfo factinfo;
22: PetscTruth equal;
23: PetscTruth TestAIJ=PETSC_FALSE,TestBAIJ=PETSC_TRUE;
24: PetscInt TestShift=0;
26: PetscInitialize(&argc,&args,(char *)0,help);
27: MPI_Comm_size(PETSC_COMM_WORLD,&size);
28: if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
29: PetscOptionsGetInt(PETSC_NULL,"-bs",&bs,PETSC_NULL);
30: PetscOptionsGetInt(PETSC_NULL,"-mbs",&mbs,PETSC_NULL);
31: PetscOptionsGetInt(PETSC_NULL,"-reorder",&reorder,PETSC_NULL);
32: PetscOptionsGetTruth(PETSC_NULL,"-testaij",&TestAIJ,PETSC_NULL);
33: PetscOptionsGetInt(PETSC_NULL,"-testShift",&TestShift,PETSC_NULL);
34: PetscOptionsGetInt(PETSC_NULL,"-displ",&displ,PETSC_NULL);
36: n = mbs*bs;
37: if (TestAIJ){ /* A is in aij format */
38: MatCreateSeqAIJ(PETSC_COMM_WORLD,n,n,nz,PETSC_NULL,&A);
39: TestBAIJ = PETSC_FALSE;
40: } else { /* A is in baij format */
41: ierr=MatCreateSeqBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,PETSC_NULL,&A);
42: TestAIJ = PETSC_FALSE;
43: }
44: MatCreateSeqSBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,PETSC_NULL,&sA);
46: /* Test MatGetOwnershipRange() */
47: MatGetOwnershipRange(A,&I,&J);
48: MatGetOwnershipRange(sA,&i,&j);
49: if (i-I || j-J){
50: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetOwnershipRange() in MatSBAIJ format\n");
51: }
53: /* Assemble matrix */
54: if (bs == 1){
55: PetscOptionsGetInt(PETSC_NULL,"-test_problem",&prob,PETSC_NULL);
56: if (prob == 1){ /* tridiagonal matrix */
57: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
58: for (i=1; i<n-1; i++) {
59: col[0] = i-1; col[1] = i; col[2] = i+1;
60: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
61: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
62: }
63: i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;
64: value[0]= 0.1; value[1]=-1; value[2]=2;
65: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
66: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
68: i = 0; col[0] = 0; col[1] = 1; col[2]=n-1;
69: value[0] = 2.0; value[1] = -1.0; value[2]=0.1;
70: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
71: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
72: } else if (prob ==2){ /* matrix for the five point stencil */
73: n1 = (int) (sqrt((PetscReal)n) + 0.001);
74: if (n1*n1 - n) SETERRQ(PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive interger!");
75: for (i=0; i<n1; i++) {
76: for (j=0; j<n1; j++) {
77: I = j + n1*i;
78: if (i>0) {
79: J = I - n1;
80: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
81: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
82: }
83: if (i<n1-1) {
84: J = I + n1;
85: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
86: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
87: }
88: if (j>0) {
89: J = I - 1;
90: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
91: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
92: }
93: if (j<n1-1) {
94: J = I + 1;
95: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
96: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
97: }
98: MatSetValues(A,1,&I,1,&I,&four,INSERT_VALUES);
99: MatSetValues(sA,1,&I,1,&I,&four,INSERT_VALUES);
100: }
101: }
102: }
103: } else { /* bs > 1 */
104: for (block=0; block<n/bs; block++){
105: /* diagonal blocks */
106: value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
107: for (i=1+block*bs; i<bs-1+block*bs; i++) {
108: col[0] = i-1; col[1] = i; col[2] = i+1;
109: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
110: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
111: }
112: i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
113: value[0]=-1.0; value[1]=4.0;
114: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
115: MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
117: i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
118: value[0]=4.0; value[1] = -1.0;
119: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
120: MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
121: }
122: /* off-diagonal blocks */
123: value[0]=-1.0;
124: for (i=0; i<(n/bs-1)*bs; i++){
125: col[0]=i+bs;
126: MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);
127: MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES);
128: col[0]=i; row=i+bs;
129: MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
130: MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);
131: }
132: }
134: if (TestShift){
135: /* set diagonals in the 0-th block as 0 for testing shift numerical factor */
136: for (i=0; i<bs; i++){
137: row = i; col[0] = i; value[0] = 0.0;
138: MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
139: MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);
140: }
141: }
143: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
144: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
145:
146: MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY);
147: MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY);
149: MatMultEqual(A,sA,5,&equal);
150: if (!equal) SETERRQ(PETSC_ERR_USER,"A != sA");
152: /* Vectors */
153: PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&rdm);
154: VecCreateSeq(PETSC_COMM_SELF,n,&x);
155: VecDuplicate(x,&b);
156: VecDuplicate(x,&y);
157: VecSetRandom(x,rdm);
159: /* Test MatReordering() on a symmetric ordering */
160: PetscMalloc(mbs*sizeof(PetscInt),&ip_ptr);
161: for (i=0; i<mbs; i++) ip_ptr[i] = i;
162: switch (reorder){
163: case 0: break;
164: case 1:
165: i = ip_ptr[2]; ip_ptr[2] = ip_ptr[mbs-3]; ip_ptr[mbs-3] = i;
166: case 2:
167: i = ip_ptr[1]; ip_ptr[1] = ip_ptr[mbs-2]; ip_ptr[mbs-2] = i;
168: case 3:
169: i = ip_ptr[0]; ip_ptr[0] = ip_ptr[mbs-1]; ip_ptr[mbs-1] = i;
170: }
172: ISCreateGeneral(PETSC_COMM_SELF,mbs,ip_ptr,&perm);
173: ISSetPermutation(perm);
175: /* initialize factinfo */
176: factinfo.shiftnz = 0.0;
177: factinfo.shiftpd = 0.0; /* false */
178: factinfo.zeropivot = 1.e-12;
179: if (TestShift == 1){
180: factinfo.shiftnz = 0.1;
181: } else if (TestShift == 2){
182: factinfo.shiftpd = 1.0; /* true*/
183: }
184:
185: /* Test MatCholeskyFactor(), MatICCFactor() */
186: /*------------------------------------------*/
187: /* Test aij matrix A */
188: if (TestAIJ){
189: if (displ>0){PetscPrintf(PETSC_COMM_SELF,"AIJ: \n");}
190: i = 0;
191: for (lvl=-1; lvl<10; lvl++){
192: if (lvl==-1) { /* Cholesky factor */
193: factinfo.fill = 5.0;
194: MatCholeskyFactorSymbolic(A,perm,&factinfo,&sC);
195: } else { /* incomplete Cholesky factor */
196: factinfo.fill = 5.0;
197: factinfo.levels = lvl;
198: MatICCFactorSymbolic(A,perm,&factinfo,&sC);
199: }
200: MatCholeskyFactorNumeric(A,&factinfo,&sC);
201: MatMult(A,x,b);
202: MatSolve(sC,b,y);
203: MatDestroy(sC);
205: /* Check the error */
206: VecAXPY(y,neg_one,x);
207: VecNorm(y,NORM_2,&norm2);
208: if (displ>0){PetscPrintf(PETSC_COMM_SELF," lvl: %d, error: %G\n", lvl,norm2);}
209: err[i++] = norm2;
210: }
211: }
212:
213: /* Test baij matrix A */
214: if (TestBAIJ){
215: if (displ>0){PetscPrintf(PETSC_COMM_SELF,"BAIJ: \n");}
216: i = 0;
217: for (lvl=-1; lvl<10; lvl++){
218: if (lvl==-1) { /* Cholesky factor */
219: factinfo.fill = 5.0;
220: MatCholeskyFactorSymbolic(A,perm,&factinfo,&sC);
221: } else { /* incomplete Cholesky factor */
222: factinfo.fill = 5.0;
223: factinfo.levels = lvl;
224: MatICCFactorSymbolic(A,perm,&factinfo,&sC);
225: }
226: MatCholeskyFactorNumeric(A,&factinfo,&sC);
228: MatMult(A,x,b);
229: MatSolve(sC,b,y);
230: MatDestroy(sC);
232: /* Check the error */
233: VecAXPY(y,neg_one,x);
234: VecNorm(y,NORM_2,&norm2);
235: if (displ>0){PetscPrintf(PETSC_COMM_SELF," lvl: %d, error: %G\n", lvl,norm2);}
236: err[i++] = norm2;
237: }
238: }
240: /* Test sbaij matrix sA */
241: if (displ>0){PetscPrintf(PETSC_COMM_SELF,"SBAIJ: \n");}
242: i = 0;
243: for (lvl=-1; lvl<10; lvl++){
244: if (lvl==-1) { /* Cholesky factor */
245: factinfo.fill = 5.0;
246: MatCholeskyFactorSymbolic(sA,perm,&factinfo,&sC);
247: } else { /* incomplete Cholesky factor */
248: factinfo.fill = 5.0;
249: factinfo.levels = lvl;
250: MatICCFactorSymbolic(sA,perm,&factinfo,&sC);
251: }
252: MatCholeskyFactorNumeric(sA,&factinfo,&sC);
253: MatMult(sA,x,b);
254: MatSolve(sC,b,y);
256: /* Test MatSolves() */
257: if (bs == 1) {
258: Vecs xx,bb;
259: VecsCreateSeq(PETSC_COMM_SELF,n,4,&xx);
260: VecsDuplicate(xx,&bb);
261: MatSolves(sC,bb,xx);
262: VecsDestroy(xx);
263: VecsDestroy(bb);
264: }
265: MatDestroy(sC);
267: /* Check the error */
268: VecAXPY(y,neg_one,x);
269: VecNorm(y,NORM_2,&norm2);
270: if (displ>0){PetscPrintf(PETSC_COMM_SELF," lvl: %d, error: %G\n", lvl,norm2); }
271: err[i] -= norm2;
272: if (err[i] > tol) SETERRQ2(PETSC_ERR_USER," level: %d, err: %G\n", lvl,err[i]);
273: }
275: ISDestroy(perm);
276: PetscFree(ip_ptr);
277: MatDestroy(A);
278: MatDestroy(sA);
279: VecDestroy(x);
280: VecDestroy(y);
281: VecDestroy(b);
282: PetscRandomDestroy(rdm);
284: PetscFinalize();
285: return 0;
286: }