1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for parallel AIJ and BAIJ formats.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat A,B,*submatA,*submatB; 9c4762a1bSJed Brown PetscInt bs=1,m=11,ov=1,i,j,k,*rows,*cols,nd=5,*idx,rstart,rend,sz,mm,nn,M,N,Mbs; 10c4762a1bSJed Brown PetscMPIInt size,rank; 11c4762a1bSJed Brown PetscScalar *vals,rval; 12c4762a1bSJed Brown IS *is1,*is2; 13c4762a1bSJed Brown PetscRandom rdm; 14c4762a1bSJed Brown Vec xx,s1,s2; 15c4762a1bSJed Brown PetscReal s1norm,s2norm,rnorm,tol = 100*PETSC_SMALL; 16c4762a1bSJed Brown PetscBool flg,test_nd0=PETSC_FALSE; 17c4762a1bSJed Brown 18*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 195f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 205f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 21c4762a1bSJed Brown 225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL)); 255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_nd0",&test_nd0,NULL)); 27c4762a1bSJed Brown 28c4762a1bSJed Brown /* Create a AIJ matrix A */ 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,m*bs,m*bs,PETSC_DECIDE,PETSC_DECIDE)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATAIJ)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqAIJSetPreallocation(A,PETSC_DEFAULT,NULL)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIAIJSetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL)); 345f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 36c4762a1bSJed Brown 37c4762a1bSJed Brown /* Create a BAIJ matrix B */ 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B,m*bs,m*bs,PETSC_DECIDE,PETSC_DECIDE)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B,MATBAIJ)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(MatSeqBAIJSetPreallocation(B,bs,PETSC_DEFAULT,NULL)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatMPIBAIJSetPreallocation(B,bs,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 45c4762a1bSJed Brown 465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rdm)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rdm)); 48c4762a1bSJed Brown 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&M,&N)); 51c4762a1bSJed Brown Mbs = M/bs; 52c4762a1bSJed Brown 535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(bs,&rows)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(bs,&cols)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(bs*bs,&vals)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(M,&idx)); 57c4762a1bSJed Brown 58c4762a1bSJed Brown /* Now set blocks of values */ 59c4762a1bSJed Brown for (i=0; i<40*bs; i++) { 605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rdm,&rval)); 61c4762a1bSJed Brown cols[0] = bs*(int)(PetscRealPart(rval)*Mbs); 625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rdm,&rval)); 63c4762a1bSJed Brown rows[0] = rstart + bs*(int)(PetscRealPart(rval)*m); 64c4762a1bSJed Brown for (j=1; j<bs; j++) { 65c4762a1bSJed Brown rows[j] = rows[j-1]+1; 66c4762a1bSJed Brown cols[j] = cols[j-1]+1; 67c4762a1bSJed Brown } 68c4762a1bSJed Brown 69c4762a1bSJed Brown for (j=0; j<bs*bs; j++) { 705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rdm,&rval)); 71c4762a1bSJed Brown vals[j] = rval; 72c4762a1bSJed Brown } 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,bs,rows,bs,cols,vals,ADD_VALUES)); 75c4762a1bSJed Brown } 76c4762a1bSJed Brown 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 81c4762a1bSJed Brown 82c4762a1bSJed Brown /* Test MatIncreaseOverlap() */ 835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nd,&is1)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nd,&is2)); 85c4762a1bSJed Brown 86dd400576SPatrick Sanan if (rank == 0 && test_nd0) nd = 0; /* test case */ 87c4762a1bSJed Brown 88c4762a1bSJed Brown for (i=0; i<nd; i++) { 895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rdm,&rval)); 90c4762a1bSJed Brown sz = (int)(PetscRealPart(rval)*m); 91c4762a1bSJed Brown for (j=0; j<sz; j++) { 925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rdm,&rval)); 93c4762a1bSJed Brown idx[j*bs] = bs*(int)(PetscRealPart(rval)*Mbs); 94c4762a1bSJed Brown for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k; 95c4762a1bSJed Brown } 965f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,sz*bs,idx,PETSC_COPY_VALUES,is1+i)); 975f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,sz*bs,idx,PETSC_COPY_VALUES,is2+i)); 98c4762a1bSJed Brown } 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatIncreaseOverlap(A,nd,is1,ov)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(MatIncreaseOverlap(B,nd,is2,ov)); 101c4762a1bSJed Brown 102c4762a1bSJed Brown for (i=0; i<nd; ++i) { 1035f80ce2aSJacob Faibussowitsch CHKERRQ(ISEqual(is1[i],is2[i],&flg)); 104c4762a1bSJed Brown 105c4762a1bSJed Brown if (!flg) { 1065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"i=%" PetscInt_FMT ", flg=%d :bs=%" PetscInt_FMT " m=%" PetscInt_FMT " ov=%" PetscInt_FMT " nd=%" PetscInt_FMT " np=%d\n",i,flg,bs,m,ov,nd,size)); 107c4762a1bSJed Brown } 108c4762a1bSJed Brown } 109c4762a1bSJed Brown 110c4762a1bSJed Brown for (i=0; i<nd; ++i) { 1115f80ce2aSJacob Faibussowitsch CHKERRQ(ISSort(is1[i])); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(ISSort(is2[i])); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown 1155f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrices(B,nd,is2,is2,MAT_INITIAL_MATRIX,&submatB)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA)); 117c4762a1bSJed Brown 118c4762a1bSJed Brown /* Test MatMult() */ 119c4762a1bSJed Brown for (i=0; i<nd; i++) { 1205f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(submatA[i],&mm,&nn)); 1215f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,mm,&xx)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(xx,&s1)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(xx,&s2)); 124c4762a1bSJed Brown for (j=0; j<3; j++) { 1255f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(xx,rdm)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(submatA[i],xx,s1)); 1275f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(submatB[i],xx,s2)); 1285f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s1,NORM_2,&s1norm)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s2,NORM_2,&s2norm)); 130c4762a1bSJed Brown rnorm = s2norm-s1norm; 131c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 1325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",rank,(double)s1norm,(double)s2norm)); 133c4762a1bSJed Brown } 134c4762a1bSJed Brown } 1355f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&xx)); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&s1)); 1375f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&s2)); 138c4762a1bSJed Brown } 139c4762a1bSJed Brown 140c4762a1bSJed Brown /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */ 1415f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrices(B,nd,is2,is2,MAT_REUSE_MATRIX,&submatB)); 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* Test MatMult() */ 145c4762a1bSJed Brown for (i=0; i<nd; i++) { 1465f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(submatA[i],&mm,&nn)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,mm,&xx)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(xx,&s1)); 1495f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(xx,&s2)); 150c4762a1bSJed Brown for (j=0; j<3; j++) { 1515f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(xx,rdm)); 1525f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(submatA[i],xx,s1)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(submatB[i],xx,s2)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s1,NORM_2,&s1norm)); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(s2,NORM_2,&s2norm)); 156c4762a1bSJed Brown rnorm = s2norm-s1norm; 157c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 1585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",rank,(double)s1norm,(double)s2norm)); 159c4762a1bSJed Brown } 160c4762a1bSJed Brown } 1615f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&xx)); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&s1)); 1635f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&s2)); 164c4762a1bSJed Brown } 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* Free allocated memory */ 167c4762a1bSJed Brown for (i=0; i<nd; ++i) { 1685f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is1[i])); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&is2[i])); 170c4762a1bSJed Brown } 1715f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroySubMatrices(nd,&submatA)); 1725f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroySubMatrices(nd,&submatB)); 173c4762a1bSJed Brown 1745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(is1)); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(is2)); 1765f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(idx)); 1775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(rows)); 1785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(cols)); 1795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(vals)); 1805f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 1815f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 1825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rdm)); 183*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 184*b122ec5aSJacob Faibussowitsch return 0; 185c4762a1bSJed Brown } 186c4762a1bSJed Brown 187c4762a1bSJed Brown /*TEST 188c4762a1bSJed Brown 189c4762a1bSJed Brown test: 190c4762a1bSJed Brown nsize: {{1 3}} 191c4762a1bSJed Brown args: -mat_block_size {{1 3 4 6 8}} -ov {{1 3}} -mat_size {{11 13}} -nd {{7}} ; done 192c4762a1bSJed Brown output_file: output/ex54.out 193c4762a1bSJed Brown 194c4762a1bSJed Brown test: 195c4762a1bSJed Brown suffix: 2 196c4762a1bSJed Brown args: -nd 2 -test_nd0 197c4762a1bSJed Brown output_file: output/ex54.out 198c4762a1bSJed Brown 199c4762a1bSJed Brown test: 200c4762a1bSJed Brown suffix: 3 201c4762a1bSJed Brown nsize: 3 202c4762a1bSJed Brown args: -nd 2 -test_nd0 203c4762a1bSJed Brown output_file: output/ex54.out 204c4762a1bSJed Brown 205c4762a1bSJed Brown TEST*/ 206