1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests the use of interface functions for MATIS matrices.\n\ 3c4762a1bSJed Brown This example tests: MatZeroRows(), MatZeroRowsLocal(), MatView(), MatDuplicate(),\n\ 4c4762a1bSJed Brown MatCopy(), MatCreateSubMatrix(), MatGetLocalSubMatrix(), MatAXPY(), MatShift()\n\ 5c4762a1bSJed Brown MatDiagonalSet(), MatTranspose() and MatPtAP(). It also tests some\n\ 6c4762a1bSJed Brown conversion routines.\n"; 7c4762a1bSJed Brown 8c4762a1bSJed Brown #include <petscmat.h> 9c4762a1bSJed Brown 10c4762a1bSJed Brown PetscErrorCode TestMatZeroRows(Mat,Mat,PetscBool,IS,PetscScalar); 11c4762a1bSJed Brown PetscErrorCode CheckMat(Mat,Mat,PetscBool,const char*); 12c4762a1bSJed Brown 13c4762a1bSJed Brown int main(int argc,char **args) 14c4762a1bSJed Brown { 15c4762a1bSJed Brown Mat A,B,A2,B2,T; 16c4762a1bSJed Brown Mat Aee,Aeo,Aoe,Aoo; 17c4762a1bSJed Brown Mat *mats; 18c4762a1bSJed Brown Vec x,y; 19c4762a1bSJed Brown MatInfo info; 20c4762a1bSJed Brown ISLocalToGlobalMapping cmap,rmap; 21c4762a1bSJed Brown IS is,is2,reven,rodd,ceven,codd; 22c4762a1bSJed Brown IS *rows,*cols; 23c4762a1bSJed Brown MatType lmtype; 24c4762a1bSJed Brown PetscScalar diag = 2.; 25c4762a1bSJed Brown PetscInt n,m,i; 26c4762a1bSJed Brown PetscInt rst,ren,cst,cen,nr,nc; 27c4762a1bSJed Brown PetscMPIInt rank,size; 28c4762a1bSJed Brown PetscBool testT,squaretest,isaij; 29c4762a1bSJed Brown PetscBool permute = PETSC_FALSE; 30c4762a1bSJed Brown PetscBool diffmap = PETSC_TRUE, symmetric = PETSC_FALSE; 31c4762a1bSJed Brown PetscErrorCode ierr; 32c4762a1bSJed Brown 33c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 34ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 35ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 36c4762a1bSJed Brown m = n = 2*size; 37c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-symmetric",&symmetric,NULL);CHKERRQ(ierr); 38c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr); 39c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr); 40c4762a1bSJed Brown if (size > 1 && m < 4) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Number of rows should be larger or equal 4 for parallel runs"); 41c4762a1bSJed Brown if (size == 1 && m < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Number of rows should be larger or equal 2 for uniprocessor runs"); 42c4762a1bSJed Brown if (n < 2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_WRONG,"Number of cols should be larger or equal 2"); 43c4762a1bSJed Brown if (symmetric) m = n = PetscMax(m,n); 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* create a MATIS matrix */ 46c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 47c4762a1bSJed Brown ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr); 48c4762a1bSJed Brown ierr = MatSetType(A,MATIS);CHKERRQ(ierr); 49c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 50fc989267SStefano Zampini /* This is not the proper setting for MATIS for finite elements, it is just used to test the routines 51fc989267SStefano Zampini Here we use a one-to-one correspondence between local row/column spaces and global row/column spaces 52fc989267SStefano Zampini Equivalent to passing NULL for the mapping */ 53c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,n,0,1,&is);CHKERRQ(ierr); 54c4762a1bSJed Brown ierr = ISLocalToGlobalMappingCreateIS(is,&cmap);CHKERRQ(ierr); 55c4762a1bSJed Brown ierr = ISDestroy(&is);CHKERRQ(ierr); 56c4762a1bSJed Brown 57c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-permmap",&permute,NULL);CHKERRQ(ierr); 58c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-diffmap",&diffmap,NULL);CHKERRQ(ierr); 59c4762a1bSJed Brown if (!symmetric && (diffmap || m != n)) { 60c4762a1bSJed Brown 61c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,m,permute ? m -1 : 0,permute ? -1 : 1,&is);CHKERRQ(ierr); 62c4762a1bSJed Brown ierr = ISLocalToGlobalMappingCreateIS(is,&rmap);CHKERRQ(ierr); 63c4762a1bSJed Brown ierr = ISDestroy(&is);CHKERRQ(ierr); 64c4762a1bSJed Brown if (m==n && !permute) squaretest = PETSC_TRUE; 65c4762a1bSJed Brown else squaretest = PETSC_FALSE; 66c4762a1bSJed Brown } else { 67c4762a1bSJed Brown ierr = PetscObjectReference((PetscObject)cmap);CHKERRQ(ierr); 68c4762a1bSJed Brown rmap = cmap; 69c4762a1bSJed Brown squaretest = PETSC_TRUE; 70c4762a1bSJed Brown } 71c4762a1bSJed Brown ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr); 72c4762a1bSJed Brown ierr = MatISStoreL2L(A,PETSC_FALSE);CHKERRQ(ierr); 73c4762a1bSJed Brown ierr = MatISSetPreallocation(A,3,NULL,0,NULL);CHKERRQ(ierr); 74c4762a1bSJed Brown for (i=0; i<m; i++) { 75c4762a1bSJed Brown PetscScalar v[3]; 76c4762a1bSJed Brown PetscInt cols[3]; 77c4762a1bSJed Brown 78c4762a1bSJed Brown cols[0] = (i-1+n)%n; 79c4762a1bSJed Brown cols[1] = i%n; 80c4762a1bSJed Brown cols[2] = (i+1)%n; 81c4762a1bSJed Brown v[0] = -1.*(symmetric ? PetscMin(i+1,cols[0]+1) : i+1); 82c4762a1bSJed Brown v[1] = 2.*(symmetric ? PetscMin(i+1,cols[1]+1) : i+1); 83c4762a1bSJed Brown v[2] = -1.*(symmetric ? PetscMin(i+1,cols[2]+1) : i+1); 84c4762a1bSJed Brown ierr = MatSetValuesLocal(A,1,&i,3,cols,v,ADD_VALUES);CHKERRQ(ierr); 85c4762a1bSJed Brown } 86c4762a1bSJed Brown if (symmetric) { 87c4762a1bSJed Brown ierr = MatSetOption(A,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 88c4762a1bSJed Brown } 89c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 91c4762a1bSJed Brown 92c4762a1bSJed Brown ierr = MatISGetLocalMat(A,&B);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = MatGetType(B,&lmtype);CHKERRQ(ierr); 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* test MatGetInfo */ 96c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatGetInfo\n");CHKERRQ(ierr); 97c4762a1bSJed Brown if (!PetscGlobalRank) { 98c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 99c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown ierr = MatGetInfo(A,MAT_LOCAL,&info);CHKERRQ(ierr); 102c4762a1bSJed Brown ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1038fffc762SJacob Faibussowitsch ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"Process %2d: %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",PetscGlobalRank,(PetscInt)info.nz_used, 104c4762a1bSJed Brown (PetscInt)info.nz_allocated,(PetscInt)info.nz_unneeded,(PetscInt)info.assemblies,(PetscInt)info.mallocs);CHKERRQ(ierr); 105c4762a1bSJed Brown ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 106c4762a1bSJed Brown ierr = MatGetInfo(A,MAT_GLOBAL_MAX,&info);CHKERRQ(ierr); 1078fffc762SJacob Faibussowitsch ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"GlobalMax : %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",(PetscInt)info.nz_used, 108c4762a1bSJed Brown (PetscInt)info.nz_allocated,(PetscInt)info.nz_unneeded,(PetscInt)info.assemblies,(PetscInt)info.mallocs);CHKERRQ(ierr); 109c4762a1bSJed Brown ierr = MatGetInfo(A,MAT_GLOBAL_SUM,&info);CHKERRQ(ierr); 1108fffc762SJacob Faibussowitsch ierr = PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"GlobalSum : %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",(PetscInt)info.nz_used, 111c4762a1bSJed Brown (PetscInt)info.nz_allocated,(PetscInt)info.nz_unneeded,(PetscInt)info.assemblies,(PetscInt)info.mallocs);CHKERRQ(ierr); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* test MatView */ 114c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatView\n");CHKERRQ(ierr); 115c4762a1bSJed Brown ierr = MatView(A,NULL);CHKERRQ(ierr); 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* Create a MPIAIJ matrix, same as A */ 118c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr); 119c4762a1bSJed Brown ierr = MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,m,n);CHKERRQ(ierr); 120c4762a1bSJed Brown ierr = MatSetType(B,MATAIJ);CHKERRQ(ierr); 121c4762a1bSJed Brown ierr = MatSetFromOptions(B);CHKERRQ(ierr); 122c4762a1bSJed Brown ierr = MatSetUp(B);CHKERRQ(ierr); 123c4762a1bSJed Brown ierr = MatSetLocalToGlobalMapping(B,rmap,cmap);CHKERRQ(ierr); 124c4762a1bSJed Brown ierr = MatMPIAIJSetPreallocation(B,3,NULL,3,NULL);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = MatMPIBAIJSetPreallocation(B,1,3,NULL,3,NULL);CHKERRQ(ierr); 126c4762a1bSJed Brown #if defined(PETSC_HAVE_HYPRE) 127c4762a1bSJed Brown ierr = MatHYPRESetPreallocation(B,3,NULL,3,NULL);CHKERRQ(ierr); 128c4762a1bSJed Brown #endif 129c4762a1bSJed Brown ierr = MatISSetPreallocation(B,3,NULL,3,NULL);CHKERRQ(ierr); 130c4762a1bSJed Brown for (i=0; i<m; i++) { 131c4762a1bSJed Brown PetscScalar v[3]; 132c4762a1bSJed Brown PetscInt cols[3]; 133c4762a1bSJed Brown 134c4762a1bSJed Brown cols[0] = (i-1+n)%n; 135c4762a1bSJed Brown cols[1] = i%n; 136c4762a1bSJed Brown cols[2] = (i+1)%n; 137c4762a1bSJed Brown v[0] = -1.*(symmetric ? PetscMin(i+1,cols[0]+1) : i+1); 138c4762a1bSJed Brown v[1] = 2.*(symmetric ? PetscMin(i+1,cols[1]+1) : i+1); 139c4762a1bSJed Brown v[2] = -1.*(symmetric ? PetscMin(i+1,cols[2]+1) : i+1); 140c4762a1bSJed Brown ierr = MatSetValuesLocal(B,1,&i,3,cols,v,ADD_VALUES);CHKERRQ(ierr); 141c4762a1bSJed Brown } 142c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 144c4762a1bSJed Brown ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr); 145c4762a1bSJed Brown ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr); 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* test CheckMat */ 148c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test CheckMat\n");CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = CheckMat(A,B,PETSC_FALSE,"CheckMat");CHKERRQ(ierr); 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* test MatDuplicate and MatAXPY */ 152c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatDuplicate and MatAXPY\n");CHKERRQ(ierr); 153c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr); 154c4762a1bSJed Brown ierr = CheckMat(A,A2,PETSC_FALSE,"MatDuplicate and MatAXPY");CHKERRQ(ierr); 155c4762a1bSJed Brown 156c4762a1bSJed Brown /* test MatConvert */ 157c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatConvert_IS_XAIJ\n");CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = MatConvert(A2,MATAIJ,MAT_INITIAL_MATRIX,&B2);CHKERRQ(ierr); 159c4762a1bSJed Brown ierr = CheckMat(B,B2,PETSC_TRUE,"MatConvert_IS_XAIJ MAT_INITIAL_MATRIX");CHKERRQ(ierr); 160c4762a1bSJed Brown ierr = MatConvert(A2,MATAIJ,MAT_REUSE_MATRIX,&B2);CHKERRQ(ierr); 161c4762a1bSJed Brown ierr = CheckMat(B,B2,PETSC_TRUE,"MatConvert_IS_XAIJ MAT_REUSE_MATRIX");CHKERRQ(ierr); 162c4762a1bSJed Brown ierr = MatConvert(A2,MATAIJ,MAT_INPLACE_MATRIX,&A2);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = CheckMat(B,A2,PETSC_TRUE,"MatConvert_IS_XAIJ MAT_INPLACE_MATRIX");CHKERRQ(ierr); 164c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 165c4762a1bSJed Brown ierr = MatDestroy(&B2);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatConvert_XAIJ_IS\n");CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = MatDuplicate(B,MAT_COPY_VALUES,&B2);CHKERRQ(ierr); 168c4762a1bSJed Brown ierr = MatConvert(B2,MATIS,MAT_INITIAL_MATRIX,&A2);CHKERRQ(ierr); 169c4762a1bSJed Brown ierr = CheckMat(A,A2,PETSC_TRUE,"MatConvert_XAIJ_IS MAT_INITIAL_MATRIX");CHKERRQ(ierr); 170c4762a1bSJed Brown ierr = MatConvert(B2,MATIS,MAT_REUSE_MATRIX,&A2);CHKERRQ(ierr); 171c4762a1bSJed Brown ierr = CheckMat(A,A2,PETSC_TRUE,"MatConvert_XAIJ_IS MAT_REUSE_MATRIX");CHKERRQ(ierr); 172c4762a1bSJed Brown ierr = MatConvert(B2,MATIS,MAT_INPLACE_MATRIX,&B2);CHKERRQ(ierr); 173c4762a1bSJed Brown ierr = CheckMat(A,B2,PETSC_TRUE,"MatConvert_XAIJ_IS MAT_INPLACE_MATRIX");CHKERRQ(ierr); 174c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 175c4762a1bSJed Brown ierr = MatDestroy(&B2);CHKERRQ(ierr); 176c4762a1bSJed Brown ierr = PetscStrcmp(lmtype,MATSEQAIJ,&isaij);CHKERRQ(ierr); 177c4762a1bSJed Brown if (size == 1 && isaij) { /* tests special code paths in MatConvert_IS_XAIJ */ 178c4762a1bSJed 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}; 179c4762a1bSJed Brown 180c4762a1bSJed Brown for (ri = 0; ri < 2; ri++) { 181c4762a1bSJed Brown PetscInt *r; 182c4762a1bSJed Brown 183c4762a1bSJed Brown r = (PetscInt*)(ri == 0 ? rr : rk); 184c4762a1bSJed Brown for (ci = 0; ci < 2; ci++) { 185c4762a1bSJed Brown PetscInt *c,rb,cb; 186c4762a1bSJed Brown 187c4762a1bSJed Brown c = (PetscInt*)(ci == 0 ? cr : ck); 188c4762a1bSJed Brown for (rb = 1; rb < 4; rb++) { 189c4762a1bSJed Brown ierr = ISCreateBlock(PETSC_COMM_SELF,rb,3,r,PETSC_COPY_VALUES,&is);CHKERRQ(ierr); 190c4762a1bSJed Brown ierr = ISLocalToGlobalMappingCreateIS(is,&rmap);CHKERRQ(ierr); 191c4762a1bSJed Brown ierr = ISDestroy(&is);CHKERRQ(ierr); 192c4762a1bSJed Brown for (cb = 1; cb < 4; cb++) { 193c4762a1bSJed Brown Mat T,lT,T2; 194c4762a1bSJed Brown char testname[256]; 195c4762a1bSJed Brown 1968fffc762SJacob Faibussowitsch ierr = PetscSNPrintf(testname,sizeof(testname),"MatConvert_IS_XAIJ special case (%" PetscInt_FMT " %" PetscInt_FMT ", bs %" PetscInt_FMT " %" PetscInt_FMT ")",ri,ci,rb,cb);CHKERRQ(ierr); 197c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test %s\n",testname);CHKERRQ(ierr); 198c4762a1bSJed Brown 199c4762a1bSJed Brown ierr = ISCreateBlock(PETSC_COMM_SELF,cb,4,c,PETSC_COPY_VALUES,&is);CHKERRQ(ierr); 200c4762a1bSJed Brown ierr = ISLocalToGlobalMappingCreateIS(is,&cmap);CHKERRQ(ierr); 201c4762a1bSJed Brown ierr = ISDestroy(&is);CHKERRQ(ierr); 202c4762a1bSJed Brown 203c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_SELF,&T);CHKERRQ(ierr); 204c4762a1bSJed Brown ierr = MatSetSizes(T,PETSC_DECIDE,PETSC_DECIDE,rb*3,cb*4);CHKERRQ(ierr); 205c4762a1bSJed Brown ierr = MatSetType(T,MATIS);CHKERRQ(ierr); 206c4762a1bSJed Brown ierr = MatSetLocalToGlobalMapping(T,rmap,cmap);CHKERRQ(ierr); 207c4762a1bSJed Brown ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr); 208c4762a1bSJed Brown ierr = MatISGetLocalMat(T,&lT);CHKERRQ(ierr); 209c4762a1bSJed Brown ierr = MatSetType(lT,MATSEQAIJ);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = MatSeqAIJSetPreallocation(lT,cb*4,NULL);CHKERRQ(ierr); 211c4762a1bSJed Brown ierr = MatSetRandom(lT,NULL);CHKERRQ(ierr); 212c4762a1bSJed Brown ierr = MatConvert(lT,lmtype,MAT_INPLACE_MATRIX,&lT);CHKERRQ(ierr); 213c4762a1bSJed Brown ierr = MatISRestoreLocalMat(T,&lT);CHKERRQ(ierr); 214c4762a1bSJed Brown ierr = MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 215c4762a1bSJed Brown ierr = MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 216c4762a1bSJed Brown 217c4762a1bSJed Brown ierr = MatConvert(T,MATAIJ,MAT_INITIAL_MATRIX,&T2);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = CheckMat(T,T2,PETSC_TRUE,"MAT_INITIAL_MATRIX");CHKERRQ(ierr); 219c4762a1bSJed Brown ierr = MatConvert(T,MATAIJ,MAT_REUSE_MATRIX,&T2);CHKERRQ(ierr); 220c4762a1bSJed Brown ierr = CheckMat(T,T2,PETSC_TRUE,"MAT_REUSE_MATRIX");CHKERRQ(ierr); 221c4762a1bSJed Brown ierr = MatDestroy(&T2);CHKERRQ(ierr); 222c4762a1bSJed Brown ierr = MatDuplicate(T,MAT_COPY_VALUES,&T2);CHKERRQ(ierr); 223c4762a1bSJed Brown ierr = MatConvert(T2,MATAIJ,MAT_INPLACE_MATRIX,&T2);CHKERRQ(ierr); 224c4762a1bSJed Brown ierr = CheckMat(T,T2,PETSC_TRUE,"MAT_INPLACE_MATRIX");CHKERRQ(ierr); 225c4762a1bSJed Brown ierr = MatDestroy(&T);CHKERRQ(ierr); 226c4762a1bSJed Brown ierr = MatDestroy(&T2);CHKERRQ(ierr); 227c4762a1bSJed Brown } 228c4762a1bSJed Brown ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr); 229c4762a1bSJed Brown } 230c4762a1bSJed Brown } 231c4762a1bSJed Brown } 232c4762a1bSJed Brown } 233c4762a1bSJed Brown 234c4762a1bSJed Brown /* test MatDiagonalScale */ 235c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatDiagonalScale\n");CHKERRQ(ierr); 236c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr); 237c4762a1bSJed Brown ierr = MatDuplicate(B,MAT_COPY_VALUES,&B2);CHKERRQ(ierr); 238c4762a1bSJed Brown ierr = MatCreateVecs(A,&x,&y);CHKERRQ(ierr); 239c4762a1bSJed Brown ierr = VecSetRandom(x,NULL);CHKERRQ(ierr); 240c4762a1bSJed Brown if (symmetric) { 241c4762a1bSJed Brown ierr = VecCopy(x,y);CHKERRQ(ierr); 242c4762a1bSJed Brown } else { 243c4762a1bSJed Brown ierr = VecSetRandom(y,NULL);CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = VecScale(y,8.);CHKERRQ(ierr); 245c4762a1bSJed Brown } 246c4762a1bSJed Brown ierr = MatDiagonalScale(A2,y,x);CHKERRQ(ierr); 247c4762a1bSJed Brown ierr = MatDiagonalScale(B2,y,x);CHKERRQ(ierr); 248c4762a1bSJed Brown ierr = CheckMat(A2,B2,PETSC_FALSE,"MatDiagonalScale");CHKERRQ(ierr); 249c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 250c4762a1bSJed Brown ierr = MatDestroy(&B2);CHKERRQ(ierr); 251c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 252c4762a1bSJed Brown ierr = VecDestroy(&y);CHKERRQ(ierr); 253c4762a1bSJed Brown 254c4762a1bSJed Brown /* test MatPtAP (A IS and B AIJ) */ 255c4762a1bSJed Brown if (isaij && m == n) { 256c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatPtAP\n");CHKERRQ(ierr); 257c4762a1bSJed Brown ierr = MatISStoreL2L(A,PETSC_TRUE);CHKERRQ(ierr); 258c4762a1bSJed Brown ierr = MatPtAP(A,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&A2);CHKERRQ(ierr); 259c4762a1bSJed Brown ierr = MatPtAP(B,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B2);CHKERRQ(ierr); 260c4762a1bSJed Brown ierr = CheckMat(A2,B2,PETSC_FALSE,"MatPtAP MAT_INITIAL_MATRIX");CHKERRQ(ierr); 261c4762a1bSJed Brown ierr = MatPtAP(A,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&A2);CHKERRQ(ierr); 262c4762a1bSJed Brown ierr = CheckMat(A2,B2,PETSC_FALSE,"MatPtAP MAT_REUSE_MATRIX");CHKERRQ(ierr); 263c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 264c4762a1bSJed Brown ierr = MatDestroy(&B2);CHKERRQ(ierr); 265c4762a1bSJed Brown } 266c4762a1bSJed Brown 267c4762a1bSJed Brown /* test MatGetLocalSubMatrix */ 268c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatGetLocalSubMatrix\n");CHKERRQ(ierr); 269c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&A2);CHKERRQ(ierr); 270c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,m/2+m%2,0,2,&reven);CHKERRQ(ierr); 271c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,m/2,1,2,&rodd);CHKERRQ(ierr); 272c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,n/2+n%2,0,2,&ceven);CHKERRQ(ierr); 273c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,n/2,1,2,&codd);CHKERRQ(ierr); 274c4762a1bSJed Brown ierr = MatGetLocalSubMatrix(A2,reven,ceven,&Aee);CHKERRQ(ierr); 275c4762a1bSJed Brown ierr = MatGetLocalSubMatrix(A2,reven,codd,&Aeo);CHKERRQ(ierr); 276c4762a1bSJed Brown ierr = MatGetLocalSubMatrix(A2,rodd,ceven,&Aoe);CHKERRQ(ierr); 277c4762a1bSJed Brown ierr = MatGetLocalSubMatrix(A2,rodd,codd,&Aoo);CHKERRQ(ierr); 278c4762a1bSJed Brown for (i=0; i<m; i++) { 279c4762a1bSJed Brown PetscInt j,je,jo,colse[3], colso[3]; 280c4762a1bSJed Brown PetscScalar ve[3], vo[3]; 281c4762a1bSJed Brown PetscScalar v[3]; 282c4762a1bSJed Brown PetscInt cols[3]; 283c4762a1bSJed Brown 284c4762a1bSJed Brown cols[0] = (i-1+n)%n; 285c4762a1bSJed Brown cols[1] = i%n; 286c4762a1bSJed Brown cols[2] = (i+1)%n; 287c4762a1bSJed Brown v[0] = -1.*(symmetric ? PetscMin(i+1,cols[0]+1) : i+1); 288c4762a1bSJed Brown v[1] = 2.*(symmetric ? PetscMin(i+1,cols[1]+1) : i+1); 289c4762a1bSJed Brown v[2] = -1.*(symmetric ? PetscMin(i+1,cols[2]+1) : i+1); 290c4762a1bSJed Brown for (j=0,je=0,jo=0;j<3;j++) { 291c4762a1bSJed Brown if (cols[j]%2) { 292c4762a1bSJed Brown vo[jo] = v[j]; 293c4762a1bSJed Brown colso[jo++] = cols[j]/2; 294c4762a1bSJed Brown } else { 295c4762a1bSJed Brown ve[je] = v[j]; 296c4762a1bSJed Brown colse[je++] = cols[j]/2; 297c4762a1bSJed Brown } 298c4762a1bSJed Brown } 299c4762a1bSJed Brown if (i%2) { 300c4762a1bSJed Brown PetscInt row = i/2; 301c4762a1bSJed Brown ierr = MatSetValuesLocal(Aoe,1,&row,je,colse,ve,ADD_VALUES);CHKERRQ(ierr); 302c4762a1bSJed Brown ierr = MatSetValuesLocal(Aoo,1,&row,jo,colso,vo,ADD_VALUES);CHKERRQ(ierr); 303c4762a1bSJed Brown } else { 304c4762a1bSJed Brown PetscInt row = i/2; 305c4762a1bSJed Brown ierr = MatSetValuesLocal(Aee,1,&row,je,colse,ve,ADD_VALUES);CHKERRQ(ierr); 306c4762a1bSJed Brown ierr = MatSetValuesLocal(Aeo,1,&row,jo,colso,vo,ADD_VALUES);CHKERRQ(ierr); 307c4762a1bSJed Brown } 308c4762a1bSJed Brown } 309c4762a1bSJed Brown ierr = MatRestoreLocalSubMatrix(A2,reven,ceven,&Aee);CHKERRQ(ierr); 310c4762a1bSJed Brown ierr = MatRestoreLocalSubMatrix(A2,reven,codd,&Aeo);CHKERRQ(ierr); 311c4762a1bSJed Brown ierr = MatRestoreLocalSubMatrix(A2,rodd,ceven,&Aoe);CHKERRQ(ierr); 312c4762a1bSJed Brown ierr = MatRestoreLocalSubMatrix(A2,rodd,codd,&Aoo);CHKERRQ(ierr); 313c4762a1bSJed Brown ierr = ISDestroy(&reven);CHKERRQ(ierr); 314c4762a1bSJed Brown ierr = ISDestroy(&ceven);CHKERRQ(ierr); 315c4762a1bSJed Brown ierr = ISDestroy(&rodd);CHKERRQ(ierr); 316c4762a1bSJed Brown ierr = ISDestroy(&codd);CHKERRQ(ierr); 317c4762a1bSJed Brown ierr = MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 318c4762a1bSJed Brown ierr = MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 319c4762a1bSJed Brown ierr = MatAXPY(A2,-1.,A,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 320c4762a1bSJed Brown ierr = CheckMat(A2,NULL,PETSC_FALSE,"MatGetLocalSubMatrix");CHKERRQ(ierr); 321c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 322c4762a1bSJed Brown 323c4762a1bSJed Brown /* test MatConvert_Nest_IS */ 324c4762a1bSJed Brown testT = PETSC_FALSE; 325c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-test_trans",&testT,NULL);CHKERRQ(ierr); 326c4762a1bSJed Brown 327c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatConvert_Nest_IS\n");CHKERRQ(ierr); 328c4762a1bSJed Brown nr = 2; 329c4762a1bSJed Brown nc = 2; 330c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-nr",&nr,NULL);CHKERRQ(ierr); 331c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-nc",&nc,NULL);CHKERRQ(ierr); 332c4762a1bSJed Brown if (testT) { 333c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&cst,&cen);CHKERRQ(ierr); 334c4762a1bSJed Brown ierr = MatGetOwnershipRangeColumn(A,&rst,&ren);CHKERRQ(ierr); 335c4762a1bSJed Brown } else { 336c4762a1bSJed Brown ierr = MatGetOwnershipRange(A,&rst,&ren);CHKERRQ(ierr); 337c4762a1bSJed Brown ierr = MatGetOwnershipRangeColumn(A,&cst,&cen);CHKERRQ(ierr); 338c4762a1bSJed Brown } 339c4762a1bSJed Brown ierr = PetscMalloc3(nr,&rows,nc,&cols,2*nr*nc,&mats);CHKERRQ(ierr); 340c4762a1bSJed Brown for (i=0;i<nr*nc;i++) { 341c4762a1bSJed Brown if (testT) { 342c4762a1bSJed Brown ierr = MatCreateTranspose(A,&mats[i]);CHKERRQ(ierr); 343c4762a1bSJed Brown ierr = MatTranspose(B,MAT_INITIAL_MATRIX,&mats[i+nr*nc]);CHKERRQ(ierr); 344c4762a1bSJed Brown } else { 345c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_COPY_VALUES,&mats[i]);CHKERRQ(ierr); 346c4762a1bSJed Brown ierr = MatDuplicate(B,MAT_COPY_VALUES,&mats[i+nr*nc]);CHKERRQ(ierr); 347c4762a1bSJed Brown } 348c4762a1bSJed Brown } 349c4762a1bSJed Brown for (i=0;i<nr;i++) { 350c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,ren-rst,i+rst,nr,&rows[i]);CHKERRQ(ierr); 351c4762a1bSJed Brown } 352c4762a1bSJed Brown for (i=0;i<nc;i++) { 353c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,cen-cst,i+cst,nc,&cols[i]);CHKERRQ(ierr); 354c4762a1bSJed Brown } 355c4762a1bSJed Brown ierr = MatCreateNest(PETSC_COMM_WORLD,nr,rows,nc,cols,mats,&A2);CHKERRQ(ierr); 356c4762a1bSJed Brown ierr = MatCreateNest(PETSC_COMM_WORLD,nr,rows,nc,cols,mats+nr*nc,&B2);CHKERRQ(ierr); 357c4762a1bSJed Brown for (i=0;i<nr;i++) { 358c4762a1bSJed Brown ierr = ISDestroy(&rows[i]);CHKERRQ(ierr); 359c4762a1bSJed Brown } 360c4762a1bSJed Brown for (i=0;i<nc;i++) { 361c4762a1bSJed Brown ierr = ISDestroy(&cols[i]);CHKERRQ(ierr); 362c4762a1bSJed Brown } 363c4762a1bSJed Brown for (i=0;i<2*nr*nc;i++) { 364c4762a1bSJed Brown ierr = MatDestroy(&mats[i]);CHKERRQ(ierr); 365c4762a1bSJed Brown } 366c4762a1bSJed Brown ierr = PetscFree3(rows,cols,mats);CHKERRQ(ierr); 367c4762a1bSJed Brown ierr = MatConvert(B2,MATAIJ,MAT_INITIAL_MATRIX,&T);CHKERRQ(ierr); 368c4762a1bSJed Brown ierr = MatDestroy(&B2);CHKERRQ(ierr); 369c4762a1bSJed Brown ierr = MatConvert(A2,MATIS,MAT_INITIAL_MATRIX,&B2);CHKERRQ(ierr); 370c4762a1bSJed Brown ierr = CheckMat(B2,T,PETSC_TRUE,"MatConvert_Nest_IS MAT_INITIAL_MATRIX");CHKERRQ(ierr); 371c4762a1bSJed Brown ierr = MatConvert(A2,MATIS,MAT_REUSE_MATRIX,&B2);CHKERRQ(ierr); 372c4762a1bSJed Brown ierr = CheckMat(B2,T,PETSC_TRUE,"MatConvert_Nest_IS MAT_REUSE_MATRIX");CHKERRQ(ierr); 373c4762a1bSJed Brown ierr = MatDestroy(&B2);CHKERRQ(ierr); 374c4762a1bSJed Brown ierr = MatConvert(A2,MATIS,MAT_INPLACE_MATRIX,&A2);CHKERRQ(ierr); 375c4762a1bSJed Brown ierr = CheckMat(A2,T,PETSC_TRUE,"MatConvert_Nest_IS MAT_INPLACE_MATRIX");CHKERRQ(ierr); 376c4762a1bSJed Brown ierr = MatDestroy(&T);CHKERRQ(ierr); 377c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 378c4762a1bSJed Brown 379c4762a1bSJed Brown /* test MatCreateSubMatrix */ 380c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatCreateSubMatrix\n");CHKERRQ(ierr); 381dd400576SPatrick Sanan if (rank == 0) { 382c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,1,1,1,&is);CHKERRQ(ierr); 383c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,2,0,1,&is2);CHKERRQ(ierr); 384c4762a1bSJed Brown } else if (rank == 1) { 385c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,1,0,1,&is);CHKERRQ(ierr); 386c4762a1bSJed Brown if (n > 3) { 387c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,1,3,1,&is2);CHKERRQ(ierr); 388c4762a1bSJed Brown } else { 389c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is2);CHKERRQ(ierr); 390c4762a1bSJed Brown } 391c4762a1bSJed Brown } else if (rank == 2 && n > 4) { 392c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is);CHKERRQ(ierr); 393c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,n-4,4,1,&is2);CHKERRQ(ierr); 394c4762a1bSJed Brown } else { 395c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is);CHKERRQ(ierr); 396c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is2);CHKERRQ(ierr); 397c4762a1bSJed Brown } 398c4762a1bSJed Brown ierr = MatCreateSubMatrix(A,is,is,MAT_INITIAL_MATRIX,&A2);CHKERRQ(ierr); 399c4762a1bSJed Brown ierr = MatCreateSubMatrix(B,is,is,MAT_INITIAL_MATRIX,&B2);CHKERRQ(ierr); 400c4762a1bSJed Brown ierr = CheckMat(A2,B2,PETSC_TRUE,"first MatCreateSubMatrix");CHKERRQ(ierr); 401c4762a1bSJed Brown 402c4762a1bSJed Brown ierr = MatCreateSubMatrix(A,is,is,MAT_REUSE_MATRIX,&A2);CHKERRQ(ierr); 403c4762a1bSJed Brown ierr = MatCreateSubMatrix(B,is,is,MAT_REUSE_MATRIX,&B2);CHKERRQ(ierr); 404c4762a1bSJed Brown ierr = CheckMat(A2,B2,PETSC_FALSE,"reuse MatCreateSubMatrix");CHKERRQ(ierr); 405c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 406c4762a1bSJed Brown ierr = MatDestroy(&B2);CHKERRQ(ierr); 407c4762a1bSJed Brown 408c4762a1bSJed Brown if (!symmetric) { 409c4762a1bSJed Brown ierr = MatCreateSubMatrix(A,is,is2,MAT_INITIAL_MATRIX,&A2);CHKERRQ(ierr); 410c4762a1bSJed Brown ierr = MatCreateSubMatrix(B,is,is2,MAT_INITIAL_MATRIX,&B2);CHKERRQ(ierr); 411c4762a1bSJed Brown ierr = MatCreateSubMatrix(A,is,is2,MAT_REUSE_MATRIX,&A2);CHKERRQ(ierr); 412c4762a1bSJed Brown ierr = MatCreateSubMatrix(B,is,is2,MAT_REUSE_MATRIX,&B2);CHKERRQ(ierr); 413c4762a1bSJed Brown ierr = CheckMat(A2,B2,PETSC_FALSE,"second MatCreateSubMatrix");CHKERRQ(ierr); 414c4762a1bSJed Brown } 415c4762a1bSJed Brown 416c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 417c4762a1bSJed Brown ierr = MatDestroy(&B2);CHKERRQ(ierr); 418c4762a1bSJed Brown ierr = ISDestroy(&is);CHKERRQ(ierr); 419c4762a1bSJed Brown ierr = ISDestroy(&is2);CHKERRQ(ierr); 420c4762a1bSJed Brown 421c4762a1bSJed Brown /* Create an IS required by MatZeroRows(): just rank zero provides the rows to be eliminated */ 422c4762a1bSJed Brown if (size > 1) { 423dd400576SPatrick Sanan if (rank == 0) { 424c4762a1bSJed Brown PetscInt st,len; 425c4762a1bSJed Brown 426c4762a1bSJed Brown st = (m+1)/2; 427c4762a1bSJed Brown len = PetscMin(m/2,PetscMax(m-(m+1)/2-1,0)); 428c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,len,st,1,&is);CHKERRQ(ierr); 429c4762a1bSJed Brown } else { 430c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,0,0,1,&is);CHKERRQ(ierr); 431c4762a1bSJed Brown } 432c4762a1bSJed Brown } else { 433c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,1,0,1,&is);CHKERRQ(ierr); 434c4762a1bSJed Brown } 435c4762a1bSJed Brown 436c4762a1bSJed Brown if (squaretest) { /* tests for square matrices only, with same maps for rows and columns */ 437c4762a1bSJed Brown /* test MatDiagonalSet */ 438c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatDiagonalSet\n");CHKERRQ(ierr); 439c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr); 440c4762a1bSJed Brown ierr = MatDuplicate(B,MAT_COPY_VALUES,&B2);CHKERRQ(ierr); 441c4762a1bSJed Brown ierr = MatCreateVecs(A,NULL,&x);CHKERRQ(ierr); 442c4762a1bSJed Brown ierr = VecSetRandom(x,NULL);CHKERRQ(ierr); 443c4762a1bSJed Brown ierr = MatDiagonalSet(A2,x,INSERT_VALUES);CHKERRQ(ierr); 444c4762a1bSJed Brown ierr = MatDiagonalSet(B2,x,INSERT_VALUES);CHKERRQ(ierr); 445c4762a1bSJed Brown ierr = CheckMat(A2,B2,PETSC_FALSE,"MatDiagonalSet");CHKERRQ(ierr); 446c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 447c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 448c4762a1bSJed Brown ierr = MatDestroy(&B2);CHKERRQ(ierr); 449c4762a1bSJed Brown 450c4762a1bSJed Brown /* test MatShift (MatShift_IS internally uses MatDiagonalSet_IS with ADD_VALUES) */ 451c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatShift\n");CHKERRQ(ierr); 452c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr); 453c4762a1bSJed Brown ierr = MatDuplicate(B,MAT_COPY_VALUES,&B2);CHKERRQ(ierr); 454c4762a1bSJed Brown ierr = MatShift(A2,2.0);CHKERRQ(ierr); 455c4762a1bSJed Brown ierr = MatShift(B2,2.0);CHKERRQ(ierr); 456c4762a1bSJed Brown ierr = CheckMat(A2,B2,PETSC_FALSE,"MatShift");CHKERRQ(ierr); 457c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 458c4762a1bSJed Brown ierr = MatDestroy(&B2);CHKERRQ(ierr); 459c4762a1bSJed Brown 460c4762a1bSJed Brown /* nonzero diag value is supported for square matrices only */ 461c4762a1bSJed Brown ierr = TestMatZeroRows(A,B,PETSC_TRUE,is,diag);CHKERRQ(ierr); 462c4762a1bSJed Brown } 463c4762a1bSJed Brown ierr = TestMatZeroRows(A,B,squaretest,is,0.0);CHKERRQ(ierr); 464c4762a1bSJed Brown ierr = ISDestroy(&is);CHKERRQ(ierr); 465c4762a1bSJed Brown 466c4762a1bSJed Brown /* test MatTranspose */ 467c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatTranspose\n");CHKERRQ(ierr); 468c4762a1bSJed Brown ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&A2);CHKERRQ(ierr); 469c4762a1bSJed Brown ierr = MatTranspose(B,MAT_INITIAL_MATRIX,&B2);CHKERRQ(ierr); 470c4762a1bSJed Brown ierr = CheckMat(A2,B2,PETSC_FALSE,"initial matrix MatTranspose");CHKERRQ(ierr); 471c4762a1bSJed Brown 472c4762a1bSJed Brown ierr = MatTranspose(A,MAT_REUSE_MATRIX,&A2);CHKERRQ(ierr); 473c4762a1bSJed Brown ierr = CheckMat(A2,B2,PETSC_FALSE,"reuse matrix (not in place) MatTranspose");CHKERRQ(ierr); 474c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 475c4762a1bSJed Brown 476c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr); 477c4762a1bSJed Brown ierr = MatTranspose(A2,MAT_INPLACE_MATRIX,&A2);CHKERRQ(ierr); 478c4762a1bSJed Brown ierr = CheckMat(A2,B2,PETSC_FALSE,"reuse matrix (in place) MatTranspose");CHKERRQ(ierr); 479c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 480c4762a1bSJed Brown 481c4762a1bSJed Brown ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&A2);CHKERRQ(ierr); 482c4762a1bSJed Brown ierr = CheckMat(A2,B2,PETSC_FALSE,"reuse matrix (different type) MatTranspose");CHKERRQ(ierr); 483c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 484c4762a1bSJed Brown ierr = MatDestroy(&B2);CHKERRQ(ierr); 485c4762a1bSJed Brown 486c4762a1bSJed Brown /* test MatISFixLocalEmpty */ 487c4762a1bSJed Brown if (isaij) { 488c4762a1bSJed Brown PetscInt r[2]; 489c4762a1bSJed Brown 490c4762a1bSJed Brown r[0] = 0; 491c4762a1bSJed Brown r[1] = PetscMin(m,n)-1; 492c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatISFixLocalEmpty\n");CHKERRQ(ierr); 493c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr); 494c4762a1bSJed Brown ierr = MatISFixLocalEmpty(A2,PETSC_TRUE);CHKERRQ(ierr); 495c4762a1bSJed Brown ierr = MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 496c4762a1bSJed Brown ierr = MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 497c4762a1bSJed Brown ierr = CheckMat(A2,B,PETSC_FALSE,"MatISFixLocalEmpty (null)");CHKERRQ(ierr); 498c4762a1bSJed Brown 499c4762a1bSJed Brown ierr = MatZeroRows(A2,2,r,0.0,NULL,NULL);CHKERRQ(ierr); 500c4762a1bSJed Brown ierr = MatViewFromOptions(A2,NULL,"-fixempty_view");CHKERRQ(ierr); 501c4762a1bSJed Brown ierr = MatDuplicate(B,MAT_COPY_VALUES,&B2);CHKERRQ(ierr); 502c4762a1bSJed Brown ierr = MatZeroRows(B2,2,r,0.0,NULL,NULL);CHKERRQ(ierr); 503c4762a1bSJed Brown ierr = MatISFixLocalEmpty(A2,PETSC_TRUE);CHKERRQ(ierr); 504c4762a1bSJed Brown ierr = MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 505c4762a1bSJed Brown ierr = MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 506c4762a1bSJed Brown ierr = MatViewFromOptions(A2,NULL,"-fixempty_view");CHKERRQ(ierr); 507c4762a1bSJed Brown ierr = CheckMat(A2,B2,PETSC_FALSE,"MatISFixLocalEmpty (rows)");CHKERRQ(ierr); 508c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 509c4762a1bSJed Brown 510c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr); 511c4762a1bSJed Brown ierr = MatZeroRows(A2,2,r,0.0,NULL,NULL);CHKERRQ(ierr); 512c4762a1bSJed Brown ierr = MatTranspose(A2,MAT_INPLACE_MATRIX,&A2);CHKERRQ(ierr); 513c4762a1bSJed Brown ierr = MatTranspose(B2,MAT_INPLACE_MATRIX,&B2);CHKERRQ(ierr); 514c4762a1bSJed Brown ierr = MatViewFromOptions(A2,NULL,"-fixempty_view");CHKERRQ(ierr); 515c4762a1bSJed Brown ierr = MatISFixLocalEmpty(A2,PETSC_TRUE);CHKERRQ(ierr); 516c4762a1bSJed Brown ierr = MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 517c4762a1bSJed Brown ierr = MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 518c4762a1bSJed Brown ierr = MatViewFromOptions(A2,NULL,"-fixempty_view");CHKERRQ(ierr); 519c4762a1bSJed Brown ierr = CheckMat(A2,B2,PETSC_FALSE,"MatISFixLocalEmpty (cols)");CHKERRQ(ierr); 520c4762a1bSJed Brown 521c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 522c4762a1bSJed Brown ierr = MatDestroy(&B2);CHKERRQ(ierr); 523c4762a1bSJed Brown 524c4762a1bSJed Brown if (squaretest) { 525c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_COPY_VALUES,&A2);CHKERRQ(ierr); 526c4762a1bSJed Brown ierr = MatDuplicate(B,MAT_COPY_VALUES,&B2);CHKERRQ(ierr); 527c4762a1bSJed Brown ierr = MatZeroRowsColumns(A2,2,r,0.0,NULL,NULL);CHKERRQ(ierr); 528c4762a1bSJed Brown ierr = MatZeroRowsColumns(B2,2,r,0.0,NULL,NULL);CHKERRQ(ierr); 529c4762a1bSJed Brown ierr = MatViewFromOptions(A2,NULL,"-fixempty_view");CHKERRQ(ierr); 530c4762a1bSJed Brown ierr = MatISFixLocalEmpty(A2,PETSC_TRUE);CHKERRQ(ierr); 531c4762a1bSJed Brown ierr = MatAssemblyBegin(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 532c4762a1bSJed Brown ierr = MatAssemblyEnd(A2,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 533c4762a1bSJed Brown ierr = MatViewFromOptions(A2,NULL,"-fixempty_view");CHKERRQ(ierr); 534c4762a1bSJed Brown ierr = CheckMat(A2,B2,PETSC_FALSE,"MatISFixLocalEmpty (rows+cols)");CHKERRQ(ierr); 535c4762a1bSJed Brown ierr = MatDestroy(&A2);CHKERRQ(ierr); 536c4762a1bSJed Brown ierr = MatDestroy(&B2);CHKERRQ(ierr); 537c4762a1bSJed Brown } 538c4762a1bSJed Brown 539c4762a1bSJed Brown } 540c4762a1bSJed Brown 541c4762a1bSJed Brown /* test MatInvertBlockDiagonal 542c4762a1bSJed Brown special cases for block-diagonal matrices */ 543c4762a1bSJed Brown if (m == n) { 544c4762a1bSJed Brown ISLocalToGlobalMapping map; 545c4762a1bSJed Brown Mat Abd,Bbd; 546c4762a1bSJed Brown IS is,bis; 547c4762a1bSJed Brown const PetscScalar *isbd,*aijbd; 548c4762a1bSJed Brown PetscScalar *vals; 549c4762a1bSJed Brown const PetscInt *sts,*idxs; 550c4762a1bSJed Brown PetscInt *idxs2,diff,perm,nl,bs,st,en,in; 551c4762a1bSJed Brown PetscBool ok; 552c4762a1bSJed Brown 553c4762a1bSJed Brown for (diff = 0; diff < 3; diff++) { 554c4762a1bSJed Brown for (perm = 0; perm < 3; perm++) { 555c4762a1bSJed Brown for (bs = 1; bs < 4; bs++) { 5568fffc762SJacob Faibussowitsch ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatInvertBlockDiagonal blockdiag %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT " %" PetscInt_FMT "\n",n,diff,perm,bs);CHKERRQ(ierr); 557c4762a1bSJed Brown ierr = PetscMalloc1(bs*bs,&vals);CHKERRQ(ierr); 558c4762a1bSJed Brown ierr = MatGetOwnershipRanges(A,&sts);CHKERRQ(ierr); 559c4762a1bSJed Brown switch (diff) { 560c4762a1bSJed Brown case 1: /* inverted layout by processes */ 561c4762a1bSJed Brown in = 1; 562c4762a1bSJed Brown st = sts[size - rank - 1]; 563c4762a1bSJed Brown en = sts[size - rank]; 564c4762a1bSJed Brown nl = en - st; 565c4762a1bSJed Brown break; 566c4762a1bSJed Brown case 2: /* round-robin layout */ 567c4762a1bSJed Brown in = size; 568c4762a1bSJed Brown st = rank; 569c4762a1bSJed Brown nl = n/size; 570c4762a1bSJed Brown if (rank < n%size) nl++; 571c4762a1bSJed Brown break; 572c4762a1bSJed Brown default: /* same layout */ 573c4762a1bSJed Brown in = 1; 574c4762a1bSJed Brown st = sts[rank]; 575c4762a1bSJed Brown en = sts[rank + 1]; 576c4762a1bSJed Brown nl = en - st; 577c4762a1bSJed Brown break; 578c4762a1bSJed Brown } 579c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,nl,st,in,&is);CHKERRQ(ierr); 580c4762a1bSJed Brown ierr = ISGetLocalSize(is,&nl);CHKERRQ(ierr); 581c4762a1bSJed Brown ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr); 582c4762a1bSJed Brown ierr = PetscMalloc1(nl,&idxs2);CHKERRQ(ierr); 583c4762a1bSJed Brown for (i=0;i<nl;i++) { 584c4762a1bSJed Brown switch (perm) { /* invert some of the indices */ 585c4762a1bSJed Brown case 2: 586c4762a1bSJed Brown idxs2[i] = rank%2 ? idxs[i] : idxs[nl-i-1]; 587c4762a1bSJed Brown break; 588c4762a1bSJed Brown case 1: 589c4762a1bSJed Brown idxs2[i] = rank%2 ? idxs[nl-i-1] : idxs[i]; 590c4762a1bSJed Brown break; 591c4762a1bSJed Brown default: 592c4762a1bSJed Brown idxs2[i] = idxs[i]; 593c4762a1bSJed Brown break; 594c4762a1bSJed Brown } 595c4762a1bSJed Brown } 596c4762a1bSJed Brown ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr); 597c4762a1bSJed Brown ierr = ISCreateBlock(PETSC_COMM_WORLD,bs,nl,idxs2,PETSC_OWN_POINTER,&bis);CHKERRQ(ierr); 598c4762a1bSJed Brown ierr = ISLocalToGlobalMappingCreateIS(bis,&map);CHKERRQ(ierr); 599fc989267SStefano Zampini ierr = MatCreateIS(PETSC_COMM_WORLD,bs,PETSC_DECIDE,PETSC_DECIDE,bs*n,bs*n,map,map,&Abd);CHKERRQ(ierr); 600c4762a1bSJed Brown ierr = ISLocalToGlobalMappingDestroy(&map);CHKERRQ(ierr); 601c4762a1bSJed Brown ierr = MatISSetPreallocation(Abd,bs,NULL,0,NULL);CHKERRQ(ierr); 602c4762a1bSJed Brown for (i=0;i<nl;i++) { 603c4762a1bSJed Brown PetscInt b1,b2; 604c4762a1bSJed Brown 605c4762a1bSJed Brown for (b1=0;b1<bs;b1++) { 606c4762a1bSJed Brown for (b2=0;b2<bs;b2++) { 607c4762a1bSJed Brown vals[b1*bs + b2] = i*bs*bs + b1*bs + b2 + 1 + (b1 == b2 ? 1.0 : 0); 608c4762a1bSJed Brown } 609c4762a1bSJed Brown } 610c4762a1bSJed Brown ierr = MatSetValuesBlockedLocal(Abd,1,&i,1,&i,vals,INSERT_VALUES);CHKERRQ(ierr); 611c4762a1bSJed Brown } 612c4762a1bSJed Brown ierr = MatAssemblyBegin(Abd,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 613c4762a1bSJed Brown ierr = MatAssemblyEnd(Abd,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 614c4762a1bSJed Brown ierr = MatConvert(Abd,MATAIJ,MAT_INITIAL_MATRIX,&Bbd);CHKERRQ(ierr); 615c4762a1bSJed Brown ierr = MatInvertBlockDiagonal(Abd,&isbd);CHKERRQ(ierr); 616c4762a1bSJed Brown ierr = MatInvertBlockDiagonal(Bbd,&aijbd);CHKERRQ(ierr); 617c4762a1bSJed Brown ierr = MatGetLocalSize(Bbd,&nl,NULL);CHKERRQ(ierr); 618c4762a1bSJed Brown ok = PETSC_TRUE; 619c4762a1bSJed Brown for (i=0;i<nl/bs;i++) { 620c4762a1bSJed Brown PetscInt b1,b2; 621c4762a1bSJed Brown 622c4762a1bSJed Brown for (b1=0;b1<bs;b1++) { 623c4762a1bSJed Brown for (b2=0;b2<bs;b2++) { 624c4762a1bSJed Brown if (PetscAbsScalar(isbd[i*bs*bs+b1*bs + b2]-aijbd[i*bs*bs+b1*bs + b2]) > PETSC_SMALL) ok = PETSC_FALSE; 625c4762a1bSJed Brown if (!ok) { 6268fffc762SJacob Faibussowitsch ierr = 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]));CHKERRQ(ierr); 627c4762a1bSJed Brown break; 628c4762a1bSJed Brown } 629c4762a1bSJed Brown } 630c4762a1bSJed Brown if (!ok) break; 631c4762a1bSJed Brown } 632c4762a1bSJed Brown if (!ok) break; 633c4762a1bSJed Brown } 634c4762a1bSJed Brown ierr = MatDestroy(&Abd);CHKERRQ(ierr); 635c4762a1bSJed Brown ierr = MatDestroy(&Bbd);CHKERRQ(ierr); 636c4762a1bSJed Brown ierr = PetscFree(vals);CHKERRQ(ierr); 637c4762a1bSJed Brown ierr = ISDestroy(&is);CHKERRQ(ierr); 638c4762a1bSJed Brown ierr = ISDestroy(&bis);CHKERRQ(ierr); 639c4762a1bSJed Brown } 640c4762a1bSJed Brown } 641c4762a1bSJed Brown } 642c4762a1bSJed Brown } 643c4762a1bSJed Brown /* free testing matrices */ 644c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 645c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 646c4762a1bSJed Brown ierr = PetscFinalize(); 647c4762a1bSJed Brown return ierr; 648c4762a1bSJed Brown } 649c4762a1bSJed Brown 650c4762a1bSJed Brown PetscErrorCode CheckMat(Mat A, Mat B, PetscBool usemult, const char* func) 651c4762a1bSJed Brown { 652c4762a1bSJed Brown Mat Bcheck; 653c4762a1bSJed Brown PetscReal error; 654c4762a1bSJed Brown PetscErrorCode ierr; 655c4762a1bSJed Brown 656c4762a1bSJed Brown PetscFunctionBeginUser; 657c4762a1bSJed Brown if (!usemult && B) { 658c4762a1bSJed Brown PetscBool hasnorm; 659c4762a1bSJed Brown 660c4762a1bSJed Brown ierr = MatHasOperation(B,MATOP_NORM,&hasnorm);CHKERRQ(ierr); 661c4762a1bSJed Brown if (!hasnorm) usemult = PETSC_TRUE; 662c4762a1bSJed Brown } 663c4762a1bSJed Brown if (!usemult) { 664c4762a1bSJed Brown if (B) { 665c4762a1bSJed Brown MatType Btype; 666c4762a1bSJed Brown 667c4762a1bSJed Brown ierr = MatGetType(B,&Btype);CHKERRQ(ierr); 668c4762a1bSJed Brown ierr = MatConvert(A,Btype,MAT_INITIAL_MATRIX,&Bcheck);CHKERRQ(ierr); 669c4762a1bSJed Brown } else { 670c4762a1bSJed Brown ierr = MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&Bcheck);CHKERRQ(ierr); 671c4762a1bSJed Brown } 672c4762a1bSJed Brown if (B) { /* if B is present, subtract it */ 673c4762a1bSJed Brown ierr = MatAXPY(Bcheck,-1.,B,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 674c4762a1bSJed Brown } 675c4762a1bSJed Brown ierr = MatNorm(Bcheck,NORM_INFINITY,&error);CHKERRQ(ierr); 676c4762a1bSJed Brown if (error > PETSC_SQRT_MACHINE_EPSILON) { 677c4762a1bSJed Brown ISLocalToGlobalMapping rl2g,cl2g; 678c4762a1bSJed Brown 679c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)Bcheck,"Assembled Bcheck");CHKERRQ(ierr); 680c4762a1bSJed Brown ierr = MatView(Bcheck,NULL);CHKERRQ(ierr); 681c4762a1bSJed Brown if (B) { 682c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)B,"Assembled AIJ");CHKERRQ(ierr); 683c4762a1bSJed Brown ierr = MatView(B,NULL);CHKERRQ(ierr); 684c4762a1bSJed Brown ierr = MatDestroy(&Bcheck);CHKERRQ(ierr); 685c4762a1bSJed Brown ierr = MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&Bcheck);CHKERRQ(ierr); 686c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)Bcheck,"Assembled IS");CHKERRQ(ierr); 687c4762a1bSJed Brown ierr = MatView(Bcheck,NULL);CHKERRQ(ierr); 688c4762a1bSJed Brown } 689c4762a1bSJed Brown ierr = MatDestroy(&Bcheck);CHKERRQ(ierr); 690c4762a1bSJed Brown ierr = PetscObjectSetName((PetscObject)A,"MatIS");CHKERRQ(ierr); 691c4762a1bSJed Brown ierr = MatView(A,NULL);CHKERRQ(ierr); 692c4762a1bSJed Brown ierr = MatGetLocalToGlobalMapping(A,&rl2g,&cl2g);CHKERRQ(ierr); 693c4762a1bSJed Brown ierr = ISLocalToGlobalMappingView(rl2g,NULL);CHKERRQ(ierr); 694c4762a1bSJed Brown ierr = ISLocalToGlobalMappingView(cl2g,NULL);CHKERRQ(ierr); 695*98921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ERROR ON %s: %g",func,(double)error); 696c4762a1bSJed Brown } 697c4762a1bSJed Brown ierr = MatDestroy(&Bcheck);CHKERRQ(ierr); 698c4762a1bSJed Brown } else { 699c4762a1bSJed Brown PetscBool ok,okt; 700c4762a1bSJed Brown 701c4762a1bSJed Brown ierr = MatMultEqual(A,B,3,&ok);CHKERRQ(ierr); 702c4762a1bSJed Brown ierr = MatMultTransposeEqual(A,B,3,&okt);CHKERRQ(ierr); 703*98921bdaSJacob Faibussowitsch if (!ok || !okt) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ERROR ON %s: mult ok ? %d, multtranspose ok ? %d",func,ok,okt); 704c4762a1bSJed Brown } 705c4762a1bSJed Brown PetscFunctionReturn(0); 706c4762a1bSJed Brown } 707c4762a1bSJed Brown 708c4762a1bSJed Brown PetscErrorCode TestMatZeroRows(Mat A, Mat Afull, PetscBool squaretest, IS is, PetscScalar diag) 709c4762a1bSJed Brown { 710c4762a1bSJed Brown Mat B,Bcheck,B2 = NULL,lB; 711c4762a1bSJed Brown Vec x = NULL, b = NULL, b2 = NULL; 712c4762a1bSJed Brown ISLocalToGlobalMapping l2gr,l2gc; 713c4762a1bSJed Brown PetscReal error; 714c4762a1bSJed Brown char diagstr[16]; 715c4762a1bSJed Brown const PetscInt *idxs; 716c4762a1bSJed Brown PetscInt rst,ren,i,n,N,d; 717c4762a1bSJed Brown PetscMPIInt rank; 718c4762a1bSJed Brown PetscBool miss,haszerorows; 719c4762a1bSJed Brown PetscErrorCode ierr; 720c4762a1bSJed Brown 721c4762a1bSJed Brown PetscFunctionBeginUser; 722c4762a1bSJed Brown if (diag == 0.) { 723c4762a1bSJed Brown ierr = PetscStrcpy(diagstr,"zero");CHKERRQ(ierr); 724c4762a1bSJed Brown } else { 725c4762a1bSJed Brown ierr = PetscStrcpy(diagstr,"nonzero");CHKERRQ(ierr); 726c4762a1bSJed Brown } 727c4762a1bSJed Brown ierr = ISView(is,NULL);CHKERRQ(ierr); 728c4762a1bSJed Brown ierr = MatGetLocalToGlobalMapping(A,&l2gr,&l2gc);CHKERRQ(ierr); 729c4762a1bSJed Brown /* tests MatDuplicate and MatCopy */ 730c4762a1bSJed Brown if (diag == 0.) { 731c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 732c4762a1bSJed Brown } else { 733c4762a1bSJed Brown ierr = MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&B);CHKERRQ(ierr); 734c4762a1bSJed Brown ierr = MatCopy(A,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 735c4762a1bSJed Brown } 736c4762a1bSJed Brown ierr = MatISGetLocalMat(B,&lB);CHKERRQ(ierr); 737c4762a1bSJed Brown ierr = MatHasOperation(lB,MATOP_ZERO_ROWS,&haszerorows);CHKERRQ(ierr); 738c4762a1bSJed Brown if (squaretest && haszerorows) { 739c4762a1bSJed Brown 740c4762a1bSJed Brown ierr = MatCreateVecs(B,&x,&b);CHKERRQ(ierr); 741c4762a1bSJed Brown ierr = MatDuplicate(B,MAT_COPY_VALUES,&B2);CHKERRQ(ierr); 742c4762a1bSJed Brown ierr = VecSetLocalToGlobalMapping(b,l2gr);CHKERRQ(ierr); 743c4762a1bSJed Brown ierr = VecSetLocalToGlobalMapping(x,l2gc);CHKERRQ(ierr); 744c4762a1bSJed Brown ierr = VecSetRandom(x,NULL);CHKERRQ(ierr); 745c4762a1bSJed Brown ierr = VecSetRandom(b,NULL);CHKERRQ(ierr); 746c4762a1bSJed Brown /* mimic b[is] = x[is] */ 747c4762a1bSJed Brown ierr = VecDuplicate(b,&b2);CHKERRQ(ierr); 748c4762a1bSJed Brown ierr = VecSetLocalToGlobalMapping(b2,l2gr);CHKERRQ(ierr); 749c4762a1bSJed Brown ierr = VecCopy(b,b2);CHKERRQ(ierr); 750c4762a1bSJed Brown ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr); 751c4762a1bSJed Brown ierr = ISGetIndices(is,&idxs);CHKERRQ(ierr); 752c4762a1bSJed Brown ierr = VecGetSize(x,&N);CHKERRQ(ierr); 753c4762a1bSJed Brown for (i=0;i<n;i++) { 754c4762a1bSJed Brown if (0 <= idxs[i] && idxs[i] < N) { 755c4762a1bSJed Brown ierr = VecSetValue(b2,idxs[i],diag,INSERT_VALUES);CHKERRQ(ierr); 756c4762a1bSJed Brown ierr = VecSetValue(x,idxs[i],1.,INSERT_VALUES);CHKERRQ(ierr); 757c4762a1bSJed Brown } 758c4762a1bSJed Brown } 759c4762a1bSJed Brown ierr = VecAssemblyBegin(b2);CHKERRQ(ierr); 760c4762a1bSJed Brown ierr = VecAssemblyEnd(b2);CHKERRQ(ierr); 761c4762a1bSJed Brown ierr = VecAssemblyBegin(x);CHKERRQ(ierr); 762c4762a1bSJed Brown ierr = VecAssemblyEnd(x);CHKERRQ(ierr); 763c4762a1bSJed Brown ierr = ISRestoreIndices(is,&idxs);CHKERRQ(ierr); 764c4762a1bSJed Brown /* test ZeroRows on MATIS */ 765c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatZeroRows (diag %s)\n",diagstr);CHKERRQ(ierr); 766c4762a1bSJed Brown ierr = MatZeroRowsIS(B,is,diag,x,b);CHKERRQ(ierr); 767c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatZeroRowsColumns (diag %s)\n",diagstr);CHKERRQ(ierr); 768c4762a1bSJed Brown ierr = MatZeroRowsColumnsIS(B2,is,diag,NULL,NULL);CHKERRQ(ierr); 769c4762a1bSJed Brown } else if (haszerorows) { 770c4762a1bSJed Brown /* test ZeroRows on MATIS */ 771c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatZeroRows (diag %s)\n",diagstr);CHKERRQ(ierr); 772c4762a1bSJed Brown ierr = MatZeroRowsIS(B,is,diag,NULL,NULL);CHKERRQ(ierr); 773c4762a1bSJed Brown b = b2 = x = NULL; 774c4762a1bSJed Brown } else { 775c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Skipping MatZeroRows (diag %s)\n",diagstr);CHKERRQ(ierr); 776c4762a1bSJed Brown b = b2 = x = NULL; 777c4762a1bSJed Brown } 778c4762a1bSJed Brown 779c4762a1bSJed Brown if (squaretest && haszerorows) { 780c4762a1bSJed Brown ierr = VecAXPY(b2,-1.,b);CHKERRQ(ierr); 781c4762a1bSJed Brown ierr = VecNorm(b2,NORM_INFINITY,&error);CHKERRQ(ierr); 782*98921bdaSJacob Faibussowitsch if (error > PETSC_SQRT_MACHINE_EPSILON) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"ERROR IN ZEROROWS ON B %g (diag %s)",(double)error,diagstr); 783c4762a1bSJed Brown } 784c4762a1bSJed Brown 785c4762a1bSJed Brown /* test MatMissingDiagonal */ 786c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Test MatMissingDiagonal\n");CHKERRQ(ierr); 787ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 788c4762a1bSJed Brown ierr = MatMissingDiagonal(B,&miss,&d);CHKERRQ(ierr); 789c4762a1bSJed Brown ierr = MatGetOwnershipRange(B,&rst,&ren);CHKERRQ(ierr); 790c4762a1bSJed Brown ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 7918fffc762SJacob Faibussowitsch ierr = 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);CHKERRQ(ierr); 792c4762a1bSJed Brown ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 793c4762a1bSJed Brown ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 794c4762a1bSJed Brown 795c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 796c4762a1bSJed Brown ierr = VecDestroy(&b);CHKERRQ(ierr); 797c4762a1bSJed Brown ierr = VecDestroy(&b2);CHKERRQ(ierr); 798c4762a1bSJed Brown 799c4762a1bSJed Brown /* check the result of ZeroRows with that from MPIAIJ routines 800c4762a1bSJed Brown assuming that MatConvert_IS_XAIJ and MatZeroRows_MPIAIJ work fine */ 801c4762a1bSJed Brown if (haszerorows) { 802c4762a1bSJed Brown ierr = MatDuplicate(Afull,MAT_COPY_VALUES,&Bcheck);CHKERRQ(ierr); 803c4762a1bSJed Brown ierr = MatSetOption(Bcheck,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 804c4762a1bSJed Brown ierr = MatZeroRowsIS(Bcheck,is,diag,NULL,NULL);CHKERRQ(ierr); 805c4762a1bSJed Brown ierr = CheckMat(B,Bcheck,PETSC_FALSE,"Zerorows");CHKERRQ(ierr); 806c4762a1bSJed Brown ierr = MatDestroy(&Bcheck);CHKERRQ(ierr); 807c4762a1bSJed Brown } 808c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 809c4762a1bSJed Brown 810c4762a1bSJed Brown if (B2) { /* test MatZeroRowsColumns */ 811c4762a1bSJed Brown ierr = MatDuplicate(Afull,MAT_COPY_VALUES,&B);CHKERRQ(ierr); 812c4762a1bSJed Brown ierr = MatSetOption(B,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 813c4762a1bSJed Brown ierr = MatZeroRowsColumnsIS(B,is,diag,NULL,NULL);CHKERRQ(ierr); 814c4762a1bSJed Brown ierr = CheckMat(B2,B,PETSC_FALSE,"MatZeroRowsColumns");CHKERRQ(ierr); 815c4762a1bSJed Brown ierr = MatDestroy(&B);CHKERRQ(ierr); 816c4762a1bSJed Brown ierr = MatDestroy(&B2);CHKERRQ(ierr); 817c4762a1bSJed Brown } 818c4762a1bSJed Brown PetscFunctionReturn(0); 819c4762a1bSJed Brown } 820c4762a1bSJed Brown 821c4762a1bSJed Brown /*TEST 822c4762a1bSJed Brown 823c4762a1bSJed Brown test: 824c4762a1bSJed Brown args: -test_trans 825c4762a1bSJed Brown 826c4762a1bSJed Brown test: 827c4762a1bSJed Brown suffix: 2 828c4762a1bSJed Brown nsize: 4 829c4762a1bSJed Brown args: -matis_convert_local_nest -nr 3 -nc 4 830c4762a1bSJed Brown 831c4762a1bSJed Brown test: 832c4762a1bSJed Brown suffix: 3 833c4762a1bSJed Brown nsize: 5 834c4762a1bSJed Brown args: -m 11 -n 10 -matis_convert_local_nest -nr 2 -nc 1 835c4762a1bSJed Brown 836c4762a1bSJed Brown test: 837c4762a1bSJed Brown suffix: 4 838c4762a1bSJed Brown nsize: 6 839c4762a1bSJed Brown args: -m 9 -n 12 -test_trans -nr 2 -nc 7 840c4762a1bSJed Brown 841c4762a1bSJed Brown test: 842c4762a1bSJed Brown suffix: 5 843c4762a1bSJed Brown nsize: 6 844c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 3 -nc 1 845c4762a1bSJed Brown 846c4762a1bSJed Brown test: 847c4762a1bSJed Brown suffix: 6 848c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap 849c4762a1bSJed Brown 850c4762a1bSJed Brown test: 851c4762a1bSJed Brown suffix: 7 852c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap -permmap 853c4762a1bSJed Brown 854c4762a1bSJed Brown test: 855c4762a1bSJed Brown suffix: 8 856c4762a1bSJed Brown args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -permmap 857c4762a1bSJed Brown 858c4762a1bSJed Brown test: 859c4762a1bSJed Brown suffix: 9 860c4762a1bSJed Brown nsize: 5 861c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap 862c4762a1bSJed Brown 863c4762a1bSJed Brown test: 864c4762a1bSJed Brown suffix: 10 865c4762a1bSJed Brown nsize: 5 866c4762a1bSJed Brown args: -m 12 -n 12 -test_trans -nr 2 -nc 3 -diffmap -permmap 867c4762a1bSJed Brown 868c20d7725SJed Brown test: 869c20d7725SJed Brown suffix: vscat_default 870c4762a1bSJed Brown nsize: 5 871c4762a1bSJed Brown args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -permmap 872c4762a1bSJed Brown output_file: output/ex23_11.out 873c4762a1bSJed Brown 874c4762a1bSJed Brown test: 875c4762a1bSJed Brown suffix: 12 876c4762a1bSJed Brown nsize: 3 877c4762a1bSJed Brown args: -m 12 -n 12 -symmetric -matis_localmat_type sbaij -test_trans -nr 2 -nc 3 878c4762a1bSJed Brown 879c4762a1bSJed Brown testset: 880c4762a1bSJed Brown output_file: output/ex23_13.out 881c4762a1bSJed Brown nsize: 3 882c4762a1bSJed Brown args: -m 12 -n 17 -test_trans -nr 2 -nc 3 -diffmap -permmap 883c4762a1bSJed Brown filter: grep -v "type:" 884c4762a1bSJed Brown test: 885c4762a1bSJed Brown suffix: baij 886c4762a1bSJed Brown args: -matis_localmat_type baij 887c4762a1bSJed Brown test: 888c4762a1bSJed Brown requires: viennacl 889c4762a1bSJed Brown suffix: viennacl 890c4762a1bSJed Brown args: -matis_localmat_type aijviennacl 891c4762a1bSJed Brown test: 892c4762a1bSJed Brown requires: cuda 893c4762a1bSJed Brown suffix: cusparse 894c4762a1bSJed Brown args: -matis_localmat_type aijcusparse 895c4762a1bSJed Brown 896c4762a1bSJed Brown TEST*/ 897