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