static char help[] = "Creates a matrix using 9 pt stencil, and uses it to test MatIncreaseOverlap (needed for additive Schwarz preconditioner). \n\ -m : problem size\n\ -x1, -x2 : no of subdomains in x and y directions\n\n"; #include PetscErrorCode FormElementStiffness(PetscReal H,PetscScalar *Ke) { Ke[0] = H/6.0; Ke[1] = -.125*H; Ke[2] = H/12.0; Ke[3] = -.125*H; Ke[4] = -.125*H; Ke[5] = H/6.0; Ke[6] = -.125*H; Ke[7] = H/12.0; Ke[8] = H/12.0; Ke[9] = -.125*H; Ke[10] = H/6.0; Ke[11] = -.125*H; Ke[12] = -.125*H; Ke[13] = H/12.0; Ke[14] = -.125*H; Ke[15] = H/6.0; return 0; } PetscErrorCode FormElementRhs(PetscReal x,PetscReal y,PetscReal H,PetscScalar *r) { r[0] = 0.; r[1] = 0.; r[2] = 0.; r[3] = 0.0; return 0; } int main(int argc,char **args) { Mat C; PetscErrorCode ierr; PetscInt i,m = 2,N,M,idx[4],Nsub1,Nsub2,ol=1,x1,x2; PetscScalar Ke[16]; PetscReal h; IS *is1,*is2,*islocal1,*islocal2; PetscBool flg; ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); N = (m+1)*(m+1); /* dimension of matrix */ M = m*m; /* number of elements */ h = 1.0/m; /* mesh width */ x1 = (m+1)/2; x2 = x1; ierr = PetscOptionsGetInt(NULL,NULL,"-x1",&x1,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-x2",&x2,NULL);CHKERRQ(ierr); /* create stiffness matrix */ ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,9,NULL,&C);CHKERRQ(ierr); /* forms the element stiffness for the Laplacian */ ierr = FormElementStiffness(h*h,Ke);CHKERRQ(ierr); for (i=0; i both index sets are same\n");CHKERRQ(ierr); if (Nsub1 != Nsub2) { ierr = PetscPrintf(PETSC_COMM_SELF,"Error: No of indes sets don't match\n");CHKERRQ(ierr); } for (i=0; i