1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for MatBAIJ format.\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=43,ov=1,i,j,k,*rows,*cols,M,nd=5,*idx,mm,nn,lsize; 10c4762a1bSJed Brown PetscScalar *vals,rval; 11c4762a1bSJed Brown IS *is1,*is2; 12c4762a1bSJed Brown PetscRandom rdm; 13c4762a1bSJed Brown Vec xx,s1,s2; 14c4762a1bSJed Brown PetscReal s1norm,s2norm,rnorm,tol = PETSC_SQRT_MACHINE_EPSILON; 15c4762a1bSJed Brown PetscBool flg; 16c4762a1bSJed Brown 17*327415f7SBarry Smith PetscFunctionBeginUser; 189566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL)); 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL)); 219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL)); 229566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL)); 23c4762a1bSJed Brown M = m*bs; 24c4762a1bSJed Brown 259566063dSJacob Faibussowitsch PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,1,NULL,&A)); 269566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 279566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,M,M,15,NULL,&B)); 289566063dSJacob Faibussowitsch PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 299566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rdm)); 309566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 31c4762a1bSJed Brown 329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs,&rows)); 339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs,&cols)); 349566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs*bs,&vals)); 359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M,&idx)); 36c4762a1bSJed Brown 37c4762a1bSJed Brown /* Now set blocks of values */ 38c4762a1bSJed Brown for (i=0; i<20*bs; i++) { 399566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 40c4762a1bSJed Brown cols[0] = bs*(int)(PetscRealPart(rval)*m); 419566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 42c4762a1bSJed Brown rows[0] = bs*(int)(PetscRealPart(rval)*m); 43c4762a1bSJed Brown for (j=1; j<bs; j++) { 44c4762a1bSJed Brown rows[j] = rows[j-1]+1; 45c4762a1bSJed Brown cols[j] = cols[j-1]+1; 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 48c4762a1bSJed Brown for (j=0; j<bs*bs; j++) { 499566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 50c4762a1bSJed Brown vals[j] = rval; 51c4762a1bSJed Brown } 529566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES)); 539566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,bs,rows,bs,cols,vals,ADD_VALUES)); 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 569566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 579566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 589566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 599566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* Test MatIncreaseOverlap() */ 629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd,&is1)); 639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd,&is2)); 64c4762a1bSJed Brown 65c4762a1bSJed Brown for (i=0; i<nd; i++) { 669566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 67c4762a1bSJed Brown lsize = (int)(PetscRealPart(rval)*m); 68c4762a1bSJed Brown for (j=0; j<lsize; j++) { 699566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 70c4762a1bSJed Brown idx[j*bs] = bs*(int)(PetscRealPart(rval)*m); 71c4762a1bSJed Brown for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k; 72c4762a1bSJed Brown } 739566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,lsize*bs,idx,PETSC_COPY_VALUES,is1+i)); 749566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,lsize*bs,idx,PETSC_COPY_VALUES,is2+i)); 75c4762a1bSJed Brown } 769566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A,nd,is1,ov)); 779566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(B,nd,is2,ov)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown for (i=0; i<nd; ++i) { 809566063dSJacob Faibussowitsch PetscCall(ISEqual(is1[i],is2[i],&flg)); 8128b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"i=%" PetscInt_FMT ", flg =%d",i,(int)flg); 82c4762a1bSJed Brown } 83c4762a1bSJed Brown 84c4762a1bSJed Brown for (i=0; i<nd; ++i) { 859566063dSJacob Faibussowitsch PetscCall(ISSort(is1[i])); 869566063dSJacob Faibussowitsch PetscCall(ISSort(is2[i])); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown 899566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA)); 909566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B,nd,is2,is2,MAT_INITIAL_MATRIX,&submatB)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* Test MatMult() */ 93c4762a1bSJed Brown for (i=0; i<nd; i++) { 949566063dSJacob Faibussowitsch PetscCall(MatGetSize(submatA[i],&mm,&nn)); 959566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,mm,&xx)); 969566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s1)); 979566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s2)); 98c4762a1bSJed Brown for (j=0; j<3; j++) { 999566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rdm)); 1009566063dSJacob Faibussowitsch PetscCall(MatMult(submatA[i],xx,s1)); 1019566063dSJacob Faibussowitsch PetscCall(MatMult(submatB[i],xx,s2)); 1029566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_2,&s1norm)); 1039566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&s2norm)); 104c4762a1bSJed Brown rnorm = s2norm-s1norm; 105c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 1069566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm)); 107c4762a1bSJed Brown } 108c4762a1bSJed Brown } 1099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 1109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 1119566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */ 1149566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA)); 1159566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B,nd,is2,is2,MAT_REUSE_MATRIX,&submatB)); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* Test MatMult() */ 118c4762a1bSJed Brown for (i=0; i<nd; i++) { 1199566063dSJacob Faibussowitsch PetscCall(MatGetSize(submatA[i],&mm,&nn)); 1209566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,mm,&xx)); 1219566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s1)); 1229566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s2)); 123c4762a1bSJed Brown for (j=0; j<3; j++) { 1249566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rdm)); 1259566063dSJacob Faibussowitsch PetscCall(MatMult(submatA[i],xx,s1)); 1269566063dSJacob Faibussowitsch PetscCall(MatMult(submatB[i],xx,s2)); 1279566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_2,&s1norm)); 1289566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&s2norm)); 129c4762a1bSJed Brown rnorm = s2norm-s1norm; 130c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 1319566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm)); 132c4762a1bSJed Brown } 133c4762a1bSJed Brown } 1349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 1359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 1369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* Free allocated memory */ 140c4762a1bSJed Brown for (i=0; i<nd; ++i) { 1419566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1[i])); 1429566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2[i])); 143c4762a1bSJed Brown } 1449566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd,&submatA)); 1459566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd,&submatB)); 1469566063dSJacob Faibussowitsch PetscCall(PetscFree(is1)); 1479566063dSJacob Faibussowitsch PetscCall(PetscFree(is2)); 1489566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 1499566063dSJacob Faibussowitsch PetscCall(PetscFree(rows)); 1509566063dSJacob Faibussowitsch PetscCall(PetscFree(cols)); 1519566063dSJacob Faibussowitsch PetscCall(PetscFree(vals)); 1529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1539566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1549566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 1559566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 156b122ec5aSJacob Faibussowitsch return 0; 157c4762a1bSJed Brown } 158c4762a1bSJed Brown 159c4762a1bSJed Brown /*TEST 160c4762a1bSJed Brown 161c4762a1bSJed Brown test: 162c4762a1bSJed Brown args: -mat_block_size {{1 2 5 7 8}} -ov {{1 3}} -mat_size {{11 13}} -nd {{7}} 163c4762a1bSJed Brown 164c4762a1bSJed Brown TEST*/ 165