1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for parallel AIJ and BAIJ formats.\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=11,ov=1,i,j,k,*rows,*cols,nd=5,*idx,rstart,rend,sz,mm,nn,M,N,Mbs; 10*c4762a1bSJed Brown PetscErrorCode ierr; 11*c4762a1bSJed Brown PetscMPIInt size,rank; 12*c4762a1bSJed Brown PetscScalar *vals,rval; 13*c4762a1bSJed Brown IS *is1,*is2; 14*c4762a1bSJed Brown PetscRandom rdm; 15*c4762a1bSJed Brown Vec xx,s1,s2; 16*c4762a1bSJed Brown PetscReal s1norm,s2norm,rnorm,tol = 100*PETSC_SMALL; 17*c4762a1bSJed Brown PetscBool flg,test_nd0=PETSC_FALSE; 18*c4762a1bSJed Brown 19*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 20*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 21*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL);CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-mat_size",&m,NULL);CHKERRQ(ierr); 25*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL);CHKERRQ(ierr); 26*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-test_nd0",&test_nd0,NULL);CHKERRQ(ierr); 28*c4762a1bSJed Brown 29*c4762a1bSJed Brown /* Create a AIJ matrix A */ 30*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = MatSetSizes(A,m*bs,m*bs,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = MatSetType(A,MATAIJ);CHKERRQ(ierr); 33*c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(A,PETSC_DEFAULT,NULL);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(A,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);CHKERRQ(ierr); 35*c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown /* Create a BAIJ matrix B */ 39*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = MatSetSizes(B,m*bs,m*bs,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = MatSetType(B,MATBAIJ);CHKERRQ(ierr); 42*c4762a1bSJed Brown ierr = MatSeqBAIJSetPreallocation(B,bs,PETSC_DEFAULT,NULL);CHKERRQ(ierr); 43*c4762a1bSJed Brown ierr = MatMPIBAIJSetPreallocation(B,bs,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);CHKERRQ(ierr); 44*c4762a1bSJed Brown ierr = MatSetFromOptions(B);CHKERRQ(ierr); 45*c4762a1bSJed Brown ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rdm);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rdm);CHKERRQ(ierr); 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 51*c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 52*c4762a1bSJed Brown Mbs = M/bs; 53*c4762a1bSJed Brown 54*c4762a1bSJed Brown ierr = PetscMalloc1(bs,&rows);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = PetscMalloc1(bs,&cols);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = PetscMalloc1(bs*bs,&vals);CHKERRQ(ierr); 57*c4762a1bSJed Brown ierr = PetscMalloc1(M,&idx);CHKERRQ(ierr); 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown /* Now set blocks of values */ 60*c4762a1bSJed Brown for (i=0; i<40*bs; i++) { 61*c4762a1bSJed Brown ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); 62*c4762a1bSJed Brown cols[0] = bs*(int)(PetscRealPart(rval)*Mbs); 63*c4762a1bSJed Brown ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); 64*c4762a1bSJed Brown rows[0] = rstart + bs*(int)(PetscRealPart(rval)*m); 65*c4762a1bSJed Brown for (j=1; j<bs; j++) { 66*c4762a1bSJed Brown rows[j] = rows[j-1]+1; 67*c4762a1bSJed Brown cols[j] = cols[j-1]+1; 68*c4762a1bSJed Brown } 69*c4762a1bSJed Brown 70*c4762a1bSJed Brown for (j=0; j<bs*bs; j++) { 71*c4762a1bSJed Brown ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); 72*c4762a1bSJed Brown vals[j] = rval; 73*c4762a1bSJed Brown } 74*c4762a1bSJed Brown ierr = MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES);CHKERRQ(ierr); 75*c4762a1bSJed Brown ierr = MatSetValues(B,bs,rows,bs,cols,vals,ADD_VALUES);CHKERRQ(ierr); 76*c4762a1bSJed Brown } 77*c4762a1bSJed Brown 78*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 81*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 82*c4762a1bSJed Brown 83*c4762a1bSJed Brown /* Test MatIncreaseOverlap() */ 84*c4762a1bSJed Brown ierr = PetscMalloc1(nd,&is1);CHKERRQ(ierr); 85*c4762a1bSJed Brown ierr = PetscMalloc1(nd,&is2);CHKERRQ(ierr); 86*c4762a1bSJed Brown 87*c4762a1bSJed Brown if (!rank && test_nd0) nd = 0; /* test case */ 88*c4762a1bSJed Brown 89*c4762a1bSJed Brown for (i=0; i<nd; i++) { 90*c4762a1bSJed Brown ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); 91*c4762a1bSJed Brown sz = (int)(PetscRealPart(rval)*m); 92*c4762a1bSJed Brown for (j=0; j<sz; j++) { 93*c4762a1bSJed Brown ierr = PetscRandomGetValue(rdm,&rval);CHKERRQ(ierr); 94*c4762a1bSJed Brown idx[j*bs] = bs*(int)(PetscRealPart(rval)*Mbs); 95*c4762a1bSJed Brown for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k; 96*c4762a1bSJed Brown } 97*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,sz*bs,idx,PETSC_COPY_VALUES,is1+i);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,sz*bs,idx,PETSC_COPY_VALUES,is2+i);CHKERRQ(ierr); 99*c4762a1bSJed Brown } 100*c4762a1bSJed Brown ierr = MatIncreaseOverlap(A,nd,is1,ov);CHKERRQ(ierr); 101*c4762a1bSJed Brown ierr = MatIncreaseOverlap(B,nd,is2,ov);CHKERRQ(ierr); 102*c4762a1bSJed Brown 103*c4762a1bSJed Brown for (i=0; i<nd; ++i) { 104*c4762a1bSJed Brown ierr = ISEqual(is1[i],is2[i],&flg);CHKERRQ(ierr); 105*c4762a1bSJed Brown 106*c4762a1bSJed Brown if (!flg) { 107*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"i=%D, flg=%d :bs=%D m=%D ov=%D nd=%D np=%D\n",i,flg,bs,m,ov,nd,size);CHKERRQ(ierr); 108*c4762a1bSJed Brown } 109*c4762a1bSJed Brown } 110*c4762a1bSJed Brown 111*c4762a1bSJed Brown for (i=0; i<nd; ++i) { 112*c4762a1bSJed Brown ierr = ISSort(is1[i]);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = ISSort(is2[i]);CHKERRQ(ierr); 114*c4762a1bSJed Brown } 115*c4762a1bSJed Brown 116*c4762a1bSJed Brown ierr = MatCreateSubMatrices(B,nd,is2,is2,MAT_INITIAL_MATRIX,&submatB);CHKERRQ(ierr); 117*c4762a1bSJed Brown ierr = MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA);CHKERRQ(ierr); 118*c4762a1bSJed Brown 119*c4762a1bSJed Brown /* Test MatMult() */ 120*c4762a1bSJed Brown for (i=0; i<nd; i++) { 121*c4762a1bSJed Brown ierr = MatGetSize(submatA[i],&mm,&nn);CHKERRQ(ierr); 122*c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,mm,&xx);CHKERRQ(ierr); 123*c4762a1bSJed Brown ierr = VecDuplicate(xx,&s1);CHKERRQ(ierr); 124*c4762a1bSJed Brown ierr = VecDuplicate(xx,&s2);CHKERRQ(ierr); 125*c4762a1bSJed Brown for (j=0; j<3; j++) { 126*c4762a1bSJed Brown ierr = VecSetRandom(xx,rdm);CHKERRQ(ierr); 127*c4762a1bSJed Brown ierr = MatMult(submatA[i],xx,s1);CHKERRQ(ierr); 128*c4762a1bSJed Brown ierr = MatMult(submatB[i],xx,s2);CHKERRQ(ierr); 129*c4762a1bSJed Brown ierr = VecNorm(s1,NORM_2,&s1norm);CHKERRQ(ierr); 130*c4762a1bSJed Brown ierr = VecNorm(s2,NORM_2,&s2norm);CHKERRQ(ierr); 131*c4762a1bSJed Brown rnorm = s2norm-s1norm; 132*c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 133*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",rank,s1norm,s2norm);CHKERRQ(ierr); 134*c4762a1bSJed Brown } 135*c4762a1bSJed Brown } 136*c4762a1bSJed Brown ierr = VecDestroy(&xx);CHKERRQ(ierr); 137*c4762a1bSJed Brown ierr = VecDestroy(&s1);CHKERRQ(ierr); 138*c4762a1bSJed Brown ierr = VecDestroy(&s2);CHKERRQ(ierr); 139*c4762a1bSJed Brown } 140*c4762a1bSJed Brown 141*c4762a1bSJed Brown /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */ 142*c4762a1bSJed Brown ierr = MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA);CHKERRQ(ierr); 143*c4762a1bSJed Brown ierr = MatCreateSubMatrices(B,nd,is2,is2,MAT_REUSE_MATRIX,&submatB);CHKERRQ(ierr); 144*c4762a1bSJed Brown 145*c4762a1bSJed Brown /* Test MatMult() */ 146*c4762a1bSJed Brown for (i=0; i<nd; i++) { 147*c4762a1bSJed Brown ierr = MatGetSize(submatA[i],&mm,&nn);CHKERRQ(ierr); 148*c4762a1bSJed Brown ierr = VecCreateSeq(PETSC_COMM_SELF,mm,&xx);CHKERRQ(ierr); 149*c4762a1bSJed Brown ierr = VecDuplicate(xx,&s1);CHKERRQ(ierr); 150*c4762a1bSJed Brown ierr = VecDuplicate(xx,&s2);CHKERRQ(ierr); 151*c4762a1bSJed Brown for (j=0; j<3; j++) { 152*c4762a1bSJed Brown ierr = VecSetRandom(xx,rdm);CHKERRQ(ierr); 153*c4762a1bSJed Brown ierr = MatMult(submatA[i],xx,s1);CHKERRQ(ierr); 154*c4762a1bSJed Brown ierr = MatMult(submatB[i],xx,s2);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = VecNorm(s1,NORM_2,&s1norm);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = VecNorm(s2,NORM_2,&s2norm);CHKERRQ(ierr); 157*c4762a1bSJed Brown rnorm = s2norm-s1norm; 158*c4762a1bSJed Brown if (rnorm<-tol || rnorm>tol) { 159*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"[%d]Error:MatMult - Norm1=%16.14e Norm2=%16.14e\n",rank,s1norm,s2norm);CHKERRQ(ierr); 160*c4762a1bSJed Brown } 161*c4762a1bSJed Brown } 162*c4762a1bSJed Brown ierr = VecDestroy(&xx);CHKERRQ(ierr); 163*c4762a1bSJed Brown ierr = VecDestroy(&s1);CHKERRQ(ierr); 164*c4762a1bSJed Brown ierr = VecDestroy(&s2);CHKERRQ(ierr); 165*c4762a1bSJed Brown } 166*c4762a1bSJed Brown 167*c4762a1bSJed Brown /* Free allocated memory */ 168*c4762a1bSJed Brown for (i=0; i<nd; ++i) { 169*c4762a1bSJed Brown ierr = ISDestroy(&is1[i]);CHKERRQ(ierr); 170*c4762a1bSJed Brown ierr = ISDestroy(&is2[i]);CHKERRQ(ierr); 171*c4762a1bSJed Brown } 172*c4762a1bSJed Brown ierr = MatDestroySubMatrices(nd,&submatA);CHKERRQ(ierr); 173*c4762a1bSJed Brown ierr = MatDestroySubMatrices(nd,&submatB);CHKERRQ(ierr); 174*c4762a1bSJed Brown 175*c4762a1bSJed Brown ierr = PetscFree(is1);CHKERRQ(ierr); 176*c4762a1bSJed Brown ierr = PetscFree(is2);CHKERRQ(ierr); 177*c4762a1bSJed Brown ierr = PetscFree(idx);CHKERRQ(ierr); 178*c4762a1bSJed Brown ierr = PetscFree(rows);CHKERRQ(ierr); 179*c4762a1bSJed Brown ierr = PetscFree(cols);CHKERRQ(ierr); 180*c4762a1bSJed Brown ierr = PetscFree(vals);CHKERRQ(ierr); 181*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 182*c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 183*c4762a1bSJed Brown ierr = PetscRandomDestroy(&rdm);CHKERRQ(ierr); 184*c4762a1bSJed Brown ierr = PetscFinalize(); 185*c4762a1bSJed Brown return ierr; 186*c4762a1bSJed Brown } 187*c4762a1bSJed Brown 188*c4762a1bSJed Brown 189*c4762a1bSJed Brown /*TEST 190*c4762a1bSJed Brown 191*c4762a1bSJed Brown test: 192*c4762a1bSJed Brown nsize: {{1 3}} 193*c4762a1bSJed Brown args: -mat_block_size {{1 3 4 6 8}} -ov {{1 3}} -mat_size {{11 13}} -nd {{7}} ; done 194*c4762a1bSJed Brown output_file: output/ex54.out 195*c4762a1bSJed Brown 196*c4762a1bSJed Brown test: 197*c4762a1bSJed Brown suffix: 2 198*c4762a1bSJed Brown args: -nd 2 -test_nd0 199*c4762a1bSJed Brown output_file: output/ex54.out 200*c4762a1bSJed Brown 201*c4762a1bSJed Brown test: 202*c4762a1bSJed Brown suffix: 3 203*c4762a1bSJed Brown nsize: 3 204*c4762a1bSJed Brown args: -nd 2 -test_nd0 205*c4762a1bSJed Brown output_file: output/ex54.out 206*c4762a1bSJed Brown 207*c4762a1bSJed Brown TEST*/ 208