Actual source code: ex75.c
1: /*$Id: ex75.c,v 1.28 2001/09/11 16:32:50 bsmith Exp $*/
3: /* Program usage: mpirun -np <procs> ex75 [-help] [all PETSc options] */
5: static char help[] = "Tests the vatious routines in MatMPISBAIJ format.n";
7: #include petscmat.h
9: int main(int argc,char **args)
10: {
11: int bs=1, mbs=16, d_nz=3, o_nz=3, prob=2;
12: Vec x,y,u,s1,s2;
13: Mat A,sA;
14: PetscRandom rctx;
15: double r1,r2,tol=1.e-10;
16: int i,j,i2,j2,I,J,ierr;
17: PetscScalar one=1.0, neg_one=-1.0, value[3], four=4.0,alpha=0.1,*vr;
18: int n,rank,size,col[3],n1,block,row;
19: int ncols,*cols,rstart,rend;
20: IS isrow;
22: PetscInitialize(&argc,&args,(char *)0,help);
23: PetscOptionsGetInt(PETSC_NULL,"-mbs",&mbs,PETSC_NULL);
24: PetscOptionsGetInt(PETSC_NULL,"-bs",&bs,PETSC_NULL);
25:
26: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
27: MPI_Comm_size(PETSC_COMM_WORLD,&size);
28:
29: n = mbs*bs;
30:
31: /* Assemble MPISBAIJ matrix sA */
32: MatCreateMPISBAIJ(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,n,n,d_nz,PETSC_NULL,o_nz,PETSC_NULL,&sA);
34: if (bs == 1){
35: if (prob == 1){ /* tridiagonal matrix */
36: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
37: for (i=1; i<n-1; i++) {
38: col[0] = i-1; col[1] = i; col[2] = i+1;
39: /* PetscTrValid(0,0,0,0); */
40: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
41: }
42: i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;
43: value[0]= 0.1; value[1]=-1; value[2]=2;
44: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
46: i = 0; col[0] = 0; col[1] = 1; col[2]=n-1;
47: value[0] = 2.0; value[1] = -1.0; value[2]=0.1;
48: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
49: }
50: else if (prob ==2){ /* matrix for the five point stencil */
51: n1 = (int) sqrt(n);
52: if (n1*n1 != n){
53: SETERRQ(PETSC_ERR_ARG_SIZ,"n must be a perfect square of n1");
54: }
55:
56: for (i=0; i<n1; i++) {
57: for (j=0; j<n1; j++) {
58: I = j + n1*i;
59: if (i>0) {J = I - n1; MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);}
60: if (i<n1-1) {J = I + n1; MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);}
61: if (j>0) {J = I - 1; MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);}
62: if (j<n1-1) {J = I + 1; MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);}
63: MatSetValues(sA,1,&I,1,&I,&four,INSERT_VALUES);
64: }
65: }
66: }
67: } /* end of if (bs == 1) */
68: else { /* bs > 1 */
69: for (block=0; block<n/bs; block++){
70: /* diagonal blocks */
71: value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
72: for (i=1+block*bs; i<bs-1+block*bs; i++) {
73: col[0] = i-1; col[1] = i; col[2] = i+1;
74: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
75: }
76: i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
77: value[0]=-1.0; value[1]=4.0;
78: MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
80: i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
81: value[0]=4.0; value[1] = -1.0;
82: MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
83: }
84: /* off-diagonal blocks */
85: value[0]=-1.0;
86: for (i=0; i<(n/bs-1)*bs; i++){
87: col[0]=i+bs;
88: MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES);
89: col[0]=i; row=i+bs;
90: MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);
91: }
92: }
93: MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY);
94: MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY);
96: /* Test MatView() */
97: /*
98: MatView(sA, PETSC_VIEWER_STDOUT_WORLD);
99: MatView(sA, PETSC_VIEWER_DRAW_WORLD);
100: */
101: /* Assemble MPIBAIJ matrix A */
102: MatCreateMPIBAIJ(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,n,n,d_nz,PETSC_NULL,o_nz,PETSC_NULL,&A);
104: if (bs == 1){
105: if (prob == 1){ /* tridiagonal matrix */
106: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
107: for (i=1; i<n-1; i++) {
108: col[0] = i-1; col[1] = i; col[2] = i+1;
109: /* PetscTrValid(0,0,0,0); */
110: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
111: }
112: i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;
113: value[0]= 0.1; value[1]=-1; value[2]=2;
114: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
116: i = 0; col[0] = 0; col[1] = 1; col[2]=n-1;
117: value[0] = 2.0; value[1] = -1.0; value[2]=0.1;
118: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
119: }
120: else if (prob ==2){ /* matrix for the five point stencil */
121: n1 = (int) sqrt(n);
122: for (i=0; i<n1; i++) {
123: for (j=0; j<n1; j++) {
124: I = j + n1*i;
125: if (i>0) {J = I - n1; MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);}
126: if (i<n1-1) {J = I + n1; MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);}
127: if (j>0) {J = I - 1; MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);}
128: if (j<n1-1) {J = I + 1; MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);}
129: MatSetValues(A,1,&I,1,&I,&four,INSERT_VALUES);
130: }
131: }
132: }
133: } /* end of if (bs == 1) */
134: else { /* bs > 1 */
135: for (block=0; block<n/bs; block++){
136: /* diagonal blocks */
137: value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
138: for (i=1+block*bs; i<bs-1+block*bs; i++) {
139: col[0] = i-1; col[1] = i; col[2] = i+1;
140: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
141: }
142: i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
143: value[0]=-1.0; value[1]=4.0;
144: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
146: i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
147: value[0]=4.0; value[1] = -1.0;
148: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
149: }
150: /* off-diagonal blocks */
151: value[0]=-1.0;
152: for (i=0; i<(n/bs-1)*bs; i++){
153: col[0]=i+bs;
154: MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);
155: col[0]=i; row=i+bs;
156: MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
157: }
158: }
159: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
160: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
162: /* Test MatGetSize(), MatGetLocalSize() */
163: MatGetSize(sA, &i,&j); MatGetSize(A, &i2,&j2);
164: i -= i2; j -= j2;
165: if (i || j) {
166: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatGetSize()n",rank);
167: PetscSynchronizedFlush(PETSC_COMM_WORLD);
168: }
169:
170: MatGetLocalSize(sA, &i,&j); MatGetLocalSize(A, &i2,&j2);
171: i2 -= i; j2 -= j;
172: if (i2 || j2) {
173: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatGetLocalSize()n",rank);
174: PetscSynchronizedFlush(PETSC_COMM_WORLD);
175: }
177: /* vectors */
178: /*--------------------*/
179: /* i is obtained from MatGetLocalSize() */
180: VecCreate(PETSC_COMM_WORLD,&x);
181: VecSetSizes(x,i,PETSC_DECIDE);
182: VecSetFromOptions(x);
183: VecDuplicate(x,&y);
184: VecDuplicate(x,&u);
185: VecDuplicate(x,&s1);
186: VecDuplicate(x,&s2);
188: PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rctx);
189: VecSetRandom(rctx,x);
190: VecSet(&one,u);
192: /* Test MatNorm() */
193: MatNorm(A,NORM_FROBENIUS,&r1);
194: MatNorm(sA,NORM_FROBENIUS,&r2);
195: r1 -= r2;
196: if (r1<-tol || r1>tol){
197: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatNorm(), A_fnorm - sA_fnorm = %16.14en",rank,r1);
198: PetscSynchronizedFlush(PETSC_COMM_WORLD);
199: }
201: /* Test MatGetOwnershipRange() */
202: MatGetOwnershipRange(sA,&rstart,&rend);
203: MatGetOwnershipRange(A,&i2,&j2);
204: i2 -= rstart; j2 -= rend;
205: if (i2 || j2) {
206: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MaGetOwnershipRange()n",rank);
207: PetscSynchronizedFlush(PETSC_COMM_WORLD);
208: }
210: /* Test MatGetRow(): can only obtain rows associated with the given processor */
211: for (i=rstart; i<rstart+1; i++) {
212: MatGetRow(sA,i,&ncols,&cols,&vr);
213: /*
214: PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"[%d] get row %d: ",rank,i);
215: for (j=0; j<ncols; j++) {
216: PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"%d %g ",cols[j],vr[j]);
217: }
218: PetscSynchronizedFPrintf(PETSC_COMM_WORLD,stdout,"n");
219: PetscSynchronizedFlush(PETSC_COMM_WORLD);
220: */
221: MatRestoreRow(sA,i,&ncols,&cols,&vr);
222: }
224: /* Test MatDiagonalScale(), MatGetDiagonal(), MatScale() */
225: MatDiagonalScale(A,x,x);
226: MatDiagonalScale(sA,x,x);
227: MatGetDiagonal(A,s1);
228: MatGetDiagonal(sA,s2);
229: VecNorm(s1,NORM_1,&r1);
230: VecNorm(s2,NORM_1,&r2);
231: r1 -= r2;
232: if (r1<-tol || r1>tol) {
233: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDiagonalScale() or MatGetDiagonal(), r1=%g n",rank,r1);
234: PetscSynchronizedFlush(PETSC_COMM_WORLD);
235: }
236:
237: MatScale(&alpha,A);
238: MatScale(&alpha,sA);
240: /* Test MatGetRowMax() */
241: MatGetRowMax(A,s1);
242: MatGetRowMax(sA,s2);
244: VecNorm(s1,NORM_1,&r1);
245: VecNorm(s2,NORM_1,&r2);
246: /* MatView(A,PETSC_VIEWER_STDOUT_WORLD);
247: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], r1: %g, r2: %gn",rank,r1,r2);
248: PetscSynchronizedFlush(PETSC_COMM_WORLD);
249: */
250: r1 -= r2;
251: if (r1<-tol || r1>tol) {
252: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetRowMax() n");
253: }
255: /* Test MatMult(), MatMultAdd() */
256: for (i=0; i<10; i++) {
257: VecSetRandom(rctx,x);
258: MatMult(A,x,s1);
259: MatMult(sA,x,s2);
260: VecNorm(s1,NORM_1,&r1);
261: VecNorm(s2,NORM_1,&r2);
262: r1 -= r2;
263: if (r1<-tol || r1>tol) {
264: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMult() or MatScale(), err=%gn",rank,r1);
265: PetscSynchronizedFlush(PETSC_COMM_WORLD);
266: }
267: }
269: for (i=0; i<10; i++) {
270: VecSetRandom(rctx,x);
271: VecSetRandom(rctx,y);
272: MatMultAdd(A,x,y,s1);
273: MatMultAdd(sA,x,y,s2);
274: VecNorm(s1,NORM_1,&r1);
275: VecNorm(s2,NORM_1,&r2);
276: r1 -= r2;
277: if (r1<-tol || r1>tol) {
278: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMultAdd(), err=%g n",rank,r1);
279: PetscSynchronizedFlush(PETSC_COMM_WORLD);
280: }
281: }
283: /* Test MatZeroRows() */
284: ISCreateStride(PETSC_COMM_SELF,2,rstart,2,&isrow);
285: /* ISView(isrow, PETSC_VIEWER_STDOUT_SELF); */
286: MatZeroRows(sA,isrow,PETSC_NULL);
287: ISDestroy(isrow);
288: /* MatView(sA, PETSC_VIEWER_STDOUT_WORLD); */
289: /* MatView(sA, PETSC_VIEWER_DRAW_WORLD); */
290:
291: VecDestroy(u);
292: VecDestroy(x);
293: VecDestroy(y);
294: VecDestroy(s1);
295: VecDestroy(s2);
296: MatDestroy(sA);
297: MatDestroy(A);
298: PetscRandomDestroy(rctx);
299:
300: PetscFinalize();
301: return 0;
302: }