1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for sequential MatSBAIJ format. Derived from ex51.c\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat A,Atrans,sA,*submatA,*submatsA; 9c4762a1bSJed Brown PetscInt bs=1,m=43,ov=1,i,j,k,*rows,*cols,M,nd=5,*idx,mm,nn; 10c4762a1bSJed Brown PetscMPIInt size; 11c4762a1bSJed Brown PetscScalar *vals,rval,one=1.0; 12c4762a1bSJed Brown IS *is1,*is2; 13c4762a1bSJed Brown PetscRandom rand; 14c4762a1bSJed Brown Vec xx,s1,s2; 15c4762a1bSJed Brown PetscReal s1norm,s2norm,rnorm,tol = 10*PETSC_SMALL; 16c4762a1bSJed Brown PetscBool flg; 17c4762a1bSJed Brown 18*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 19*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL)); 20*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL)); 21*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL)); 22*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL)); 23c4762a1bSJed Brown 24c4762a1bSJed Brown /* create a SeqBAIJ matrix A */ 25c4762a1bSJed Brown M = m*bs; 26*9566063dSJacob Faibussowitsch PetscCall(MatCreateSeqBAIJ(PETSC_COMM_SELF,bs,M,M,1,NULL,&A)); 27*9566063dSJacob Faibussowitsch PetscCall(MatSetOption(A, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 28*9566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&rand)); 29*9566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 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 random values */ 37c4762a1bSJed Brown /* first, set diagonal blocks as zero */ 38c4762a1bSJed Brown for (j=0; j<bs*bs; j++) vals[j] = 0.0; 39c4762a1bSJed Brown for (i=0; i<m; i++) { 40c4762a1bSJed Brown cols[0] = i*bs; rows[0] = i*bs; 41c4762a1bSJed Brown for (j=1; j<bs; j++) { 42c4762a1bSJed Brown rows[j] = rows[j-1]+1; 43c4762a1bSJed Brown cols[j] = cols[j-1]+1; 44c4762a1bSJed Brown } 45*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES)); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown /* second, add random blocks */ 48c4762a1bSJed Brown for (i=0; i<20*bs; i++) { 49*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 50c4762a1bSJed Brown cols[0] = bs*(int)(PetscRealPart(rval)*m); 51*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 52c4762a1bSJed Brown rows[0] = bs*(int)(PetscRealPart(rval)*m); 53c4762a1bSJed Brown for (j=1; j<bs; j++) { 54c4762a1bSJed Brown rows[j] = rows[j-1]+1; 55c4762a1bSJed Brown cols[j] = cols[j-1]+1; 56c4762a1bSJed Brown } 57c4762a1bSJed Brown 58c4762a1bSJed Brown for (j=0; j<bs*bs; j++) { 59*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 60c4762a1bSJed Brown vals[j] = rval; 61c4762a1bSJed Brown } 62*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES)); 63c4762a1bSJed Brown } 64c4762a1bSJed Brown 65*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 66*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* make A a symmetric matrix: A <- A^T + A */ 69*9566063dSJacob Faibussowitsch PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX, &Atrans)); 70*9566063dSJacob Faibussowitsch PetscCall(MatAXPY(A,one,Atrans,DIFFERENT_NONZERO_PATTERN)); 71*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Atrans)); 72*9566063dSJacob Faibussowitsch PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX, &Atrans)); 73*9566063dSJacob Faibussowitsch PetscCall(MatEqual(A, Atrans, &flg)); 7428b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"A+A^T is non-symmetric"); 75*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Atrans)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* create a SeqSBAIJ matrix sA (= A) */ 78*9566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE)); 79*9566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&sA)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* Test sA==A through MatMult() */ 82c4762a1bSJed Brown for (i=0; i<nd; i++) { 83*9566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&mm,&nn)); 84*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,mm,&xx)); 85*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s1)); 86*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s2)); 87c4762a1bSJed Brown for (j=0; j<3; j++) { 88*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rand)); 89*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,xx,s1)); 90*9566063dSJacob Faibussowitsch PetscCall(MatMult(sA,xx,s2)); 91*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_2,&s1norm)); 92*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&s2norm)); 93c4762a1bSJed Brown rnorm = s2norm-s1norm; 94c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 95*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm)); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown } 98*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 99*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 100*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 101c4762a1bSJed Brown } 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* Test MatIncreaseOverlap() */ 104*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd,&is1)); 105*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nd,&is2)); 106c4762a1bSJed Brown 107c4762a1bSJed Brown for (i=0; i<nd; i++) { 108*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 109c4762a1bSJed Brown size = (int)(PetscRealPart(rval)*m); 110c4762a1bSJed Brown for (j=0; j<size; j++) { 111*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 112c4762a1bSJed Brown idx[j*bs] = bs*(int)(PetscRealPart(rval)*m); 113c4762a1bSJed Brown for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k; 114c4762a1bSJed Brown } 115*9566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,size*bs,idx,PETSC_COPY_VALUES,is1+i)); 116*9566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,size*bs,idx,PETSC_COPY_VALUES,is2+i)); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown /* for debugging */ 119c4762a1bSJed Brown /* 120*9566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_SELF)); 121*9566063dSJacob Faibussowitsch PetscCall(MatView(sA,PETSC_VIEWER_STDOUT_SELF)); 122c4762a1bSJed Brown */ 123c4762a1bSJed Brown 124*9566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(A,nd,is1,ov)); 125*9566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(sA,nd,is2,ov)); 126c4762a1bSJed Brown 127c4762a1bSJed Brown for (i=0; i<nd; ++i) { 128*9566063dSJacob Faibussowitsch PetscCall(ISSort(is1[i])); 129*9566063dSJacob Faibussowitsch PetscCall(ISSort(is2[i])); 130c4762a1bSJed Brown } 131c4762a1bSJed Brown 132c4762a1bSJed Brown for (i=0; i<nd; ++i) { 133*9566063dSJacob Faibussowitsch PetscCall(ISEqual(is1[i],is2[i],&flg)); 13428b400f6SJacob Faibussowitsch PetscCheck(flg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"i=%" PetscInt_FMT ", is1 != is2",i); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown 137*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA)); 138*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(sA,nd,is2,is2,MAT_INITIAL_MATRIX,&submatsA)); 139c4762a1bSJed Brown 140c4762a1bSJed Brown /* Test MatMult() */ 141c4762a1bSJed Brown for (i=0; i<nd; i++) { 142*9566063dSJacob Faibussowitsch PetscCall(MatGetSize(submatA[i],&mm,&nn)); 143*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,mm,&xx)); 144*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s1)); 145*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s2)); 146c4762a1bSJed Brown for (j=0; j<3; j++) { 147*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rand)); 148*9566063dSJacob Faibussowitsch PetscCall(MatMult(submatA[i],xx,s1)); 149*9566063dSJacob Faibussowitsch PetscCall(MatMult(submatsA[i],xx,s2)); 150*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_2,&s1norm)); 151*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&s2norm)); 152c4762a1bSJed Brown rnorm = s2norm-s1norm; 153c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 154*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm)); 155c4762a1bSJed Brown } 156c4762a1bSJed Brown } 157*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 158*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 159*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 162c4762a1bSJed Brown /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */ 163*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA)); 164*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(sA,nd,is2,is2,MAT_REUSE_MATRIX,&submatsA)); 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* Test MatMult() */ 167c4762a1bSJed Brown for (i=0; i<nd; i++) { 168*9566063dSJacob Faibussowitsch PetscCall(MatGetSize(submatA[i],&mm,&nn)); 169*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,mm,&xx)); 170*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s1)); 171*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx,&s2)); 172c4762a1bSJed Brown for (j=0; j<3; j++) { 173*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(xx,rand)); 174*9566063dSJacob Faibussowitsch PetscCall(MatMult(submatA[i],xx,s1)); 175*9566063dSJacob Faibussowitsch PetscCall(MatMult(submatsA[i],xx,s2)); 176*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s1,NORM_2,&s1norm)); 177*9566063dSJacob Faibussowitsch PetscCall(VecNorm(s2,NORM_2,&s2norm)); 178c4762a1bSJed Brown rnorm = s2norm-s1norm; 179c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 180*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",(double)s1norm,(double)s2norm)); 181c4762a1bSJed Brown } 182c4762a1bSJed Brown } 183*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 184*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s1)); 185*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&s2)); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* Free allocated memory */ 189c4762a1bSJed Brown for (i=0; i<nd; ++i) { 190*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is1[i])); 191*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2[i])); 192c4762a1bSJed Brown } 193*9566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd,&submatA)); 194*9566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(nd,&submatsA)); 195c4762a1bSJed Brown 196*9566063dSJacob Faibussowitsch PetscCall(PetscFree(is1)); 197*9566063dSJacob Faibussowitsch PetscCall(PetscFree(is2)); 198*9566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 199*9566063dSJacob Faibussowitsch PetscCall(PetscFree(rows)); 200*9566063dSJacob Faibussowitsch PetscCall(PetscFree(cols)); 201*9566063dSJacob Faibussowitsch PetscCall(PetscFree(vals)); 202*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 203*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sA)); 204*9566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 205*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 206b122ec5aSJacob Faibussowitsch return 0; 207c4762a1bSJed Brown } 208c4762a1bSJed Brown 209c4762a1bSJed Brown /*TEST 210c4762a1bSJed Brown 211c4762a1bSJed Brown test: 212c4762a1bSJed Brown args: -ov 2 213c4762a1bSJed Brown 214c4762a1bSJed Brown TEST*/ 215