static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for parallel MatSBAIJ format.\n"; /* Example of usage: mpiexec -n 2 ./ex92 -nd 2 -ov 3 -mat_block_size 2 -view_id 0 -test_overlap -test_submat */ #include int main(int argc,char **args) { Mat A,Atrans,sA,*submatA,*submatsA; PetscErrorCode ierr; PetscMPIInt size,rank; PetscInt bs=1,mbs=10,ov=1,i,j,k,*rows,*cols,nd=2,*idx,rstart,rend,sz,M,N,Mbs; PetscScalar *vals,rval,one=1.0; IS *is1,*is2; PetscRandom rand; PetscBool flg,TestOverlap,TestSubMat,TestAllcols,test_sorted=PETSC_FALSE; PetscInt vid = -1; #if defined(PETSC_USE_LOG) PetscLogStage stages[2]; #endif ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL)); CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_mbs",&mbs,NULL)); CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL)); CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL)); CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-view_id",&vid,NULL)); CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_overlap", &TestOverlap)); CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_submat", &TestSubMat)); CHKERRQ(PetscOptionsHasName(NULL,NULL, "-test_allcols", &TestAllcols)); CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_sorted",&test_sorted,NULL)); CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); CHKERRQ(MatSetSizes(A,mbs*bs,mbs*bs,PETSC_DECIDE,PETSC_DECIDE)); CHKERRQ(MatSetType(A,MATBAIJ)); CHKERRQ(MatSeqBAIJSetPreallocation(A,bs,PETSC_DEFAULT,NULL)); CHKERRQ(MatMPIBAIJSetPreallocation(A,bs,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL)); CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); CHKERRQ(PetscRandomSetFromOptions(rand)); CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend)); CHKERRQ(MatGetSize(A,&M,&N)); Mbs = M/bs; CHKERRQ(PetscMalloc1(bs,&rows)); CHKERRQ(PetscMalloc1(bs,&cols)); CHKERRQ(PetscMalloc1(bs*bs,&vals)); CHKERRQ(PetscMalloc1(M,&idx)); /* Now set blocks of values */ for (j=0; j= 0 && vid < size) { CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"A:\n")); CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"sA:\n")); CHKERRQ(MatView(sA,PETSC_VIEWER_STDOUT_WORLD)); } /* Test sA==A through MatMult() */ CHKERRQ(MatMultEqual(A,sA,10,&flg)); PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Error in MatConvert(): A != sA"); /* Test MatIncreaseOverlap() */ CHKERRQ(PetscMalloc1(nd,&is1)); CHKERRQ(PetscMalloc1(nd,&is2)); for (i=0; i