1c4762a1bSJed Brown 2*d0dbe9f7SStefano Zampini static char help[] = "Tests the use of interface functions for MATIS matrices and conversion routines.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown PetscErrorCode TestMatZeroRows(Mat,Mat,PetscBool,IS,PetscScalar); 7c4762a1bSJed Brown PetscErrorCode CheckMat(Mat,Mat,PetscBool,const char*); 8c4762a1bSJed Brown 9c4762a1bSJed Brown int main(int argc,char **args) 10c4762a1bSJed Brown { 11c4762a1bSJed Brown Mat A,B,A2,B2,T; 12c4762a1bSJed Brown Mat Aee,Aeo,Aoe,Aoo; 13*d0dbe9f7SStefano Zampini Mat *mats,*Asub,*Bsub; 14c4762a1bSJed Brown Vec x,y; 15c4762a1bSJed Brown MatInfo info; 16c4762a1bSJed Brown ISLocalToGlobalMapping cmap,rmap; 17c4762a1bSJed Brown IS is,is2,reven,rodd,ceven,codd; 18c4762a1bSJed Brown IS *rows,*cols; 19*d0dbe9f7SStefano Zampini IS irow[2], icol[2]; 20*d0dbe9f7SStefano Zampini PetscLayout rlayout, clayout; 21*d0dbe9f7SStefano Zampini const PetscInt *rrange, *crange; 22c4762a1bSJed Brown MatType lmtype; 23c4762a1bSJed Brown PetscScalar diag = 2.; 24e432b41dSStefano Zampini PetscInt n,m,i,lm,ln; 25c4762a1bSJed Brown PetscInt rst,ren,cst,cen,nr,nc; 26*d0dbe9f7SStefano Zampini PetscMPIInt rank,size,lrank,rrank; 27c4762a1bSJed Brown PetscBool testT,squaretest,isaij; 28e432b41dSStefano Zampini PetscBool permute = PETSC_FALSE, negmap = PETSC_FALSE, repmap = PETSC_FALSE; 29e432b41dSStefano Zampini PetscBool diffmap = PETSC_TRUE, symmetric = PETSC_FALSE, issymmetric; 30c4762a1bSJed Brown 319566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 34c4762a1bSJed Brown m = n = 2*size; 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-symmetric",&symmetric,NULL)); 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 379566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 389566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-negmap",&negmap,NULL)); 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-repmap",&repmap,NULL)); 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-permmap",&permute,NULL)); 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-diffmap",&diffmap,NULL)); 42e00437b9SBarry Smith PetscCheck(size == 1 || m >= 4,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Number of rows should be larger or equal 4 for parallel runs"); 43e00437b9SBarry Smith PetscCheck(size != 1 || m >= 2,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Number of rows should be larger or equal 2 for uniprocessor runs"); 44e00437b9SBarry Smith PetscCheck(n >= 2,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Number of cols should be larger or equal 2"); 45c4762a1bSJed Brown 46c4762a1bSJed Brown /* create a MATIS matrix */ 479566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 489566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n)); 499566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATIS)); 509566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 51e432b41dSStefano Zampini if (!negmap && !repmap) { 52fc989267SStefano Zampini /* This is not the proper setting for MATIS for finite elements, it is just used to test the routines 53fc989267SStefano Zampini Here we use a one-to-one correspondence between local row/column spaces and global row/column spaces 54fc989267SStefano Zampini Equivalent to passing NULL for the mapping */ 559566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n,0,1,&is)); 56e432b41dSStefano Zampini } else if (negmap && !repmap) { /* non repeated but with negative indices */ 579566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n+2,-2,1,&is)); 58e432b41dSStefano Zampini } else if (!negmap && repmap) { /* non negative but repeated indices */ 59e432b41dSStefano Zampini IS isl[2]; 60e432b41dSStefano Zampini 619566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n,0,1,&isl[0])); 629566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n,n-1,-1,&isl[1])); 639566063dSJacob Faibussowitsch PetscCall(ISConcatenate(PETSC_COMM_WORLD,2,isl,&is)); 649566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl[0])); 659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl[1])); 66e432b41dSStefano Zampini } else { /* negative and repeated indices */ 67e432b41dSStefano Zampini IS isl[2]; 68e432b41dSStefano Zampini 699566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n+1,-1,1,&isl[0])); 709566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n+1,n-1,-1,&isl[1])); 719566063dSJacob Faibussowitsch PetscCall(ISConcatenate(PETSC_COMM_WORLD,2,isl,&is)); 729566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl[0])); 739566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl[1])); 74e432b41dSStefano Zampini } 759566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(is,&cmap)); 769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 77c4762a1bSJed Brown 78e432b41dSStefano Zampini if (m != n || diffmap) { 799566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,m,permute ? m-1 : 0,permute ? -1 : 1,&is)); 809566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(is,&rmap)); 819566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 82c4762a1bSJed Brown } else { 839566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)cmap)); 84c4762a1bSJed Brown rmap = cmap; 85c4762a1bSJed Brown } 86e432b41dSStefano Zampini 879566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(A,rmap,cmap)); 889566063dSJacob Faibussowitsch PetscCall(MatISStoreL2L(A,PETSC_FALSE)); 899566063dSJacob Faibussowitsch PetscCall(MatISSetPreallocation(A,3,NULL,3,NULL)); 909566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,(PetscBool)!(repmap || negmap))); /* I do not want to precompute the pattern */ 919566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(rmap,&lm)); 929566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(cmap,&ln)); 93e432b41dSStefano Zampini for (i=0; i<lm; i++) { 94c4762a1bSJed Brown PetscScalar v[3]; 95c4762a1bSJed Brown PetscInt cols[3]; 96c4762a1bSJed Brown 97c4762a1bSJed Brown cols[0] = (i-1+n)%n; 98c4762a1bSJed Brown cols[1] = i%n; 99c4762a1bSJed Brown cols[2] = (i+1)%n; 100c4762a1bSJed Brown v[0] = -1.*(symmetric ? PetscMin(i+1,cols[0]+1) : i+1); 101c4762a1bSJed Brown v[1] = 2.*(symmetric ? PetscMin(i+1,cols[1]+1) : i+1); 102c4762a1bSJed Brown v[2] = -1.*(symmetric ? PetscMin(i+1,cols[2]+1) : i+1); 1039566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(cmap,IS_GTOLM_MASK,3,cols,NULL,cols)); 1049566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A,1,&i,3,cols,v,ADD_VALUES)); 105c4762a1bSJed Brown } 1069566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1079566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 108c4762a1bSJed Brown 109e432b41dSStefano Zampini /* activate tests for square matrices with same maps only */ 1109566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(A,&squaretest)); 111e432b41dSStefano Zampini if (squaretest && rmap != cmap) { 112e432b41dSStefano Zampini PetscInt nr, nc; 113e432b41dSStefano Zampini 1149566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(rmap,&nr)); 1159566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(cmap,&nc)); 116e432b41dSStefano Zampini if (nr != nc) squaretest = PETSC_FALSE; 117e432b41dSStefano Zampini else { 118e432b41dSStefano Zampini const PetscInt *idxs1,*idxs2; 119e432b41dSStefano Zampini 1209566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(rmap,&idxs1)); 1219566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(cmap,&idxs2)); 1229566063dSJacob Faibussowitsch PetscCall(PetscArraycmp(idxs1,idxs2,nr,&squaretest)); 1239566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(rmap,&idxs1)); 1249566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(cmap,&idxs2)); 125e432b41dSStefano Zampini } 1261c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE,&squaretest,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A))); 127e432b41dSStefano Zampini } 128e432b41dSStefano Zampini 129e432b41dSStefano Zampini /* test MatISGetLocalMat */ 1309566063dSJacob Faibussowitsch PetscCall(MatISGetLocalMat(A,&B)); 1319566063dSJacob Faibussowitsch PetscCall(MatGetType(B,&lmtype)); 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* test MatGetInfo */ 1349566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatGetInfo\n")); 1359566063dSJacob Faibussowitsch PetscCall(MatGetInfo(A,MAT_LOCAL,&info)); 1369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 137d0609cedSBarry Smith PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"Process %2d: %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",PetscGlobalRank,(PetscInt)info.nz_used, 138d0609cedSBarry Smith (PetscInt)info.nz_allocated,(PetscInt)info.nz_unneeded,(PetscInt)info.assemblies,(PetscInt)info.mallocs)); 1399566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 1409566063dSJacob Faibussowitsch PetscCall(MatGetInfo(A,MAT_GLOBAL_MAX,&info)); 141d0609cedSBarry Smith PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"GlobalMax : %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",(PetscInt)info.nz_used, 142d0609cedSBarry Smith (PetscInt)info.nz_allocated,(PetscInt)info.nz_unneeded,(PetscInt)info.assemblies,(PetscInt)info.mallocs)); 1439566063dSJacob Faibussowitsch PetscCall(MatGetInfo(A,MAT_GLOBAL_SUM,&info)); 144d0609cedSBarry Smith PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"GlobalSum : %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",(PetscInt)info.nz_used, 145d0609cedSBarry Smith (PetscInt)info.nz_allocated,(PetscInt)info.nz_unneeded,(PetscInt)info.assemblies,(PetscInt)info.mallocs)); 146c4762a1bSJed Brown 147e432b41dSStefano Zampini /* test MatIsSymmetric */ 1489566063dSJacob Faibussowitsch PetscCall(MatIsSymmetric(A,0.0,&issymmetric)); 1499566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatIsSymmetric: %d\n",issymmetric)); 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* Create a MPIAIJ matrix, same as A */ 1529566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); 1539566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n)); 1549566063dSJacob Faibussowitsch PetscCall(MatSetType(B,MATAIJ)); 1559566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 1569566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 1579566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(B,rmap,cmap)); 1589566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(B,3,NULL,3,NULL)); 1599566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(B,1,3,NULL,3,NULL)); 160c4762a1bSJed Brown #if defined(PETSC_HAVE_HYPRE) 1619566063dSJacob Faibussowitsch PetscCall(MatHYPRESetPreallocation(B,3,NULL,3,NULL)); 162c4762a1bSJed Brown #endif 1639566063dSJacob Faibussowitsch PetscCall(MatISSetPreallocation(B,3,NULL,3,NULL)); 1649566063dSJacob Faibussowitsch PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,(PetscBool)!(repmap || negmap))); /* I do not want to precompute the pattern */ 165e432b41dSStefano Zampini for (i=0; i<lm; i++) { 166c4762a1bSJed Brown PetscScalar v[3]; 167c4762a1bSJed Brown PetscInt cols[3]; 168c4762a1bSJed Brown 169c4762a1bSJed Brown cols[0] = (i-1+n)%n; 170c4762a1bSJed Brown cols[1] = i%n; 171c4762a1bSJed Brown cols[2] = (i+1)%n; 172c4762a1bSJed Brown v[0] = -1.*(symmetric ? PetscMin(i+1,cols[0]+1) : i+1); 173c4762a1bSJed Brown v[1] = 2.*(symmetric ? PetscMin(i+1,cols[1]+1) : i+1); 174c4762a1bSJed Brown v[2] = -1.*(symmetric ? PetscMin(i+1,cols[2]+1) : i+1); 1759566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(cmap,IS_GTOLM_MASK,3,cols,NULL,cols)); 1769566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(B,1,&i,3,cols,v,ADD_VALUES)); 177c4762a1bSJed Brown } 1789566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1799566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 180e432b41dSStefano Zampini 181e432b41dSStefano Zampini /* test MatView */ 1829566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatView\n")); 1839566063dSJacob Faibussowitsch PetscCall(MatView(A,NULL)); 1849566063dSJacob Faibussowitsch PetscCall(MatView(B,NULL)); 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* test CheckMat */ 1879566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test CheckMat\n")); 1889566063dSJacob Faibussowitsch PetscCall(CheckMat(A,B,PETSC_FALSE,"CheckMat")); 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* test MatDuplicate and MatAXPY */ 1919566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatDuplicate and MatAXPY\n")); 1929566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 1939566063dSJacob Faibussowitsch PetscCall(CheckMat(A,A2,PETSC_FALSE,"MatDuplicate and MatAXPY")); 194c4762a1bSJed Brown 195c4762a1bSJed Brown /* test MatConvert */ 1969566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatConvert_IS_XAIJ\n")); 1979566063dSJacob Faibussowitsch PetscCall(MatConvert(A2,MATAIJ,MAT_INITIAL_MATRIX,&B2)); 1989566063dSJacob Faibussowitsch PetscCall(CheckMat(B,B2,PETSC_TRUE,"MatConvert_IS_XAIJ MAT_INITIAL_MATRIX")); 1999566063dSJacob Faibussowitsch PetscCall(MatConvert(A2,MATAIJ,MAT_REUSE_MATRIX,&B2)); 2009566063dSJacob Faibussowitsch PetscCall(CheckMat(B,B2,PETSC_TRUE,"MatConvert_IS_XAIJ MAT_REUSE_MATRIX")); 2019566063dSJacob Faibussowitsch PetscCall(MatConvert(A2,MATAIJ,MAT_INPLACE_MATRIX,&A2)); 2029566063dSJacob Faibussowitsch PetscCall(CheckMat(B,A2,PETSC_TRUE,"MatConvert_IS_XAIJ MAT_INPLACE_MATRIX")); 2039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 2049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 2059566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatConvert_XAIJ_IS\n")); 2069566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&B2)); 2079566063dSJacob Faibussowitsch PetscCall(MatConvert(B2,MATIS,MAT_INITIAL_MATRIX,&A2)); 2089566063dSJacob Faibussowitsch PetscCall(CheckMat(A,A2,PETSC_TRUE,"MatConvert_XAIJ_IS MAT_INITIAL_MATRIX")); 2099566063dSJacob Faibussowitsch PetscCall(MatConvert(B2,MATIS,MAT_REUSE_MATRIX,&A2)); 2109566063dSJacob Faibussowitsch PetscCall(CheckMat(A,A2,PETSC_TRUE,"MatConvert_XAIJ_IS MAT_REUSE_MATRIX")); 2119566063dSJacob Faibussowitsch PetscCall(MatConvert(B2,MATIS,MAT_INPLACE_MATRIX,&B2)); 2129566063dSJacob Faibussowitsch PetscCall(CheckMat(A,B2,PETSC_TRUE,"MatConvert_XAIJ_IS MAT_INPLACE_MATRIX")); 2139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 2149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 2159566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(lmtype,MATSEQAIJ,&isaij)); 216c4762a1bSJed Brown if (size == 1 && isaij) { /* tests special code paths in MatConvert_IS_XAIJ */ 217c4762a1bSJed Brown PetscInt ri, ci, rr[3] = {0,1,0}, cr[4] = {1,2,0,1}, rk[3] = {0,2,1}, ck[4] = {1,0,3,2}; 218e432b41dSStefano Zampini ISLocalToGlobalMapping tcmap,trmap; 219c4762a1bSJed Brown 220c4762a1bSJed Brown for (ri = 0; ri < 2; ri++) { 221c4762a1bSJed Brown PetscInt *r; 222c4762a1bSJed Brown 223c4762a1bSJed Brown r = (PetscInt*)(ri == 0 ? rr : rk); 224c4762a1bSJed Brown for (ci = 0; ci < 2; ci++) { 225c4762a1bSJed Brown PetscInt *c,rb,cb; 226c4762a1bSJed Brown 227c4762a1bSJed Brown c = (PetscInt*)(ci == 0 ? cr : ck); 228c4762a1bSJed Brown for (rb = 1; rb < 4; rb++) { 2299566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF,rb,3,r,PETSC_COPY_VALUES,&is)); 2309566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(is,&trmap)); 2319566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 232c4762a1bSJed Brown for (cb = 1; cb < 4; cb++) { 233c4762a1bSJed Brown Mat T,lT,T2; 234c4762a1bSJed Brown char testname[256]; 235c4762a1bSJed Brown 2369566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(testname,sizeof(testname),"MatConvert_IS_XAIJ special case (%" PetscInt_FMT " %" PetscInt_FMT ", bs %" PetscInt_FMT " %" PetscInt_FMT ")",ri,ci,rb,cb)); 2379566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test %s\n",testname)); 238c4762a1bSJed Brown 2399566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF,cb,4,c,PETSC_COPY_VALUES,&is)); 2409566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(is,&tcmap)); 2419566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 242c4762a1bSJed Brown 2439566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&T)); 2449566063dSJacob Faibussowitsch PetscCall(MatSetSizes(T,PETSC_DECIDE,PETSC_DECIDE,rb*3,cb*4)); 2459566063dSJacob Faibussowitsch PetscCall(MatSetType(T,MATIS)); 2469566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(T,trmap,tcmap)); 2479566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&tcmap)); 2489566063dSJacob Faibussowitsch PetscCall(MatISGetLocalMat(T,&lT)); 2499566063dSJacob Faibussowitsch PetscCall(MatSetType(lT,MATSEQAIJ)); 2509566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(lT,cb*4,NULL)); 2519566063dSJacob Faibussowitsch PetscCall(MatSetRandom(lT,NULL)); 2529566063dSJacob Faibussowitsch PetscCall(MatConvert(lT,lmtype,MAT_INPLACE_MATRIX,&lT)); 2539566063dSJacob Faibussowitsch PetscCall(MatISRestoreLocalMat(T,&lT)); 2549566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY)); 2559566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY)); 256c4762a1bSJed Brown 2579566063dSJacob Faibussowitsch PetscCall(MatConvert(T,MATAIJ,MAT_INITIAL_MATRIX,&T2)); 2589566063dSJacob Faibussowitsch PetscCall(CheckMat(T,T2,PETSC_TRUE,"MAT_INITIAL_MATRIX")); 2599566063dSJacob Faibussowitsch PetscCall(MatConvert(T,MATAIJ,MAT_REUSE_MATRIX,&T2)); 2609566063dSJacob Faibussowitsch PetscCall(CheckMat(T,T2,PETSC_TRUE,"MAT_REUSE_MATRIX")); 2619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T2)); 2629566063dSJacob Faibussowitsch PetscCall(MatDuplicate(T,MAT_COPY_VALUES,&T2)); 2639566063dSJacob Faibussowitsch PetscCall(MatConvert(T2,MATAIJ,MAT_INPLACE_MATRIX,&T2)); 2649566063dSJacob Faibussowitsch PetscCall(CheckMat(T,T2,PETSC_TRUE,"MAT_INPLACE_MATRIX")); 2659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T)); 2669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T2)); 267c4762a1bSJed Brown } 2689566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&trmap)); 269c4762a1bSJed Brown } 270c4762a1bSJed Brown } 271c4762a1bSJed Brown } 272c4762a1bSJed Brown } 273c4762a1bSJed Brown 274c4762a1bSJed Brown /* test MatDiagonalScale */ 2759566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatDiagonalScale\n")); 2769566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 2779566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&B2)); 2789566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&x,&y)); 2799566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,NULL)); 280e432b41dSStefano Zampini if (issymmetric) { 2819566063dSJacob Faibussowitsch PetscCall(VecCopy(x,y)); 282c4762a1bSJed Brown } else { 2839566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y,NULL)); 2849566063dSJacob Faibussowitsch PetscCall(VecScale(y,8.)); 285c4762a1bSJed Brown } 2869566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A2,y,x)); 2879566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(B2,y,x)); 2889566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatDiagonalScale")); 2899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 2909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 2919566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 293c4762a1bSJed Brown 294c4762a1bSJed Brown /* test MatPtAP (A IS and B AIJ) */ 295c4762a1bSJed Brown if (isaij && m == n) { 2969566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatPtAP\n")); 2979566063dSJacob Faibussowitsch PetscCall(MatISStoreL2L(A,PETSC_TRUE)); 2989566063dSJacob Faibussowitsch PetscCall(MatPtAP(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&A2)); 2999566063dSJacob Faibussowitsch PetscCall(MatPtAP(B,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B2)); 3009566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatPtAP MAT_INITIAL_MATRIX")); 3019566063dSJacob Faibussowitsch PetscCall(MatPtAP(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&A2)); 3029566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatPtAP MAT_REUSE_MATRIX")); 3039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 3049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 305c4762a1bSJed Brown } 306c4762a1bSJed Brown 307c4762a1bSJed Brown /* test MatGetLocalSubMatrix */ 3089566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatGetLocalSubMatrix\n")); 3099566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&A2)); 3109566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,lm/2+lm%2,0,2,&reven)); 3119566063dSJacob Faibussowitsch PetscCall(ISComplement(reven,0,lm,&rodd)); 3129566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,ln/2+ln%2,0,2,&ceven)); 3139566063dSJacob Faibussowitsch PetscCall(ISComplement(ceven,0,ln,&codd)); 3149566063dSJacob Faibussowitsch PetscCall(MatGetLocalSubMatrix(A2,reven,ceven,&Aee)); 3159566063dSJacob Faibussowitsch PetscCall(MatGetLocalSubMatrix(A2,reven,codd,&Aeo)); 3169566063dSJacob Faibussowitsch PetscCall(MatGetLocalSubMatrix(A2,rodd,ceven,&Aoe)); 3179566063dSJacob Faibussowitsch PetscCall(MatGetLocalSubMatrix(A2,rodd,codd,&Aoo)); 318e432b41dSStefano Zampini for (i=0; i<lm; i++) { 319c4762a1bSJed Brown PetscInt j,je,jo,colse[3], colso[3]; 320c4762a1bSJed Brown PetscScalar ve[3], vo[3]; 321c4762a1bSJed Brown PetscScalar v[3]; 322c4762a1bSJed Brown PetscInt cols[3]; 323e432b41dSStefano Zampini PetscInt row = i/2; 324c4762a1bSJed Brown 325c4762a1bSJed Brown cols[0] = (i-1+n)%n; 326c4762a1bSJed Brown cols[1] = i%n; 327c4762a1bSJed Brown cols[2] = (i+1)%n; 328c4762a1bSJed Brown v[0] = -1.*(symmetric ? PetscMin(i+1,cols[0]+1) : i+1); 329c4762a1bSJed Brown v[1] = 2.*(symmetric ? PetscMin(i+1,cols[1]+1) : i+1); 330c4762a1bSJed Brown v[2] = -1.*(symmetric ? PetscMin(i+1,cols[2]+1) : i+1); 3319566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(cmap,IS_GTOLM_MASK,3,cols,NULL,cols)); 332c4762a1bSJed Brown for (j=0,je=0,jo=0;j<3;j++) { 333c4762a1bSJed Brown if (cols[j]%2) { 334c4762a1bSJed Brown vo[jo] = v[j]; 335c4762a1bSJed Brown colso[jo++] = cols[j]/2; 336c4762a1bSJed Brown } else { 337c4762a1bSJed Brown ve[je] = v[j]; 338c4762a1bSJed Brown colse[je++] = cols[j]/2; 339c4762a1bSJed Brown } 340c4762a1bSJed Brown } 341c4762a1bSJed Brown if (i%2) { 3429566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(Aoe,1,&row,je,colse,ve,ADD_VALUES)); 3439566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(Aoo,1,&row,jo,colso,vo,ADD_VALUES)); 344c4762a1bSJed Brown } else { 3459566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(Aee,1,&row,je,colse,ve,ADD_VALUES)); 3469566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(Aeo,1,&row,jo,colso,vo,ADD_VALUES)); 347c4762a1bSJed Brown } 348c4762a1bSJed Brown } 3499566063dSJacob Faibussowitsch PetscCall(MatRestoreLocalSubMatrix(A2,reven,ceven,&Aee)); 3509566063dSJacob Faibussowitsch PetscCall(MatRestoreLocalSubMatrix(A2,reven,codd,&Aeo)); 3519566063dSJacob Faibussowitsch PetscCall(MatRestoreLocalSubMatrix(A2,rodd,ceven,&Aoe)); 3529566063dSJacob Faibussowitsch PetscCall(MatRestoreLocalSubMatrix(A2,rodd,codd,&Aoo)); 3539566063dSJacob Faibussowitsch PetscCall(ISDestroy(&reven)); 3549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ceven)); 3559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&rodd)); 3569566063dSJacob Faibussowitsch PetscCall(ISDestroy(&codd)); 3579566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY)); 3589566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY)); 3599566063dSJacob Faibussowitsch PetscCall(MatAXPY(A2,-1.,A,SAME_NONZERO_PATTERN)); 3609566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,NULL,PETSC_FALSE,"MatGetLocalSubMatrix")); 3619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 362c4762a1bSJed Brown 363c4762a1bSJed Brown /* test MatConvert_Nest_IS */ 364c4762a1bSJed Brown testT = PETSC_FALSE; 3659566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_trans",&testT,NULL)); 366c4762a1bSJed Brown 3679566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatConvert_Nest_IS\n")); 368c4762a1bSJed Brown nr = 2; 369c4762a1bSJed Brown nc = 2; 3709566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nr",&nr,NULL)); 3719566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nc",&nc,NULL)); 372c4762a1bSJed Brown if (testT) { 3739566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&cst,&cen)); 3749566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(A,&rst,&ren)); 375c4762a1bSJed Brown } else { 3769566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rst,&ren)); 3779566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(A,&cst,&cen)); 378c4762a1bSJed Brown } 3799566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(nr,&rows,nc,&cols,2*nr*nc,&mats)); 380c4762a1bSJed Brown for (i=0;i<nr*nc;i++) { 381c4762a1bSJed Brown if (testT) { 3829566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(A,&mats[i])); 3839566063dSJacob Faibussowitsch PetscCall(MatTranspose(B,MAT_INITIAL_MATRIX,&mats[i+nr*nc])); 384c4762a1bSJed Brown } else { 3859566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&mats[i])); 3869566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&mats[i+nr*nc])); 387c4762a1bSJed Brown } 388c4762a1bSJed Brown } 389c4762a1bSJed Brown for (i=0;i<nr;i++) { 3909566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,ren-rst,i+rst,nr,&rows[i])); 391c4762a1bSJed Brown } 392c4762a1bSJed Brown for (i=0;i<nc;i++) { 3939566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,cen-cst,i+cst,nc,&cols[i])); 394c4762a1bSJed Brown } 3959566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD,nr,rows,nc,cols,mats,&A2)); 3969566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD,nr,rows,nc,cols,mats+nr*nc,&B2)); 397c4762a1bSJed Brown for (i=0;i<nr;i++) { 3989566063dSJacob Faibussowitsch PetscCall(ISDestroy(&rows[i])); 399c4762a1bSJed Brown } 400c4762a1bSJed Brown for (i=0;i<nc;i++) { 4019566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cols[i])); 402c4762a1bSJed Brown } 403c4762a1bSJed Brown for (i=0;i<2*nr*nc;i++) { 4049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mats[i])); 405c4762a1bSJed Brown } 4069566063dSJacob Faibussowitsch PetscCall(PetscFree3(rows,cols,mats)); 4079566063dSJacob Faibussowitsch PetscCall(MatConvert(B2,MATAIJ,MAT_INITIAL_MATRIX,&T)); 4089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 4099566063dSJacob Faibussowitsch PetscCall(MatConvert(A2,MATIS,MAT_INITIAL_MATRIX,&B2)); 4109566063dSJacob Faibussowitsch PetscCall(CheckMat(B2,T,PETSC_TRUE,"MatConvert_Nest_IS MAT_INITIAL_MATRIX")); 4119566063dSJacob Faibussowitsch PetscCall(MatConvert(A2,MATIS,MAT_REUSE_MATRIX,&B2)); 4129566063dSJacob Faibussowitsch PetscCall(CheckMat(B2,T,PETSC_TRUE,"MatConvert_Nest_IS MAT_REUSE_MATRIX")); 4139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 4149566063dSJacob Faibussowitsch PetscCall(MatConvert(A2,MATIS,MAT_INPLACE_MATRIX,&A2)); 4159566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,T,PETSC_TRUE,"MatConvert_Nest_IS MAT_INPLACE_MATRIX")); 4169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T)); 4179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 418c4762a1bSJed Brown 419c4762a1bSJed Brown /* test MatCreateSubMatrix */ 4209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatCreateSubMatrix\n")); 421dd400576SPatrick Sanan if (rank == 0) { 4229566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,1,1,1,&is)); 4239566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,2,0,1,&is2)); 424c4762a1bSJed Brown } else if (rank == 1) { 4259566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,1,0,1,&is)); 426c4762a1bSJed Brown if (n > 3) { 4279566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,1,3,1,&is2)); 428c4762a1bSJed Brown } else { 4299566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is2)); 430c4762a1bSJed Brown } 431c4762a1bSJed Brown } else if (rank == 2 && n > 4) { 4329566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is)); 4339566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n-4,4,1,&is2)); 434c4762a1bSJed Brown } else { 4359566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is)); 4369566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is2)); 437c4762a1bSJed Brown } 4389566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&A2)); 4399566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(B,is,is,MAT_INITIAL_MATRIX,&B2)); 4409566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_TRUE,"first MatCreateSubMatrix")); 441c4762a1bSJed Brown 4429566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A,is,is,MAT_REUSE_MATRIX,&A2)); 4439566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(B,is,is,MAT_REUSE_MATRIX,&B2)); 4449566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"reuse MatCreateSubMatrix")); 4459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 4469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 447c4762a1bSJed Brown 448e432b41dSStefano Zampini if (!issymmetric) { 4499566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A,is,is2,MAT_INITIAL_MATRIX,&A2)); 4509566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(B,is,is2,MAT_INITIAL_MATRIX,&B2)); 4519566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A,is,is2,MAT_REUSE_MATRIX,&A2)); 4529566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(B,is,is2,MAT_REUSE_MATRIX,&B2)); 4539566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"second MatCreateSubMatrix")); 454c4762a1bSJed Brown } 455c4762a1bSJed Brown 4569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 4579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 4589566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 4599566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2)); 460c4762a1bSJed Brown 461*d0dbe9f7SStefano Zampini /* test MatCreateSubMatrices */ 462*d0dbe9f7SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatCreateSubMatrices\n")); 463*d0dbe9f7SStefano Zampini PetscCall(MatGetLayouts(A,&rlayout,&clayout)); 464*d0dbe9f7SStefano Zampini PetscCall(PetscLayoutGetRanges(rlayout,&rrange)); 465*d0dbe9f7SStefano Zampini PetscCall(PetscLayoutGetRanges(clayout,&crange)); 466*d0dbe9f7SStefano Zampini lrank = (size + rank - 1)%size; 467*d0dbe9f7SStefano Zampini rrank = (rank + 1)%size; 468*d0dbe9f7SStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF,(rrange[lrank+1] - rrange[lrank]),rrange[lrank],1,&irow[0])); 469*d0dbe9f7SStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF,(crange[rrank+1] - crange[rrank]),crange[rrank],1,&icol[0])); 470*d0dbe9f7SStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF,(rrange[rrank+1] - rrange[rrank]),rrange[rrank],1,&irow[1])); 471*d0dbe9f7SStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF,(crange[lrank+1] - crange[lrank]),crange[lrank],1,&icol[1])); 472*d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(A,2,irow,icol,MAT_INITIAL_MATRIX,&Asub)); 473*d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(B,2,irow,icol,MAT_INITIAL_MATRIX,&Bsub)); 474*d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[0],Bsub[0],PETSC_FALSE,"MatCreateSubMatrices[0]")); 475*d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[1],Bsub[1],PETSC_FALSE,"MatCreateSubMatrices[1]")); 476*d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(A,2,irow,icol,MAT_REUSE_MATRIX,&Asub)); 477*d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(B,2,irow,icol,MAT_REUSE_MATRIX,&Bsub)); 478*d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[0],Bsub[0],PETSC_FALSE,"MatCreateSubMatrices[0]")); 479*d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[1],Bsub[1],PETSC_FALSE,"MatCreateSubMatrices[1]")); 480*d0dbe9f7SStefano Zampini PetscCall(MatDestroySubMatrices(2,&Asub)); 481*d0dbe9f7SStefano Zampini PetscCall(MatDestroySubMatrices(2,&Bsub)); 482*d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&irow[0])); 483*d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&irow[1])); 484*d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&icol[0])); 485*d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&icol[1])); 486*d0dbe9f7SStefano Zampini 487c4762a1bSJed Brown /* Create an IS required by MatZeroRows(): just rank zero provides the rows to be eliminated */ 488c4762a1bSJed Brown if (size > 1) { 489dd400576SPatrick Sanan if (rank == 0) { 490c4762a1bSJed Brown PetscInt st,len; 491c4762a1bSJed Brown 492c4762a1bSJed Brown st = (m+1)/2; 493c4762a1bSJed Brown len = PetscMin(m/2,PetscMax(m-(m+1)/2-1,0)); 4949566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,len,st,1,&is)); 495c4762a1bSJed Brown } else { 4969566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is)); 497c4762a1bSJed Brown } 498c4762a1bSJed Brown } else { 4999566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,1,0,1,&is)); 500c4762a1bSJed Brown } 501c4762a1bSJed Brown 502c4762a1bSJed Brown if (squaretest) { /* tests for square matrices only, with same maps for rows and columns */ 503*d0dbe9f7SStefano Zampini PetscInt *idx0, *idx1, n0, n1; 504*d0dbe9f7SStefano Zampini IS Ais[2], Bis[2]; 505*d0dbe9f7SStefano Zampini 506c4762a1bSJed Brown /* test MatDiagonalSet */ 5079566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatDiagonalSet\n")); 5089566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 5099566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&B2)); 5109566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,NULL,&x)); 5119566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,NULL)); 5129566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(A2,x,INSERT_VALUES)); 5139566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(B2,x,INSERT_VALUES)); 5149566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatDiagonalSet")); 5159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 5169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 5179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 518c4762a1bSJed Brown 519c4762a1bSJed Brown /* test MatShift (MatShift_IS internally uses MatDiagonalSet_IS with ADD_VALUES) */ 5209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatShift\n")); 5219566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 5229566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&B2)); 5239566063dSJacob Faibussowitsch PetscCall(MatShift(A2,2.0)); 5249566063dSJacob Faibussowitsch PetscCall(MatShift(B2,2.0)); 5259566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatShift")); 5269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 5279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 528c4762a1bSJed Brown 529c4762a1bSJed Brown /* nonzero diag value is supported for square matrices only */ 5309566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows(A,B,PETSC_TRUE,is,diag)); 531*d0dbe9f7SStefano Zampini 532*d0dbe9f7SStefano Zampini /* test MatIncreaseOverlap */ 533*d0dbe9f7SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatIncreaseOverlap\n")); 534*d0dbe9f7SStefano Zampini PetscCall(MatGetOwnershipRange(A,&rst,&ren)); 535*d0dbe9f7SStefano Zampini n0 = (ren - rst)/2; 536*d0dbe9f7SStefano Zampini n1 = (ren - rst)/3; 537*d0dbe9f7SStefano Zampini PetscCall(PetscMalloc1(n0,&idx0)); 538*d0dbe9f7SStefano Zampini PetscCall(PetscMalloc1(n1,&idx1)); 539*d0dbe9f7SStefano Zampini for (i = 0; i < n0; i++) idx0[i] = ren - i - 1; 540*d0dbe9f7SStefano Zampini for (i = 0; i < n1; i++) idx1[i] = rst + i; 541*d0dbe9f7SStefano Zampini PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,n0,idx0,PETSC_OWN_POINTER,&Ais[0])); 542*d0dbe9f7SStefano Zampini PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,n1,idx1,PETSC_OWN_POINTER,&Ais[1])); 543*d0dbe9f7SStefano Zampini PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,n0,idx0,PETSC_COPY_VALUES,&Bis[0])); 544*d0dbe9f7SStefano Zampini PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,n1,idx1,PETSC_COPY_VALUES,&Bis[1])); 545*d0dbe9f7SStefano Zampini PetscCall(MatIncreaseOverlap(A,2,Ais,3)); 546*d0dbe9f7SStefano Zampini PetscCall(MatIncreaseOverlap(B,2,Bis,3)); 547*d0dbe9f7SStefano Zampini /* Non deterministic output! */ 548*d0dbe9f7SStefano Zampini PetscCall(ISSort(Ais[0])); 549*d0dbe9f7SStefano Zampini PetscCall(ISSort(Ais[1])); 550*d0dbe9f7SStefano Zampini PetscCall(ISSort(Bis[0])); 551*d0dbe9f7SStefano Zampini PetscCall(ISSort(Bis[1])); 552*d0dbe9f7SStefano Zampini PetscCall(ISView(Ais[0],NULL)); 553*d0dbe9f7SStefano Zampini PetscCall(ISView(Bis[0],NULL)); 554*d0dbe9f7SStefano Zampini PetscCall(ISView(Ais[1],NULL)); 555*d0dbe9f7SStefano Zampini PetscCall(ISView(Bis[1],NULL)); 556*d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(A,2,Ais,Ais,MAT_INITIAL_MATRIX,&Asub)); 557*d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(B,2,Bis,Ais,MAT_INITIAL_MATRIX,&Bsub)); 558*d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[0],Bsub[0],PETSC_FALSE,"MatIncreaseOverlap[0]")); 559*d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[1],Bsub[1],PETSC_FALSE,"MatIncreaseOverlap[1]")); 560*d0dbe9f7SStefano Zampini PetscCall(MatDestroySubMatrices(2,&Asub)); 561*d0dbe9f7SStefano Zampini PetscCall(MatDestroySubMatrices(2,&Bsub)); 562*d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&Ais[0])); 563*d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&Ais[1])); 564*d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&Bis[0])); 565*d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&Bis[1])); 566*d0dbe9f7SStefano Zampini 567c4762a1bSJed Brown } 5689566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows(A,B,squaretest,is,0.0)); 5699566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 570c4762a1bSJed Brown 571c4762a1bSJed Brown /* test MatTranspose */ 5729566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatTranspose\n")); 5739566063dSJacob Faibussowitsch PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX,&A2)); 5749566063dSJacob Faibussowitsch PetscCall(MatTranspose(B,MAT_INITIAL_MATRIX,&B2)); 5759566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"initial matrix MatTranspose")); 576c4762a1bSJed Brown 5779566063dSJacob Faibussowitsch PetscCall(MatTranspose(A,MAT_REUSE_MATRIX,&A2)); 5789566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"reuse matrix (not in place) MatTranspose")); 5799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 580c4762a1bSJed Brown 5819566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 5829566063dSJacob Faibussowitsch PetscCall(MatTranspose(A2,MAT_INPLACE_MATRIX,&A2)); 5839566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"reuse matrix (in place) MatTranspose")); 5849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 585c4762a1bSJed Brown 5869566063dSJacob Faibussowitsch PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX,&A2)); 5879566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"reuse matrix (different type) MatTranspose")); 5889566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 5899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 590c4762a1bSJed Brown 591c4762a1bSJed Brown /* test MatISFixLocalEmpty */ 592c4762a1bSJed Brown if (isaij) { 593c4762a1bSJed Brown PetscInt r[2]; 594c4762a1bSJed Brown 595c4762a1bSJed Brown r[0] = 0; 596c4762a1bSJed Brown r[1] = PetscMin(m,n)-1; 5979566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatISFixLocalEmpty\n")); 5989566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 599e432b41dSStefano Zampini 6009566063dSJacob Faibussowitsch PetscCall(MatISFixLocalEmpty(A2,PETSC_TRUE)); 6019566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY)); 6029566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY)); 6039566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B,PETSC_FALSE,"MatISFixLocalEmpty (null)")); 604c4762a1bSJed Brown 6059566063dSJacob Faibussowitsch PetscCall(MatZeroRows(A2,2,r,0.0,NULL,NULL)); 6069566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2,NULL,"-fixempty_view")); 6079566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&B2)); 6089566063dSJacob Faibussowitsch PetscCall(MatZeroRows(B2,2,r,0.0,NULL,NULL)); 6099566063dSJacob Faibussowitsch PetscCall(MatISFixLocalEmpty(A2,PETSC_TRUE)); 6109566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY)); 6119566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY)); 6129566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2,NULL,"-fixempty_view")); 6139566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatISFixLocalEmpty (rows)")); 6149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 615c4762a1bSJed Brown 6169566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 6179566063dSJacob Faibussowitsch PetscCall(MatZeroRows(A2,2,r,0.0,NULL,NULL)); 6189566063dSJacob Faibussowitsch PetscCall(MatTranspose(A2,MAT_INPLACE_MATRIX,&A2)); 6199566063dSJacob Faibussowitsch PetscCall(MatTranspose(B2,MAT_INPLACE_MATRIX,&B2)); 6209566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2,NULL,"-fixempty_view")); 6219566063dSJacob Faibussowitsch PetscCall(MatISFixLocalEmpty(A2,PETSC_TRUE)); 6229566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY)); 6239566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY)); 6249566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2,NULL,"-fixempty_view")); 6259566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatISFixLocalEmpty (cols)")); 626c4762a1bSJed Brown 6279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 6289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 629c4762a1bSJed Brown 630c4762a1bSJed Brown if (squaretest) { 6319566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 6329566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&B2)); 6339566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumns(A2,2,r,0.0,NULL,NULL)); 6349566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumns(B2,2,r,0.0,NULL,NULL)); 6359566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2,NULL,"-fixempty_view")); 6369566063dSJacob Faibussowitsch PetscCall(MatISFixLocalEmpty(A2,PETSC_TRUE)); 6379566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY)); 6389566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY)); 6399566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2,NULL,"-fixempty_view")); 6409566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatISFixLocalEmpty (rows+cols)")); 6419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 6429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 643c4762a1bSJed Brown } 644c4762a1bSJed Brown 645c4762a1bSJed Brown } 646c4762a1bSJed Brown 647c4762a1bSJed Brown /* test MatInvertBlockDiagonal 648c4762a1bSJed Brown special cases for block-diagonal matrices */ 649c4762a1bSJed Brown if (m == n) { 650c4762a1bSJed Brown ISLocalToGlobalMapping map; 651c4762a1bSJed Brown Mat Abd,Bbd; 652c4762a1bSJed Brown IS is,bis; 653c4762a1bSJed Brown const PetscScalar *isbd,*aijbd; 654c4762a1bSJed Brown PetscScalar *vals; 655c4762a1bSJed Brown const PetscInt *sts,*idxs; 656c4762a1bSJed Brown PetscInt *idxs2,diff,perm,nl,bs,st,en,in; 657c4762a1bSJed Brown PetscBool ok; 658c4762a1bSJed Brown 659c4762a1bSJed Brown for (diff = 0; diff < 3; diff++) { 660c4762a1bSJed Brown for (perm = 0; perm < 3; perm++) { 661c4762a1bSJed Brown for (bs = 1; bs < 4; bs++) { 6629566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatInvertBlockDiagonal blockdiag %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",n,diff,perm,bs)); 6639566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs*bs,&vals)); 6649566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(A,&sts)); 665c4762a1bSJed Brown switch (diff) { 666c4762a1bSJed Brown case 1: /* inverted layout by processes */ 667c4762a1bSJed Brown in = 1; 668c4762a1bSJed Brown st = sts[size - rank - 1]; 669c4762a1bSJed Brown en = sts[size - rank]; 670c4762a1bSJed Brown nl = en - st; 671c4762a1bSJed Brown break; 672c4762a1bSJed Brown case 2: /* round-robin layout */ 673c4762a1bSJed Brown in = size; 674c4762a1bSJed Brown st = rank; 675c4762a1bSJed Brown nl = n/size; 676c4762a1bSJed Brown if (rank < n%size) nl++; 677c4762a1bSJed Brown break; 678c4762a1bSJed Brown default: /* same layout */ 679c4762a1bSJed Brown in = 1; 680c4762a1bSJed Brown st = sts[rank]; 681c4762a1bSJed Brown en = sts[rank + 1]; 682c4762a1bSJed Brown nl = en - st; 683c4762a1bSJed Brown break; 684c4762a1bSJed Brown } 6859566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,nl,st,in,&is)); 6869566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is,&nl)); 6879566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is,&idxs)); 6889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nl,&idxs2)); 689c4762a1bSJed Brown for (i=0;i<nl;i++) { 690c4762a1bSJed Brown switch (perm) { /* invert some of the indices */ 691c4762a1bSJed Brown case 2: 692c4762a1bSJed Brown idxs2[i] = rank%2 ? idxs[i] : idxs[nl-i-1]; 693c4762a1bSJed Brown break; 694c4762a1bSJed Brown case 1: 695c4762a1bSJed Brown idxs2[i] = rank%2 ? idxs[nl-i-1] : idxs[i]; 696c4762a1bSJed Brown break; 697c4762a1bSJed Brown default: 698c4762a1bSJed Brown idxs2[i] = idxs[i]; 699c4762a1bSJed Brown break; 700c4762a1bSJed Brown } 701c4762a1bSJed Brown } 7029566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is,&idxs)); 7039566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_WORLD,bs,nl,idxs2,PETSC_OWN_POINTER,&bis)); 7049566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(bis,&map)); 7059566063dSJacob Faibussowitsch PetscCall(MatCreateIS(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,bs*n,bs*n,map,map,&Abd)); 7069566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&map)); 7079566063dSJacob Faibussowitsch PetscCall(MatISSetPreallocation(Abd,bs,NULL,0,NULL)); 708c4762a1bSJed Brown for (i=0;i<nl;i++) { 709c4762a1bSJed Brown PetscInt b1,b2; 710c4762a1bSJed Brown 711c4762a1bSJed Brown for (b1=0;b1<bs;b1++) { 712c4762a1bSJed Brown for (b2=0;b2<bs;b2++) { 713c4762a1bSJed Brown vals[b1*bs + b2] = i*bs*bs + b1*bs + b2 + 1 + (b1 == b2 ? 1.0 : 0); 714c4762a1bSJed Brown } 715c4762a1bSJed Brown } 7169566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedLocal(Abd,1,&i,1,&i,vals,INSERT_VALUES)); 717c4762a1bSJed Brown } 7189566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Abd,MAT_FINAL_ASSEMBLY)); 7199566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Abd,MAT_FINAL_ASSEMBLY)); 7209566063dSJacob Faibussowitsch PetscCall(MatConvert(Abd,MATAIJ,MAT_INITIAL_MATRIX,&Bbd)); 7219566063dSJacob Faibussowitsch PetscCall(MatInvertBlockDiagonal(Abd,&isbd)); 7229566063dSJacob Faibussowitsch PetscCall(MatInvertBlockDiagonal(Bbd,&aijbd)); 7239566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Bbd,&nl,NULL)); 724c4762a1bSJed Brown ok = PETSC_TRUE; 725c4762a1bSJed Brown for (i=0;i<nl/bs;i++) { 726c4762a1bSJed Brown PetscInt b1,b2; 727c4762a1bSJed Brown 728c4762a1bSJed Brown for (b1=0;b1<bs;b1++) { 729c4762a1bSJed Brown for (b2=0;b2<bs;b2++) { 730c4762a1bSJed Brown if (PetscAbsScalar(isbd[i*bs*bs+b1*bs + b2]-aijbd[i*bs*bs+b1*bs + b2]) > PETSC_SMALL) ok = PETSC_FALSE; 731c4762a1bSJed Brown if (!ok) { 7329566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"[%d] ERROR block %" PetscInt_FMT ", entry %" PetscInt_FMT " %" PetscInt_FMT ": %g %g\n",rank,i,b1,b2,(double)PetscAbsScalar(isbd[i*bs*bs+b1*bs + b2]),(double)PetscAbsScalar(aijbd[i*bs*bs+b1*bs + b2]))); 733c4762a1bSJed Brown break; 734c4762a1bSJed Brown } 735c4762a1bSJed Brown } 736c4762a1bSJed Brown if (!ok) break; 737c4762a1bSJed Brown } 738c4762a1bSJed Brown if (!ok) break; 739c4762a1bSJed Brown } 7409566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Abd)); 7419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bbd)); 7429566063dSJacob Faibussowitsch PetscCall(PetscFree(vals)); 7439566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 7449566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bis)); 745c4762a1bSJed Brown } 746c4762a1bSJed Brown } 747c4762a1bSJed Brown } 748c4762a1bSJed Brown } 749*d0dbe9f7SStefano Zampini 750*d0dbe9f7SStefano Zampini /* test MatGetDiagonalBlock */ 751*d0dbe9f7SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatGetDiagonalBlock\n")); 752*d0dbe9f7SStefano Zampini PetscCall(MatGetDiagonalBlock(A,&A2)); 753*d0dbe9f7SStefano Zampini PetscCall(MatGetDiagonalBlock(B,&B2)); 754*d0dbe9f7SStefano Zampini PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatGetDiagonalBlock")); 755*d0dbe9f7SStefano Zampini PetscCall(MatScale(A,2.0)); 756*d0dbe9f7SStefano Zampini PetscCall(MatScale(B,2.0)); 757*d0dbe9f7SStefano Zampini PetscCall(MatGetDiagonalBlock(A,&A2)); 758*d0dbe9f7SStefano Zampini PetscCall(MatGetDiagonalBlock(B,&B2)); 759*d0dbe9f7SStefano Zampini PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatGetDiagonalBlock")); 760*d0dbe9f7SStefano Zampini 761c4762a1bSJed Brown /* free testing matrices */ 7629566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&cmap)); 7639566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&rmap)); 7649566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 7659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 7669566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 767b122ec5aSJacob Faibussowitsch return 0; 768c4762a1bSJed Brown } 769c4762a1bSJed Brown 770c4762a1bSJed Brown PetscErrorCode CheckMat(Mat A, Mat B, PetscBool usemult, const char* func) 771c4762a1bSJed Brown { 772c4762a1bSJed Brown Mat Bcheck; 773c4762a1bSJed Brown PetscReal error; 774c4762a1bSJed Brown 775c4762a1bSJed Brown PetscFunctionBeginUser; 776c4762a1bSJed Brown if (!usemult && B) { 777c4762a1bSJed Brown PetscBool hasnorm; 778c4762a1bSJed Brown 7799566063dSJacob Faibussowitsch PetscCall(MatHasOperation(B,MATOP_NORM,&hasnorm)); 780c4762a1bSJed Brown if (!hasnorm) usemult = PETSC_TRUE; 781c4762a1bSJed Brown } 782c4762a1bSJed Brown if (!usemult) { 783c4762a1bSJed Brown if (B) { 784c4762a1bSJed Brown MatType Btype; 785c4762a1bSJed Brown 7869566063dSJacob Faibussowitsch PetscCall(MatGetType(B,&Btype)); 7879566063dSJacob Faibussowitsch PetscCall(MatConvert(A,Btype,MAT_INITIAL_MATRIX,&Bcheck)); 788c4762a1bSJed Brown } else { 7899566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&Bcheck)); 790c4762a1bSJed Brown } 791c4762a1bSJed Brown if (B) { /* if B is present, subtract it */ 7929566063dSJacob Faibussowitsch PetscCall(MatAXPY(Bcheck,-1.,B,DIFFERENT_NONZERO_PATTERN)); 793c4762a1bSJed Brown } 7949566063dSJacob Faibussowitsch PetscCall(MatNorm(Bcheck,NORM_INFINITY,&error)); 795c4762a1bSJed Brown if (error > PETSC_SQRT_MACHINE_EPSILON) { 796c4762a1bSJed Brown ISLocalToGlobalMapping rl2g,cl2g; 797c4762a1bSJed Brown 798*d0dbe9f7SStefano Zampini PetscCall(PetscObjectSetName((PetscObject)Bcheck,"Bcheck")); 7999566063dSJacob Faibussowitsch PetscCall(MatView(Bcheck,NULL)); 800c4762a1bSJed Brown if (B) { 801*d0dbe9f7SStefano Zampini PetscCall(PetscObjectSetName((PetscObject)B,"B")); 8029566063dSJacob Faibussowitsch PetscCall(MatView(B,NULL)); 8039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bcheck)); 8049566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&Bcheck)); 805*d0dbe9f7SStefano Zampini PetscCall(PetscObjectSetName((PetscObject)Bcheck,"Assembled A")); 8069566063dSJacob Faibussowitsch PetscCall(MatView(Bcheck,NULL)); 807c4762a1bSJed Brown } 8089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bcheck)); 809*d0dbe9f7SStefano Zampini PetscCall(PetscObjectSetName((PetscObject)A,"A")); 8109566063dSJacob Faibussowitsch PetscCall(MatView(A,NULL)); 8119566063dSJacob Faibussowitsch PetscCall(MatGetLocalToGlobalMapping(A,&rl2g,&cl2g)); 812*d0dbe9f7SStefano Zampini if (rl2g) PetscCall(ISLocalToGlobalMappingView(rl2g,NULL)); 813*d0dbe9f7SStefano Zampini if (cl2g) PetscCall(ISLocalToGlobalMappingView(cl2g,NULL)); 814*d0dbe9f7SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"ERROR ON %s: %g",func,(double)error); 815c4762a1bSJed Brown } 8169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bcheck)); 817c4762a1bSJed Brown } else { 818c4762a1bSJed Brown PetscBool ok,okt; 819c4762a1bSJed Brown 8209566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A,B,3,&ok)); 8219566063dSJacob Faibussowitsch PetscCall(MatMultTransposeEqual(A,B,3,&okt)); 822e00437b9SBarry Smith PetscCheck(ok && okt,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ERROR ON %s: mult ok ? %d, multtranspose ok ? %d",func,ok,okt); 823c4762a1bSJed Brown } 824c4762a1bSJed Brown PetscFunctionReturn(0); 825c4762a1bSJed Brown } 826c4762a1bSJed Brown 827c4762a1bSJed Brown PetscErrorCode TestMatZeroRows(Mat A, Mat Afull, PetscBool squaretest, IS is, PetscScalar diag) 828c4762a1bSJed Brown { 829c4762a1bSJed Brown Mat B,Bcheck,B2 = NULL,lB; 830c4762a1bSJed Brown Vec x = NULL, b = NULL, b2 = NULL; 831c4762a1bSJed Brown ISLocalToGlobalMapping l2gr,l2gc; 832c4762a1bSJed Brown PetscReal error; 833c4762a1bSJed Brown char diagstr[16]; 834c4762a1bSJed Brown const PetscInt *idxs; 835c4762a1bSJed Brown PetscInt rst,ren,i,n,N,d; 836c4762a1bSJed Brown PetscMPIInt rank; 837c4762a1bSJed Brown PetscBool miss,haszerorows; 838c4762a1bSJed Brown 839c4762a1bSJed Brown PetscFunctionBeginUser; 840c4762a1bSJed Brown if (diag == 0.) { 8419566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(diagstr,"zero")); 842c4762a1bSJed Brown } else { 8439566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(diagstr,"nonzero")); 844c4762a1bSJed Brown } 8459566063dSJacob Faibussowitsch PetscCall(ISView(is,NULL)); 8469566063dSJacob Faibussowitsch PetscCall(MatGetLocalToGlobalMapping(A,&l2gr,&l2gc)); 847c4762a1bSJed Brown /* tests MatDuplicate and MatCopy */ 848c4762a1bSJed Brown if (diag == 0.) { 8499566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&B)); 850c4762a1bSJed Brown } else { 8519566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&B)); 8529566063dSJacob Faibussowitsch PetscCall(MatCopy(A,B,SAME_NONZERO_PATTERN)); 853c4762a1bSJed Brown } 8549566063dSJacob Faibussowitsch PetscCall(MatISGetLocalMat(B,&lB)); 8559566063dSJacob Faibussowitsch PetscCall(MatHasOperation(lB,MATOP_ZERO_ROWS,&haszerorows)); 856c4762a1bSJed Brown if (squaretest && haszerorows) { 857c4762a1bSJed Brown 8589566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B,&x,&b)); 8599566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&B2)); 8609566063dSJacob Faibussowitsch PetscCall(VecSetLocalToGlobalMapping(b,l2gr)); 8619566063dSJacob Faibussowitsch PetscCall(VecSetLocalToGlobalMapping(x,l2gc)); 8629566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,NULL)); 8639566063dSJacob Faibussowitsch PetscCall(VecSetRandom(b,NULL)); 864c4762a1bSJed Brown /* mimic b[is] = x[is] */ 8659566063dSJacob Faibussowitsch PetscCall(VecDuplicate(b,&b2)); 8669566063dSJacob Faibussowitsch PetscCall(VecSetLocalToGlobalMapping(b2,l2gr)); 8679566063dSJacob Faibussowitsch PetscCall(VecCopy(b,b2)); 8689566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is,&n)); 8699566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is,&idxs)); 8709566063dSJacob Faibussowitsch PetscCall(VecGetSize(x,&N)); 871c4762a1bSJed Brown for (i=0;i<n;i++) { 872c4762a1bSJed Brown if (0 <= idxs[i] && idxs[i] < N) { 8739566063dSJacob Faibussowitsch PetscCall(VecSetValue(b2,idxs[i],diag,INSERT_VALUES)); 8749566063dSJacob Faibussowitsch PetscCall(VecSetValue(x,idxs[i],1.,INSERT_VALUES)); 875c4762a1bSJed Brown } 876c4762a1bSJed Brown } 8779566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(b2)); 8789566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(b2)); 8799566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x)); 8809566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x)); 8819566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is,&idxs)); 882c4762a1bSJed Brown /* test ZeroRows on MATIS */ 8839566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatZeroRows (diag %s)\n",diagstr)); 8849566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(B,is,diag,x,b)); 8859566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatZeroRowsColumns (diag %s)\n",diagstr)); 8869566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumnsIS(B2,is,diag,NULL,NULL)); 887c4762a1bSJed Brown } else if (haszerorows) { 888c4762a1bSJed Brown /* test ZeroRows on MATIS */ 8899566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatZeroRows (diag %s)\n",diagstr)); 8909566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(B,is,diag,NULL,NULL)); 891c4762a1bSJed Brown b = b2 = x = NULL; 892c4762a1bSJed Brown } else { 8939566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Skipping MatZeroRows (diag %s)\n",diagstr)); 894c4762a1bSJed Brown b = b2 = x = NULL; 895c4762a1bSJed Brown } 896c4762a1bSJed Brown 897c4762a1bSJed Brown if (squaretest && haszerorows) { 8989566063dSJacob Faibussowitsch PetscCall(VecAXPY(b2,-1.,b)); 8999566063dSJacob Faibussowitsch PetscCall(VecNorm(b2,NORM_INFINITY,&error)); 900e00437b9SBarry Smith PetscCheck(error <= PETSC_SQRT_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ERROR IN ZEROROWS ON B %g (diag %s)",(double)error,diagstr); 901c4762a1bSJed Brown } 902c4762a1bSJed Brown 903c4762a1bSJed Brown /* test MatMissingDiagonal */ 9049566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatMissingDiagonal\n")); 9059566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 9069566063dSJacob Faibussowitsch PetscCall(MatMissingDiagonal(B,&miss,&d)); 9079566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(B,&rst,&ren)); 9089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 9099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD, "[%d] [%" PetscInt_FMT ",%" PetscInt_FMT ") Missing %d, row %" PetscInt_FMT " (diag %s)\n",rank,rst,ren,(int)miss,d,diagstr)); 9109566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 9119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 912c4762a1bSJed Brown 9139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 9149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 9159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b2)); 916c4762a1bSJed Brown 917c4762a1bSJed Brown /* check the result of ZeroRows with that from MPIAIJ routines 918c4762a1bSJed Brown assuming that MatConvert_IS_XAIJ and MatZeroRows_MPIAIJ work fine */ 919c4762a1bSJed Brown if (haszerorows) { 9209566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Afull,MAT_COPY_VALUES,&Bcheck)); 9219566063dSJacob Faibussowitsch PetscCall(MatSetOption(Bcheck,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 9229566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(Bcheck,is,diag,NULL,NULL)); 9239566063dSJacob Faibussowitsch PetscCall(CheckMat(B,Bcheck,PETSC_FALSE,"Zerorows")); 9249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bcheck)); 925c4762a1bSJed Brown } 9269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 927c4762a1bSJed Brown 928c4762a1bSJed Brown if (B2) { /* test MatZeroRowsColumns */ 9299566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Afull,MAT_COPY_VALUES,&B)); 9309566063dSJacob Faibussowitsch PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 9319566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumnsIS(B,is,diag,NULL,NULL)); 9329566063dSJacob Faibussowitsch PetscCall(CheckMat(B2,B,PETSC_FALSE,"MatZeroRowsColumns")); 9339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 9349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 935c4762a1bSJed Brown } 936c4762a1bSJed Brown PetscFunctionReturn(0); 937c4762a1bSJed Brown } 938c4762a1bSJed Brown 939c4762a1bSJed Brown /*TEST 940c4762a1bSJed Brown 941c4762a1bSJed Brown test: 942c4762a1bSJed Brown args: -test_trans 943c4762a1bSJed Brown 944c4762a1bSJed Brown test: 945c4762a1bSJed Brown suffix: 2 946c4762a1bSJed Brown nsize: 4 947c4762a1bSJed Brown args: -matis_convert_local_nest -nr 3 -nc 4 948c4762a1bSJed Brown 949c4762a1bSJed Brown test: 950c4762a1bSJed Brown suffix: 3 951c4762a1bSJed Brown nsize: 5 952c4762a1bSJed Brown args: -m 11 -n 10 -matis_convert_local_nest -nr 2 -nc 1 953c4762a1bSJed Brown 954c4762a1bSJed Brown test: 955c4762a1bSJed Brown suffix: 4 956c4762a1bSJed Brown nsize: 6 957c4762a1bSJed Brown args: -m 9 -n 12 -test_trans -nr 2 -nc 7 958c4762a1bSJed Brown 959c4762a1bSJed Brown test: 960c4762a1bSJed Brown suffix: 5 961c4762a1bSJed Brown nsize: 6 962c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 3 -nc 1 963c4762a1bSJed Brown 964c4762a1bSJed Brown test: 965c4762a1bSJed Brown suffix: 6 966c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap 967c4762a1bSJed Brown 968c4762a1bSJed Brown test: 969c4762a1bSJed Brown suffix: 7 970c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap -permmap 971c4762a1bSJed Brown 972c4762a1bSJed Brown test: 973c4762a1bSJed Brown suffix: 8 974c4762a1bSJed Brown args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -permmap 975c4762a1bSJed Brown 976c4762a1bSJed Brown test: 977c4762a1bSJed Brown suffix: 9 978c4762a1bSJed Brown nsize: 5 979c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap 980c4762a1bSJed Brown 981c4762a1bSJed Brown test: 982c4762a1bSJed Brown suffix: 10 983c4762a1bSJed Brown nsize: 5 984c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap -permmap 985c4762a1bSJed Brown 986c20d7725SJed Brown test: 987c20d7725SJed Brown suffix: vscat_default 988c4762a1bSJed Brown nsize: 5 989c4762a1bSJed Brown args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -permmap 990c4762a1bSJed Brown output_file: output/ex23_11.out 991c4762a1bSJed Brown 992c4762a1bSJed Brown test: 993c4762a1bSJed Brown suffix: 12 994c4762a1bSJed Brown nsize: 3 995c4762a1bSJed Brown args: -m 12 -n 12 -symmetric -matis_localmat_type sbaij -test_trans -nr 2 -nc 3 996c4762a1bSJed Brown 997c4762a1bSJed Brown testset: 998c4762a1bSJed Brown output_file: output/ex23_13.out 999c4762a1bSJed Brown nsize: 3 1000c4762a1bSJed Brown args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -diffmap -permmap 1001c4762a1bSJed Brown filter: grep -v "type:" 1002c4762a1bSJed Brown test: 1003c4762a1bSJed Brown suffix: baij 1004c4762a1bSJed Brown args: -matis_localmat_type baij 1005c4762a1bSJed Brown test: 1006c4762a1bSJed Brown requires: viennacl 1007c4762a1bSJed Brown suffix: viennacl 1008c4762a1bSJed Brown args: -matis_localmat_type aijviennacl 1009c4762a1bSJed Brown test: 1010c4762a1bSJed Brown requires: cuda 1011c4762a1bSJed Brown suffix: cusparse 1012c4762a1bSJed Brown args: -matis_localmat_type aijcusparse 1013c4762a1bSJed Brown 1014e432b41dSStefano Zampini test: 1015e432b41dSStefano Zampini suffix: negrep 1016e432b41dSStefano Zampini nsize: {{1 3}separate output} 1017e432b41dSStefano Zampini args: -m {{5 7}separate output} -n {{5 7}separate output} -test_trans -nr 2 -nc 3 -negmap {{0 1}separate output} -repmap {{0 1}separate output} -permmap -diffmap {{0 1}separate output} 1018e432b41dSStefano Zampini 1019c4762a1bSJed Brown TEST*/ 1020