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