1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatIncreaseOverlap(), MatCreateSubMatrices() for parallel MatSBAIJ format.\n"; 3c4762a1bSJed Brown /* Example of usage: 4c4762a1bSJed Brown mpiexec -n 2 ./ex92 -nd 2 -ov 3 -mat_block_size 2 -view_id 0 -test_overlap -test_submat 5c4762a1bSJed Brown */ 6c4762a1bSJed Brown #include <petscmat.h> 7c4762a1bSJed Brown 8c4762a1bSJed Brown int main(int argc,char **args) 9c4762a1bSJed Brown { 10c4762a1bSJed Brown Mat A,Atrans,sA,*submatA,*submatsA; 11c4762a1bSJed Brown PetscErrorCode ierr; 12c4762a1bSJed Brown PetscMPIInt size,rank; 13c4762a1bSJed Brown PetscInt bs=1,mbs=10,ov=1,i,j,k,*rows,*cols,nd=2,*idx,rstart,rend,sz,M,N,Mbs; 14c4762a1bSJed Brown PetscScalar *vals,rval,one=1.0; 15c4762a1bSJed Brown IS *is1,*is2; 16c4762a1bSJed Brown PetscRandom rand; 17c4762a1bSJed Brown PetscBool flg,TestOverlap,TestSubMat,TestAllcols,test_sorted=PETSC_FALSE; 18c4762a1bSJed Brown PetscInt vid = -1; 19c4762a1bSJed Brown #if defined(PETSC_USE_LOG) 20c4762a1bSJed Brown PetscLogStage stages[2]; 21c4762a1bSJed Brown #endif 22c4762a1bSJed Brown 23c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 24ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 25ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 26c4762a1bSJed Brown 27c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-mat_block_size",&bs,NULL);CHKERRQ(ierr); 28c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-mat_mbs",&mbs,NULL);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-ov",&ov,NULL);CHKERRQ(ierr); 30c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-nd",&nd,NULL);CHKERRQ(ierr); 31c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-view_id",&vid,NULL);CHKERRQ(ierr); 32c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-test_overlap", &TestOverlap);CHKERRQ(ierr); 33c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-test_submat", &TestSubMat);CHKERRQ(ierr); 34c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-test_allcols", &TestAllcols);CHKERRQ(ierr); 35c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-test_sorted",&test_sorted,NULL);CHKERRQ(ierr); 36c4762a1bSJed Brown 37c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = MatSetSizes(A,mbs*bs,mbs*bs,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 39c4762a1bSJed Brown ierr = MatSetType(A,MATBAIJ);CHKERRQ(ierr); 40c4762a1bSJed Brown ierr = MatSeqBAIJSetPreallocation(A,bs,PETSC_DEFAULT,NULL);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = MatMPIBAIJSetPreallocation(A,bs,PETSC_DEFAULT,NULL,PETSC_DEFAULT,NULL);CHKERRQ(ierr); 42c4762a1bSJed Brown 43c4762a1bSJed Brown ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr); 44c4762a1bSJed Brown ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr); 45c4762a1bSJed Brown 46c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr); 48c4762a1bSJed Brown Mbs = M/bs; 49c4762a1bSJed Brown 50c4762a1bSJed Brown ierr = PetscMalloc1(bs,&rows);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = PetscMalloc1(bs,&cols);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = PetscMalloc1(bs*bs,&vals);CHKERRQ(ierr); 53c4762a1bSJed Brown ierr = PetscMalloc1(M,&idx);CHKERRQ(ierr); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Now set blocks of values */ 56c4762a1bSJed Brown for (j=0; j<bs*bs; j++) vals[j] = 0.0; 57c4762a1bSJed Brown for (i=0; i<Mbs; i++) { 58c4762a1bSJed Brown cols[0] = i*bs; rows[0] = i*bs; 59c4762a1bSJed Brown for (j=1; j<bs; j++) { 60c4762a1bSJed Brown rows[j] = rows[j-1]+1; 61c4762a1bSJed Brown cols[j] = cols[j-1]+1; 62c4762a1bSJed Brown } 63c4762a1bSJed Brown ierr = MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES);CHKERRQ(ierr); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown /* second, add random blocks */ 66c4762a1bSJed Brown for (i=0; i<20*bs; i++) { 67c4762a1bSJed Brown ierr = PetscRandomGetValue(rand,&rval);CHKERRQ(ierr); 68c4762a1bSJed Brown cols[0] = bs*(PetscInt)(PetscRealPart(rval)*Mbs); 69c4762a1bSJed Brown ierr = PetscRandomGetValue(rand,&rval);CHKERRQ(ierr); 70c4762a1bSJed Brown rows[0] = rstart + bs*(PetscInt)(PetscRealPart(rval)*mbs); 71c4762a1bSJed Brown for (j=1; j<bs; j++) { 72c4762a1bSJed Brown rows[j] = rows[j-1]+1; 73c4762a1bSJed Brown cols[j] = cols[j-1]+1; 74c4762a1bSJed Brown } 75c4762a1bSJed Brown 76c4762a1bSJed Brown for (j=0; j<bs*bs; j++) { 77c4762a1bSJed Brown ierr = PetscRandomGetValue(rand,&rval);CHKERRQ(ierr); 78c4762a1bSJed Brown vals[j] = rval; 79c4762a1bSJed Brown } 80c4762a1bSJed Brown ierr = MatSetValues(A,bs,rows,bs,cols,vals,ADD_VALUES);CHKERRQ(ierr); 81c4762a1bSJed Brown } 82c4762a1bSJed Brown 83c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 84c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 85c4762a1bSJed Brown 86c4762a1bSJed Brown /* make A a symmetric matrix: A <- A^T + A */ 87c4762a1bSJed Brown ierr = MatTranspose(A,MAT_INITIAL_MATRIX, &Atrans);CHKERRQ(ierr); 88c4762a1bSJed Brown ierr = MatAXPY(A,one,Atrans,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 89c4762a1bSJed Brown ierr = MatDestroy(&Atrans);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = MatTranspose(A,MAT_INITIAL_MATRIX, &Atrans);CHKERRQ(ierr); 91c4762a1bSJed Brown ierr = MatEqual(A, Atrans, &flg);CHKERRQ(ierr); 92c4762a1bSJed Brown if (flg) { 93c4762a1bSJed Brown ierr = MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 94c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"A+A^T is non-symmetric"); 95c4762a1bSJed Brown ierr = MatDestroy(&Atrans);CHKERRQ(ierr); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* create a SeqSBAIJ matrix sA (= A) */ 98c4762a1bSJed Brown ierr = MatConvert(A,MATSBAIJ,MAT_INITIAL_MATRIX,&sA);CHKERRQ(ierr); 99c4762a1bSJed Brown if (vid >= 0 && vid < size) { 100c4762a1bSJed Brown if (!rank) printf("A: \n"); 101c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 102c4762a1bSJed Brown if (!rank) printf("sA: \n"); 103c4762a1bSJed Brown ierr = MatView(sA,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106c4762a1bSJed Brown /* Test sA==A through MatMult() */ 107c4762a1bSJed Brown ierr = MatMultEqual(A,sA,10,&flg);CHKERRQ(ierr); 108c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Error in MatConvert(): A != sA"); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* Test MatIncreaseOverlap() */ 111c4762a1bSJed Brown ierr = PetscMalloc1(nd,&is1);CHKERRQ(ierr); 112c4762a1bSJed Brown ierr = PetscMalloc1(nd,&is2);CHKERRQ(ierr); 113c4762a1bSJed Brown 114c4762a1bSJed Brown for (i=0; i<nd; i++) { 115c4762a1bSJed Brown if (!TestAllcols) { 116c4762a1bSJed Brown ierr = PetscRandomGetValue(rand,&rval);CHKERRQ(ierr); 117c4762a1bSJed Brown sz = (PetscInt)((0.5+0.2*PetscRealPart(rval))*mbs); /* 0.5*mbs < sz < 0.7*mbs */ 118c4762a1bSJed Brown 119c4762a1bSJed Brown for (j=0; j<sz; j++) { 120c4762a1bSJed Brown ierr = PetscRandomGetValue(rand,&rval);CHKERRQ(ierr); 121c4762a1bSJed Brown idx[j*bs] = bs*(PetscInt)(PetscRealPart(rval)*Mbs); 122c4762a1bSJed Brown for (k=1; k<bs; k++) idx[j*bs+k] = idx[j*bs]+k; 123c4762a1bSJed Brown } 124c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,sz*bs,idx,PETSC_COPY_VALUES,is1+i);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,sz*bs,idx,PETSC_COPY_VALUES,is2+i);CHKERRQ(ierr); 126c4762a1bSJed Brown if (rank == vid) { 127c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," [%d] IS sz[%d]: %d\n",rank,i,sz);CHKERRQ(ierr); 128c4762a1bSJed Brown ierr = ISView(is2[i],PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 129c4762a1bSJed Brown } 130*a5b23f4aSJose E. Roman } else { /* Test all rows and columns */ 131c4762a1bSJed Brown sz = M; 132c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,sz,0,1,is1+i);CHKERRQ(ierr); 133c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,sz,0,1,is2+i);CHKERRQ(ierr); 134c4762a1bSJed Brown 135c4762a1bSJed Brown if (rank == vid) { 136c4762a1bSJed Brown PetscBool colflag; 137c4762a1bSJed Brown ierr = ISIdentity(is2[i],&colflag);CHKERRQ(ierr); 138c4762a1bSJed Brown printf("[%d] is2[%d], colflag %d\n",rank,(int)i,(int)colflag); 139c4762a1bSJed Brown ierr = ISView(is2[i],PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown } 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144c4762a1bSJed Brown ierr = PetscLogStageRegister("MatOv_SBAIJ",&stages[0]);CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = PetscLogStageRegister("MatOv_BAIJ",&stages[1]);CHKERRQ(ierr); 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* Test MatIncreaseOverlap */ 148c4762a1bSJed Brown if (TestOverlap) { 149c4762a1bSJed Brown ierr = PetscLogStagePush(stages[0]);CHKERRQ(ierr); 150c4762a1bSJed Brown ierr = MatIncreaseOverlap(sA,nd,is2,ov);CHKERRQ(ierr); 151c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 152c4762a1bSJed Brown 153c4762a1bSJed Brown ierr = PetscLogStagePush(stages[1]);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = MatIncreaseOverlap(A,nd,is1,ov);CHKERRQ(ierr); 155c4762a1bSJed Brown ierr = PetscLogStagePop();CHKERRQ(ierr); 156c4762a1bSJed Brown 157c4762a1bSJed Brown if (rank == vid) { 158c4762a1bSJed Brown printf("\n[%d] IS from BAIJ:\n",rank); 159c4762a1bSJed Brown ierr = ISView(is1[0],PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 160c4762a1bSJed Brown printf("\n[%d] IS from SBAIJ:\n",rank); 161c4762a1bSJed Brown ierr = ISView(is2[0],PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 162c4762a1bSJed Brown } 163c4762a1bSJed Brown 164c4762a1bSJed Brown for (i=0; i<nd; ++i) { 165c4762a1bSJed Brown ierr = ISEqual(is1[i],is2[i],&flg);CHKERRQ(ierr); 166c4762a1bSJed Brown if (!flg) { 167c4762a1bSJed Brown if (!rank) { 168c4762a1bSJed Brown ierr = ISSort(is1[i]);CHKERRQ(ierr); 169c4762a1bSJed Brown /* ISView(is1[i],PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); */ 170c4762a1bSJed Brown ierr = ISSort(is2[i]);CHKERRQ(ierr); 171c4762a1bSJed Brown /* ISView(is2[i],PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); */ 172c4762a1bSJed Brown } 173d8185827SBarry Smith SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"i=%D, is1 != is2",i); 174c4762a1bSJed Brown } 175c4762a1bSJed Brown } 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* Test MatCreateSubmatrices */ 179c4762a1bSJed Brown if (TestSubMat) { 180c4762a1bSJed Brown if (test_sorted) { 181c4762a1bSJed Brown for (i = 0; i < nd; ++i) { 182c4762a1bSJed Brown ierr = ISSort(is1[i]);CHKERRQ(ierr); 183c4762a1bSJed Brown } 184c4762a1bSJed Brown } 185c4762a1bSJed Brown ierr = MatCreateSubMatrices(A,nd,is1,is1,MAT_INITIAL_MATRIX,&submatA);CHKERRQ(ierr); 186c4762a1bSJed Brown ierr = MatCreateSubMatrices(sA,nd,is1,is1,MAT_INITIAL_MATRIX,&submatsA);CHKERRQ(ierr); 187c4762a1bSJed Brown 188c4762a1bSJed Brown ierr = MatMultEqual(A,sA,10,&flg);CHKERRQ(ierr); 189c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"A != sA"); 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* Now test MatCreateSubmatrices with MAT_REUSE_MATRIX option */ 192c4762a1bSJed Brown ierr = MatCreateSubMatrices(A,nd,is1,is1,MAT_REUSE_MATRIX,&submatA);CHKERRQ(ierr); 193c4762a1bSJed Brown ierr = MatCreateSubMatrices(sA,nd,is1,is1,MAT_REUSE_MATRIX,&submatsA);CHKERRQ(ierr); 194c4762a1bSJed Brown ierr = MatMultEqual(A,sA,10,&flg);CHKERRQ(ierr); 195c4762a1bSJed Brown if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"MatCreateSubmatrices(): A != sA"); 196c4762a1bSJed Brown 197c4762a1bSJed Brown ierr = MatDestroySubMatrices(nd,&submatA);CHKERRQ(ierr); 198c4762a1bSJed Brown ierr = MatDestroySubMatrices(nd,&submatsA);CHKERRQ(ierr); 199c4762a1bSJed Brown } 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* Free allocated memory */ 202c4762a1bSJed Brown for (i=0; i<nd; ++i) { 203c4762a1bSJed Brown ierr = ISDestroy(&is1[i]);CHKERRQ(ierr); 204c4762a1bSJed Brown ierr = ISDestroy(&is2[i]);CHKERRQ(ierr); 205c4762a1bSJed Brown } 206c4762a1bSJed Brown ierr = PetscFree(is1);CHKERRQ(ierr); 207c4762a1bSJed Brown ierr = PetscFree(is2);CHKERRQ(ierr); 208c4762a1bSJed Brown ierr = PetscFree(idx);CHKERRQ(ierr); 209c4762a1bSJed Brown ierr = PetscFree(rows);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = PetscFree(cols);CHKERRQ(ierr); 211c4762a1bSJed Brown ierr = PetscFree(vals);CHKERRQ(ierr); 212c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 213c4762a1bSJed Brown ierr = MatDestroy(&sA);CHKERRQ(ierr); 214c4762a1bSJed Brown ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = PetscFinalize(); 216c4762a1bSJed Brown return ierr; 217c4762a1bSJed Brown } 218c4762a1bSJed Brown 219c4762a1bSJed Brown /*TEST 220c4762a1bSJed Brown 221c4762a1bSJed Brown test: 222c4762a1bSJed Brown args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_submat 223c4762a1bSJed Brown output_file: output/ex92_1.out 224c4762a1bSJed Brown 225c4762a1bSJed Brown test: 226c4762a1bSJed Brown suffix: 2 227c4762a1bSJed Brown nsize: {{3 4}} 228c4762a1bSJed Brown args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_submat 229c4762a1bSJed Brown output_file: output/ex92_1.out 230c4762a1bSJed Brown 231c4762a1bSJed Brown test: 232c4762a1bSJed Brown suffix: 3 233c4762a1bSJed Brown nsize: {{3 4}} 234c4762a1bSJed Brown args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_allcols 235c4762a1bSJed Brown output_file: output/ex92_1.out 236c4762a1bSJed Brown 237c4762a1bSJed Brown test: 238c4762a1bSJed Brown suffix: 3_sorted 239c4762a1bSJed Brown nsize: {{3 4}} 240c4762a1bSJed Brown args: -ov {{1 3}} -mat_block_size {{2 8}} -test_overlap -test_allcols -test_sorted 241c4762a1bSJed Brown output_file: output/ex92_1.out 242c4762a1bSJed Brown 243c4762a1bSJed Brown test: 244c4762a1bSJed Brown suffix: 4 245c4762a1bSJed Brown nsize: {{3 4}} 246c4762a1bSJed Brown args: -ov {{1 3}} -mat_block_size {{2 8}} -test_submat -test_allcols 247c4762a1bSJed Brown output_file: output/ex92_1.out 248c4762a1bSJed Brown 249c4762a1bSJed Brown TEST*/ 250