static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for MatBAIJ format.\n"; #include int main(int argc,char **args) { Mat A,B,*submatA,*submatB; PetscInt bs=1,m=43,ov=1,i,j,k,*rows,*cols,M,nd=5,*idx,mm,nn,lsize; PetscScalar *vals,rval; IS *is1,*is2; PetscRandom rdm; Vec xx,s1,s2; PetscReal s1norm,s2norm,rnorm,tol = PETSC_SQRT_MACHINE_EPSILON; PetscBool flg; PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL)); PetscCall(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL)); M = m*bs; PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,1,NULL,&A)); PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,M,M,15,NULL,&B)); PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rdm)); PetscCall(PetscRandomSetFromOptions(rdm)); PetscCall(PetscMalloc1(bs,&rows)); PetscCall(PetscMalloc1(bs,&cols)); PetscCall(PetscMalloc1(bs*bs,&vals)); PetscCall(PetscMalloc1(M,&idx)); /* Now set blocks of values */ for (i=0; i<20*bs; i++) { PetscCall(PetscRandomGetValue(rdm,&rval)); cols[0] = bs*(int)(PetscRealPart(rval)*m); PetscCall(PetscRandomGetValue(rdm,&rval)); rows[0] = bs*(int)(PetscRealPart(rval)*m); for (j=1; jtol) { PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm)); } } PetscCall(VecDestroy(&xx)); PetscCall(VecDestroy(&s1)); PetscCall(VecDestroy(&s2)); } /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */ PetscCall(MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA)); PetscCall(MatCreateSubMatrices(B,nd,is2,is2,MAT_REUSE_MATRIX,&submatB)); /* Test MatMult() */ for (i=0; itol) { PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm)); } } PetscCall(VecDestroy(&xx)); PetscCall(VecDestroy(&s1)); PetscCall(VecDestroy(&s2)); } /* Free allocated memory */ for (i=0; i