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; 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; PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 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)); PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_nd0",&test_nd0,NULL)); /* Create a AIJ matrix A */ PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); PetscCall(MatSetSizes(A,m*bs,m*bs,PETSC_DECIDE,PETSC_DECIDE)); PetscCall(MatSetType(A,MATAIJ)); PetscCall(MatSeqAIJSetPreallocation(A,PETSC_DEFAULT,NULL)); PetscCall(MatMPIAIJSetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL)); PetscCall(MatSetFromOptions(A)); PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); /* Create a BAIJ matrix B */ PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); PetscCall(MatSetSizes(B,m*bs,m*bs,PETSC_DECIDE,PETSC_DECIDE)); PetscCall(MatSetType(B,MATBAIJ)); PetscCall(MatSeqBAIJSetPreallocation(B,bs,PETSC_DEFAULT,NULL)); PetscCall(MatMPIBAIJSetPreallocation(B,bs,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL)); PetscCall(MatSetFromOptions(B)); PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rdm)); PetscCall(PetscRandomSetFromOptions(rdm)); PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); PetscCall(MatGetSize(A,&M,&N)); Mbs = M/bs; 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<40*bs; i++) { PetscCall(PetscRandomGetValue(rdm,&rval)); cols[0] = bs*(int)(PetscRealPart(rval)*Mbs); PetscCall(PetscRandomGetValue(rdm,&rval)); rows[0] = rstart + bs*(int)(PetscRealPart(rval)*m); for (j=1; jtol) { PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",rank,(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,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",rank,(double)s1norm,(double)s2norm)); } } PetscCall(VecDestroy(&xx)); PetscCall(VecDestroy(&s1)); PetscCall(VecDestroy(&s2)); } /* Free allocated memory */ for (i=0; i