Actual source code: ex77.c
1: /*$Id: ex77.c,v 1.12 2001/08/07 21:30:08 bsmith Exp $*/
3: static char help[] = "Tests the various sequential routines in MatSBAIJ format. Same as ex74.c except diagonal entries of the matrices are zeros.n";
5: #include petscmat.h
7: /* extern int MatReorderingSeqSBAIJ(Mat,IS); */
9: int main(int argc,char **args)
10: {
11: Vec x,y,b,s1,s2;
12: Mat A; /* linear system matrix */
13: Mat sA,sC; /* symmetric part of the matrices */
15: int n,mbs=16,bs=1,nz=3,prob=2;
16: PetscScalar neg_one = -1.0,value[3],alpha=0.1;
17: int ierr,i,j,col[3],size,row,I,J,n1,*ip_ptr;
18: IS ip, isrow, iscol;
19: PetscRandom rand;
20: PetscTruth reorder=PETSC_FALSE,getrow=PETSC_FALSE;
21: MatInfo minfo1,minfo2;
22: PetscScalar *vr1,*vr2,*vr1_wk,*vr2_wk;
23: int *cols1,*cols2;
24: PetscReal norm1,norm2,tol=1.e-10;
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);
32: n = mbs*bs;
33: ierr=MatCreateSeqBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,PETSC_NULL, &A);
34: ierr=MatCreateSeqSBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,PETSC_NULL, &sA);
36: /* Test MatGetOwnershipRange() */
37: MatGetOwnershipRange(A,&I,&J);
38: MatGetOwnershipRange(sA,&i,&j);
39: if (i-I || j-J){
40: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetOwnershipRange() in MatSBAIJ formatn");
41: }
43: /* Assemble matrix */
44: if (bs == 1){
45: PetscOptionsGetInt(PETSC_NULL,"-test_problem",&prob,PETSC_NULL);
46: if (prob == 1){ /* tridiagonal matrix */
47: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
48: for (i=1; i<n-1; i++) {
49: col[0] = i-1; col[1] = i; col[2] = i+1;
50: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
51: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
52: }
53: i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;
54: value[0]= 0.1; value[1]=-1; value[2]=2;
55: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
56: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
58: i = 0; col[0] = 0; col[1] = 1; col[2]=n-1;
59: value[0] = 2.0; value[1] = -1.0; value[2]=0.1;
60: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
61: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
62: }
63: else if (prob ==2){ /* matrix for the five point stencil */
64: n1 = (int) (sqrt((PetscReal)n) + 0.001);
65: if (n1*n1 - n) SETERRQ(PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive interger!");
66: for (i=0; i<n1; i++) {
67: for (j=0; j<n1; j++) {
68: I = j + n1*i;
69: if (i>0) {
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 (i<n1-1) {
75: J = I + n1;
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>0) {
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: if (j<n1-1) {
85: J = I + 1;
86: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
87: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
88: }
89: /*
90: MatSetValues(A,1,&I,1,&I,&four,INSERT_VALUES);
91: MatSetValues(sA,1,&I,1,&I,&four,INSERT_VALUES);
92: */
93: }
94: }
95: }
96: }
97: else { /* bs > 1 */
98: #ifdef DIAGB
99: for (block=0; block<n/bs; block++){
100: /* diagonal blocks */
101: value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
102: for (i=1+block*bs; i<bs-1+block*bs; i++) {
103: col[0] = i-1; col[1] = i; col[2] = i+1;
104: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
105: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
106: }
107: i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
108: value[0]=-1.0; value[1]=4.0;
109: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
110: MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
112: i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
113: value[0]=4.0; value[1] = -1.0;
114: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
115: MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
116: }
117: #endif
118: /* off-diagonal blocks */
119: value[0]=-1.0;
120: for (i=0; i<(n/bs-1)*bs; i++){
121: col[0]=i+bs;
122: MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);
123: MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES);
124: col[0]=i; row=i+bs;
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, VIEWER_DRAW_WORLD);
133: MatView(A, 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, VIEWER_DRAW_WORLD);
139: MatView(sA, VIEWER_STDOUT_WORLD);
140: */
142: /* Test MatNorm() */
143: MatNorm(A,NORM_FROBENIUS,&norm1);
144: MatNorm(sA,NORM_FROBENIUS,&norm2);
145: norm1 -= norm2;
146: if (norm1<-tol || norm1>tol){
147: PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm(), fnorm1-fnorm2=%16.14en",norm1);
148: }
149: MatNorm(A,NORM_INFINITY,&norm1);
150: MatNorm(sA,NORM_INFINITY,&norm2);
151: norm1 -= norm2;
152: if (norm1<-tol || norm1>tol){
153: PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm(), inf_norm1-inf_norm2=%16.14en",norm1);
154: }
155:
156: /* Test MatGetInfo(), MatGetSize(), MatGetBlockSize() */
157: MatGetInfo(A,MAT_LOCAL,&minfo1);
158: MatGetInfo(sA,MAT_LOCAL,&minfo2);
159: /*
160: printf("matrix nonzeros (BAIJ format) = %d, allocated nonzeros= %dn", (int)minfo1.nz_used,(int)minfo1.nz_allocated);
161: printf("matrix nonzeros(SBAIJ format) = %d, allocated nonzeros= %dn", (int)minfo2.nz_used,(int)minfo2.nz_allocated);
162: */
163: i = (int) (minfo1.nz_used - minfo2.nz_used);
164: j = (int) (minfo1.nz_allocated - minfo2.nz_allocated);
165: if (i<0 || j<0) {
166: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetInfo()n");
167: }
169: MatGetSize(A,&I,&J);
170: MatGetSize(sA,&i,&j);
171: if (i-I || j-J) {
172: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetSize()n");
173: }
174:
175: MatGetBlockSize(A, &I);
176: MatGetBlockSize(sA, &i);
177: if (i-I){
178: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetBlockSize()n");
179: }
181: /* Test MatGetRow() */
182: if (getrow){
183: row = n/2;
184: PetscMalloc(n*sizeof(PetscScalar),&vr1);
185: vr1_wk = vr1;
186: PetscMalloc(n*sizeof(PetscScalar),&vr2);
187: vr2_wk = vr2;
188: MatGetRow(A,row,&J,&cols1,&vr1);
189: vr1_wk += J-1;
190: MatGetRow(sA,row,&j,&cols2,&vr2);
191: vr2_wk += j-1;
192: VecCreateSeq(PETSC_COMM_SELF,j,&x);
193:
194: for (i=j-1; i>-1; i--){
195: VecSetValue(x,i,*vr2_wk - *vr1_wk,INSERT_VALUES);
196: vr2_wk--; vr1_wk--;
197: }
198: VecNorm(x,NORM_1,&norm2);
199: if (norm2<-tol || norm2>tol) {
200: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetRow()n");
201: }
202: VecDestroy(x);
203: MatRestoreRow(A,row,&J,&cols1,&vr1);
204: MatRestoreRow(sA,row,&j,&cols2,&vr2);
205: PetscFree(vr1);
206: PetscFree(vr2);
208: /* Test GetSubMatrix() */
209: /* get a submatrix consisting of every next block row and column of the original matrix */
210: /* for symm. matrix, iscol=isrow. */
211: PetscMalloc(n*sizeof(IS),&isrow);
212: PetscMalloc(n*sizeof(int),&ip_ptr);
213: j = 0;
214: for (n1=0; n1<mbs; n1 += 2){ /* n1: block row */
215: for (i=0; i<bs; i++) ip_ptr[j++] = n1*bs + i;
216: }
217: ISCreateGeneral(PETSC_COMM_SELF, j, ip_ptr, &isrow);
218: /* ISView(isrow, VIEWER_STDOUT_SELF); */
219:
220: MatGetSubMatrix(sA,isrow,isrow,PETSC_DECIDE,MAT_INITIAL_MATRIX,&sC);
221: ISDestroy(isrow);
222: PetscFree(ip_ptr);
223: printf("sA =n");
224: MatView(sA,PETSC_VIEWER_STDOUT_WORLD);
225: printf("submatrix of sA =n");
226: MatView(sC,PETSC_VIEWER_STDOUT_WORLD);
227: MatDestroy(sC);
228: }
230: /* Test MatDiagonalScale(), MatGetDiagonal(), MatScale() */
231: PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&rand);
232: VecCreateSeq(PETSC_COMM_SELF,n,&x);
233: VecDuplicate(x,&s1);
234: VecDuplicate(x,&s2);
235: VecDuplicate(x,&y);
236: VecDuplicate(x,&b);
238: VecSetRandom(rand,x);
240: MatDiagonalScale(A,x,x);
241: MatDiagonalScale(sA,x,x);
243: MatGetDiagonal(A,s1);
244: MatGetDiagonal(sA,s2);
245: VecNorm(s1,NORM_1,&norm1);
246: VecNorm(s2,NORM_1,&norm2);
247: norm1 -= norm2;
248: if (norm1<-tol || norm1>tol) {
249: PetscPrintf(PETSC_COMM_SELF,"Error:MatGetDiagonal() n");
250: }
252: MatScale(&alpha,A);
253: MatScale(&alpha,sA);
255: /* Test MatMult(), MatMultAdd() */
256: for (i=0; i<40; i++) {
257: VecSetRandom(rand,x);
258: MatMult(A,x,s1);
259: MatMult(sA,x,s2);
260: VecNorm(s1,NORM_1,&norm1);
261: VecNorm(s2,NORM_1,&norm2);
262: norm1 -= norm2;
263: if (norm1<-tol || norm1>tol) {
264: PetscPrintf(PETSC_COMM_SELF,"Error: MatMult(), MatDiagonalScale() or MatScale()n");
265: }
266: }
268: for (i=0; i<40; i++) {
269: VecSetRandom(rand,x);
270: VecSetRandom(rand,y);
271: MatMultAdd(A,x,y,s1);
272: MatMultAdd(sA,x,y,s2);
273: VecNorm(s1,NORM_1,&norm1);
274: VecNorm(s2,NORM_1,&norm2);
275: norm1 -= norm2;
276: if (norm1<-tol || norm1>tol) {
277: PetscPrintf(PETSC_COMM_SELF,"Error:MatMultAdd(), MatDiagonalScale() or MatScale() n");
278: }
279: }
281: /* Test MatReordering() */
282: MatGetOrdering(A,MATORDERING_NATURAL,&isrow,&iscol);
283: ip = isrow;
285: if (reorder){
286: ISGetIndices(ip,&ip_ptr);
287: i = ip_ptr[1]; ip_ptr[1] = ip_ptr[mbs-2]; ip_ptr[mbs-2] = i;
288: i = ip_ptr[0]; ip_ptr[0] = ip_ptr[mbs-1]; ip_ptr[mbs-1] = i;
289: ierr= ISRestoreIndices(ip,&ip_ptr);
291: MatReorderingSeqSBAIJ(sA, ip);
292: /* ISView(ip, VIEWER_STDOUT_SELF);
293: MatView(sA,VIEWER_DRAW_SELF); */
294: }
295:
296: ISDestroy(iscol);
297: /* ISDestroy(isrow);*/
299: ISDestroy(isrow);
300: MatDestroy(A);
301: MatDestroy(sA);
302: VecDestroy(x);
303: VecDestroy(y);
304: VecDestroy(s1);
305: VecDestroy(s2);
306: VecDestroy(b);
307: PetscRandomDestroy(rand);
309: PetscFinalize();
310: return 0;
311: }