Actual source code: ex75.c
2: /* Program usage: mpirun -np <procs> ex75 [-help] [all PETSc options] */
4: static char help[] = "Tests the vatious routines in MatMPISBAIJ format.\n";
6: #include petscmat.h
10: int main(int argc,char **args)
11: {
12: Vec x,y,u,s1,s2;
13: Mat A,sA,sB;
14: PetscRandom rctx;
15: PetscReal r1,r2,tol=1.e-10;
16: PetscScalar one=1.0, neg_one=-1.0, value[3], four=4.0,alpha=0.1;
17: const PetscScalar *vr;
18: PetscInt n,col[3],n1,block,row,i,j,i2,j2,I,J,ncols,rstart,rend,bs=1, mbs=16, d_nz=3, o_nz=3, prob=2;
19: PetscErrorCode ierr;
20: PetscMPIInt size,rank;
21: const PetscInt *cols;
22: PetscTruth flg;
23: MatType type;
25: PetscInitialize(&argc,&args,(char *)0,help);
26: PetscOptionsGetInt(PETSC_NULL,"-mbs",&mbs,PETSC_NULL);
27: PetscOptionsGetInt(PETSC_NULL,"-bs",&bs,PETSC_NULL);
28:
29: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
30: MPI_Comm_size(PETSC_COMM_WORLD,&size);
31:
32: n = mbs*bs;
33:
34: /* Assemble MPISBAIJ matrix sA */
35: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,&sA);
36: MatSetType(sA,MATSBAIJ);
37: /* -mat_type <seqsbaij_derived type>, e.g., mpisbaijspooles, sbaijmumps */
38: MatSetFromOptions(sA);
39: MatGetType(sA,&type);
40: /* printf(" mattype: %s\n",type); */
41: MatMPISBAIJSetPreallocation(sA,bs,d_nz,PETSC_NULL,o_nz,PETSC_NULL);
43: if (bs == 1){
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: /* PetscTrValid(0,0,0,0); */
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(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(sA,1,&i,3,col,value,INSERT_VALUES);
58: }
59: else if (prob ==2){ /* matrix for the five point stencil */
60: n1 = (int) sqrt((double)n);
61: if (n1*n1 != n){
62: SETERRQ(PETSC_ERR_ARG_SIZ,"n must be a perfect square of n1");
63: }
64:
65: for (i=0; i<n1; i++) {
66: for (j=0; j<n1; j++) {
67: I = j + n1*i;
68: if (i>0) {J = I - n1; MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);}
69: if (i<n1-1) {J = I + n1; MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);}
70: if (j>0) {J = I - 1; MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);}
71: if (j<n1-1) {J = I + 1; MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);}
72: MatSetValues(sA,1,&I,1,&I,&four,INSERT_VALUES);
73: }
74: }
75: }
76: } /* end of if (bs == 1) */
77: else { /* bs > 1 */
78: for (block=0; block<n/bs; block++){
79: /* diagonal blocks */
80: value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
81: for (i=1+block*bs; i<bs-1+block*bs; i++) {
82: col[0] = i-1; col[1] = i; col[2] = i+1;
83: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
84: }
85: i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
86: value[0]=-1.0; value[1]=4.0;
87: MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
89: i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
90: value[0]=4.0; value[1] = -1.0;
91: MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
92: }
93: /* off-diagonal blocks */
94: value[0]=-1.0;
95: for (i=0; i<(n/bs-1)*bs; i++){
96: col[0]=i+bs;
97: MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES);
98: col[0]=i; row=i+bs;
99: MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);
100: }
101: }
102: MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY);
103: MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY);
105: /* Test MatView() */
106: /*
107: MatView(sA, PETSC_VIEWER_STDOUT_WORLD);
108: MatView(sA, PETSC_VIEWER_DRAW_WORLD);
109: */
110: /* Assemble MPIBAIJ matrix A */
111: MatCreateMPIBAIJ(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,n,n,d_nz,PETSC_NULL,o_nz,PETSC_NULL,&A);
113: if (bs == 1){
114: if (prob == 1){ /* tridiagonal matrix */
115: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
116: for (i=1; i<n-1; i++) {
117: col[0] = i-1; col[1] = i; col[2] = i+1;
118: /* PetscTrValid(0,0,0,0); */
119: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
120: }
121: i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;
122: value[0]= 0.1; value[1]=-1; value[2]=2;
123: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
125: i = 0; col[0] = 0; col[1] = 1; col[2]=n-1;
126: value[0] = 2.0; value[1] = -1.0; value[2]=0.1;
127: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
128: }
129: else if (prob ==2){ /* matrix for the five point stencil */
130: n1 = (int) sqrt((double)n);
131: for (i=0; i<n1; i++) {
132: for (j=0; j<n1; j++) {
133: I = j + n1*i;
134: if (i>0) {J = I - n1; MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);}
135: if (i<n1-1) {J = I + n1; MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);}
136: if (j>0) {J = I - 1; MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);}
137: if (j<n1-1) {J = I + 1; MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);}
138: MatSetValues(A,1,&I,1,&I,&four,INSERT_VALUES);
139: }
140: }
141: }
142: } /* end of if (bs == 1) */
143: else { /* bs > 1 */
144: for (block=0; block<n/bs; block++){
145: /* diagonal blocks */
146: value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
147: for (i=1+block*bs; i<bs-1+block*bs; i++) {
148: col[0] = i-1; col[1] = i; col[2] = i+1;
149: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
150: }
151: i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
152: value[0]=-1.0; value[1]=4.0;
153: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
155: i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
156: value[0]=4.0; value[1] = -1.0;
157: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
158: }
159: /* off-diagonal blocks */
160: value[0]=-1.0;
161: for (i=0; i<(n/bs-1)*bs; i++){
162: col[0]=i+bs;
163: MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);
164: col[0]=i; row=i+bs;
165: MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
166: }
167: }
168: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
169: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
171: /* Test MatGetSize(), MatGetLocalSize() */
172: MatGetSize(sA, &i,&j); MatGetSize(A, &i2,&j2);
173: i -= i2; j -= j2;
174: if (i || j) {
175: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatGetSize()\n",rank);
176: PetscSynchronizedFlush(PETSC_COMM_WORLD);
177: }
178:
179: MatGetLocalSize(sA, &i,&j); MatGetLocalSize(A, &i2,&j2);
180: i2 -= i; j2 -= j;
181: if (i2 || j2) {
182: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatGetLocalSize()\n",rank);
183: PetscSynchronizedFlush(PETSC_COMM_WORLD);
184: }
186: /* vectors */
187: /*--------------------*/
188: /* i is obtained from MatGetLocalSize() */
189: VecCreate(PETSC_COMM_WORLD,&x);
190: VecSetSizes(x,i,PETSC_DECIDE);
191: VecSetFromOptions(x);
192: VecDuplicate(x,&y);
193: VecDuplicate(x,&u);
194: VecDuplicate(x,&s1);
195: VecDuplicate(x,&s2);
197: PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rctx);
198: VecSetRandom(rctx,x);
199: VecSet(&one,u);
201: /* Test MatNorm() */
202: MatNorm(A,NORM_FROBENIUS,&r1);
203: MatNorm(sA,NORM_FROBENIUS,&r2);
204: r1 -= r2;
205: if (r1<-tol || r1>tol){
206: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatNorm(), A_fnorm - sA_fnorm = %16.14e\n",rank,r1);
207: PetscSynchronizedFlush(PETSC_COMM_WORLD);
208: }
210: /* Test MatGetOwnershipRange() */
211: MatGetOwnershipRange(sA,&rstart,&rend);
212: MatGetOwnershipRange(A,&i2,&j2);
213: i2 -= rstart; j2 -= rend;
214: if (i2 || j2) {
215: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MaGetOwnershipRange()\n",rank);
216: PetscSynchronizedFlush(PETSC_COMM_WORLD);
217: }
219: /* Test MatGetRow(): can only obtain rows associated with the given processor */
220: for (i=rstart; i<rstart+1; i++) {
221: MatGetRow(sA,i,&ncols,&cols,&vr);
222: MatRestoreRow(sA,i,&ncols,&cols,&vr);
223: }
225: /* Test MatDiagonalScale(), MatGetDiagonal(), MatScale() */
226: MatDiagonalScale(A,x,x);
227: MatDiagonalScale(sA,x,x);
228: MatGetDiagonal(A,s1);
229: MatGetDiagonal(sA,s2);
230: VecNorm(s1,NORM_1,&r1);
231: VecNorm(s2,NORM_1,&r2);
232: r1 -= r2;
233: if (r1<-tol || r1>tol) {
234: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDiagonalScale() or MatGetDiagonal(), r1=%g \n",rank,r1);
235: PetscSynchronizedFlush(PETSC_COMM_WORLD);
236: }
237:
238: MatScale(&alpha,A);
239: MatScale(&alpha,sA);
241: /* Test MatGetRowMax() */
242: MatGetRowMax(A,s1);
243: MatGetRowMax(sA,s2);
245: VecNorm(s1,NORM_1,&r1);
246: VecNorm(s2,NORM_1,&r2);
247: r1 -= r2;
248: if (r1<-tol || r1>tol) {
249: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetRowMax() \n");
250: }
252: /* Test MatMult(), MatMultAdd() */
253: for (i=0; i<10; i++) {
254: VecSetRandom(rctx,x);
255: MatMult(A,x,s1);
256: MatMult(sA,x,s2);
257: VecNorm(s1,NORM_1,&r1);
258: VecNorm(s2,NORM_1,&r2);
259: r1 -= r2;
260: if (r1<-tol || r1>tol) {
261: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMult() or MatScale(), err=%g\n",rank,r1);
262: PetscSynchronizedFlush(PETSC_COMM_WORLD);
263: }
264: }
266: for (i=0; i<10; i++) {
267: VecSetRandom(rctx,x);
268: VecSetRandom(rctx,y);
269: MatMultAdd(A,x,y,s1);
270: MatMultAdd(sA,x,y,s2);
271: VecNorm(s1,NORM_1,&r1);
272: VecNorm(s2,NORM_1,&r2);
273: r1 -= r2;
274: if (r1<-tol || r1>tol) {
275: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMultAdd(), err=%g \n",rank,r1);
276: PetscSynchronizedFlush(PETSC_COMM_WORLD);
277: }
278: }
280: /* Test MatMultTranspose(), MatMultTransposeAdd() */
281: for (i=0; i<10; i++) {
282: VecSetRandom(rctx,x);
283: MatMultTranspose(A,x,s1);
284: MatMultTranspose(sA,x,s2);
285: VecNorm(s1,NORM_1,&r1);
286: VecNorm(s2,NORM_1,&r2);
287: r1 -= r2;
288: if (r1<-tol || r1>tol) {
289: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMult() or MatScale(), err=%g\n",rank,r1);
290: PetscSynchronizedFlush(PETSC_COMM_WORLD);
291: }
292: }
293: for (i=0; i<10; i++) {
294: VecSetRandom(rctx,x);
295: VecSetRandom(rctx,y);
296: MatMultTransposeAdd(A,x,y,s1);
297: MatMultTransposeAdd(sA,x,y,s2);
298: VecNorm(s1,NORM_1,&r1);
299: VecNorm(s2,NORM_1,&r2);
300: r1 -= r2;
301: if (r1<-tol || r1>tol) {
302: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatMultAdd(), err=%g \n",rank,r1);
303: PetscSynchronizedFlush(PETSC_COMM_WORLD);
304: }
305: }
306: /* MatView(sA, PETSC_VIEWER_STDOUT_WORLD); */
307: /* MatView(sA, PETSC_VIEWER_DRAW_WORLD); */
309: /* Test MatDuplicate() */
310: MatDuplicate(sA,MAT_COPY_VALUES,&sB);
311: MatEqual(sA,sB,&flg);
312: if (!flg){
313: PetscPrintf(PETSC_COMM_WORLD," Error in MatDuplicate(), sA != sB \n");
314: }
315: for (i=0; i<5; i++) {
316: VecSetRandom(rctx,x);
317: MatMult(A,x,s1);
318: MatMult(sB,x,s2);
319: VecNorm(s1,NORM_1,&r1);
320: VecNorm(s2,NORM_1,&r2);
321: r1 -= r2;
322: if (r1<-tol || r1>tol) {
323: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDuplicate() or MatMult(), err=%g\n",rank,r1);
324: PetscSynchronizedFlush(PETSC_COMM_WORLD);
325: }
326: }
328: for (i=0; i<5; i++) {
329: VecSetRandom(rctx,x);
330: VecSetRandom(rctx,y);
331: MatMultAdd(A,x,y,s1);
332: MatMultAdd(sB,x,y,s2);
333: VecNorm(s1,NORM_1,&r1);
334: VecNorm(s2,NORM_1,&r2);
335: r1 -= r2;
336: if (r1<-tol || r1>tol) {
337: PetscSynchronizedPrintf(PETSC_COMM_WORLD,"[%d], Error: MatDuplicate() or MatMultAdd(), err=%g \n",rank,r1);
338: PetscSynchronizedFlush(PETSC_COMM_WORLD);
339: }
340: }
341: MatDestroy(sB);
342:
343: VecDestroy(u);
344: VecDestroy(x);
345: VecDestroy(y);
346: VecDestroy(s1);
347: VecDestroy(s2);
348: MatDestroy(sA);
349: MatDestroy(A);
350: PetscRandomDestroy(rctx);
351:
352: PetscFinalize();
353: return 0;
354: }