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; ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetBool(NULL,NULL,"-test_nd0",&test_nd0,NULL);CHKERRQ(ierr); /* Create a AIJ matrix A */ ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); ierr = MatSetSizes(A,m*bs,m*bs,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); ierr = MatSeqAIJSetPreallocation(A,PETSC_DEFAULT,NULL);CHKERRQ(ierr); ierr = MatMPIAIJSetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);CHKERRQ(ierr); ierr = MatSetFromOptions(A);CHKERRQ(ierr); ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); /* Create a BAIJ matrix B */ ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); ierr = MatSetSizes(B,m*bs,m*bs,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); ierr = MatSetType(B,MATBAIJ);CHKERRQ(ierr); ierr = MatSeqBAIJSetPreallocation(B,bs,PETSC_DEFAULT,NULL);CHKERRQ(ierr); ierr = MatMPIBAIJSetPreallocation(B,bs,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);CHKERRQ(ierr); ierr = MatSetFromOptions(B);CHKERRQ(ierr); ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rdm);CHKERRQ(ierr); ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr); ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); Mbs = M/bs; ierr = PetscMalloc1(bs,&rows);CHKERRQ(ierr); ierr = PetscMalloc1(bs,&cols);CHKERRQ(ierr); ierr = PetscMalloc1(bs*bs,&vals);CHKERRQ(ierr); ierr = PetscMalloc1(M,&idx);CHKERRQ(ierr); /* Now set blocks of values */ for (i=0; i<40*bs; i++) { ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); cols[0] = bs*(int)(PetscRealPart(rval)*Mbs); ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); rows[0] = rstart + bs*(int)(PetscRealPart(rval)*m); for (j=1; jtol) { ierr = PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",rank,s1norm,s2norm);CHKERRQ(ierr); } } ierr = VecDestroy(&xx);CHKERRQ(ierr); ierr = VecDestroy(&s1);CHKERRQ(ierr); ierr = VecDestroy(&s2);CHKERRQ(ierr); } /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */ ierr = MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA);CHKERRQ(ierr); ierr = MatCreateSubMatrices(B,nd,is2,is2,MAT_REUSE_MATRIX,&submatB);CHKERRQ(ierr); /* Test MatMult() */ for (i=0; itol) { ierr = PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",rank,s1norm,s2norm);CHKERRQ(ierr); } } ierr = VecDestroy(&xx);CHKERRQ(ierr); ierr = VecDestroy(&s1);CHKERRQ(ierr); ierr = VecDestroy(&s2);CHKERRQ(ierr); } /* Free allocated memory */ for (i=0; i