Actual source code: ex91.c
2: static char help[] = "Tests MatIncreaseOverlap(), MatGetSubMatrices() for sequential MatSBAIJ format. Derived from ex51.c\n";
4: #include petscmat.h
8: int main(int argc,char **args)
9: {
10: Mat A,Atrans,sA,*submatA,*submatsA;
11: PetscInt bs=1,m=43,ov=1,i,j,k,*rows,*cols,M,nd=5,*idx,mm,nn;
13: PetscMPIInt size;
14: PetscScalar *vals,rval,one=1.0;
15: IS *is1,*is2;
16: PetscRandom rand;
17: Vec xx,s1,s2;
18: PetscReal s1norm,s2norm,rnorm,tol = 1.e-10;
19: PetscTruth flg;
21: PetscInitialize(&argc,&args,(char *)0,help);
22:
24: PetscOptionsGetInt(PETSC_NULL,"-mat_block_size",&bs,PETSC_NULL);
25: PetscOptionsGetInt(PETSC_NULL,"-mat_size",&m,PETSC_NULL);
26: PetscOptionsGetInt(PETSC_NULL,"-ov",&ov,PETSC_NULL);
27: PetscOptionsGetInt(PETSC_NULL,"-nd",&nd,PETSC_NULL);
29: /* create a SeqBAIJ matrix A */
30: M = m*bs;
31: MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,1,PETSC_NULL,&A);
32: PetscRandomCreate(PETSC_COMM_SELF,RANDOM_DEFAULT,&rand);
34: PetscMalloc(bs*sizeof(PetscInt),&rows);
35: PetscMalloc(bs*sizeof(PetscInt),&cols);
36: PetscMalloc(bs*bs*sizeof(PetscScalar),&vals);
37: PetscMalloc(M*sizeof(PetscScalar),&idx);
38:
39: /* Now set blocks of random values */
40: /* first, set diagonal blocks as zero */
41: for (j=0; j<bs*bs; j++) vals[j] = 0.0;
42: for (i=0; i<m; i++){
43: cols[0] = i*bs; rows[0] = i*bs;
44: for (j=1; j<bs; j++) {
45: rows[j] = rows[j-1]+1;
46: cols[j] = cols[j-1]+1;
47: }
48: MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES);
49: }
50: /* second, add random blocks */
51: for (i=0; i<20*bs; i++) {
52: PetscRandomGetValue(rand,&rval);
53: cols[0] = bs*(int)(PetscRealPart(rval)*m);
54: PetscRandomGetValue(rand,&rval);
55: rows[0] = bs*(int)(PetscRealPart(rval)*m);
56: for (j=1; j<bs; j++) {
57: rows[j] = rows[j-1]+1;
58: cols[j] = cols[j-1]+1;
59: }
61: for (j=0; j<bs*bs; j++) {
62: PetscRandomGetValue(rand,&rval);
63: vals[j] = rval;
64: }
65: MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES);
66: }
68: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
69: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
71: /* make A a symmetric matrix: A <- A^T + A */
72: MatTranspose(A, &Atrans);
73: MatAXPY(A,one,Atrans,DIFFERENT_NONZERO_PATTERN);
74: MatDestroy(Atrans);
75: MatTranspose(A, &Atrans);
76: MatEqual(A, Atrans, &flg);
77: if (!flg) {
78: SETERRQ(1,"A+A^T is non-symmetric");
79: }
80: MatDestroy(Atrans);
82: /* create a SeqSBAIJ matrix sA (= A) */
83: MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&sA);
84:
85: /* Test sA==A through MatMult() */
86: for (i=0; i<nd; i++) {
87: MatGetSize(A,&mm,&nn);
88: VecCreateSeq(PETSC_COMM_SELF,mm,&xx);
89: VecDuplicate(xx,&s1);
90: VecDuplicate(xx,&s2);
91: for (j=0; j<3; j++) {
92: VecSetRandom(xx,rand);
93: MatMult(A,xx,s1);
94: MatMult(sA,xx,s2);
95: VecNorm(s1,NORM_2,&s1norm);
96: VecNorm(s2,NORM_2,&s2norm);
97: rnorm = s2norm-s1norm;
98: if (rnorm<-tol || rnorm>tol) {
99: PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",s1norm,s2norm);
100: }
101: }
102: VecDestroy(xx);
103: VecDestroy(s1);
104: VecDestroy(s2);
105: }
107: /* Test MatIncreaseOverlap() */
108: PetscMalloc(nd*sizeof(IS **),&is1);
109: PetscMalloc(nd*sizeof(IS **),&is2);
111:
112: for (i=0; i<nd; i++) {
113: PetscRandomGetValue(rand,&rval);
114: size = (int)(PetscRealPart(rval)*m);
115: for (j=0; j<size; j++) {
116: PetscRandomGetValue(rand,&rval);
117: idx[j*bs] = bs*(int)(PetscRealPart(rval)*m);
118: for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k;
119: }
120: ISCreateGeneral(PETSC_COMM_SELF,size*bs,idx,is1+i);
121: ISCreateGeneral(PETSC_COMM_SELF,size*bs,idx,is2+i);
122: }
123: /* for debugging */
124: /*
125: MatView(A,PETSC_VIEWER_STDOUT_SELF);
126: MatView(sA,PETSC_VIEWER_STDOUT_SELF);
127: */
129: MatIncreaseOverlap(A,nd,is1,ov);
130: MatIncreaseOverlap(sA,nd,is2,ov);
132: for (i=0; i<nd; ++i) {
133: ISSort(is1[i]);
134: ISSort(is2[i]);
135: }
137: for (i=0; i<nd; ++i) {
138: ISEqual(is1[i],is2[i],&flg);
139: if (!flg){
140: /* ISView(is1[i],PETSC_VIEWER_STDOUT_SELF);
141: ISView(is2[i],PETSC_VIEWER_STDOUT_SELF); */
142: SETERRQ1(1,"i=%d, is1 != is2",i);
143: }
144: }
145:
146: MatGetSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA);
147: MatGetSubMatrices(sA,nd,is2,is2,MAT_INITIAL_MATRIX,&submatsA);
149: /* Test MatMult() */
150: for (i=0; i<nd; i++) {
151: MatGetSize(submatA[i],&mm,&nn);
152: VecCreateSeq(PETSC_COMM_SELF,mm,&xx);
153: VecDuplicate(xx,&s1);
154: VecDuplicate(xx,&s2);
155: for (j=0; j<3; j++) {
156: VecSetRandom(xx,rand);
157: MatMult(submatA[i],xx,s1);
158: MatMult(submatsA[i],xx,s2);
159: VecNorm(s1,NORM_2,&s1norm);
160: VecNorm(s2,NORM_2,&s2norm);
161: rnorm = s2norm-s1norm;
162: if (rnorm<-tol || rnorm>tol) {
163: PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",s1norm,s2norm);
164: }
165: }
166: VecDestroy(xx);
167: VecDestroy(s1);
168: VecDestroy(s2);
169: }
171: /* Now test MatGetSubmatrices with MAT_REUSE_MATRIX option */
172: MatGetSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA);
173: MatGetSubMatrices(sA,nd,is2,is2,MAT_REUSE_MATRIX,&submatsA);
174:
175: /* Test MatMult() */
176: for (i=0; i<nd; i++) {
177: MatGetSize(submatA[i],&mm,&nn);
178: VecCreateSeq(PETSC_COMM_SELF,mm,&xx);
179: VecDuplicate(xx,&s1);
180: VecDuplicate(xx,&s2);
181: for (j=0; j<3; j++) {
182: VecSetRandom(xx,rand);
183: MatMult(submatA[i],xx,s1);
184: MatMult(submatsA[i],xx,s2);
185: VecNorm(s1,NORM_2,&s1norm);
186: VecNorm(s2,NORM_2,&s2norm);
187: rnorm = s2norm-s1norm;
188: if (rnorm<-tol || rnorm>tol) {
189: PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",s1norm,s2norm);
190: }
191: }
192: VecDestroy(xx);
193: VecDestroy(s1);
194: VecDestroy(s2);
195: }
196:
197: /* Free allocated memory */
198: for (i=0; i<nd; ++i) {
199: ISDestroy(is1[i]);
200: ISDestroy(is2[i]);
201:
202: MatDestroy(submatA[i]);
203: MatDestroy(submatsA[i]);
204:
205: }
206:
207: PetscFree(submatA);
208: PetscFree(submatsA);
209:
210: PetscFree(is1);
211: PetscFree(is2);
212: PetscFree(idx);
213: PetscFree(rows);
214: PetscFree(cols);
215: PetscFree(vals);
216: MatDestroy(A);
217: MatDestroy(sA);
218: PetscRandomDestroy(rand);
219: PetscFinalize();
220: return 0;
221: }