1 2 static char help[] = "Creates a matrix using 9 pt stencil, and uses it to test MatIncreaseOverlap (needed for additive Schwarz preconditioner). \n\ 3 -m <size> : problem size\n\ 4 -x1, -x2 <size> : no of subdomains in x and y directions\n\n"; 5 6 #include <petscksp.h> 7 8 PetscErrorCode FormElementStiffness(PetscReal H,PetscScalar *Ke) 9 { 10 Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H; 11 Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0; 12 Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H; 13 Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0; 14 return 0; 15 } 16 PetscErrorCode FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r) 17 { 18 r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0; 19 return 0; 20 } 21 22 int main(int argc,char **args) 23 { 24 Mat C; 25 PetscErrorCode ierr; 26 PetscInt i,m = 2,N,M,idx[4],Nsub1,Nsub2,ol=1,x1,x2; 27 PetscScalar Ke[16]; 28 PetscReal h; 29 IS *is1,*is2,*islocal1,*islocal2; 30 PetscBool flg; 31 32 ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 33 ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 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 ierr = PetscOptionsGetInt(NULL,NULL,"-x1",&x1,NULL);CHKERRQ(ierr); 40 ierr = PetscOptionsGetInt(NULL,NULL,"-x2",&x2,NULL);CHKERRQ(ierr); 41 /* create stiffness matrix */ 42 ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,9,NULL,&C);CHKERRQ(ierr); 43 44 /* forms the element stiffness for the Laplacian */ 45 ierr = FormElementStiffness(h*h,Ke);CHKERRQ(ierr); 46 for (i=0; i<M; i++) { 47 /* node numbers for the four corners of element */ 48 idx[0] = (m+1)*(i/m) + (i % m); 49 idx[1] = idx[0]+1; idx[2] = idx[1] + m + 1; idx[3] = idx[2] - 1; 50 ierr = MatSetValues(C,4,idx,4,idx,Ke,ADD_VALUES);CHKERRQ(ierr); 51 } 52 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 53 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 54 55 for (ol=0; ol<m+2; ++ol) { 56 57 ierr = PCASMCreateSubdomains2D(m+1,m+1,x1,x2,1,0,&Nsub1,&is1,&islocal1);CHKERRQ(ierr); 58 ierr = MatIncreaseOverlap(C,Nsub1,is1,ol);CHKERRQ(ierr); 59 ierr = PCASMCreateSubdomains2D(m+1,m+1,x1,x2,1,ol,&Nsub2,&is2,&islocal2);CHKERRQ(ierr); 60 61 ierr = PetscPrintf(PETSC_COMM_SELF,"flg == 1 => both index sets are same\n");CHKERRQ(ierr); 62 if (Nsub1 != Nsub2) { 63 ierr = PetscPrintf(PETSC_COMM_SELF,"Error: No of indes sets don't match\n");CHKERRQ(ierr); 64 } 65 66 for (i=0; i<Nsub1; ++i) { 67 ierr = ISEqual(is1[i],is2[i],&flg);CHKERRQ(ierr); 68 ierr = PetscPrintf(PETSC_COMM_SELF,"i = %D,flg = %d \n",i,(int)flg);CHKERRQ(ierr); 69 70 } 71 for (i=0; i<Nsub1; ++i) {ierr = ISDestroy(&is1[i]);CHKERRQ(ierr);} 72 for (i=0; i<Nsub2; ++i) {ierr = ISDestroy(&is2[i]);CHKERRQ(ierr);} 73 for (i=0; i<Nsub1; ++i) {ierr = ISDestroy(&islocal1[i]);CHKERRQ(ierr);} 74 for (i=0; i<Nsub2; ++i) {ierr = ISDestroy(&islocal2[i]);CHKERRQ(ierr);} 75 76 ierr = PetscFree(is1);CHKERRQ(ierr); 77 ierr = PetscFree(is2);CHKERRQ(ierr); 78 ierr = PetscFree(islocal1);CHKERRQ(ierr); 79 ierr = PetscFree(islocal2);CHKERRQ(ierr); 80 } 81 ierr = MatDestroy(&C);CHKERRQ(ierr); 82 ierr = PetscFinalize(); 83 return ierr; 84 } 85 86 /*TEST 87 88 test: 89 args: -m 7 90 91 TEST*/ 92