Actual source code: ex51.c
1: /*$Id: ex51.c,v 1.22 2001/08/07 21:30:08 bsmith Exp $*/
3: static char help[] = "Tests MatIncreaseOverlap(), MatGetSubMatrices() for MatBAIJ format.n";
5: #include petscmat.h
7: int main(int argc,char **args)
8: {
9: Mat A,B,*submatA,*submatB;
10: int bs=1,m=43,ov=1,i,j,k,*rows,*cols,ierr,M,nd=5,*idx,size,mm,nn;
11: PetscScalar *vals,rval;
12: IS *is1,*is2;
13: PetscRandom rand;
14: Vec xx,s1,s2;
15: PetscReal s1norm,s2norm,rnorm,tol = 1.e-10;
16: PetscTruth flg;
18: PetscInitialize(&argc,&args,(char *)0,help);
19:
21: PetscOptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);
22: PetscOptionsGetInt(PETSC_NULL,"-mat_size",&m,PETSC_NULL);
23: PetscOptionsGetInt(PETSC_NULL,"-ov",&ov,PETSC_NULL);
24: PetscOptionsGetInt(PETSC_NULL,"-nd",&nd,PETSC_NULL);
25: M = m*bs;
27: MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,1,PETSC_NULL,&A);
28: MatCreateSeqAIJ(PETSC_COMM_SELF,M,M,15,PETSC_NULL,&B);
29: PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&rand);
31: PetscMalloc(bs*sizeof(int),&rows);
32: PetscMalloc(bs*sizeof(int),&cols);
33: PetscMalloc(bs*bs*sizeof(PetscScalar),&vals);
34: PetscMalloc(M*sizeof(PetscScalar),&idx);
35:
36: /* Now set blocks of values */
37: for (i=0; i<20*bs; i++) {
38: PetscRandomGetValue(rand,&rval);
39: cols[0] = bs*(int)(PetscRealPart(rval)*m);
40: PetscRandomGetValue(rand,&rval);
41: rows[0] = bs*(int)(PetscRealPart(rval)*m);
42: for (j=1; j<bs; j++) {
43: rows[j] = rows[j-1]+1;
44: cols[j] = cols[j-1]+1;
45: }
47: for (j=0; j<bs*bs; j++) {
48: PetscRandomGetValue(rand,&rval);
49: vals[j] = rval;
50: }
51: MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES);
52: MatSetValues(B,bs,rows,bs,cols,vals,ADD_VALUES);
53: }
55: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
56: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
57: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
58: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
60: /* Test MatIncreaseOverlap() */
61: PetscMalloc(nd*sizeof(IS **),&is1);
62: PetscMalloc(nd*sizeof(IS **),&is2);
64:
65: for (i=0; i<nd; i++) {
66: PetscRandomGetValue(rand,&rval);
67: size = (int)(PetscRealPart(rval)*m);
68: for (j=0; j<size; j++) {
69: PetscRandomGetValue(rand,&rval);
70: idx[j*bs] = bs*(int)(PetscRealPart(rval)*m);
71: for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k;
72: }
73: ISCreateGeneral(PETSC_COMM_SELF,size*bs,idx,is1+i);
74: ISCreateGeneral(PETSC_COMM_SELF,size*bs,idx,is2+i);
75: }
76: MatIncreaseOverlap(A,nd,is1,ov);
77: MatIncreaseOverlap(B,nd,is2,ov);
79: for (i=0; i<nd; ++i) {
80: ISEqual(is1[i],is2[i],&flg);
81: PetscPrintf(PETSC_COMM_SELF,"i=%d, flg =%dn",i,flg);
82: }
84: for (i=0; i<nd; ++i) {
85: ISSort(is1[i]);
86: ISSort(is2[i]);
87: }
88:
89: MatGetSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA);
90: MatGetSubMatrices(B,nd,is2,is2,MAT_INITIAL_MATRIX,&submatB);
92: /* Test MatMult() */
93: for (i=0; i<nd; i++) {
94: MatGetSize(submatA[i],&mm,&nn);
95: VecCreateSeq(PETSC_COMM_SELF,mm,&xx);
96: VecDuplicate(xx,&s1);
97: VecDuplicate(xx,&s2);
98: for (j=0; j<3; j++) {
99: VecSetRandom(rand,xx);
100: MatMult(submatA[i],xx,s1);
101: MatMult(submatB[i],xx,s2);
102: VecNorm(s1,NORM_2,&s1norm);
103: VecNorm(s2,NORM_2,&s2norm);
104: rnorm = s2norm-s1norm;
105: if (rnorm<-tol || rnorm>tol) {
106: PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14en",s1norm,s2norm);
107: }
108: }
109: VecDestroy(xx);
110: VecDestroy(s1);
111: VecDestroy(s2);
112: }
113: /* Now test MatGetSubmatrices with MAT_REUSE_MATRIX option */
114: MatGetSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA);
115: MatGetSubMatrices(B,nd,is2,is2,MAT_REUSE_MATRIX,&submatB);
116:
117: /* Test MatMult() */
118: for (i=0; i<nd; i++) {
119: MatGetSize(submatA[i],&mm,&nn);
120: VecCreateSeq(PETSC_COMM_SELF,mm,&xx);
121: VecDuplicate(xx,&s1);
122: VecDuplicate(xx,&s2);
123: for (j=0; j<3; j++) {
124: VecSetRandom(rand,xx);
125: MatMult(submatA[i],xx,s1);
126: MatMult(submatB[i],xx,s2);
127: VecNorm(s1,NORM_2,&s1norm);
128: VecNorm(s2,NORM_2,&s2norm);
129: rnorm = s2norm-s1norm;
130: if (rnorm<-tol || rnorm>tol) {
131: PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14en",s1norm,s2norm);
132: }
133: }
134: VecDestroy(xx);
135: VecDestroy(s1);
136: VecDestroy(s2);
137: }
138:
139: /* Free allocated memory */
140: for (i=0; i<nd; ++i) {
141: ISDestroy(is1[i]);
142: ISDestroy(is2[i]);
143: MatDestroy(submatA[i]);
144: MatDestroy(submatB[i]);
145: }
146: PetscFree(is1);
147: PetscFree(is2);
148: PetscFree(idx);
149: PetscFree(rows);
150: PetscFree(cols);
151: PetscFree(vals);
152: MatDestroy(A);
153: MatDestroy(B);
154: PetscFree(submatA);
155: PetscFree(submatB);
156: PetscRandomDestroy(rand);
157: PetscFinalize();
158: return 0;
159: }