Actual source code: ex6.c
1: /*$Id: ex6.c,v 1.19 2001/08/07 21:30:37 bsmith Exp $*/
3: static char help[] = "Creates a matrix using 9 pt stensil, and uses it to test MatIncreaseOverlap (needed for aditive schwarts preconditioner. n
4: -m <size> : problem sizen
5: -x1, -x2 <size> : no of subdomains in x and y directionsnn";
7: #include petscsles.h
9: int FormElementStiffness(PetscReal H,PetscScalar *Ke)
10: {
11: Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H;
12: Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0;
13: Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H;
14: Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0;
15: return 0;
16: }
17: int FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r)
18: {
19: r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0;
20: return 0;
21: }
23: int main(int argc,char **args)
24: {
25: Mat C;
26: int i,m = 2,N,M,ierr,idx[4],Nsub1,Nsub2,ol=1,x1,x2;
27: PetscScalar Ke[16];
28: PetscReal x,y,h;
29: IS *is1,*is2;
30: PetscTruth flg;
32: PetscInitialize(&argc,&args,(char *)0,help);
33: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
34: N = (m+1)*(m+1); /* dimension of matrix */
35: M = m*m; /* number of elements */
36: h = 1.0/m; /* mesh width */
37: x1= (m+1)/2;
38: x2= x1;
39: PetscOptionsGetInt(PETSC_NULL,"-x1",&x1,PETSC_NULL);
40: PetscOptionsGetInt(PETSC_NULL,"-x2",&x2,PETSC_NULL);
41: /* create stiffness matrix */
42: MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,9,PETSC_NULL,&C);
44: /* forms the element stiffness for the Laplacian */
45: FormElementStiffness(h*h,Ke);
46: for (i=0; i<M; i++) {
47: /* location of lower left corner of element */
48: x = h*(i % m); y = h*(i/m);
49: /* node numbers for the four corners of element */
50: idx[0] = (m+1)*(i/m) + (i % m);
51: idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1;
52: MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);
53: }
54: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
55: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
58: for (ol=0; ol<m+2; ++ol) {
60: PCASMCreateSubdomains2D(m+1,m+1,x1,x2,1,0,&Nsub1,&is1);
61: MatIncreaseOverlap(C,Nsub1,is1,ol);
62: PCASMCreateSubdomains2D(m+1,m+1,x1,x2,1,ol,&Nsub2,&is2);
63:
64: PetscPrintf(PETSC_COMM_SELF,"flg == 1 => both index sets are samen");
65: if(Nsub1 != Nsub2){
66: PetscPrintf(PETSC_COMM_SELF,"Error: No of indes sets don't matchn");
67: }
68:
69: for (i=0; i<Nsub1; ++i) {
70: ISEqual(is1[i],is2[i],&flg);
71: PetscPrintf(PETSC_COMM_SELF,"i = %d,flg = %d n",i,flg);
72:
73: }
74: for (i=0; i<Nsub1; ++i) {ISDestroy(is1[i]);}
75: for (i=0; i<Nsub2; ++i) {ISDestroy(is2[i]);}
76:
78: PetscFree(is1);
79: PetscFree(is2);
80: }
81: MatDestroy(C);
82: PetscFinalize();
83: return 0;
84: }