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