1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests converting a matrix to another format with MatConvert().\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown /* Usage: mpiexec -n <np> ex55 -verbose <0 or 1> */ 6c4762a1bSJed Brown 7c4762a1bSJed Brown int main(int argc,char **args) 8c4762a1bSJed Brown { 9c4762a1bSJed Brown Mat C,A,B,D; 10c4762a1bSJed Brown PetscInt i,j,ntypes,bs,mbs,m,block,d_nz=6, o_nz=3,col[3],row,verbose=0; 11c4762a1bSJed Brown PetscMPIInt size,rank; 12c4762a1bSJed Brown MatType type[9]; 13c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 14c4762a1bSJed Brown PetscViewer fd; 15c4762a1bSJed Brown PetscBool equal,flg_loadmat,flg,issymmetric; 16c4762a1bSJed Brown PetscScalar value[3]; 17c4762a1bSJed Brown 189566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-verbose",&verbose,NULL)); 209566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg_loadmat)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 23c4762a1bSJed Brown 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-testseqaij",&flg)); 25c4762a1bSJed Brown if (flg) { 26c4762a1bSJed Brown if (size == 1) { 27c4762a1bSJed Brown type[0] = MATSEQAIJ; 28c4762a1bSJed Brown } else { 29c4762a1bSJed Brown type[0] = MATMPIAIJ; 30c4762a1bSJed Brown } 31c4762a1bSJed Brown } else { 32c4762a1bSJed Brown type[0] = MATAIJ; 33c4762a1bSJed Brown } 34c4762a1bSJed Brown if (size == 1) { 35c4762a1bSJed Brown ntypes = 3; 36c4762a1bSJed Brown type[1] = MATSEQBAIJ; 37c4762a1bSJed Brown type[2] = MATSEQSBAIJ; 38c4762a1bSJed Brown } else { 39c4762a1bSJed Brown ntypes = 3; 40c4762a1bSJed Brown type[1] = MATMPIBAIJ; 41c4762a1bSJed Brown type[2] = MATMPISBAIJ; 42c4762a1bSJed Brown } 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* input matrix C */ 45c4762a1bSJed Brown if (flg_loadmat) { 46c4762a1bSJed Brown /* Open binary file. */ 479566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* Load a baij matrix, then destroy the viewer. */ 509566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 51c4762a1bSJed Brown if (size == 1) { 529566063dSJacob Faibussowitsch PetscCall(MatSetType(C,MATSEQBAIJ)); 53c4762a1bSJed Brown } else { 549566063dSJacob Faibussowitsch PetscCall(MatSetType(C,MATMPIBAIJ)); 55c4762a1bSJed Brown } 569566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 579566063dSJacob Faibussowitsch PetscCall(MatLoad(C,fd)); 589566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&fd)); 59c4762a1bSJed Brown } else { /* Create a baij mat with bs>1 */ 60c4762a1bSJed Brown bs = 2; mbs=8; 619566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL)); 629566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL)); 63*08401ef6SPierre Jolivet PetscCheck(bs > 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG," bs must be >1 in this case"); 64c4762a1bSJed Brown m = mbs*bs; 659566063dSJacob Faibussowitsch PetscCall(MatCreateBAIJ(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,m,m,d_nz,NULL,o_nz,NULL,&C)); 66c4762a1bSJed Brown for (block=0; block<mbs; block++) { 67c4762a1bSJed Brown /* diagonal blocks */ 68c4762a1bSJed Brown value[0] = -1.0; value[1] = 4.0; value[2] = -1.0; 69c4762a1bSJed Brown for (i=1+block*bs; i<bs-1+block*bs; i++) { 70c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 719566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&i,3,col,value,INSERT_VALUES)); 72c4762a1bSJed Brown } 73c4762a1bSJed Brown i = bs - 1+block*bs; col[0] = bs - 2+block*bs; col[1] = bs - 1+block*bs; 74c4762a1bSJed Brown value[0]=-1.0; value[1]=4.0; 759566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&i,2,col,value,INSERT_VALUES)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown i = 0+block*bs; col[0] = 0+block*bs; col[1] = 1+block*bs; 78c4762a1bSJed Brown value[0]=4.0; value[1] = -1.0; 799566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&i,2,col,value,INSERT_VALUES)); 80c4762a1bSJed Brown } 81c4762a1bSJed Brown /* off-diagonal blocks */ 82c4762a1bSJed Brown value[0]=-1.0; 83c4762a1bSJed Brown for (i=0; i<(mbs-1)*bs; i++) { 84c4762a1bSJed Brown col[0]=i+bs; 859566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&i,1,col,value,INSERT_VALUES)); 86c4762a1bSJed Brown col[0]=i; row=i+bs; 879566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&row,1,col,value,INSERT_VALUES)); 88c4762a1bSJed Brown } 899566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 909566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 91c4762a1bSJed Brown } 92c4762a1bSJed Brown { 93c4762a1bSJed Brown /* Check the symmetry of C because it will be converted to a sbaij matrix */ 94c4762a1bSJed Brown Mat Ctrans; 959566063dSJacob Faibussowitsch PetscCall(MatTranspose(C,MAT_INITIAL_MATRIX,&Ctrans)); 969566063dSJacob Faibussowitsch PetscCall(MatEqual(C,Ctrans,&flg)); 97c4762a1bSJed Brown /* 98c4762a1bSJed Brown { 999566063dSJacob Faibussowitsch PetscCall(MatAXPY(C,1.,Ctrans,DIFFERENT_NONZERO_PATTERN)); 100c4762a1bSJed Brown flg = PETSC_TRUE; 101c4762a1bSJed Brown } 102c4762a1bSJed Brown */ 1039566063dSJacob Faibussowitsch PetscCall(MatSetOption(C,MAT_SYMMETRIC,flg)); 1049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Ctrans)); 105c4762a1bSJed Brown } 1069566063dSJacob Faibussowitsch PetscCall(MatIsSymmetric(C,0.0,&issymmetric)); 1079566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(C,NULL,"-view_mat")); 108c4762a1bSJed Brown 109c4762a1bSJed Brown /* convert C to other formats */ 110c4762a1bSJed Brown for (i=0; i<ntypes; i++) { 111c4762a1bSJed Brown PetscBool ismpisbaij,isseqsbaij; 112c4762a1bSJed Brown 1139566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(type[i],MATMPISBAIJ,&ismpisbaij)); 1149566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(type[i],MATMPISBAIJ,&isseqsbaij)); 115c4762a1bSJed Brown if (!issymmetric && (ismpisbaij || isseqsbaij)) continue; 1169566063dSJacob Faibussowitsch PetscCall(MatConvert(C,type[i],MAT_INITIAL_MATRIX,&A)); 1179566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A,C,10,&equal)); 11828b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in conversion from BAIJ to %s",type[i]); 119c4762a1bSJed Brown for (j=i+1; j<ntypes; j++) { 1209566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(type[j],MATMPISBAIJ,&ismpisbaij)); 1219566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(type[j],MATMPISBAIJ,&isseqsbaij)); 122c4762a1bSJed Brown if (!issymmetric && (ismpisbaij || isseqsbaij)) continue; 123c4762a1bSJed Brown if (verbose>0) { 1249566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," \n[%d] test conversion between %s and %s\n",rank,type[i],type[j])); 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127dd400576SPatrick Sanan if (rank == 0 && verbose) printf("Convert %s A to %s B\n",type[i],type[j]); 1289566063dSJacob Faibussowitsch PetscCall(MatConvert(A,type[j],MAT_INITIAL_MATRIX,&B)); 129c4762a1bSJed Brown /* 130c4762a1bSJed Brown if (j == 2) { 1319566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," A: %s\n",type[i])); 1329566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 1339566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," B: %s\n",type[j])); 1349566063dSJacob Faibussowitsch PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 135c4762a1bSJed Brown } 136c4762a1bSJed Brown */ 1379566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A,B,10,&equal)); 13828b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in conversion from %s to %s",type[i],type[j]); 139c4762a1bSJed Brown 140c4762a1bSJed Brown if (size == 1 || j != 2) { /* Matconvert from mpisbaij mat to other formats are not supported */ 141dd400576SPatrick Sanan if (rank == 0 && verbose) printf("Convert %s B to %s D\n",type[j],type[i]); 1429566063dSJacob Faibussowitsch PetscCall(MatConvert(B,type[i],MAT_INITIAL_MATRIX,&D)); 1439566063dSJacob Faibussowitsch PetscCall(MatMultEqual(B,D,10,&equal)); 14428b400f6SJacob Faibussowitsch PetscCheck(equal,PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Error in conversion from %s to %s",type[j],type[i]); 145c4762a1bSJed Brown 1469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 147c4762a1bSJed Brown } 1489566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&D)); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* Test in-place convert */ 153c4762a1bSJed Brown if (size == 1) { /* size > 1 is not working yet! */ 154c4762a1bSJed Brown j = (i+1)%ntypes; 1559566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(type[j],MATMPISBAIJ,&ismpisbaij)); 1569566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(type[j],MATMPISBAIJ,&isseqsbaij)); 157c4762a1bSJed Brown if (!issymmetric && (ismpisbaij || isseqsbaij)) continue; 158c4762a1bSJed Brown /* printf("[%d] i: %d, j: %d\n",rank,i,j); */ 1599566063dSJacob Faibussowitsch PetscCall(MatConvert(A,type[j],MAT_INPLACE_MATRIX,&A)); 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 1629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 163c4762a1bSJed Brown } 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* test BAIJ to MATIS */ 166c4762a1bSJed Brown if (size > 1) { 167c4762a1bSJed Brown MatType ctype; 168c4762a1bSJed Brown 1699566063dSJacob Faibussowitsch PetscCall(MatGetType(C,&ctype)); 1709566063dSJacob Faibussowitsch PetscCall(MatConvert(C,MATIS,MAT_INITIAL_MATRIX,&A)); 1719566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A,C,10,&equal)); 1729566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A,NULL,"-view_conv")); 173c4762a1bSJed Brown if (!equal) { 174c4762a1bSJed Brown Mat C2; 175c4762a1bSJed Brown 1769566063dSJacob Faibussowitsch PetscCall(MatConvert(A,ctype,MAT_INITIAL_MATRIX,&C2)); 1779566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(C2,NULL,"-view_conv_assembled")); 1789566063dSJacob Faibussowitsch PetscCall(MatAXPY(C2,-1.,C,DIFFERENT_NONZERO_PATTERN)); 1799566063dSJacob Faibussowitsch PetscCall(MatChop(C2,PETSC_SMALL)); 1809566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(C2,NULL,"-view_err")); 1819566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C2)); 182c4762a1bSJed Brown SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in conversion from BAIJ to MATIS"); 183c4762a1bSJed Brown } 1849566063dSJacob Faibussowitsch PetscCall(MatConvert(C,MATIS,MAT_REUSE_MATRIX,&A)); 1859566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A,C,10,&equal)); 1869566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A,NULL,"-view_conv")); 187c4762a1bSJed Brown if (!equal) { 188c4762a1bSJed Brown Mat C2; 189c4762a1bSJed Brown 1909566063dSJacob Faibussowitsch PetscCall(MatConvert(A,ctype,MAT_INITIAL_MATRIX,&C2)); 1919566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(C2,NULL,"-view_conv_assembled")); 1929566063dSJacob Faibussowitsch PetscCall(MatAXPY(C2,-1.,C,DIFFERENT_NONZERO_PATTERN)); 1939566063dSJacob Faibussowitsch PetscCall(MatChop(C2,PETSC_SMALL)); 1949566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(C2,NULL,"-view_err")); 1959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C2)); 196c4762a1bSJed Brown SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in conversion from BAIJ to MATIS"); 197c4762a1bSJed Brown } 1989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1999566063dSJacob Faibussowitsch PetscCall(MatDuplicate(C,MAT_COPY_VALUES,&A)); 2009566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATIS,MAT_INPLACE_MATRIX,&A)); 2019566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A,NULL,"-view_conv")); 2029566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A,C,10,&equal)); 203c4762a1bSJed Brown if (!equal) { 204c4762a1bSJed Brown Mat C2; 205c4762a1bSJed Brown 2069566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A,NULL,"-view_conv")); 2079566063dSJacob Faibussowitsch PetscCall(MatConvert(A,ctype,MAT_INITIAL_MATRIX,&C2)); 2089566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(C2,NULL,"-view_conv_assembled")); 2099566063dSJacob Faibussowitsch PetscCall(MatAXPY(C2,-1.,C,DIFFERENT_NONZERO_PATTERN)); 2109566063dSJacob Faibussowitsch PetscCall(MatChop(C2,PETSC_SMALL)); 2119566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(C2,NULL,"-view_err")); 2129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C2)); 213c4762a1bSJed Brown SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Error in conversion from BAIJ to MATIS"); 214c4762a1bSJed Brown } 2159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 216c4762a1bSJed Brown } 2179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 218c4762a1bSJed Brown 2199566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 220b122ec5aSJacob Faibussowitsch return 0; 221c4762a1bSJed Brown } 222c4762a1bSJed Brown 223c4762a1bSJed Brown /*TEST 224c4762a1bSJed Brown 225c4762a1bSJed Brown test: 226c4762a1bSJed Brown 227c4762a1bSJed Brown test: 228c4762a1bSJed Brown suffix: 2 229c4762a1bSJed Brown nsize: 3 230c4762a1bSJed Brown 231c4762a1bSJed Brown testset: 232c4762a1bSJed Brown requires: parmetis 233c4762a1bSJed Brown output_file: output/ex55_1.out 234c4762a1bSJed Brown nsize: 3 235c4762a1bSJed Brown args: -mat_is_disassemble_l2g_type nd -mat_partitioning_type parmetis 236c4762a1bSJed Brown test: 237c4762a1bSJed Brown suffix: matis_baij_parmetis_nd 238c4762a1bSJed Brown test: 239c4762a1bSJed Brown suffix: matis_aij_parmetis_nd 240c4762a1bSJed Brown args: -testseqaij 241c4762a1bSJed Brown test: 242dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 243c4762a1bSJed Brown suffix: matis_poisson1_parmetis_nd 244c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/poisson1 245c4762a1bSJed Brown 246c4762a1bSJed Brown testset: 247dfd57a17SPierre Jolivet requires: ptscotch defined(PETSC_HAVE_SCOTCH_PARMETIS_V3_NODEND) 248c4762a1bSJed Brown output_file: output/ex55_1.out 249c4762a1bSJed Brown nsize: 4 250c4762a1bSJed Brown args: -mat_is_disassemble_l2g_type nd -mat_partitioning_type ptscotch 251c4762a1bSJed Brown test: 252c4762a1bSJed Brown suffix: matis_baij_ptscotch_nd 253c4762a1bSJed Brown test: 254c4762a1bSJed Brown suffix: matis_aij_ptscotch_nd 255c4762a1bSJed Brown args: -testseqaij 256c4762a1bSJed Brown test: 257dfd57a17SPierre Jolivet requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES) 258c4762a1bSJed Brown suffix: matis_poisson1_ptscotch_nd 259c4762a1bSJed Brown args: -f ${DATAFILESPATH}/matrices/poisson1 260c4762a1bSJed Brown 261c4762a1bSJed Brown TEST*/ 262