Actual source code: ex74.c
2: static char help[] = "Tests the various sequential routines in MatSBAIJ format.\n";
4: #include petscmat.h
8: int main(int argc,char **args)
9: {
10: PetscMPIInt size;
11: PetscErrorCode ierr;
12: Vec x,y,b,s1,s2;
13: Mat A; /* linear system matrix */
14: Mat sA,sB,sC; /* symmetric part of the matrices */
15: PetscInt n,mbs=16,bs=1,nz=3,prob=1,i,j,col[3],lf,block, row,I,J,n1,inc;
16: PetscReal norm1,norm2,tol=1.e-10;
17: PetscScalar neg_one = -1.0,four=4.0,value[3],alpha=0.1;
18: IS perm, iscol;
19: PetscRandom rdm;
20: PetscTruth doIcc=PETSC_TRUE;
21: MatInfo minfo1,minfo2;
22: MatFactorInfo factinfo;
23: MatType type;
25: PetscInitialize(&argc,&args,(char *)0,help);
26: MPI_Comm_size(PETSC_COMM_WORLD,&size);
27: if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
28: PetscOptionsGetInt(PETSC_NULL,"-bs",&bs,PETSC_NULL);
29: PetscOptionsGetInt(PETSC_NULL,"-mbs",&mbs,PETSC_NULL);
31: n = mbs*bs;
32: ierr=MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,n,n,nz,PETSC_NULL, &A);
33: MatCreate(PETSC_COMM_SELF,n,n,PETSC_NULL,PETSC_NULL,&sA);
34: MatSetType(sA,MATSBAIJ);
35: /* -mat_type <seqsbaij_derived type>, e.g., seqsbaijspooles, sbaijmumps */
36: MatSetFromOptions(sA);
37: MatGetType(sA,&type);
38: PetscTypeCompare((PetscObject)sA,MATSEQSBAIJ,&doIcc);
39: /* printf(" mattype: %s\n",type); */
40: MatSeqSBAIJSetPreallocation(sA,bs,nz,PETSC_NULL);
42: /* Test MatGetOwnershipRange() */
43: MatGetOwnershipRange(A,&I,&J);
44: MatGetOwnershipRange(sA,&i,&j);
45: if (i-I || j-J){
46: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetOwnershipRange() in MatSBAIJ format\n");
47: }
49: /* Assemble matrix */
50: if (bs == 1){
51: PetscOptionsGetInt(PETSC_NULL,"-test_problem",&prob,PETSC_NULL);
52: if (prob == 1){ /* tridiagonal matrix */
53: value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
54: for (i=1; i<n-1; i++) {
55: col[0] = i-1; col[1] = i; col[2] = i+1;
56: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
57: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
58: }
59: i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;
60: value[0]= 0.1; value[1]=-1; value[2]=2;
61: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
62: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
64: i = 0;
65: col[0] = n-1; col[1] = 1; col[2]=0;
66: value[0] = 0.1; value[1] = -1.0; value[2]=2;
67: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
68: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
69: }
70: else if (prob ==2){ /* matrix for the five point stencil */
71: n1 = (int) (sqrt((PetscReal)n) + 0.001);
72: if (n1*n1 - n) SETERRQ(PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive interger!");
73: for (i=0; i<n1; i++) {
74: for (j=0; j<n1; j++) {
75: I = j + n1*i;
76: if (i>0) {
77: J = I - n1;
78: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
79: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
80: }
81: if (i<n1-1) {
82: J = I + n1;
83: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
84: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
85: }
86: if (j>0) {
87: J = I - 1;
88: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
89: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
90: }
91: if (j<n1-1) {
92: J = I + 1;
93: MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
94: MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
95: }
96: MatSetValues(A,1,&I,1,&I,&four,INSERT_VALUES);
97: MatSetValues(sA,1,&I,1,&I,&four,INSERT_VALUES);
98: }
99: }
100: }
101: }
102: else { /* bs > 1 */
103: for (block=0; block<n/bs; block++){
104: /* diagonal blocks */
105: value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
106: for (i=1+block*bs; i<bs-1+block*bs; i++) {
107: col[0] = i-1; col[1] = i; col[2] = i+1;
108: MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
109: MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
110: }
111: i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
112: value[0]=-1.0; value[1]=4.0;
113: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
114: MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
116: i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
117: value[0]=4.0; value[1] = -1.0;
118: MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
119: MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
120: }
121: /* off-diagonal blocks */
122: value[0]=-1.0;
123: for (i=0; i<(n/bs-1)*bs; i++){
124: col[0]=i+bs;
125: MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);
126: MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES);
127: col[0]=i; row=i+bs;
128: MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
129: MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);
130: }
131: }
132: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
133: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
134: /* PetscPrintf(PETSC_COMM_SELF,"\n The Matrix: \n");
135: MatView(A, PETSC_VIEWER_DRAW_WORLD);
136: MatView(A, PETSC_VIEWER_STDOUT_WORLD); */
138: MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY);
139: MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY);
140: /* PetscPrintf(PETSC_COMM_SELF,"\n Symmetric Part of Matrix: \n");
141: MatView(sA, PETSC_VIEWER_DRAW_WORLD);
142: MatView(sA, PETSC_VIEWER_STDOUT_WORLD);
143: */
145: /* Test MatNorm(), MatDuplicate() */
146: MatNorm(A,NORM_FROBENIUS,&norm1);
147: MatDuplicate(sA,MAT_COPY_VALUES,&sB);
148: MatNorm(sB,NORM_FROBENIUS,&norm2);
149: norm1 -= norm2;
150: if (norm1<-tol || norm1>tol){
151: PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm(), fnorm1-fnorm2=%16.14e\n",norm1);
152: }
153: MatNorm(A,NORM_INFINITY,&norm1);
154: MatNorm(sB,NORM_INFINITY,&norm2);
155: norm1 -= norm2;
156: if (norm1<-tol || norm1>tol){
157: PetscPrintf(PETSC_COMM_SELF,"Error: MatNorm(), inf_norm1-inf_norm2=%16.14e\n",norm1);
158: }
160: /* Test MatGetInfo(), MatGetSize(), MatGetBlockSize() */
161: MatGetInfo(A,MAT_LOCAL,&minfo1);
162: MatGetInfo(sB,MAT_LOCAL,&minfo2);
163: /*
164: printf("matrix nonzeros (BAIJ format) = %d, allocated nonzeros= %d\n", (int)minfo1.nz_used,(int)minfo1.nz_allocated);
165: printf("matrix nonzeros(SBAIJ format) = %d, allocated nonzeros= %d\n", (int)minfo2.nz_used,(int)minfo2.nz_allocated);
166: */
167: i = (int) (minfo1.nz_used - minfo2.nz_used);
168: j = (int) (minfo2.nz_allocated - minfo2.nz_used);
169: if (i<0 || j<0) {
170: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetInfo()\n");
171: }
173: MatGetSize(A,&I,&J);
174: MatGetSize(sB,&i,&j);
175: if (i-I || j-J) {
176: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetSize()\n");
177: }
178:
179: MatGetBlockSize(A, &I);
180: MatGetBlockSize(sB, &i);
181: if (i-I){
182: PetscPrintf(PETSC_COMM_SELF,"Error: MatGetBlockSize()\n");
183: }
185: /* Test MatDiagonalScale(), MatGetDiagonal(), MatScale() */
186: PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&rdm);
187: VecCreateSeq(PETSC_COMM_SELF,n,&x);
188: VecDuplicate(x,&s1);
189: VecDuplicate(x,&s2);
190: VecDuplicate(x,&y);
191: VecDuplicate(x,&b);
192: VecSetRandom(rdm,x);
194: MatDiagonalScale(A,x,x);
195: MatDiagonalScale(sB,x,x);
197: MatGetDiagonal(A,s1);
198: MatGetDiagonal(sB,s2);
199:
200: VecAXPY(&neg_one,s1,s2);
201: VecNorm(s2,NORM_1,&norm1);
202: if ( norm1>tol) {
203: PetscPrintf(PETSC_COMM_SELF,"Error:MatGetDiagonal(), ||s1-s2||=%g\n",norm1);
204: }
205: /*
206: VecNorm(s1,NORM_1,&norm1);
207: VecNorm(s2,NORM_1,&norm2);
208: norm1 -= norm2;
209: if (norm1<-tol || norm1>tol) {
210: PetscPrintf(PETSC_COMM_SELF,"Error:MatGetDiagonal() \n");
211: }
212: */
214: MatScale(&alpha,A);
215: MatScale(&alpha,sB);
217: /* Test MatGetRowMax() */
218: MatGetRowMax(A,s1);
219: MatGetRowMax(sB,s2);
220: VecNorm(s1,NORM_1,&norm1);
221: VecNorm(s2,NORM_1,&norm2);
222: norm1 -= norm2;
223: if (norm1<-tol || norm1>tol) {
224: PetscPrintf(PETSC_COMM_SELF,"Error:MatGetRowMax() \n");
225: }
227: /* Test MatMult() */
228: for (i=0; i<40; i++) {
229: VecSetRandom(rdm,x);
230: MatMult(A,x,s1);
231: MatMult(sB,x,s2);
232: VecNorm(s1,NORM_1,&norm1);
233: VecNorm(s2,NORM_1,&norm2);
234: norm1 -= norm2;
235: if (norm1<-tol || norm1>tol) {
236: PetscPrintf(PETSC_COMM_SELF,"Error: MatMult(), norm1-norm2: %g\n",norm1);
237: }
238: }
240: /* MatMultAdd() */
241: for (i=0; i<40; i++) {
242: VecSetRandom(rdm,x);
243: VecSetRandom(rdm,y);
244: MatMultAdd(A,x,y,s1);
245: MatMultAdd(sB,x,y,s2);
246: VecNorm(s1,NORM_1,&norm1);
247: VecNorm(s2,NORM_1,&norm2);
248: norm1 -= norm2;
249: if (norm1<-tol || norm1>tol) {
250: PetscPrintf(PETSC_COMM_SELF,"Error:MatMultAdd(), norm1-norm2: %g\n",norm1);
251: }
252: }
254: /* Test MatCholeskyFactor(), MatICCFactor() with natural ordering */
255: MatGetOrdering(A,MATORDERING_NATURAL,&perm,&iscol);
256: ISDestroy(iscol);
257: norm1 = tol;
258: inc = bs;
260: /* initialize factinfo */
261: PetscMemzero(&factinfo,sizeof(MatFactorInfo));
263: for (lf=-1; lf<10; lf += inc){
264: if (lf==-1) { /* Cholesky factor of sB (duplicate sA) */
265: factinfo.fill = 5.0;
266: MatCholeskyFactorSymbolic(sB,perm,&factinfo,&sC);
267: } else if (!doIcc){
268: break;
269: } else { /* incomplete Cholesky factor */
270: factinfo.fill = 5.0;
271: factinfo.levels = lf;
272: MatICCFactorSymbolic(sB,perm,&factinfo,&sC);
273: }
274: MatCholeskyFactorNumeric(sB,&sC);
275: /* MatView(sC, PETSC_VIEWER_DRAW_WORLD); */
277: /* test MatGetDiagonal on numeric factor */
278: /*
279: if (lf == -1) {
280: MatGetDiagonal(sC,s1);
281: printf(" in ex74.c, diag: \n");
282: VecView(s1,PETSC_VIEWER_STDOUT_SELF);
283: }
284: */
286: MatMult(sB,x,b);
287: MatSolve(sC,b,y);
288: MatDestroy(sC);
289:
290: /* Check the error */
291: VecAXPY(&neg_one,x,y);
292: VecNorm(y,NORM_2,&norm2);
293: /* printf("lf: %d, error: %g\n", lf,norm2); */
294: if (10*norm1 < norm2 && lf-inc != -1){
295: PetscPrintf(PETSC_COMM_SELF,"lf=%D, %D, Norm of error=%g, %g\n",lf-inc,lf,norm1,norm2);
296: }
297: norm1 = norm2;
298: if (norm2 < tol && lf != -1) break;
299: }
301: ISDestroy(perm);
303: MatDestroy(A);
304: MatDestroy(sB);
305: MatDestroy(sA);
306: VecDestroy(x);
307: VecDestroy(y);
308: VecDestroy(s1);
309: VecDestroy(s2);
310: VecDestroy(b);
311: PetscRandomDestroy(rdm);
313: PetscFinalize();
314: return 0;
315: }