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