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