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