Actual source code: ex76.c

  2: static char help[] = "Tests matrix permutation for factorization and solve on matrix with MatSBAIJ format. Modified from ex74.c\n";

 4:  #include petscmat.h

  8: int main(int argc,char **args)
  9: {
 10:   Vec            x,y,b;
 11:   Mat            A;           /* linear system matrix */
 12:   Mat            sA,sC;       /* symmetric part of the matrices */
 13:   PetscInt       n,mbs=16,bs=1,nz=3,prob=1,i,j,col[3],block, row,I,J,n1,*ip_ptr,lf;
 15:   PetscMPIInt    size;
 16:   PetscReal      norm1,norm2,tol=1.e-10;
 17:   PetscScalar    neg_one = -1.0,four=4.0,value[3];
 18:   IS             perm;
 19:   PetscRandom    rdm;
 20:   PetscTruth     reorder=PETSC_TRUE;
 21:   MatFactorInfo  factinfo;

 23:   PetscInitialize(&argc,&args,(char *)0,help);
 24:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 25:   if (size != 1) SETERRQ(1,"This is a uniprocessor example only!");
 26:   PetscOptionsGetInt(PETSC_NULL,"-bs",&bs,PETSC_NULL);
 27:   PetscOptionsGetInt(PETSC_NULL,"-mbs",&mbs,PETSC_NULL);

 29:   n = mbs*bs;
 30:   ierr=MatCreateSeqBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,PETSC_NULL, &A);
 31:   ierr=MatCreateSeqSBAIJ(PETSC_COMM_WORLD,bs,n,n,nz,PETSC_NULL, &sA);

 33:   /* Test MatGetOwnershipRange() */
 34:   MatGetOwnershipRange(A,&I,&J);
 35:   MatGetOwnershipRange(sA,&i,&j);
 36:   if (i-I || j-J){
 37:     PetscPrintf(PETSC_COMM_SELF,"Error: MatGetOwnershipRange() in MatSBAIJ format\n");
 38:   }

 40:   /* Assemble matrix */
 41:   if (bs == 1){
 42:     PetscOptionsGetInt(PETSC_NULL,"-test_problem",&prob,PETSC_NULL);
 43:     if (prob == 1){ /* tridiagonal matrix */
 44:       value[0] = -1.0; value[1] = 2.0; value[2] = -1.0;
 45:       for (i=1; i<n-1; i++) {
 46:         col[0] = i-1; col[1] = i; col[2] = i+1;
 47:         MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 48:         MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
 49:       }
 50:       i = n - 1; col[0]=0; col[1] = n - 2; col[2] = n - 1;
 51:       value[0]= 0.1; value[1]=-1; value[2]=2;
 52:       MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 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(A,1,&i,3,col,value,INSERT_VALUES);
 58:       MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
 59:     }
 60:     else if (prob ==2){ /* matrix for the five point stencil */
 61:       n1 = (int) (sqrt((PetscReal)n) + 0.001);
 62:       if (n1*n1 - n) SETERRQ(PETSC_ERR_ARG_WRONG,"sqrt(n) must be a positive interger!");
 63:       for (i=0; i<n1; i++) {
 64:         for (j=0; j<n1; j++) {
 65:           I = j + n1*i;
 66:           if (i>0)   {
 67:             J = I - n1;
 68:             MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
 69:             MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
 70:           }
 71:           if (i<n1-1) {
 72:             J = I + n1;
 73:             MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
 74:             MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
 75:           }
 76:           if (j>0)   {
 77:             J = I - 1;
 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 (j<n1-1) {
 82:             J = I + 1;
 83:             MatSetValues(A,1,&I,1,&J,&neg_one,INSERT_VALUES);
 84:             MatSetValues(sA,1,&I,1,&J,&neg_one,INSERT_VALUES);
 85:           }
 86:           MatSetValues(A,1,&I,1,&I,&four,INSERT_VALUES);
 87:           MatSetValues(sA,1,&I,1,&I,&four,INSERT_VALUES);
 88:         }
 89:       }
 90:     }
 91:   }
 92:   else { /* bs > 1 */
 93:     for (block=0; block<n/bs; block++){
 94:       /* diagonal blocks */
 95:       value[0] = -1.0; value[1] = 4.0; value[2] = -1.0;
 96:       for (i=1+block*bs; i<bs-1+block*bs; i++) {
 97:         col[0] = i-1; col[1] = i; col[2] = i+1;
 98:         MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 99:         MatSetValues(sA,1,&i,3,col,value,INSERT_VALUES);
100:       }
101:       i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs;
102:       value[0]=-1.0; value[1]=4.0;
103:       MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
104:       MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);

106:       i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs;
107:       value[0]=4.0; value[1] = -1.0;
108:       MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
109:       MatSetValues(sA,1,&i,2,col,value,INSERT_VALUES);
110:     }
111:     /* off-diagonal blocks */
112:     value[0]=-1.0;
113:     for (i=0; i<(n/bs-1)*bs; i++){
114:       col[0]=i+bs;
115:       MatSetValues(A,1,&i,1,col,value,INSERT_VALUES);
116:       MatSetValues(sA,1,&i,1,col,value,INSERT_VALUES);
117:       col[0]=i; row=i+bs;
118:       MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
119:       MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);
120:     }
121:     if (bs == 2){
122:       /* insert a value to off-diag blocks */
123:       row = 2; col[0] = 5; value[0] = 0.01;
124:       MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
125:       MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);
126:       row = 0; col[0] = 3; value[0] = 0.01;
127:       MatSetValues(A,1,&row,1,col,value,INSERT_VALUES);
128:       MatSetValues(sA,1,&row,1,col,value,INSERT_VALUES);
129:     }
130:   }
131:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
132:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
133:   /* PetscPrintf(PETSC_COMM_SELF,"\n The Matrix: \n");
134:   MatView(A, PETSC_VIEWER_DRAW_WORLD);
135:   MatView(A, PETSC_VIEWER_STDOUT_WORLD); */
136:   MatAssemblyBegin(sA,MAT_FINAL_ASSEMBLY);
137:   MatAssemblyEnd(sA,MAT_FINAL_ASSEMBLY);
138:   /* PetscPrintf(PETSC_COMM_SELF,"\n Symmetric Part of Matrix: \n"); */
139:   /* MatView(sA, PETSC_VIEWER_DRAW_WORLD); */
140:   /* MatView(sA, PETSC_VIEWER_STDOUT_WORLD); */

142:   /* Vectors */
143:   PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&rdm);
144:   VecCreateSeq(PETSC_COMM_SELF,n,&x);
145:   VecDuplicate(x,&b);
146:   VecDuplicate(x,&y);
147:   VecSetRandom(rdm,x);

149:   /* Test MatReordering() */
150:   PetscMalloc(mbs*sizeof(PetscInt),&ip_ptr);

152:   for (i=0; i<mbs; i++) ip_ptr[i] = i;
153:   if(reorder){
154:     i = ip_ptr[1]; ip_ptr[1] = ip_ptr[mbs-2]; ip_ptr[mbs-2] = i;
155:     /* i = ip_ptr[0]; ip_ptr[0] = ip_ptr[mbs-1]; ip_ptr[mbs-1] = i; */
156:     /* i = ip_ptr[2]; ip_ptr[2] = ip_ptr[mbs-3]; ip_ptr[mbs-3] = i; */
157:   }

159:   ISCreateGeneral(PETSC_COMM_SELF,mbs,ip_ptr,&perm);
160:   ISSetPermutation(perm);
161: 
162:   /* Test MatCholeskyFactor(), MatICCFactor() */
163:   norm1 = tol;
164:   for (lf=-1; lf<10*bs; lf += bs){
165:     if (lf==-1) {  /* Cholesky factor */
166:       factinfo.fill = 5.0;
167:       MatCholeskyFactorSymbolic(sA,perm,&factinfo,&sC);
168:     } else {       /* incomplete Cholesky factor */
169:       factinfo.fill   = 5.0;
170:       factinfo.levels = lf;
171:       MatICCFactorSymbolic(sA,perm,&factinfo,&sC);
172:     }
173:     MatCholeskyFactorNumeric(sA,&sC);
174:     /* MatView(sC, PETSC_VIEWER_DRAW_WORLD);  */ /* view factored matrix */
175:     /* MatView(sC, PETSC_VIEWER_STDOUT_WORLD); */
176: 
177:     MatMult(sA,x,b);
178:     MatSolve(sC,b,y);
179:     if (bs == 1) {
180:       Vecs xx,bb;
181:       VecsCreateSeq(PETSC_COMM_SELF,n,4,&xx);
182:       VecsDuplicate(xx,&bb);
183:       MatSolves(sC,bb,xx);
184:       VecsDestroy(xx);
185:       VecsDestroy(bb);
186:     }
187:     MatDestroy(sC);

189:     /* Check the error */
190:     VecAXPY(&neg_one,x,y);
191:     VecNorm(y,NORM_2,&norm2);
192:     /* printf("lf: %d, error: %g\n", lf,norm2); */
193:     if (10*norm1 < norm2 && lf-bs != -1){
194:       PetscPrintf(PETSC_COMM_SELF,"lf=%D, %D, Norm of error=%g, %g\n",lf-bs,lf,norm1,norm2);
195:     }
196:     norm1 = norm2;
197:     if (norm2 < tol && lf != -1) break;
198:   }

200:   ISDestroy(perm);
201:   PetscFree(ip_ptr);
202:   MatDestroy(A);
203:   MatDestroy(sA);
204:   VecDestroy(x);
205:   VecDestroy(y);
206:   VecDestroy(b);
207:   PetscRandomDestroy(rdm);

209:   PetscFinalize();
210:   return 0;
211: }