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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 18*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL)); 19*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL)); 20*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL)); 21*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL)); 22c4762a1bSJed Brown M = m*bs; 23c4762a1bSJed Brown 24*9566063dSJacob Faibussowitsch PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,1,NULL,&A)); 25*9566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 26*9566063dSJacob Faibussowitsch PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF,M,M,15,NULL,&B)); 27*9566063dSJacob Faibussowitsch PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 28*9566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rdm)); 29*9566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rdm)); 30c4762a1bSJed Brown 31*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs,&rows)); 32*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs,&cols)); 33*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs*bs,&vals)); 34*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(M,&idx)); 35c4762a1bSJed Brown 36c4762a1bSJed Brown /* Now set blocks of values */ 37c4762a1bSJed Brown for (i=0; i<20*bs; i++) { 38*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 39c4762a1bSJed Brown cols[0] = bs*(int)(PetscRealPart(rval)*m); 40*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 41c4762a1bSJed Brown rows[0] = bs*(int)(PetscRealPart(rval)*m); 42c4762a1bSJed Brown for (j=1; j<bs; j++) { 43c4762a1bSJed Brown rows[j] = rows[j-1]+1; 44c4762a1bSJed Brown cols[j] = cols[j-1]+1; 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47c4762a1bSJed Brown for (j=0; j<bs*bs; j++) { 48*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 49c4762a1bSJed Brown vals[j] = rval; 50c4762a1bSJed Brown } 51*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES)); 52*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,bs,rows,bs,cols,vals,ADD_VALUES)); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown 55*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 56*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 57*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 58*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Test MatIncreaseOverlap() */ 61*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd,&is1)); 62*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd,&is2)); 63c4762a1bSJed Brown 64c4762a1bSJed Brown for (i=0; i<nd; i++) { 65*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 66c4762a1bSJed Brown lsize = (int)(PetscRealPart(rval)*m); 67c4762a1bSJed Brown for (j=0; j<lsize; j++) { 68*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rdm,&rval)); 69c4762a1bSJed Brown idx[j*bs] = bs*(int)(PetscRealPart(rval)*m); 70c4762a1bSJed Brown for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k; 71c4762a1bSJed Brown } 72*9566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,lsize*bs,idx,PETSC_COPY_VALUES,is1+i)); 73*9566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,lsize*bs,idx,PETSC_COPY_VALUES,is2+i)); 74c4762a1bSJed Brown } 75*9566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A,nd,is1,ov)); 76*9566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(B,nd,is2,ov)); 77c4762a1bSJed Brown 78c4762a1bSJed Brown for (i=0; i<nd; ++i) { 79*9566063dSJacob Faibussowitsch PetscCall(ISEqual(is1[i],is2[i],&flg)); 8028b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"i=%" PetscInt_FMT ", flg =%d",i,(int)flg); 81c4762a1bSJed Brown } 82c4762a1bSJed Brown 83c4762a1bSJed Brown for (i=0; i<nd; ++i) { 84*9566063dSJacob Faibussowitsch PetscCall(ISSort(is1[i])); 85*9566063dSJacob Faibussowitsch PetscCall(ISSort(is2[i])); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 88*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA)); 89*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B,nd,is2,is2,MAT_INITIAL_MATRIX,&submatB)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* Test MatMult() */ 92c4762a1bSJed Brown for (i=0; i<nd; i++) { 93*9566063dSJacob Faibussowitsch PetscCall(MatGetSize(submatA[i],&mm,&nn)); 94*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,mm,&xx)); 95*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s1)); 96*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s2)); 97c4762a1bSJed Brown for (j=0; j<3; j++) { 98*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rdm)); 99*9566063dSJacob Faibussowitsch PetscCall(MatMult(submatA[i],xx,s1)); 100*9566063dSJacob Faibussowitsch PetscCall(MatMult(submatB[i],xx,s2)); 101*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_2,&s1norm)); 102*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&s2norm)); 103c4762a1bSJed Brown rnorm = s2norm-s1norm; 104c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 105*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm)); 106c4762a1bSJed Brown } 107c4762a1bSJed Brown } 108*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 109*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 110*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */ 113*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA)); 114*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(B,nd,is2,is2,MAT_REUSE_MATRIX,&submatB)); 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* Test MatMult() */ 117c4762a1bSJed Brown for (i=0; i<nd; i++) { 118*9566063dSJacob Faibussowitsch PetscCall(MatGetSize(submatA[i],&mm,&nn)); 119*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,mm,&xx)); 120*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s1)); 121*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s2)); 122c4762a1bSJed Brown for (j=0; j<3; j++) { 123*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rdm)); 124*9566063dSJacob Faibussowitsch PetscCall(MatMult(submatA[i],xx,s1)); 125*9566063dSJacob Faibussowitsch PetscCall(MatMult(submatB[i],xx,s2)); 126*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_2,&s1norm)); 127*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&s2norm)); 128c4762a1bSJed Brown rnorm = s2norm-s1norm; 129c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 130*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm)); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown } 133*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 134*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 135*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* Free allocated memory */ 139c4762a1bSJed Brown for (i=0; i<nd; ++i) { 140*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1[i])); 141*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2[i])); 142c4762a1bSJed Brown } 143*9566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd,&submatA)); 144*9566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd,&submatB)); 145*9566063dSJacob Faibussowitsch PetscCall(PetscFree(is1)); 146*9566063dSJacob Faibussowitsch PetscCall(PetscFree(is2)); 147*9566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 148*9566063dSJacob Faibussowitsch PetscCall(PetscFree(rows)); 149*9566063dSJacob Faibussowitsch PetscCall(PetscFree(cols)); 150*9566063dSJacob Faibussowitsch PetscCall(PetscFree(vals)); 151*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 152*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 153*9566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rdm)); 154*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 155b122ec5aSJacob Faibussowitsch return 0; 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158c4762a1bSJed Brown /*TEST 159c4762a1bSJed Brown 160c4762a1bSJed Brown test: 161c4762a1bSJed Brown args: -mat_block_size {{1 2 5 7 8}} -ov {{1 3}} -mat_size {{11 13}} -nd {{7}} 162c4762a1bSJed Brown 163c4762a1bSJed Brown TEST*/ 164