1c4762a1bSJed Brown 2d0dbe9f7SStefano 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; 13d0dbe9f7SStefano 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; 19d0dbe9f7SStefano Zampini IS irow[2], icol[2]; 20d0dbe9f7SStefano Zampini PetscLayout rlayout, clayout; 21d0dbe9f7SStefano 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; 26d0dbe9f7SStefano 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 31*327415f7SBarry Smith PetscFunctionBeginUser; 329566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 349566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 35c4762a1bSJed Brown m = n = 2*size; 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-symmetric",&symmetric,NULL)); 379566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 389566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-negmap",&negmap,NULL)); 409566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-repmap",&repmap,NULL)); 419566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-permmap",&permute,NULL)); 429566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-diffmap",&diffmap,NULL)); 43e00437b9SBarry 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"); 44e00437b9SBarry 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"); 45e00437b9SBarry Smith PetscCheck(n >= 2,PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Number of cols should be larger or equal 2"); 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* create a MATIS matrix */ 489566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 499566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n)); 509566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATIS)); 519566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 52e432b41dSStefano Zampini if (!negmap && !repmap) { 53fc989267SStefano Zampini /* This is not the proper setting for MATIS for finite elements, it is just used to test the routines 54fc989267SStefano Zampini Here we use a one-to-one correspondence between local row/column spaces and global row/column spaces 55fc989267SStefano Zampini Equivalent to passing NULL for the mapping */ 569566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n,0,1,&is)); 57e432b41dSStefano Zampini } else if (negmap && !repmap) { /* non repeated but with negative indices */ 589566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n+2,-2,1,&is)); 59e432b41dSStefano Zampini } else if (!negmap && repmap) { /* non negative but repeated indices */ 60e432b41dSStefano Zampini IS isl[2]; 61e432b41dSStefano Zampini 629566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n,0,1,&isl[0])); 639566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n,n-1,-1,&isl[1])); 649566063dSJacob Faibussowitsch PetscCall(ISConcatenate(PETSC_COMM_WORLD,2,isl,&is)); 659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl[0])); 669566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl[1])); 67e432b41dSStefano Zampini } else { /* negative and repeated indices */ 68e432b41dSStefano Zampini IS isl[2]; 69e432b41dSStefano Zampini 709566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n+1,-1,1,&isl[0])); 719566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n+1,n-1,-1,&isl[1])); 729566063dSJacob Faibussowitsch PetscCall(ISConcatenate(PETSC_COMM_WORLD,2,isl,&is)); 739566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl[0])); 749566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl[1])); 75e432b41dSStefano Zampini } 769566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(is,&cmap)); 779566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 78c4762a1bSJed Brown 79e432b41dSStefano Zampini if (m != n || diffmap) { 809566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,m,permute ? m-1 : 0,permute ? -1 : 1,&is)); 819566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(is,&rmap)); 829566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 83c4762a1bSJed Brown } else { 849566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)cmap)); 85c4762a1bSJed Brown rmap = cmap; 86c4762a1bSJed Brown } 87e432b41dSStefano Zampini 889566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(A,rmap,cmap)); 899566063dSJacob Faibussowitsch PetscCall(MatISStoreL2L(A,PETSC_FALSE)); 909566063dSJacob Faibussowitsch PetscCall(MatISSetPreallocation(A,3,NULL,3,NULL)); 919566063dSJacob Faibussowitsch PetscCall(MatSetOption(A,MAT_NEW_NONZERO_ALLOCATION_ERR,(PetscBool)!(repmap || negmap))); /* I do not want to precompute the pattern */ 929566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(rmap,&lm)); 939566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(cmap,&ln)); 94e432b41dSStefano Zampini for (i=0; i<lm; i++) { 95c4762a1bSJed Brown PetscScalar v[3]; 96c4762a1bSJed Brown PetscInt cols[3]; 97c4762a1bSJed Brown 98c4762a1bSJed Brown cols[0] = (i-1+n)%n; 99c4762a1bSJed Brown cols[1] = i%n; 100c4762a1bSJed Brown cols[2] = (i+1)%n; 101c4762a1bSJed Brown v[0] = -1.*(symmetric ? PetscMin(i+1,cols[0]+1) : i+1); 102c4762a1bSJed Brown v[1] = 2.*(symmetric ? PetscMin(i+1,cols[1]+1) : i+1); 103c4762a1bSJed Brown v[2] = -1.*(symmetric ? PetscMin(i+1,cols[2]+1) : i+1); 1049566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(cmap,IS_GTOLM_MASK,3,cols,NULL,cols)); 1059566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(A,1,&i,3,cols,v,ADD_VALUES)); 106c4762a1bSJed Brown } 1079566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1089566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 109c4762a1bSJed Brown 110e432b41dSStefano Zampini /* activate tests for square matrices with same maps only */ 1119566063dSJacob Faibussowitsch PetscCall(MatHasCongruentLayouts(A,&squaretest)); 112e432b41dSStefano Zampini if (squaretest && rmap != cmap) { 113e432b41dSStefano Zampini PetscInt nr, nc; 114e432b41dSStefano Zampini 1159566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(rmap,&nr)); 1169566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetSize(cmap,&nc)); 117e432b41dSStefano Zampini if (nr != nc) squaretest = PETSC_FALSE; 118e432b41dSStefano Zampini else { 119e432b41dSStefano Zampini const PetscInt *idxs1,*idxs2; 120e432b41dSStefano Zampini 1219566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(rmap,&idxs1)); 1229566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingGetIndices(cmap,&idxs2)); 1239566063dSJacob Faibussowitsch PetscCall(PetscArraycmp(idxs1,idxs2,nr,&squaretest)); 1249566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(rmap,&idxs1)); 1259566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingRestoreIndices(cmap,&idxs2)); 126e432b41dSStefano Zampini } 1271c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(MPI_IN_PLACE,&squaretest,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A))); 128e432b41dSStefano Zampini } 129e432b41dSStefano Zampini 130e432b41dSStefano Zampini /* test MatISGetLocalMat */ 1319566063dSJacob Faibussowitsch PetscCall(MatISGetLocalMat(A,&B)); 1329566063dSJacob Faibussowitsch PetscCall(MatGetType(B,&lmtype)); 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* test MatGetInfo */ 1359566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatGetInfo\n")); 1369566063dSJacob Faibussowitsch PetscCall(MatGetInfo(A,MAT_LOCAL,&info)); 1379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 138d0609cedSBarry 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, 139d0609cedSBarry Smith (PetscInt)info.nz_allocated,(PetscInt)info.nz_unneeded,(PetscInt)info.assemblies,(PetscInt)info.mallocs)); 1409566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 1419566063dSJacob Faibussowitsch PetscCall(MatGetInfo(A,MAT_GLOBAL_MAX,&info)); 142d0609cedSBarry Smith PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"GlobalMax : %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",(PetscInt)info.nz_used, 143d0609cedSBarry Smith (PetscInt)info.nz_allocated,(PetscInt)info.nz_unneeded,(PetscInt)info.assemblies,(PetscInt)info.mallocs)); 1449566063dSJacob Faibussowitsch PetscCall(MatGetInfo(A,MAT_GLOBAL_SUM,&info)); 145d0609cedSBarry Smith PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"GlobalSum : %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",(PetscInt)info.nz_used, 146d0609cedSBarry Smith (PetscInt)info.nz_allocated,(PetscInt)info.nz_unneeded,(PetscInt)info.assemblies,(PetscInt)info.mallocs)); 147c4762a1bSJed Brown 148e432b41dSStefano Zampini /* test MatIsSymmetric */ 1499566063dSJacob Faibussowitsch PetscCall(MatIsSymmetric(A,0.0,&issymmetric)); 1509566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatIsSymmetric: %d\n",issymmetric)); 151c4762a1bSJed Brown 152c4762a1bSJed Brown /* Create a MPIAIJ matrix, same as A */ 1539566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); 1549566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n)); 1559566063dSJacob Faibussowitsch PetscCall(MatSetType(B,MATAIJ)); 1569566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 1579566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 1589566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(B,rmap,cmap)); 1599566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(B,3,NULL,3,NULL)); 1609566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(B,1,3,NULL,3,NULL)); 161c4762a1bSJed Brown #if defined(PETSC_HAVE_HYPRE) 1629566063dSJacob Faibussowitsch PetscCall(MatHYPRESetPreallocation(B,3,NULL,3,NULL)); 163c4762a1bSJed Brown #endif 1649566063dSJacob Faibussowitsch PetscCall(MatISSetPreallocation(B,3,NULL,3,NULL)); 1659566063dSJacob Faibussowitsch PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,(PetscBool)!(repmap || negmap))); /* I do not want to precompute the pattern */ 166e432b41dSStefano Zampini for (i=0; i<lm; i++) { 167c4762a1bSJed Brown PetscScalar v[3]; 168c4762a1bSJed Brown PetscInt cols[3]; 169c4762a1bSJed Brown 170c4762a1bSJed Brown cols[0] = (i-1+n)%n; 171c4762a1bSJed Brown cols[1] = i%n; 172c4762a1bSJed Brown cols[2] = (i+1)%n; 173c4762a1bSJed Brown v[0] = -1.*(symmetric ? PetscMin(i+1,cols[0]+1) : i+1); 174c4762a1bSJed Brown v[1] = 2.*(symmetric ? PetscMin(i+1,cols[1]+1) : i+1); 175c4762a1bSJed Brown v[2] = -1.*(symmetric ? PetscMin(i+1,cols[2]+1) : i+1); 1769566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(cmap,IS_GTOLM_MASK,3,cols,NULL,cols)); 1779566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(B,1,&i,3,cols,v,ADD_VALUES)); 178c4762a1bSJed Brown } 1799566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1809566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 181e432b41dSStefano Zampini 182e432b41dSStefano Zampini /* test MatView */ 1839566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatView\n")); 1849566063dSJacob Faibussowitsch PetscCall(MatView(A,NULL)); 1859566063dSJacob Faibussowitsch PetscCall(MatView(B,NULL)); 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* test CheckMat */ 1889566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test CheckMat\n")); 1899566063dSJacob Faibussowitsch PetscCall(CheckMat(A,B,PETSC_FALSE,"CheckMat")); 190c4762a1bSJed Brown 191c4762a1bSJed Brown /* test MatDuplicate and MatAXPY */ 1929566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatDuplicate and MatAXPY\n")); 1939566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 1949566063dSJacob Faibussowitsch PetscCall(CheckMat(A,A2,PETSC_FALSE,"MatDuplicate and MatAXPY")); 195c4762a1bSJed Brown 196c4762a1bSJed Brown /* test MatConvert */ 1979566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatConvert_IS_XAIJ\n")); 1989566063dSJacob Faibussowitsch PetscCall(MatConvert(A2,MATAIJ,MAT_INITIAL_MATRIX,&B2)); 1999566063dSJacob Faibussowitsch PetscCall(CheckMat(B,B2,PETSC_TRUE,"MatConvert_IS_XAIJ MAT_INITIAL_MATRIX")); 2009566063dSJacob Faibussowitsch PetscCall(MatConvert(A2,MATAIJ,MAT_REUSE_MATRIX,&B2)); 2019566063dSJacob Faibussowitsch PetscCall(CheckMat(B,B2,PETSC_TRUE,"MatConvert_IS_XAIJ MAT_REUSE_MATRIX")); 2029566063dSJacob Faibussowitsch PetscCall(MatConvert(A2,MATAIJ,MAT_INPLACE_MATRIX,&A2)); 2039566063dSJacob Faibussowitsch PetscCall(CheckMat(B,A2,PETSC_TRUE,"MatConvert_IS_XAIJ MAT_INPLACE_MATRIX")); 2049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 2059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 2069566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatConvert_XAIJ_IS\n")); 2079566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&B2)); 2089566063dSJacob Faibussowitsch PetscCall(MatConvert(B2,MATIS,MAT_INITIAL_MATRIX,&A2)); 2099566063dSJacob Faibussowitsch PetscCall(CheckMat(A,A2,PETSC_TRUE,"MatConvert_XAIJ_IS MAT_INITIAL_MATRIX")); 2109566063dSJacob Faibussowitsch PetscCall(MatConvert(B2,MATIS,MAT_REUSE_MATRIX,&A2)); 2119566063dSJacob Faibussowitsch PetscCall(CheckMat(A,A2,PETSC_TRUE,"MatConvert_XAIJ_IS MAT_REUSE_MATRIX")); 2129566063dSJacob Faibussowitsch PetscCall(MatConvert(B2,MATIS,MAT_INPLACE_MATRIX,&B2)); 2139566063dSJacob Faibussowitsch PetscCall(CheckMat(A,B2,PETSC_TRUE,"MatConvert_XAIJ_IS MAT_INPLACE_MATRIX")); 2149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 2159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 2169566063dSJacob Faibussowitsch PetscCall(PetscStrcmp(lmtype,MATSEQAIJ,&isaij)); 217c4762a1bSJed Brown if (size == 1 && isaij) { /* tests special code paths in MatConvert_IS_XAIJ */ 218c4762a1bSJed 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}; 219e432b41dSStefano Zampini ISLocalToGlobalMapping tcmap,trmap; 220c4762a1bSJed Brown 221c4762a1bSJed Brown for (ri = 0; ri < 2; ri++) { 222c4762a1bSJed Brown PetscInt *r; 223c4762a1bSJed Brown 224c4762a1bSJed Brown r = (PetscInt*)(ri == 0 ? rr : rk); 225c4762a1bSJed Brown for (ci = 0; ci < 2; ci++) { 226c4762a1bSJed Brown PetscInt *c,rb,cb; 227c4762a1bSJed Brown 228c4762a1bSJed Brown c = (PetscInt*)(ci == 0 ? cr : ck); 229c4762a1bSJed Brown for (rb = 1; rb < 4; rb++) { 2309566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF,rb,3,r,PETSC_COPY_VALUES,&is)); 2319566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(is,&trmap)); 2329566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 233c4762a1bSJed Brown for (cb = 1; cb < 4; cb++) { 234c4762a1bSJed Brown Mat T,lT,T2; 235c4762a1bSJed Brown char testname[256]; 236c4762a1bSJed Brown 2379566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(testname,sizeof(testname),"MatConvert_IS_XAIJ special case (%" PetscInt_FMT " %" PetscInt_FMT ", bs %" PetscInt_FMT " %" PetscInt_FMT ")",ri,ci,rb,cb)); 2389566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test %s\n",testname)); 239c4762a1bSJed Brown 2409566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_SELF,cb,4,c,PETSC_COPY_VALUES,&is)); 2419566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(is,&tcmap)); 2429566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 243c4762a1bSJed Brown 2449566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_SELF,&T)); 2459566063dSJacob Faibussowitsch PetscCall(MatSetSizes(T,PETSC_DECIDE,PETSC_DECIDE,rb*3,cb*4)); 2469566063dSJacob Faibussowitsch PetscCall(MatSetType(T,MATIS)); 2479566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(T,trmap,tcmap)); 2489566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&tcmap)); 2499566063dSJacob Faibussowitsch PetscCall(MatISGetLocalMat(T,&lT)); 2509566063dSJacob Faibussowitsch PetscCall(MatSetType(lT,MATSEQAIJ)); 2519566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(lT,cb*4,NULL)); 2529566063dSJacob Faibussowitsch PetscCall(MatSetRandom(lT,NULL)); 2539566063dSJacob Faibussowitsch PetscCall(MatConvert(lT,lmtype,MAT_INPLACE_MATRIX,&lT)); 2549566063dSJacob Faibussowitsch PetscCall(MatISRestoreLocalMat(T,&lT)); 2559566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY)); 2569566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY)); 257c4762a1bSJed Brown 2589566063dSJacob Faibussowitsch PetscCall(MatConvert(T,MATAIJ,MAT_INITIAL_MATRIX,&T2)); 2599566063dSJacob Faibussowitsch PetscCall(CheckMat(T,T2,PETSC_TRUE,"MAT_INITIAL_MATRIX")); 2609566063dSJacob Faibussowitsch PetscCall(MatConvert(T,MATAIJ,MAT_REUSE_MATRIX,&T2)); 2619566063dSJacob Faibussowitsch PetscCall(CheckMat(T,T2,PETSC_TRUE,"MAT_REUSE_MATRIX")); 2629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T2)); 2639566063dSJacob Faibussowitsch PetscCall(MatDuplicate(T,MAT_COPY_VALUES,&T2)); 2649566063dSJacob Faibussowitsch PetscCall(MatConvert(T2,MATAIJ,MAT_INPLACE_MATRIX,&T2)); 2659566063dSJacob Faibussowitsch PetscCall(CheckMat(T,T2,PETSC_TRUE,"MAT_INPLACE_MATRIX")); 2669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T)); 2679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T2)); 268c4762a1bSJed Brown } 2699566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&trmap)); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown } 272c4762a1bSJed Brown } 273c4762a1bSJed Brown } 274c4762a1bSJed Brown 275c4762a1bSJed Brown /* test MatDiagonalScale */ 2769566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatDiagonalScale\n")); 2779566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 2789566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&B2)); 2799566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&x,&y)); 2809566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,NULL)); 281e432b41dSStefano Zampini if (issymmetric) { 2829566063dSJacob Faibussowitsch PetscCall(VecCopy(x,y)); 283c4762a1bSJed Brown } else { 2849566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y,NULL)); 2859566063dSJacob Faibussowitsch PetscCall(VecScale(y,8.)); 286c4762a1bSJed Brown } 2879566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A2,y,x)); 2889566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(B2,y,x)); 2899566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatDiagonalScale")); 2909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 2919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 2929566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 294c4762a1bSJed Brown 295c4762a1bSJed Brown /* test MatPtAP (A IS and B AIJ) */ 296c4762a1bSJed Brown if (isaij && m == n) { 2979566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatPtAP\n")); 2989566063dSJacob Faibussowitsch PetscCall(MatISStoreL2L(A,PETSC_TRUE)); 2999566063dSJacob Faibussowitsch PetscCall(MatPtAP(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&A2)); 3009566063dSJacob Faibussowitsch PetscCall(MatPtAP(B,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B2)); 3019566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatPtAP MAT_INITIAL_MATRIX")); 3029566063dSJacob Faibussowitsch PetscCall(MatPtAP(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&A2)); 3039566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatPtAP MAT_REUSE_MATRIX")); 3049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 3059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 306c4762a1bSJed Brown } 307c4762a1bSJed Brown 308c4762a1bSJed Brown /* test MatGetLocalSubMatrix */ 3099566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatGetLocalSubMatrix\n")); 3109566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&A2)); 3119566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,lm/2+lm%2,0,2,&reven)); 3129566063dSJacob Faibussowitsch PetscCall(ISComplement(reven,0,lm,&rodd)); 3139566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,ln/2+ln%2,0,2,&ceven)); 3149566063dSJacob Faibussowitsch PetscCall(ISComplement(ceven,0,ln,&codd)); 3159566063dSJacob Faibussowitsch PetscCall(MatGetLocalSubMatrix(A2,reven,ceven,&Aee)); 3169566063dSJacob Faibussowitsch PetscCall(MatGetLocalSubMatrix(A2,reven,codd,&Aeo)); 3179566063dSJacob Faibussowitsch PetscCall(MatGetLocalSubMatrix(A2,rodd,ceven,&Aoe)); 3189566063dSJacob Faibussowitsch PetscCall(MatGetLocalSubMatrix(A2,rodd,codd,&Aoo)); 319e432b41dSStefano Zampini for (i=0; i<lm; i++) { 320c4762a1bSJed Brown PetscInt j,je,jo,colse[3], colso[3]; 321c4762a1bSJed Brown PetscScalar ve[3], vo[3]; 322c4762a1bSJed Brown PetscScalar v[3]; 323c4762a1bSJed Brown PetscInt cols[3]; 324e432b41dSStefano Zampini PetscInt row = i/2; 325c4762a1bSJed Brown 326c4762a1bSJed Brown cols[0] = (i-1+n)%n; 327c4762a1bSJed Brown cols[1] = i%n; 328c4762a1bSJed Brown cols[2] = (i+1)%n; 329c4762a1bSJed Brown v[0] = -1.*(symmetric ? PetscMin(i+1,cols[0]+1) : i+1); 330c4762a1bSJed Brown v[1] = 2.*(symmetric ? PetscMin(i+1,cols[1]+1) : i+1); 331c4762a1bSJed Brown v[2] = -1.*(symmetric ? PetscMin(i+1,cols[2]+1) : i+1); 3329566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(cmap,IS_GTOLM_MASK,3,cols,NULL,cols)); 333c4762a1bSJed Brown for (j=0,je=0,jo=0;j<3;j++) { 334c4762a1bSJed Brown if (cols[j]%2) { 335c4762a1bSJed Brown vo[jo] = v[j]; 336c4762a1bSJed Brown colso[jo++] = cols[j]/2; 337c4762a1bSJed Brown } else { 338c4762a1bSJed Brown ve[je] = v[j]; 339c4762a1bSJed Brown colse[je++] = cols[j]/2; 340c4762a1bSJed Brown } 341c4762a1bSJed Brown } 342c4762a1bSJed Brown if (i%2) { 3439566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(Aoe,1,&row,je,colse,ve,ADD_VALUES)); 3449566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(Aoo,1,&row,jo,colso,vo,ADD_VALUES)); 345c4762a1bSJed Brown } else { 3469566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(Aee,1,&row,je,colse,ve,ADD_VALUES)); 3479566063dSJacob Faibussowitsch PetscCall(MatSetValuesLocal(Aeo,1,&row,jo,colso,vo,ADD_VALUES)); 348c4762a1bSJed Brown } 349c4762a1bSJed Brown } 3509566063dSJacob Faibussowitsch PetscCall(MatRestoreLocalSubMatrix(A2,reven,ceven,&Aee)); 3519566063dSJacob Faibussowitsch PetscCall(MatRestoreLocalSubMatrix(A2,reven,codd,&Aeo)); 3529566063dSJacob Faibussowitsch PetscCall(MatRestoreLocalSubMatrix(A2,rodd,ceven,&Aoe)); 3539566063dSJacob Faibussowitsch PetscCall(MatRestoreLocalSubMatrix(A2,rodd,codd,&Aoo)); 3549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&reven)); 3559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ceven)); 3569566063dSJacob Faibussowitsch PetscCall(ISDestroy(&rodd)); 3579566063dSJacob Faibussowitsch PetscCall(ISDestroy(&codd)); 3589566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY)); 3599566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY)); 3609566063dSJacob Faibussowitsch PetscCall(MatAXPY(A2,-1.,A,SAME_NONZERO_PATTERN)); 3619566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,NULL,PETSC_FALSE,"MatGetLocalSubMatrix")); 3629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 363c4762a1bSJed Brown 364c4762a1bSJed Brown /* test MatConvert_Nest_IS */ 365c4762a1bSJed Brown testT = PETSC_FALSE; 3669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_trans",&testT,NULL)); 367c4762a1bSJed Brown 3689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatConvert_Nest_IS\n")); 369c4762a1bSJed Brown nr = 2; 370c4762a1bSJed Brown nc = 2; 3719566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nr",&nr,NULL)); 3729566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-nc",&nc,NULL)); 373c4762a1bSJed Brown if (testT) { 3749566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&cst,&cen)); 3759566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(A,&rst,&ren)); 376c4762a1bSJed Brown } else { 3779566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rst,&ren)); 3789566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(A,&cst,&cen)); 379c4762a1bSJed Brown } 3809566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(nr,&rows,nc,&cols,2*nr*nc,&mats)); 381c4762a1bSJed Brown for (i=0;i<nr*nc;i++) { 382c4762a1bSJed Brown if (testT) { 3839566063dSJacob Faibussowitsch PetscCall(MatCreateTranspose(A,&mats[i])); 3849566063dSJacob Faibussowitsch PetscCall(MatTranspose(B,MAT_INITIAL_MATRIX,&mats[i+nr*nc])); 385c4762a1bSJed Brown } else { 3869566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&mats[i])); 3879566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&mats[i+nr*nc])); 388c4762a1bSJed Brown } 389c4762a1bSJed Brown } 390c4762a1bSJed Brown for (i=0;i<nr;i++) { 3919566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,ren-rst,i+rst,nr,&rows[i])); 392c4762a1bSJed Brown } 393c4762a1bSJed Brown for (i=0;i<nc;i++) { 3949566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,cen-cst,i+cst,nc,&cols[i])); 395c4762a1bSJed Brown } 3969566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD,nr,rows,nc,cols,mats,&A2)); 3979566063dSJacob Faibussowitsch PetscCall(MatCreateNest(PETSC_COMM_WORLD,nr,rows,nc,cols,mats+nr*nc,&B2)); 398c4762a1bSJed Brown for (i=0;i<nr;i++) { 3999566063dSJacob Faibussowitsch PetscCall(ISDestroy(&rows[i])); 400c4762a1bSJed Brown } 401c4762a1bSJed Brown for (i=0;i<nc;i++) { 4029566063dSJacob Faibussowitsch PetscCall(ISDestroy(&cols[i])); 403c4762a1bSJed Brown } 404c4762a1bSJed Brown for (i=0;i<2*nr*nc;i++) { 4059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mats[i])); 406c4762a1bSJed Brown } 4079566063dSJacob Faibussowitsch PetscCall(PetscFree3(rows,cols,mats)); 4089566063dSJacob Faibussowitsch PetscCall(MatConvert(B2,MATAIJ,MAT_INITIAL_MATRIX,&T)); 4099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 4109566063dSJacob Faibussowitsch PetscCall(MatConvert(A2,MATIS,MAT_INITIAL_MATRIX,&B2)); 4119566063dSJacob Faibussowitsch PetscCall(CheckMat(B2,T,PETSC_TRUE,"MatConvert_Nest_IS MAT_INITIAL_MATRIX")); 4129566063dSJacob Faibussowitsch PetscCall(MatConvert(A2,MATIS,MAT_REUSE_MATRIX,&B2)); 4139566063dSJacob Faibussowitsch PetscCall(CheckMat(B2,T,PETSC_TRUE,"MatConvert_Nest_IS MAT_REUSE_MATRIX")); 4149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 4159566063dSJacob Faibussowitsch PetscCall(MatConvert(A2,MATIS,MAT_INPLACE_MATRIX,&A2)); 4169566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,T,PETSC_TRUE,"MatConvert_Nest_IS MAT_INPLACE_MATRIX")); 4179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&T)); 4189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 419c4762a1bSJed Brown 420c4762a1bSJed Brown /* test MatCreateSubMatrix */ 4219566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatCreateSubMatrix\n")); 422dd400576SPatrick Sanan if (rank == 0) { 4239566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,1,1,1,&is)); 4249566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,2,0,1,&is2)); 425c4762a1bSJed Brown } else if (rank == 1) { 4269566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,1,0,1,&is)); 427c4762a1bSJed Brown if (n > 3) { 4289566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,1,3,1,&is2)); 429c4762a1bSJed Brown } else { 4309566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is2)); 431c4762a1bSJed Brown } 432c4762a1bSJed Brown } else if (rank == 2 && n > 4) { 4339566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is)); 4349566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,n-4,4,1,&is2)); 435c4762a1bSJed Brown } else { 4369566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is)); 4379566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is2)); 438c4762a1bSJed Brown } 4399566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&A2)); 4409566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(B,is,is,MAT_INITIAL_MATRIX,&B2)); 4419566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_TRUE,"first MatCreateSubMatrix")); 442c4762a1bSJed Brown 4439566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A,is,is,MAT_REUSE_MATRIX,&A2)); 4449566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(B,is,is,MAT_REUSE_MATRIX,&B2)); 4459566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"reuse MatCreateSubMatrix")); 4469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 4479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 448c4762a1bSJed Brown 449e432b41dSStefano Zampini if (!issymmetric) { 4509566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A,is,is2,MAT_INITIAL_MATRIX,&A2)); 4519566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(B,is,is2,MAT_INITIAL_MATRIX,&B2)); 4529566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(A,is,is2,MAT_REUSE_MATRIX,&A2)); 4539566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(B,is,is2,MAT_REUSE_MATRIX,&B2)); 4549566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"second MatCreateSubMatrix")); 455c4762a1bSJed Brown } 456c4762a1bSJed Brown 4579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 4589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 4599566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 4609566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is2)); 461c4762a1bSJed Brown 462d0dbe9f7SStefano Zampini /* test MatCreateSubMatrices */ 463d0dbe9f7SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatCreateSubMatrices\n")); 464d0dbe9f7SStefano Zampini PetscCall(MatGetLayouts(A,&rlayout,&clayout)); 465d0dbe9f7SStefano Zampini PetscCall(PetscLayoutGetRanges(rlayout,&rrange)); 466d0dbe9f7SStefano Zampini PetscCall(PetscLayoutGetRanges(clayout,&crange)); 467d0dbe9f7SStefano Zampini lrank = (size + rank - 1)%size; 468d0dbe9f7SStefano Zampini rrank = (rank + 1)%size; 469d0dbe9f7SStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF,(rrange[lrank+1] - rrange[lrank]),rrange[lrank],1,&irow[0])); 470d0dbe9f7SStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF,(crange[rrank+1] - crange[rrank]),crange[rrank],1,&icol[0])); 471d0dbe9f7SStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF,(rrange[rrank+1] - rrange[rrank]),rrange[rrank],1,&irow[1])); 472d0dbe9f7SStefano Zampini PetscCall(ISCreateStride(PETSC_COMM_SELF,(crange[lrank+1] - crange[lrank]),crange[lrank],1,&icol[1])); 473d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(A,2,irow,icol,MAT_INITIAL_MATRIX,&Asub)); 474d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(B,2,irow,icol,MAT_INITIAL_MATRIX,&Bsub)); 475d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[0],Bsub[0],PETSC_FALSE,"MatCreateSubMatrices[0]")); 476d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[1],Bsub[1],PETSC_FALSE,"MatCreateSubMatrices[1]")); 477d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(A,2,irow,icol,MAT_REUSE_MATRIX,&Asub)); 478d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(B,2,irow,icol,MAT_REUSE_MATRIX,&Bsub)); 479d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[0],Bsub[0],PETSC_FALSE,"MatCreateSubMatrices[0]")); 480d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[1],Bsub[1],PETSC_FALSE,"MatCreateSubMatrices[1]")); 481d0dbe9f7SStefano Zampini PetscCall(MatDestroySubMatrices(2,&Asub)); 482d0dbe9f7SStefano Zampini PetscCall(MatDestroySubMatrices(2,&Bsub)); 483d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&irow[0])); 484d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&irow[1])); 485d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&icol[0])); 486d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&icol[1])); 487d0dbe9f7SStefano Zampini 488c4762a1bSJed Brown /* Create an IS required by MatZeroRows(): just rank zero provides the rows to be eliminated */ 489c4762a1bSJed Brown if (size > 1) { 490dd400576SPatrick Sanan if (rank == 0) { 491c4762a1bSJed Brown PetscInt st,len; 492c4762a1bSJed Brown 493c4762a1bSJed Brown st = (m+1)/2; 494c4762a1bSJed Brown len = PetscMin(m/2,PetscMax(m-(m+1)/2-1,0)); 4959566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,len,st,1,&is)); 496c4762a1bSJed Brown } else { 4979566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is)); 498c4762a1bSJed Brown } 499c4762a1bSJed Brown } else { 5009566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,1,0,1,&is)); 501c4762a1bSJed Brown } 502c4762a1bSJed Brown 503c4762a1bSJed Brown if (squaretest) { /* tests for square matrices only, with same maps for rows and columns */ 504d0dbe9f7SStefano Zampini PetscInt *idx0, *idx1, n0, n1; 505d0dbe9f7SStefano Zampini IS Ais[2], Bis[2]; 506d0dbe9f7SStefano Zampini 507c4762a1bSJed Brown /* test MatDiagonalSet */ 5089566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatDiagonalSet\n")); 5099566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 5109566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&B2)); 5119566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,NULL,&x)); 5129566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,NULL)); 5139566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(A2,x,INSERT_VALUES)); 5149566063dSJacob Faibussowitsch PetscCall(MatDiagonalSet(B2,x,INSERT_VALUES)); 5159566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatDiagonalSet")); 5169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 5179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 5189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 519c4762a1bSJed Brown 520c4762a1bSJed Brown /* test MatShift (MatShift_IS internally uses MatDiagonalSet_IS with ADD_VALUES) */ 5219566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatShift\n")); 5229566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 5239566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&B2)); 5249566063dSJacob Faibussowitsch PetscCall(MatShift(A2,2.0)); 5259566063dSJacob Faibussowitsch PetscCall(MatShift(B2,2.0)); 5269566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatShift")); 5279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 5289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 529c4762a1bSJed Brown 530c4762a1bSJed Brown /* nonzero diag value is supported for square matrices only */ 5319566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows(A,B,PETSC_TRUE,is,diag)); 532d0dbe9f7SStefano Zampini 533d0dbe9f7SStefano Zampini /* test MatIncreaseOverlap */ 534d0dbe9f7SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatIncreaseOverlap\n")); 535d0dbe9f7SStefano Zampini PetscCall(MatGetOwnershipRange(A,&rst,&ren)); 536d0dbe9f7SStefano Zampini n0 = (ren - rst)/2; 537d0dbe9f7SStefano Zampini n1 = (ren - rst)/3; 538d0dbe9f7SStefano Zampini PetscCall(PetscMalloc1(n0,&idx0)); 539d0dbe9f7SStefano Zampini PetscCall(PetscMalloc1(n1,&idx1)); 540d0dbe9f7SStefano Zampini for (i = 0; i < n0; i++) idx0[i] = ren - i - 1; 541d0dbe9f7SStefano Zampini for (i = 0; i < n1; i++) idx1[i] = rst + i; 542d0dbe9f7SStefano Zampini PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,n0,idx0,PETSC_OWN_POINTER,&Ais[0])); 543d0dbe9f7SStefano Zampini PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,n1,idx1,PETSC_OWN_POINTER,&Ais[1])); 544d0dbe9f7SStefano Zampini PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,n0,idx0,PETSC_COPY_VALUES,&Bis[0])); 545d0dbe9f7SStefano Zampini PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,n1,idx1,PETSC_COPY_VALUES,&Bis[1])); 546d0dbe9f7SStefano Zampini PetscCall(MatIncreaseOverlap(A,2,Ais,3)); 547d0dbe9f7SStefano Zampini PetscCall(MatIncreaseOverlap(B,2,Bis,3)); 548d0dbe9f7SStefano Zampini /* Non deterministic output! */ 549d0dbe9f7SStefano Zampini PetscCall(ISSort(Ais[0])); 550d0dbe9f7SStefano Zampini PetscCall(ISSort(Ais[1])); 551d0dbe9f7SStefano Zampini PetscCall(ISSort(Bis[0])); 552d0dbe9f7SStefano Zampini PetscCall(ISSort(Bis[1])); 553d0dbe9f7SStefano Zampini PetscCall(ISView(Ais[0],NULL)); 554d0dbe9f7SStefano Zampini PetscCall(ISView(Bis[0],NULL)); 555d0dbe9f7SStefano Zampini PetscCall(ISView(Ais[1],NULL)); 556d0dbe9f7SStefano Zampini PetscCall(ISView(Bis[1],NULL)); 557d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(A,2,Ais,Ais,MAT_INITIAL_MATRIX,&Asub)); 558d0dbe9f7SStefano Zampini PetscCall(MatCreateSubMatrices(B,2,Bis,Ais,MAT_INITIAL_MATRIX,&Bsub)); 559d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[0],Bsub[0],PETSC_FALSE,"MatIncreaseOverlap[0]")); 560d0dbe9f7SStefano Zampini PetscCall(CheckMat(Asub[1],Bsub[1],PETSC_FALSE,"MatIncreaseOverlap[1]")); 561d0dbe9f7SStefano Zampini PetscCall(MatDestroySubMatrices(2,&Asub)); 562d0dbe9f7SStefano Zampini PetscCall(MatDestroySubMatrices(2,&Bsub)); 563d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&Ais[0])); 564d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&Ais[1])); 565d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&Bis[0])); 566d0dbe9f7SStefano Zampini PetscCall(ISDestroy(&Bis[1])); 567d0dbe9f7SStefano Zampini 568c4762a1bSJed Brown } 5699566063dSJacob Faibussowitsch PetscCall(TestMatZeroRows(A,B,squaretest,is,0.0)); 5709566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 571c4762a1bSJed Brown 572c4762a1bSJed Brown /* test MatTranspose */ 5739566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatTranspose\n")); 5749566063dSJacob Faibussowitsch PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX,&A2)); 5759566063dSJacob Faibussowitsch PetscCall(MatTranspose(B,MAT_INITIAL_MATRIX,&B2)); 5769566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"initial matrix MatTranspose")); 577c4762a1bSJed Brown 5789566063dSJacob Faibussowitsch PetscCall(MatTranspose(A,MAT_REUSE_MATRIX,&A2)); 5799566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"reuse matrix (not in place) MatTranspose")); 5809566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 581c4762a1bSJed Brown 5829566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 5839566063dSJacob Faibussowitsch PetscCall(MatTranspose(A2,MAT_INPLACE_MATRIX,&A2)); 5849566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"reuse matrix (in place) MatTranspose")); 5859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 586c4762a1bSJed Brown 5879566063dSJacob Faibussowitsch PetscCall(MatTranspose(A,MAT_INITIAL_MATRIX,&A2)); 5889566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"reuse matrix (different type) MatTranspose")); 5899566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 5909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 591c4762a1bSJed Brown 592c4762a1bSJed Brown /* test MatISFixLocalEmpty */ 593c4762a1bSJed Brown if (isaij) { 594c4762a1bSJed Brown PetscInt r[2]; 595c4762a1bSJed Brown 596c4762a1bSJed Brown r[0] = 0; 597c4762a1bSJed Brown r[1] = PetscMin(m,n)-1; 5989566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatISFixLocalEmpty\n")); 5999566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 600e432b41dSStefano Zampini 6019566063dSJacob Faibussowitsch PetscCall(MatISFixLocalEmpty(A2,PETSC_TRUE)); 6029566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY)); 6039566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY)); 6049566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B,PETSC_FALSE,"MatISFixLocalEmpty (null)")); 605c4762a1bSJed Brown 6069566063dSJacob Faibussowitsch PetscCall(MatZeroRows(A2,2,r,0.0,NULL,NULL)); 6079566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2,NULL,"-fixempty_view")); 6089566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&B2)); 6099566063dSJacob Faibussowitsch PetscCall(MatZeroRows(B2,2,r,0.0,NULL,NULL)); 6109566063dSJacob Faibussowitsch PetscCall(MatISFixLocalEmpty(A2,PETSC_TRUE)); 6119566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY)); 6129566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY)); 6139566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2,NULL,"-fixempty_view")); 6149566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatISFixLocalEmpty (rows)")); 6159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 616c4762a1bSJed Brown 6179566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 6189566063dSJacob Faibussowitsch PetscCall(MatZeroRows(A2,2,r,0.0,NULL,NULL)); 6199566063dSJacob Faibussowitsch PetscCall(MatTranspose(A2,MAT_INPLACE_MATRIX,&A2)); 6209566063dSJacob Faibussowitsch PetscCall(MatTranspose(B2,MAT_INPLACE_MATRIX,&B2)); 6219566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2,NULL,"-fixempty_view")); 6229566063dSJacob Faibussowitsch PetscCall(MatISFixLocalEmpty(A2,PETSC_TRUE)); 6239566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY)); 6249566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY)); 6259566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2,NULL,"-fixempty_view")); 6269566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatISFixLocalEmpty (cols)")); 627c4762a1bSJed Brown 6289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 6299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 630c4762a1bSJed Brown 631c4762a1bSJed Brown if (squaretest) { 6329566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&A2)); 6339566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&B2)); 6349566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumns(A2,2,r,0.0,NULL,NULL)); 6359566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumns(B2,2,r,0.0,NULL,NULL)); 6369566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2,NULL,"-fixempty_view")); 6379566063dSJacob Faibussowitsch PetscCall(MatISFixLocalEmpty(A2,PETSC_TRUE)); 6389566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY)); 6399566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY)); 6409566063dSJacob Faibussowitsch PetscCall(MatViewFromOptions(A2,NULL,"-fixempty_view")); 6419566063dSJacob Faibussowitsch PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatISFixLocalEmpty (rows+cols)")); 6429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A2)); 6439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 644c4762a1bSJed Brown } 645c4762a1bSJed Brown 646c4762a1bSJed Brown } 647c4762a1bSJed Brown 648c4762a1bSJed Brown /* test MatInvertBlockDiagonal 649c4762a1bSJed Brown special cases for block-diagonal matrices */ 650c4762a1bSJed Brown if (m == n) { 651c4762a1bSJed Brown ISLocalToGlobalMapping map; 652c4762a1bSJed Brown Mat Abd,Bbd; 653c4762a1bSJed Brown IS is,bis; 654c4762a1bSJed Brown const PetscScalar *isbd,*aijbd; 655c4762a1bSJed Brown PetscScalar *vals; 656c4762a1bSJed Brown const PetscInt *sts,*idxs; 657c4762a1bSJed Brown PetscInt *idxs2,diff,perm,nl,bs,st,en,in; 658c4762a1bSJed Brown PetscBool ok; 659c4762a1bSJed Brown 660c4762a1bSJed Brown for (diff = 0; diff < 3; diff++) { 661c4762a1bSJed Brown for (perm = 0; perm < 3; perm++) { 662c4762a1bSJed Brown for (bs = 1; bs < 4; bs++) { 6639566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatInvertBlockDiagonal blockdiag %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",n,diff,perm,bs)); 6649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(bs*bs,&vals)); 6659566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRanges(A,&sts)); 666c4762a1bSJed Brown switch (diff) { 667c4762a1bSJed Brown case 1: /* inverted layout by processes */ 668c4762a1bSJed Brown in = 1; 669c4762a1bSJed Brown st = sts[size - rank - 1]; 670c4762a1bSJed Brown en = sts[size - rank]; 671c4762a1bSJed Brown nl = en - st; 672c4762a1bSJed Brown break; 673c4762a1bSJed Brown case 2: /* round-robin layout */ 674c4762a1bSJed Brown in = size; 675c4762a1bSJed Brown st = rank; 676c4762a1bSJed Brown nl = n/size; 677c4762a1bSJed Brown if (rank < n%size) nl++; 678c4762a1bSJed Brown break; 679c4762a1bSJed Brown default: /* same layout */ 680c4762a1bSJed Brown in = 1; 681c4762a1bSJed Brown st = sts[rank]; 682c4762a1bSJed Brown en = sts[rank + 1]; 683c4762a1bSJed Brown nl = en - st; 684c4762a1bSJed Brown break; 685c4762a1bSJed Brown } 6869566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,nl,st,in,&is)); 6879566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is,&nl)); 6889566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is,&idxs)); 6899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nl,&idxs2)); 690c4762a1bSJed Brown for (i=0;i<nl;i++) { 691c4762a1bSJed Brown switch (perm) { /* invert some of the indices */ 692c4762a1bSJed Brown case 2: 693c4762a1bSJed Brown idxs2[i] = rank%2 ? idxs[i] : idxs[nl-i-1]; 694c4762a1bSJed Brown break; 695c4762a1bSJed Brown case 1: 696c4762a1bSJed Brown idxs2[i] = rank%2 ? idxs[nl-i-1] : idxs[i]; 697c4762a1bSJed Brown break; 698c4762a1bSJed Brown default: 699c4762a1bSJed Brown idxs2[i] = idxs[i]; 700c4762a1bSJed Brown break; 701c4762a1bSJed Brown } 702c4762a1bSJed Brown } 7039566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is,&idxs)); 7049566063dSJacob Faibussowitsch PetscCall(ISCreateBlock(PETSC_COMM_WORLD,bs,nl,idxs2,PETSC_OWN_POINTER,&bis)); 7059566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(bis,&map)); 7069566063dSJacob Faibussowitsch PetscCall(MatCreateIS(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,bs*n,bs*n,map,map,&Abd)); 7079566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&map)); 7089566063dSJacob Faibussowitsch PetscCall(MatISSetPreallocation(Abd,bs,NULL,0,NULL)); 709c4762a1bSJed Brown for (i=0;i<nl;i++) { 710c4762a1bSJed Brown PetscInt b1,b2; 711c4762a1bSJed Brown 712c4762a1bSJed Brown for (b1=0;b1<bs;b1++) { 713c4762a1bSJed Brown for (b2=0;b2<bs;b2++) { 714c4762a1bSJed Brown vals[b1*bs + b2] = i*bs*bs + b1*bs + b2 + 1 + (b1 == b2 ? 1.0 : 0); 715c4762a1bSJed Brown } 716c4762a1bSJed Brown } 7179566063dSJacob Faibussowitsch PetscCall(MatSetValuesBlockedLocal(Abd,1,&i,1,&i,vals,INSERT_VALUES)); 718c4762a1bSJed Brown } 7199566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Abd,MAT_FINAL_ASSEMBLY)); 7209566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Abd,MAT_FINAL_ASSEMBLY)); 7219566063dSJacob Faibussowitsch PetscCall(MatConvert(Abd,MATAIJ,MAT_INITIAL_MATRIX,&Bbd)); 7229566063dSJacob Faibussowitsch PetscCall(MatInvertBlockDiagonal(Abd,&isbd)); 7239566063dSJacob Faibussowitsch PetscCall(MatInvertBlockDiagonal(Bbd,&aijbd)); 7249566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(Bbd,&nl,NULL)); 725c4762a1bSJed Brown ok = PETSC_TRUE; 726c4762a1bSJed Brown for (i=0;i<nl/bs;i++) { 727c4762a1bSJed Brown PetscInt b1,b2; 728c4762a1bSJed Brown 729c4762a1bSJed Brown for (b1=0;b1<bs;b1++) { 730c4762a1bSJed Brown for (b2=0;b2<bs;b2++) { 731c4762a1bSJed Brown if (PetscAbsScalar(isbd[i*bs*bs+b1*bs + b2]-aijbd[i*bs*bs+b1*bs + b2]) > PETSC_SMALL) ok = PETSC_FALSE; 732c4762a1bSJed Brown if (!ok) { 7339566063dSJacob 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]))); 734c4762a1bSJed Brown break; 735c4762a1bSJed Brown } 736c4762a1bSJed Brown } 737c4762a1bSJed Brown if (!ok) break; 738c4762a1bSJed Brown } 739c4762a1bSJed Brown if (!ok) break; 740c4762a1bSJed Brown } 7419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Abd)); 7429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bbd)); 7439566063dSJacob Faibussowitsch PetscCall(PetscFree(vals)); 7449566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 7459566063dSJacob Faibussowitsch PetscCall(ISDestroy(&bis)); 746c4762a1bSJed Brown } 747c4762a1bSJed Brown } 748c4762a1bSJed Brown } 749c4762a1bSJed Brown } 750d0dbe9f7SStefano Zampini 751d0dbe9f7SStefano Zampini /* test MatGetDiagonalBlock */ 752d0dbe9f7SStefano Zampini PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatGetDiagonalBlock\n")); 753d0dbe9f7SStefano Zampini PetscCall(MatGetDiagonalBlock(A,&A2)); 754d0dbe9f7SStefano Zampini PetscCall(MatGetDiagonalBlock(B,&B2)); 755d0dbe9f7SStefano Zampini PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatGetDiagonalBlock")); 756d0dbe9f7SStefano Zampini PetscCall(MatScale(A,2.0)); 757d0dbe9f7SStefano Zampini PetscCall(MatScale(B,2.0)); 758d0dbe9f7SStefano Zampini PetscCall(MatGetDiagonalBlock(A,&A2)); 759d0dbe9f7SStefano Zampini PetscCall(MatGetDiagonalBlock(B,&B2)); 760d0dbe9f7SStefano Zampini PetscCall(CheckMat(A2,B2,PETSC_FALSE,"MatGetDiagonalBlock")); 761d0dbe9f7SStefano Zampini 762c4762a1bSJed Brown /* free testing matrices */ 7639566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&cmap)); 7649566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&rmap)); 7659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 7669566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 7679566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 768b122ec5aSJacob Faibussowitsch return 0; 769c4762a1bSJed Brown } 770c4762a1bSJed Brown 771c4762a1bSJed Brown PetscErrorCode CheckMat(Mat A, Mat B, PetscBool usemult, const char* func) 772c4762a1bSJed Brown { 773c4762a1bSJed Brown Mat Bcheck; 774c4762a1bSJed Brown PetscReal error; 775c4762a1bSJed Brown 776c4762a1bSJed Brown PetscFunctionBeginUser; 777c4762a1bSJed Brown if (!usemult && B) { 778c4762a1bSJed Brown PetscBool hasnorm; 779c4762a1bSJed Brown 7809566063dSJacob Faibussowitsch PetscCall(MatHasOperation(B,MATOP_NORM,&hasnorm)); 781c4762a1bSJed Brown if (!hasnorm) usemult = PETSC_TRUE; 782c4762a1bSJed Brown } 783c4762a1bSJed Brown if (!usemult) { 784c4762a1bSJed Brown if (B) { 785c4762a1bSJed Brown MatType Btype; 786c4762a1bSJed Brown 7879566063dSJacob Faibussowitsch PetscCall(MatGetType(B,&Btype)); 7889566063dSJacob Faibussowitsch PetscCall(MatConvert(A,Btype,MAT_INITIAL_MATRIX,&Bcheck)); 789c4762a1bSJed Brown } else { 7909566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&Bcheck)); 791c4762a1bSJed Brown } 792c4762a1bSJed Brown if (B) { /* if B is present, subtract it */ 7939566063dSJacob Faibussowitsch PetscCall(MatAXPY(Bcheck,-1.,B,DIFFERENT_NONZERO_PATTERN)); 794c4762a1bSJed Brown } 7959566063dSJacob Faibussowitsch PetscCall(MatNorm(Bcheck,NORM_INFINITY,&error)); 796c4762a1bSJed Brown if (error > PETSC_SQRT_MACHINE_EPSILON) { 797c4762a1bSJed Brown ISLocalToGlobalMapping rl2g,cl2g; 798c4762a1bSJed Brown 799d0dbe9f7SStefano Zampini PetscCall(PetscObjectSetName((PetscObject)Bcheck,"Bcheck")); 8009566063dSJacob Faibussowitsch PetscCall(MatView(Bcheck,NULL)); 801c4762a1bSJed Brown if (B) { 802d0dbe9f7SStefano Zampini PetscCall(PetscObjectSetName((PetscObject)B,"B")); 8039566063dSJacob Faibussowitsch PetscCall(MatView(B,NULL)); 8049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bcheck)); 8059566063dSJacob Faibussowitsch PetscCall(MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&Bcheck)); 806d0dbe9f7SStefano Zampini PetscCall(PetscObjectSetName((PetscObject)Bcheck,"Assembled A")); 8079566063dSJacob Faibussowitsch PetscCall(MatView(Bcheck,NULL)); 808c4762a1bSJed Brown } 8099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bcheck)); 810d0dbe9f7SStefano Zampini PetscCall(PetscObjectSetName((PetscObject)A,"A")); 8119566063dSJacob Faibussowitsch PetscCall(MatView(A,NULL)); 8129566063dSJacob Faibussowitsch PetscCall(MatGetLocalToGlobalMapping(A,&rl2g,&cl2g)); 813d0dbe9f7SStefano Zampini if (rl2g) PetscCall(ISLocalToGlobalMappingView(rl2g,NULL)); 814d0dbe9f7SStefano Zampini if (cl2g) PetscCall(ISLocalToGlobalMappingView(cl2g,NULL)); 815d0dbe9f7SStefano Zampini SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"ERROR ON %s: %g",func,(double)error); 816c4762a1bSJed Brown } 8179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bcheck)); 818c4762a1bSJed Brown } else { 819c4762a1bSJed Brown PetscBool ok,okt; 820c4762a1bSJed Brown 8219566063dSJacob Faibussowitsch PetscCall(MatMultEqual(A,B,3,&ok)); 8229566063dSJacob Faibussowitsch PetscCall(MatMultTransposeEqual(A,B,3,&okt)); 823e00437b9SBarry Smith PetscCheck(ok && okt,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ERROR ON %s: mult ok ? %d, multtranspose ok ? %d",func,ok,okt); 824c4762a1bSJed Brown } 825c4762a1bSJed Brown PetscFunctionReturn(0); 826c4762a1bSJed Brown } 827c4762a1bSJed Brown 828c4762a1bSJed Brown PetscErrorCode TestMatZeroRows(Mat A, Mat Afull, PetscBool squaretest, IS is, PetscScalar diag) 829c4762a1bSJed Brown { 830c4762a1bSJed Brown Mat B,Bcheck,B2 = NULL,lB; 831c4762a1bSJed Brown Vec x = NULL, b = NULL, b2 = NULL; 832c4762a1bSJed Brown ISLocalToGlobalMapping l2gr,l2gc; 833c4762a1bSJed Brown PetscReal error; 834c4762a1bSJed Brown char diagstr[16]; 835c4762a1bSJed Brown const PetscInt *idxs; 836c4762a1bSJed Brown PetscInt rst,ren,i,n,N,d; 837c4762a1bSJed Brown PetscMPIInt rank; 838c4762a1bSJed Brown PetscBool miss,haszerorows; 839c4762a1bSJed Brown 840c4762a1bSJed Brown PetscFunctionBeginUser; 841c4762a1bSJed Brown if (diag == 0.) { 8429566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(diagstr,"zero")); 843c4762a1bSJed Brown } else { 8449566063dSJacob Faibussowitsch PetscCall(PetscStrcpy(diagstr,"nonzero")); 845c4762a1bSJed Brown } 8469566063dSJacob Faibussowitsch PetscCall(ISView(is,NULL)); 8479566063dSJacob Faibussowitsch PetscCall(MatGetLocalToGlobalMapping(A,&l2gr,&l2gc)); 848c4762a1bSJed Brown /* tests MatDuplicate and MatCopy */ 849c4762a1bSJed Brown if (diag == 0.) { 8509566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&B)); 851c4762a1bSJed Brown } else { 8529566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&B)); 8539566063dSJacob Faibussowitsch PetscCall(MatCopy(A,B,SAME_NONZERO_PATTERN)); 854c4762a1bSJed Brown } 8559566063dSJacob Faibussowitsch PetscCall(MatISGetLocalMat(B,&lB)); 8569566063dSJacob Faibussowitsch PetscCall(MatHasOperation(lB,MATOP_ZERO_ROWS,&haszerorows)); 857c4762a1bSJed Brown if (squaretest && haszerorows) { 858c4762a1bSJed Brown 8599566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(B,&x,&b)); 8609566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_COPY_VALUES,&B2)); 8619566063dSJacob Faibussowitsch PetscCall(VecSetLocalToGlobalMapping(b,l2gr)); 8629566063dSJacob Faibussowitsch PetscCall(VecSetLocalToGlobalMapping(x,l2gc)); 8639566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,NULL)); 8649566063dSJacob Faibussowitsch PetscCall(VecSetRandom(b,NULL)); 865c4762a1bSJed Brown /* mimic b[is] = x[is] */ 8669566063dSJacob Faibussowitsch PetscCall(VecDuplicate(b,&b2)); 8679566063dSJacob Faibussowitsch PetscCall(VecSetLocalToGlobalMapping(b2,l2gr)); 8689566063dSJacob Faibussowitsch PetscCall(VecCopy(b,b2)); 8699566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is,&n)); 8709566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is,&idxs)); 8719566063dSJacob Faibussowitsch PetscCall(VecGetSize(x,&N)); 872c4762a1bSJed Brown for (i=0;i<n;i++) { 873c4762a1bSJed Brown if (0 <= idxs[i] && idxs[i] < N) { 8749566063dSJacob Faibussowitsch PetscCall(VecSetValue(b2,idxs[i],diag,INSERT_VALUES)); 8759566063dSJacob Faibussowitsch PetscCall(VecSetValue(x,idxs[i],1.,INSERT_VALUES)); 876c4762a1bSJed Brown } 877c4762a1bSJed Brown } 8789566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(b2)); 8799566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(b2)); 8809566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(x)); 8819566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(x)); 8829566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is,&idxs)); 883c4762a1bSJed Brown /* test ZeroRows on MATIS */ 8849566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatZeroRows (diag %s)\n",diagstr)); 8859566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(B,is,diag,x,b)); 8869566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatZeroRowsColumns (diag %s)\n",diagstr)); 8879566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumnsIS(B2,is,diag,NULL,NULL)); 888c4762a1bSJed Brown } else if (haszerorows) { 889c4762a1bSJed Brown /* test ZeroRows on MATIS */ 8909566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatZeroRows (diag %s)\n",diagstr)); 8919566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(B,is,diag,NULL,NULL)); 892c4762a1bSJed Brown b = b2 = x = NULL; 893c4762a1bSJed Brown } else { 8949566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Skipping MatZeroRows (diag %s)\n",diagstr)); 895c4762a1bSJed Brown b = b2 = x = NULL; 896c4762a1bSJed Brown } 897c4762a1bSJed Brown 898c4762a1bSJed Brown if (squaretest && haszerorows) { 8999566063dSJacob Faibussowitsch PetscCall(VecAXPY(b2,-1.,b)); 9009566063dSJacob Faibussowitsch PetscCall(VecNorm(b2,NORM_INFINITY,&error)); 901e00437b9SBarry Smith PetscCheck(error <= PETSC_SQRT_MACHINE_EPSILON,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ERROR IN ZEROROWS ON B %g (diag %s)",(double)error,diagstr); 902c4762a1bSJed Brown } 903c4762a1bSJed Brown 904c4762a1bSJed Brown /* test MatMissingDiagonal */ 9059566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test MatMissingDiagonal\n")); 9069566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 9079566063dSJacob Faibussowitsch PetscCall(MatMissingDiagonal(B,&miss,&d)); 9089566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(B,&rst,&ren)); 9099566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 9109566063dSJacob 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)); 9119566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 9129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD)); 913c4762a1bSJed Brown 9149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 9159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 9169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b2)); 917c4762a1bSJed Brown 918c4762a1bSJed Brown /* check the result of ZeroRows with that from MPIAIJ routines 919c4762a1bSJed Brown assuming that MatConvert_IS_XAIJ and MatZeroRows_MPIAIJ work fine */ 920c4762a1bSJed Brown if (haszerorows) { 9219566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Afull,MAT_COPY_VALUES,&Bcheck)); 9229566063dSJacob Faibussowitsch PetscCall(MatSetOption(Bcheck,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 9239566063dSJacob Faibussowitsch PetscCall(MatZeroRowsIS(Bcheck,is,diag,NULL,NULL)); 9249566063dSJacob Faibussowitsch PetscCall(CheckMat(B,Bcheck,PETSC_FALSE,"Zerorows")); 9259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Bcheck)); 926c4762a1bSJed Brown } 9279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 928c4762a1bSJed Brown 929c4762a1bSJed Brown if (B2) { /* test MatZeroRowsColumns */ 9309566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Afull,MAT_COPY_VALUES,&B)); 9319566063dSJacob Faibussowitsch PetscCall(MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 9329566063dSJacob Faibussowitsch PetscCall(MatZeroRowsColumnsIS(B,is,diag,NULL,NULL)); 9339566063dSJacob Faibussowitsch PetscCall(CheckMat(B2,B,PETSC_FALSE,"MatZeroRowsColumns")); 9349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 9359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B2)); 936c4762a1bSJed Brown } 937c4762a1bSJed Brown PetscFunctionReturn(0); 938c4762a1bSJed Brown } 939c4762a1bSJed Brown 940c4762a1bSJed Brown /*TEST 941c4762a1bSJed Brown 942c4762a1bSJed Brown test: 943c4762a1bSJed Brown args: -test_trans 944c4762a1bSJed Brown 945c4762a1bSJed Brown test: 946c4762a1bSJed Brown suffix: 2 947c4762a1bSJed Brown nsize: 4 948c4762a1bSJed Brown args: -matis_convert_local_nest -nr 3 -nc 4 949c4762a1bSJed Brown 950c4762a1bSJed Brown test: 951c4762a1bSJed Brown suffix: 3 952c4762a1bSJed Brown nsize: 5 953c4762a1bSJed Brown args: -m 11 -n 10 -matis_convert_local_nest -nr 2 -nc 1 954c4762a1bSJed Brown 955c4762a1bSJed Brown test: 956c4762a1bSJed Brown suffix: 4 957c4762a1bSJed Brown nsize: 6 958c4762a1bSJed Brown args: -m 9 -n 12 -test_trans -nr 2 -nc 7 959c4762a1bSJed Brown 960c4762a1bSJed Brown test: 961c4762a1bSJed Brown suffix: 5 962c4762a1bSJed Brown nsize: 6 963c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 3 -nc 1 964c4762a1bSJed Brown 965c4762a1bSJed Brown test: 966c4762a1bSJed Brown suffix: 6 967c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap 968c4762a1bSJed Brown 969c4762a1bSJed Brown test: 970c4762a1bSJed Brown suffix: 7 971c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap -permmap 972c4762a1bSJed Brown 973c4762a1bSJed Brown test: 974c4762a1bSJed Brown suffix: 8 975c4762a1bSJed Brown args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -permmap 976c4762a1bSJed Brown 977c4762a1bSJed Brown test: 978c4762a1bSJed Brown suffix: 9 979c4762a1bSJed Brown nsize: 5 980c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap 981c4762a1bSJed Brown 982c4762a1bSJed Brown test: 983c4762a1bSJed Brown suffix: 10 984c4762a1bSJed Brown nsize: 5 985c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap -permmap 986c4762a1bSJed Brown 987c20d7725SJed Brown test: 988c20d7725SJed Brown suffix: vscat_default 989c4762a1bSJed Brown nsize: 5 990c4762a1bSJed Brown args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -permmap 991c4762a1bSJed Brown output_file: output/ex23_11.out 992c4762a1bSJed Brown 993c4762a1bSJed Brown test: 994c4762a1bSJed Brown suffix: 12 995c4762a1bSJed Brown nsize: 3 996c4762a1bSJed Brown args: -m 12 -n 12 -symmetric -matis_localmat_type sbaij -test_trans -nr 2 -nc 3 997c4762a1bSJed Brown 998c4762a1bSJed Brown testset: 999c4762a1bSJed Brown output_file: output/ex23_13.out 1000c4762a1bSJed Brown nsize: 3 1001c4762a1bSJed Brown args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -diffmap -permmap 1002c4762a1bSJed Brown filter: grep -v "type:" 1003c4762a1bSJed Brown test: 1004c4762a1bSJed Brown suffix: baij 1005c4762a1bSJed Brown args: -matis_localmat_type baij 1006c4762a1bSJed Brown test: 1007c4762a1bSJed Brown requires: viennacl 1008c4762a1bSJed Brown suffix: viennacl 1009c4762a1bSJed Brown args: -matis_localmat_type aijviennacl 1010c4762a1bSJed Brown test: 1011c4762a1bSJed Brown requires: cuda 1012c4762a1bSJed Brown suffix: cusparse 1013c4762a1bSJed Brown args: -matis_localmat_type aijcusparse 1014c4762a1bSJed Brown 1015e432b41dSStefano Zampini test: 1016e432b41dSStefano Zampini suffix: negrep 1017e432b41dSStefano Zampini nsize: {{1 3}separate output} 1018e432b41dSStefano 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} 1019e432b41dSStefano Zampini 1020c4762a1bSJed Brown TEST*/ 1021