Actual source code: ex74.c
1: /*$Id: ex74.c,v 1.47 2001/08/07 21:30:08 bsmith Exp $*/
3: static char help[] = "Tests the various sequential routines in MatSBAIJ format.n";
5: #include petscmat.h
7: int main(int argc,char **args)
8: {
9: Vec x,y,b,s1,s2;
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,inc;
14: int lf; /* level of fill for icc */
15: int *cols1,*cols2;
16: PetscReal norm1,norm2,tol=1.e-10,fill;
17: PetscScalar neg_one = -1.0,four=4.0,value[3],alpha=0.1;
18: PetscScalar *vr1,*vr2,*vr1_wk,*vr2_wk;
19: IS perm, isrow, iscol;
20: PetscRandom rand;
21: PetscTruth getrow=PETSC_FALSE;
22: MatInfo minfo1,minfo2;
24: PetscInitialize(&argc,&args,(char *)0,help);
25: MPI_Comm_size(PETSC_COMM_WORLD,&size);
26: if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
27: PetscOptionsGetInt(PETSC_NULL,"-bs",&bs,PETSC_NULL);
28: PetscOptionsGetInt(PETSC_NULL,"-mbs",&mbs,PETSC_NULL);
30: n = mbs*bs;
31: ierr=MatCreateSeqBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,PETSC_NULL, &A);
32: ierr=MatCreateSeqSBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,PETSC_NULL, &sA);
34: /* Test MatGetOwnershipRange() */
35: MatGetOwnershipRange(A,&I,&J);
36: MatGetOwnershipRange(sA,&i,&j);
37: if (i-I || j-J){
38: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetOwnershipRange() in MatSBAIJ formatn");
39: }
41: /* Assemble matrix */
42: if (bs == 1){
43: PetscOptionsGetInt(PETSC_NULL,"-test_problem",&prob,PETSC_NULL);
44: if (prob == 1){ /* tridiagonal matrix */
45: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
46: for (i=1; i<n-1; i++) {
47: col[0] = i-1; col[1] = i; col[2] = i+1;
48: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
49: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
50: }
51: i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;
52: value[0]= 0.1; value[1]=-1; value[2]=2;
53: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
54: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
56: i = 0; col[0] = 0; col[1] = 1; col[2]=n-1;
57: value[0] = 2.0; value[1] = -1.0; value[2]=0.1;
58: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
59: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
60: }
61: else if (prob ==2){ /* matrix for the five point stencil */
62: n1 = (int) (sqrt((PetscReal)n) + 0.001);
63: if (n1*n1 - n) SETERRQ(PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive interger!");
64: for (i=0; i<n1; i++) {
65: for (j=0; j<n1; j++) {
66: I = j + n1*i;
67: if (i>0) {
68: J = I - n1;
69: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
70: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
71: }
72: if (i<n1-1) {
73: J = I + n1;
74: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
75: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
76: }
77: if (j>0) {
78: J = I - 1;
79: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
80: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
81: }
82: if (j<n1-1) {
83: J = I + 1;
84: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
85: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
86: }
87: MatSetValues(A,1,&I,1,&I,&four,INSERT_VALUES);
88: MatSetValues(sA,1,&I,1,&I,&four,INSERT_VALUES);
89: }
90: }
91: }
92: }
93: else { /* bs > 1 */
94: for (block=0; block<n/bs; block++){
95: /* diagonal blocks */
96: value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
97: for (i=1+block*bs; i<bs-1+block*bs; i++) {
98: col[0] = i-1; col[1] = i; col[2] = i+1;
99: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
100: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
101: }
102: i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
103: value[0]=-1.0; value[1]=4.0;
104: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
105: MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
107: i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
108: value[0]=4.0; value[1] = -1.0;
109: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
110: MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
111: }
112: /* off-diagonal blocks */
113: value[0]=-1.0;
114: for (i=0; i<(n/bs-1)*bs; i++){
115: col[0]=i+bs;
116: MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);
117: MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES);
118: col[0]=i; row=i+bs;
119: MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
120: MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);
121: }
122: }
123: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
124: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
125: /* PetscPrintf(PETSC_COMM_SELF,"n The Matrix: n");
126: MatView(A, PETSC_VIEWER_DRAW_WORLD);
127: MatView(A, PETSC_VIEWER_STDOUT_WORLD); */
129: MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY);
130: MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY);
131: /* PetscPrintf(PETSC_COMM_SELF,"n Symmetric Part of Matrix: n");
132: MatView(sA, PETSC_VIEWER_DRAW_WORLD);
133: MatView(sA, PETSC_VIEWER_STDOUT_WORLD);
134: */
136: /* Test MatNorm(), MatDuplicate() */
137: MatNorm(A,NORM_FROBENIUS,&norm1);
138: MatDuplicate(sA,MAT_COPY_VALUES,&sC);
139: MatNorm(sC,NORM_FROBENIUS,&norm2);
140: MatDestroy(sC);
141: norm1 -= norm2;
142: if (norm1<-tol || norm1>tol){
143: PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm(), fnorm1-fnorm2=%16.14en",norm1);
144: }
145: MatNorm(A,NORM_INFINITY,&norm1);
146: MatNorm(sA,NORM_INFINITY,&norm2);
147: norm1 -= norm2;
148: if (norm1<-tol || norm1>tol){
149: PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm(), inf_norm1-inf_norm2=%16.14en",norm1);
150: }
151:
152: /* Test MatGetInfo(), MatGetSize(), MatGetBlockSize() */
153: MatGetInfo(A,MAT_LOCAL,&minfo1);
154: MatGetInfo(sA,MAT_LOCAL,&minfo2);
155: /*
156: printf("matrix nonzeros (BAIJ format) = %d, allocated nonzeros= %dn", (int)minfo1.nz_used,(int)minfo1.nz_allocated);
157: printf("matrix nonzeros(SBAIJ format) = %d, allocated nonzeros= %dn", (int)minfo2.nz_used,(int)minfo2.nz_allocated);
158: */
159: i = (int) (minfo1.nz_used - minfo2.nz_used);
160: j = (int) (minfo2.nz_allocated - minfo2.nz_used);
161: if (i<0 || j<0) {
162: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetInfo()n");
163: }
165: MatGetSize(A,&I,&J);
166: MatGetSize(sA,&i,&j);
167: if (i-I || j-J) {
168: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetSize()n");
169: }
170:
171: MatGetBlockSize(A, &I);
172: MatGetBlockSize(sA, &i);
173: if (i-I){
174: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetBlockSize()n");
175: }
177: /* Test MatGetRow() */
178: if (getrow){
179: row = n/2;
180: PetscMalloc(n*sizeof(PetscScalar),&vr1);
181: vr1_wk = vr1;
182: PetscMalloc(n*sizeof(PetscScalar),&vr2);
183: vr2_wk = vr2;
184: MatGetRow(A,row,&J,&cols1,&vr1);
185: vr1_wk += J-1;
186: MatGetRow(sA,row,&j,&cols2,&vr2);
187: vr2_wk += j-1;
188: VecCreateSeq(PETSC_COMM_SELF,j,&x);
189:
190: for (i=j-1; i>-1; i--){
191: VecSetValue(x,i,*vr2_wk - *vr1_wk,INSERT_VALUES);
192: vr2_wk--; vr1_wk--;
193: }
194: VecNorm(x,NORM_1,&norm2);
195: if (norm2<-tol || norm2>tol) {
196: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetRow()n");
197: }
198: VecDestroy(x);
199: MatRestoreRow(A,row,&J,&cols1,&vr1);
200: MatRestoreRow(sA,row,&j,&cols2,&vr2);
201: PetscFree(vr1);
202: PetscFree(vr2);
204: /* Test GetSubMatrix() */
205: /* get a submatrix consisting of every next block row and column of the original matrix */
206: /* for symm. matrix, iscol=isrow. */
207: PetscMalloc(n*sizeof(IS),&isrow);
208: PetscMalloc(n*sizeof(int),&ip_ptr);
209: j = 0;
210: for (n1=0; n1<mbs; n1 += 2){ /* n1: block row */
211: for (i=0; i<bs; i++) ip_ptr[j++] = n1*bs + i;
212: }
213: ISCreateGeneral(PETSC_COMM_SELF, j, ip_ptr, &isrow);
214: /* ISView(isrow, PETSC_VIEWER_STDOUT_SELF); */
215:
216: MatGetSubMatrix(sA,isrow,isrow,PETSC_DECIDE,MAT_INITIAL_MATRIX,&sC);
217: ISDestroy(isrow);
218: PetscFree(ip_ptr);
219: printf("sA =n");
220: MatView(sA,PETSC_VIEWER_STDOUT_WORLD);
221: printf("submatrix of sA =n");
222: MatView(sC,PETSC_VIEWER_STDOUT_WORLD);
223: MatDestroy(sC);
224: }
226: /* Test MatDiagonalScale(), MatGetDiagonal(), MatScale() */
227: PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&rand);
228: VecCreateSeq(PETSC_COMM_SELF,n,&x);
229: VecDuplicate(x,&s1);
230: VecDuplicate(x,&s2);
231: VecDuplicate(x,&y);
232: VecDuplicate(x,&b);
234: VecSetRandom(rand,x);
235: MatDiagonalScale(A,x,x);
236: MatDiagonalScale(sA,x,x);
238: MatGetDiagonal(A,s1);
239: MatGetDiagonal(sA,s2);
240: VecNorm(s1,NORM_1,&norm1);
241: VecNorm(s2,NORM_1,&norm2);
242: norm1 -= norm2;
243: if (norm1<-tol || norm1>tol) {
244: PetscPrintf(PETSC_COMM_SELF,"Error:MatGetDiagonal() n");
245: }
247: MatScale(&alpha,A);
248: MatScale(&alpha,sA);
250: /* Test MatGetRowMax() */
251: MatGetRowMax(A,s1);
252: MatGetRowMax(sA,s2);
253: VecNorm(s1,NORM_1,&norm1);
254: VecNorm(s2,NORM_1,&norm2);
255: norm1 -= norm2;
256: if (norm1<-tol || norm1>tol) {
257: PetscPrintf(PETSC_COMM_SELF,"Error:MatGetRowMax() n");
258: }
260: /* Test MatMult(), MatMultAdd() */
261: for (i=0; i<40; i++) {
262: VecSetRandom(rand,x);
263: MatMult(A,x,s1);
264: MatMult(sA,x,s2);
265: VecNorm(s1,NORM_1,&norm1);
266: VecNorm(s2,NORM_1,&norm2);
267: norm1 -= norm2;
268: if (norm1<-tol || norm1>tol) {
269: PetscPrintf(PETSC_COMM_SELF,"Error: MatMult(), MatDiagonalScale() or MatScale()n");
270: }
271: }
273: for (i=0; i<40; i++) {
274: VecSetRandom(rand,x);
275: VecSetRandom(rand,y);
276: MatMultAdd(A,x,y,s1);
277: MatMultAdd(sA,x,y,s2);
278: VecNorm(s1,NORM_1,&norm1);
279: VecNorm(s2,NORM_1,&norm2);
280: norm1 -= norm2;
281: if (norm1<-tol || norm1>tol) {
282: PetscPrintf(PETSC_COMM_SELF,"Error:MatMultAdd(), MatDiagonalScale() or MatScale() n");
283: }
284: }
286: /* Test MatCholeskyFactor(), MatICCFactor() with natural ordering */
287: MatGetOrdering(A,MATORDERING_NATURAL,&perm,&iscol);
288: ISDestroy(iscol);
289: norm1 = tol;
290: inc = bs;
291: for (lf=-1; lf<10; lf += inc){
292: if (lf==-1) { /* Cholesky factor */
293: fill = 5.0;
294: MatCholeskyFactorSymbolic(sA,perm,fill,&sC);
295: } else { /* incomplete Cholesky factor */
296: fill = 5.0;
297: MatICCFactorSymbolic(sA,perm,fill,lf,&sC);
298: }
299: MatCholeskyFactorNumeric(sA,&sC);
300: /* MatView(sC, PETSC_VIEWER_DRAW_WORLD); */
301:
302: MatMult(sA,x,b);
303: MatSolve(sC,b,y);
304: MatDestroy(sC);
305:
306: /* Check the error */
307: VecAXPY(&neg_one,x,y);
308: VecNorm(y,NORM_2,&norm2);
309: /* printf("lf: %d, error: %gn", lf,norm2); */
310: if (10*norm1 < norm2 && lf-inc != -1){
311: PetscPrintf(PETSC_COMM_SELF,"lf=%d, %d, Norm of error=%g, %gn",lf-inc,lf,norm1,norm2);
312: }
313: norm1 = norm2;
314: if (norm2 < tol && lf != -1) break;
315: }
317: ISDestroy(perm);
319: MatDestroy(A);
320: MatDestroy(sA);
321: VecDestroy(x);
322: VecDestroy(y);
323: VecDestroy(s1);
324: VecDestroy(s2);
325: VecDestroy(b);
326: PetscRandomDestroy(rand);
328: PetscFinalize();
329: return 0;
330: }