static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for parallel AIJ and BAIJ formats.\n"; #include int main(int argc,char **args) { Mat A,B,*submatA,*submatB; PetscInt bs=1,m=11,ov=1,i,j,k,*rows,*cols,nd=5,*idx,rstart,rend,sz,mm,nn,M,N,Mbs; PetscErrorCode ierr; PetscMPIInt size,rank; PetscScalar *vals,rval; IS *is1,*is2; PetscRandom rdm; Vec xx,s1,s2; PetscReal s1norm,s2norm,rnorm,tol = 100*PETSC_SMALL; PetscBool flg,test_nd0=PETSC_FALSE; 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_size",&m,NULL)); CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL)); CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL)); CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_nd0",&test_nd0,NULL)); /* Create a AIJ matrix A */ CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); CHKERRQ(MatSetSizes(A,m*bs,m*bs,PETSC_DECIDE,PETSC_DECIDE)); CHKERRQ(MatSetType(A,MATAIJ)); CHKERRQ(MatSeqAIJSetPreallocation(A,PETSC_DEFAULT,NULL)); CHKERRQ(MatMPIAIJSetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL)); CHKERRQ(MatSetFromOptions(A)); CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); /* Create a BAIJ matrix B */ CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); CHKERRQ(MatSetSizes(B,m*bs,m*bs,PETSC_DECIDE,PETSC_DECIDE)); CHKERRQ(MatSetType(B,MATBAIJ)); CHKERRQ(MatSeqBAIJSetPreallocation(B,bs,PETSC_DEFAULT,NULL)); CHKERRQ(MatMPIBAIJSetPreallocation(B,bs,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL)); CHKERRQ(MatSetFromOptions(B)); CHKERRQ(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rdm)); CHKERRQ(PetscRandomSetFromOptions(rdm)); 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 (i=0; i<40*bs; i++) { CHKERRQ(PetscRandomGetValue(rdm,&rval)); cols[0] = bs*(int)(PetscRealPart(rval)*Mbs); CHKERRQ(PetscRandomGetValue(rdm,&rval)); rows[0] = rstart + bs*(int)(PetscRealPart(rval)*m); for (j=1; jtol) { CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",rank,(double)s1norm,(double)s2norm)); } } CHKERRQ(VecDestroy(&xx)); CHKERRQ(VecDestroy(&s1)); CHKERRQ(VecDestroy(&s2)); } /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */ CHKERRQ(MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA)); CHKERRQ(MatCreateSubMatrices(B,nd,is2,is2,MAT_REUSE_MATRIX,&submatB)); /* Test MatMult() */ for (i=0; itol) { CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",rank,(double)s1norm,(double)s2norm)); } } CHKERRQ(VecDestroy(&xx)); CHKERRQ(VecDestroy(&s1)); CHKERRQ(VecDestroy(&s2)); } /* Free allocated memory */ for (i=0; i