1 2 /* 3 Defines the basic matrix operations for the ADJ adjacency list matrix data-structure. 4 */ 5 #include <../src/mat/impls/adj/mpi/mpiadj.h> /*I "petscmat.h" I*/ 6 #include <petscsf.h> 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "MatView_MPIAdj_ASCII" 10 PetscErrorCode MatView_MPIAdj_ASCII(Mat A,PetscViewer viewer) 11 { 12 Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 13 PetscErrorCode ierr; 14 PetscInt i,j,m = A->rmap->n; 15 const char *name; 16 PetscViewerFormat format; 17 18 PetscFunctionBegin; 19 ierr = PetscObjectGetName((PetscObject)A,&name);CHKERRQ(ierr); 20 ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 21 if (format == PETSC_VIEWER_ASCII_INFO) { 22 PetscFunctionReturn(0); 23 } else if (format == PETSC_VIEWER_ASCII_MATLAB) { 24 SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MATLAB format not supported"); 25 } else { 26 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 27 ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 28 for (i=0; i<m; i++) { 29 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"row %D:",i+A->rmap->rstart);CHKERRQ(ierr); 30 for (j=a->i[i]; j<a->i[i+1]; j++) { 31 ierr = PetscViewerASCIISynchronizedPrintf(viewer," %D ",a->j[j]);CHKERRQ(ierr); 32 } 33 ierr = PetscViewerASCIISynchronizedPrintf(viewer,"\n");CHKERRQ(ierr); 34 } 35 ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 36 ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 37 ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 38 } 39 PetscFunctionReturn(0); 40 } 41 42 #undef __FUNCT__ 43 #define __FUNCT__ "MatView_MPIAdj" 44 PetscErrorCode MatView_MPIAdj(Mat A,PetscViewer viewer) 45 { 46 PetscErrorCode ierr; 47 PetscBool iascii; 48 49 PetscFunctionBegin; 50 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 51 if (iascii) { 52 ierr = MatView_MPIAdj_ASCII(A,viewer);CHKERRQ(ierr); 53 } 54 PetscFunctionReturn(0); 55 } 56 57 #undef __FUNCT__ 58 #define __FUNCT__ "MatDestroy_MPIAdj" 59 PetscErrorCode MatDestroy_MPIAdj(Mat mat) 60 { 61 Mat_MPIAdj *a = (Mat_MPIAdj*)mat->data; 62 PetscErrorCode ierr; 63 64 PetscFunctionBegin; 65 #if defined(PETSC_USE_LOG) 66 PetscLogObjectState((PetscObject)mat,"Rows=%D, Cols=%D, NZ=%D",mat->rmap->n,mat->cmap->n,a->nz); 67 #endif 68 ierr = PetscFree(a->diag);CHKERRQ(ierr); 69 if (a->freeaij) { 70 if (a->freeaijwithfree) { 71 if (a->i) free(a->i); 72 if (a->j) free(a->j); 73 } else { 74 ierr = PetscFree(a->i);CHKERRQ(ierr); 75 ierr = PetscFree(a->j);CHKERRQ(ierr); 76 ierr = PetscFree(a->values);CHKERRQ(ierr); 77 } 78 } 79 ierr = PetscFree(a->rowvalues);CHKERRQ(ierr); 80 ierr = PetscFree(mat->data);CHKERRQ(ierr); 81 ierr = PetscObjectChangeTypeName((PetscObject)mat,0);CHKERRQ(ierr); 82 ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjSetPreallocation_C",NULL);CHKERRQ(ierr); 83 ierr = PetscObjectComposeFunction((PetscObject)mat,"MatMPIAdjCreateNonemptySubcommMat_C",NULL);CHKERRQ(ierr); 84 PetscFunctionReturn(0); 85 } 86 87 #undef __FUNCT__ 88 #define __FUNCT__ "MatSetOption_MPIAdj" 89 PetscErrorCode MatSetOption_MPIAdj(Mat A,MatOption op,PetscBool flg) 90 { 91 Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 92 PetscErrorCode ierr; 93 94 PetscFunctionBegin; 95 switch (op) { 96 case MAT_SYMMETRIC: 97 case MAT_STRUCTURALLY_SYMMETRIC: 98 case MAT_HERMITIAN: 99 a->symmetric = flg; 100 break; 101 case MAT_SYMMETRY_ETERNAL: 102 break; 103 default: 104 ierr = PetscInfo1(A,"Option %s ignored\n",MatOptions[op]);CHKERRQ(ierr); 105 break; 106 } 107 PetscFunctionReturn(0); 108 } 109 110 111 /* 112 Adds diagonal pointers to sparse matrix structure. 113 */ 114 115 #undef __FUNCT__ 116 #define __FUNCT__ "MatMarkDiagonal_MPIAdj" 117 PetscErrorCode MatMarkDiagonal_MPIAdj(Mat A) 118 { 119 Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 120 PetscErrorCode ierr; 121 PetscInt i,j,m = A->rmap->n; 122 123 PetscFunctionBegin; 124 ierr = PetscMalloc1(m,&a->diag);CHKERRQ(ierr); 125 ierr = PetscLogObjectMemory((PetscObject)A,m*sizeof(PetscInt));CHKERRQ(ierr); 126 for (i=0; i<A->rmap->n; i++) { 127 for (j=a->i[i]; j<a->i[i+1]; j++) { 128 if (a->j[j] == i) { 129 a->diag[i] = j; 130 break; 131 } 132 } 133 } 134 PetscFunctionReturn(0); 135 } 136 137 #undef __FUNCT__ 138 #define __FUNCT__ "MatGetRow_MPIAdj" 139 PetscErrorCode MatGetRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 140 { 141 Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 142 PetscErrorCode ierr; 143 144 PetscFunctionBegin; 145 row -= A->rmap->rstart; 146 if (row < 0 || row >= A->rmap->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Row out of range"); 147 *nz = a->i[row+1] - a->i[row]; 148 if (v) { 149 PetscInt j; 150 if (a->rowvalues_alloc < *nz) { 151 ierr = PetscFree(a->rowvalues);CHKERRQ(ierr); 152 a->rowvalues_alloc = PetscMax(a->rowvalues_alloc*2, *nz); 153 ierr = PetscMalloc1(a->rowvalues_alloc,&a->rowvalues);CHKERRQ(ierr); 154 } 155 for (j=0; j<*nz; j++){ 156 a->rowvalues[j] = a->values ? a->values[a->i[row]+j]:1.0; 157 } 158 *v = (*nz) ? a->rowvalues : NULL; 159 } 160 if (idx) *idx = (*nz) ? a->j + a->i[row] : NULL; 161 PetscFunctionReturn(0); 162 } 163 164 #undef __FUNCT__ 165 #define __FUNCT__ "MatRestoreRow_MPIAdj" 166 PetscErrorCode MatRestoreRow_MPIAdj(Mat A,PetscInt row,PetscInt *nz,PetscInt **idx,PetscScalar **v) 167 { 168 PetscFunctionBegin; 169 PetscFunctionReturn(0); 170 } 171 172 #undef __FUNCT__ 173 #define __FUNCT__ "MatEqual_MPIAdj" 174 PetscErrorCode MatEqual_MPIAdj(Mat A,Mat B,PetscBool * flg) 175 { 176 Mat_MPIAdj *a = (Mat_MPIAdj*)A->data,*b = (Mat_MPIAdj*)B->data; 177 PetscErrorCode ierr; 178 PetscBool flag; 179 180 PetscFunctionBegin; 181 /* If the matrix dimensions are not equal,or no of nonzeros */ 182 if ((A->rmap->n != B->rmap->n) ||(a->nz != b->nz)) { 183 flag = PETSC_FALSE; 184 } 185 186 /* if the a->i are the same */ 187 ierr = PetscMemcmp(a->i,b->i,(A->rmap->n+1)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 188 189 /* if a->j are the same */ 190 ierr = PetscMemcmp(a->j,b->j,(a->nz)*sizeof(PetscInt),&flag);CHKERRQ(ierr); 191 192 ierr = MPI_Allreduce(&flag,flg,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)A));CHKERRQ(ierr); 193 PetscFunctionReturn(0); 194 } 195 196 #undef __FUNCT__ 197 #define __FUNCT__ "MatGetRowIJ_MPIAdj" 198 PetscErrorCode MatGetRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 199 { 200 PetscInt i; 201 Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 202 PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 203 204 PetscFunctionBegin; 205 *m = A->rmap->n; 206 *ia = a->i; 207 *ja = a->j; 208 *done = PETSC_TRUE; 209 if (oshift) { 210 for (i=0; i<(*ia)[*m]; i++) { 211 (*ja)[i]++; 212 } 213 for (i=0; i<=(*m); i++) (*ia)[i]++; 214 } 215 PetscFunctionReturn(0); 216 } 217 218 #undef __FUNCT__ 219 #define __FUNCT__ "MatRestoreRowIJ_MPIAdj" 220 PetscErrorCode MatRestoreRowIJ_MPIAdj(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool blockcompressed,PetscInt *m,const PetscInt *inia[],const PetscInt *inja[],PetscBool *done) 221 { 222 PetscInt i; 223 Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 224 PetscInt **ia = (PetscInt**)inia,**ja = (PetscInt**)inja; 225 226 PetscFunctionBegin; 227 if (ia && a->i != *ia) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ia passed back is not one obtained with MatGetRowIJ()"); 228 if (ja && a->j != *ja) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ja passed back is not one obtained with MatGetRowIJ()"); 229 if (oshift) { 230 for (i=0; i<=(*m); i++) (*ia)[i]--; 231 for (i=0; i<(*ia)[*m]; i++) { 232 (*ja)[i]--; 233 } 234 } 235 PetscFunctionReturn(0); 236 } 237 238 #undef __FUNCT__ 239 #define __FUNCT__ "MatConvertFrom_MPIAdj" 240 PetscErrorCode MatConvertFrom_MPIAdj(Mat A,MatType type,MatReuse reuse,Mat *newmat) 241 { 242 Mat B; 243 PetscErrorCode ierr; 244 PetscInt i,m,N,nzeros = 0,*ia,*ja,len,rstart,cnt,j,*a; 245 const PetscInt *rj; 246 const PetscScalar *ra; 247 MPI_Comm comm; 248 249 PetscFunctionBegin; 250 ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr); 251 ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 252 ierr = MatGetOwnershipRange(A,&rstart,NULL);CHKERRQ(ierr); 253 254 /* count the number of nonzeros per row */ 255 for (i=0; i<m; i++) { 256 ierr = MatGetRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr); 257 for (j=0; j<len; j++) { 258 if (rj[j] == i+rstart) {len--; break;} /* don't count diagonal */ 259 } 260 nzeros += len; 261 ierr = MatRestoreRow(A,i+rstart,&len,&rj,NULL);CHKERRQ(ierr); 262 } 263 264 /* malloc space for nonzeros */ 265 ierr = PetscMalloc1(nzeros+1,&a);CHKERRQ(ierr); 266 ierr = PetscMalloc1(N+1,&ia);CHKERRQ(ierr); 267 ierr = PetscMalloc1(nzeros+1,&ja);CHKERRQ(ierr); 268 269 nzeros = 0; 270 ia[0] = 0; 271 for (i=0; i<m; i++) { 272 ierr = MatGetRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 273 cnt = 0; 274 for (j=0; j<len; j++) { 275 if (rj[j] != i+rstart) { /* if not diagonal */ 276 a[nzeros+cnt] = (PetscInt) PetscAbsScalar(ra[j]); 277 ja[nzeros+cnt++] = rj[j]; 278 } 279 } 280 ierr = MatRestoreRow(A,i+rstart,&len,&rj,&ra);CHKERRQ(ierr); 281 nzeros += cnt; 282 ia[i+1] = nzeros; 283 } 284 285 ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 286 ierr = MatCreate(comm,&B);CHKERRQ(ierr); 287 ierr = MatSetSizes(B,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 288 ierr = MatSetType(B,type);CHKERRQ(ierr); 289 ierr = MatMPIAdjSetPreallocation(B,ia,ja,a);CHKERRQ(ierr); 290 291 if (reuse == MAT_REUSE_MATRIX) { 292 ierr = MatHeaderReplace(A,B);CHKERRQ(ierr); 293 } else { 294 *newmat = B; 295 } 296 PetscFunctionReturn(0); 297 } 298 299 /* -------------------------------------------------------------------*/ 300 static struct _MatOps MatOps_Values = {0, 301 MatGetRow_MPIAdj, 302 MatRestoreRow_MPIAdj, 303 0, 304 /* 4*/ 0, 305 0, 306 0, 307 0, 308 0, 309 0, 310 /*10*/ 0, 311 0, 312 0, 313 0, 314 0, 315 /*15*/ 0, 316 MatEqual_MPIAdj, 317 0, 318 0, 319 0, 320 /*20*/ 0, 321 0, 322 MatSetOption_MPIAdj, 323 0, 324 /*24*/ 0, 325 0, 326 0, 327 0, 328 0, 329 /*29*/ 0, 330 0, 331 0, 332 0, 333 0, 334 /*34*/ 0, 335 0, 336 0, 337 0, 338 0, 339 /*39*/ 0, 340 MatGetSubMatrices_MPIAdj, 341 0, 342 0, 343 0, 344 /*44*/ 0, 345 0, 346 MatShift_Basic, 347 0, 348 0, 349 /*49*/ 0, 350 MatGetRowIJ_MPIAdj, 351 MatRestoreRowIJ_MPIAdj, 352 0, 353 0, 354 /*54*/ 0, 355 0, 356 0, 357 0, 358 0, 359 /*59*/ 0, 360 MatDestroy_MPIAdj, 361 MatView_MPIAdj, 362 MatConvertFrom_MPIAdj, 363 0, 364 /*64*/ 0, 365 0, 366 0, 367 0, 368 0, 369 /*69*/ 0, 370 0, 371 0, 372 0, 373 0, 374 /*74*/ 0, 375 0, 376 0, 377 0, 378 0, 379 /*79*/ 0, 380 0, 381 0, 382 0, 383 0, 384 /*84*/ 0, 385 0, 386 0, 387 0, 388 0, 389 /*89*/ 0, 390 0, 391 0, 392 0, 393 0, 394 /*94*/ 0, 395 0, 396 0, 397 0, 398 0, 399 /*99*/ 0, 400 0, 401 0, 402 0, 403 0, 404 /*104*/ 0, 405 0, 406 0, 407 0, 408 0, 409 /*109*/ 0, 410 0, 411 0, 412 0, 413 0, 414 /*114*/ 0, 415 0, 416 0, 417 0, 418 0, 419 /*119*/ 0, 420 0, 421 0, 422 0, 423 0, 424 /*124*/ 0, 425 0, 426 0, 427 0, 428 MatGetSubMatricesMPI_MPIAdj, 429 /*129*/ 0, 430 0, 431 0, 432 0, 433 0, 434 /*134*/ 0, 435 0, 436 0, 437 0, 438 0, 439 /*139*/ 0, 440 0, 441 0 442 }; 443 444 #undef __FUNCT__ 445 #define __FUNCT__ "MatMPIAdjSetPreallocation_MPIAdj" 446 static PetscErrorCode MatMPIAdjSetPreallocation_MPIAdj(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 447 { 448 Mat_MPIAdj *b = (Mat_MPIAdj*)B->data; 449 PetscErrorCode ierr; 450 #if defined(PETSC_USE_DEBUG) 451 PetscInt ii; 452 #endif 453 454 PetscFunctionBegin; 455 ierr = PetscLayoutSetUp(B->rmap);CHKERRQ(ierr); 456 ierr = PetscLayoutSetUp(B->cmap);CHKERRQ(ierr); 457 458 #if defined(PETSC_USE_DEBUG) 459 if (i[0] != 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"First i[] index must be zero, instead it is %D\n",i[0]); 460 for (ii=1; ii<B->rmap->n; ii++) { 461 if (i[ii] < 0 || i[ii] < i[ii-1]) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i[%D]=%D index is out of range: i[%D]=%D",ii,i[ii],ii-1,i[ii-1]); 462 } 463 for (ii=0; ii<i[B->rmap->n]; ii++) { 464 if (j[ii] < 0 || j[ii] >= B->cmap->N) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Column index %D out of range %D\n",ii,j[ii]); 465 } 466 #endif 467 B->preallocated = PETSC_TRUE; 468 469 b->j = j; 470 b->i = i; 471 b->values = values; 472 473 b->nz = i[B->rmap->n]; 474 b->diag = 0; 475 b->symmetric = PETSC_FALSE; 476 b->freeaij = PETSC_TRUE; 477 478 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 479 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 480 PetscFunctionReturn(0); 481 } 482 483 static PetscErrorCode MatGetSubMatrix_MPIAdj_data(Mat,IS,IS, PetscInt **,PetscInt **,PetscInt **); 484 static PetscErrorCode MatGetSubMatrices_MPIAdj_Private(Mat,PetscInt,const IS[],const IS[],PetscBool,MatReuse,Mat **); 485 486 #undef __FUNCT__ 487 #define __FUNCT__ "MatGetSubMatricesMPI_MPIAdj" 488 PetscErrorCode MatGetSubMatricesMPI_MPIAdj(Mat mat,PetscInt n, const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 489 { 490 PetscErrorCode ierr; 491 /*get sub-matrices across a sub communicator */ 492 PetscFunctionBegin; 493 ierr = MatGetSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_TRUE,scall,submat);CHKERRQ(ierr); 494 PetscFunctionReturn(0); 495 } 496 497 498 #undef __FUNCT__ 499 #define __FUNCT__ "MatGetSubMatrices_MPIAdj" 500 PetscErrorCode MatGetSubMatrices_MPIAdj(Mat mat,PetscInt n,const IS irow[],const IS icol[],MatReuse scall,Mat *submat[]) 501 { 502 PetscErrorCode ierr; 503 504 PetscFunctionBegin; 505 /*get sub-matrices based on PETSC_COMM_SELF */ 506 ierr = MatGetSubMatrices_MPIAdj_Private(mat,n,irow,icol,PETSC_FALSE,scall,submat);CHKERRQ(ierr); 507 PetscFunctionReturn(0); 508 } 509 510 #undef __FUNCT__ 511 #define __FUNCT__ "MatGetSubMatrices_MPIAdj_Private" 512 static PetscErrorCode MatGetSubMatrices_MPIAdj_Private(Mat mat,PetscInt n,const IS irow[],const IS icol[],PetscBool subcomm,MatReuse scall,Mat *submat[]) 513 { 514 PetscInt i,irow_n,icol_n,*sxadj,*sadjncy,*svalues; 515 PetscInt *indices,nindx,j,k,loc; 516 PetscMPIInt issame; 517 const PetscInt *irow_indices,*icol_indices; 518 MPI_Comm scomm_row,scomm_col,scomm_mat; 519 PetscErrorCode ierr; 520 521 PetscFunctionBegin; 522 nindx = 0; 523 /* 524 * Estimate a maximum number for allocating memory 525 */ 526 for(i=0; i<n; i++){ 527 ierr = ISGetLocalSize(irow[i],&irow_n);CHKERRQ(ierr); 528 ierr = ISGetLocalSize(icol[i],&icol_n);CHKERRQ(ierr); 529 nindx = nindx>(irow_n+icol_n)? nindx:(irow_n+icol_n); 530 } 531 ierr = PetscCalloc1(nindx,&indices);CHKERRQ(ierr); 532 /* construct a submat */ 533 for(i=0; i<n; i++){ 534 /*comms */ 535 if(subcomm){ 536 ierr = PetscObjectGetComm((PetscObject)irow[i],&scomm_row);CHKERRQ(ierr); 537 ierr = PetscObjectGetComm((PetscObject)icol[i],&scomm_col);CHKERRQ(ierr); 538 ierr = MPI_Comm_compare(scomm_row,scomm_col,&issame);CHKERRQ(ierr); 539 if(issame != MPI_IDENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"row index set must have the same comm as the col index set\n"); 540 ierr = MPI_Comm_compare(scomm_row,PETSC_COMM_SELF,&issame);CHKERRQ(ierr); 541 if(issame == MPI_IDENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP," can not use PETSC_COMM_SELF as comm when extracting a parallel submatrix\n"); 542 }else{ 543 scomm_row = PETSC_COMM_SELF; 544 } 545 /*get sub-matrix data*/ 546 ierr = MatGetSubMatrix_MPIAdj_data(mat,irow[i],icol[i],&sxadj,&sadjncy,&svalues);CHKERRQ(ierr); 547 ierr = ISGetLocalSize(irow[i],&irow_n);CHKERRQ(ierr); 548 ierr = ISGetLocalSize(icol[i],&icol_n);CHKERRQ(ierr); 549 ierr = ISGetIndices(irow[i],&irow_indices);CHKERRQ(ierr); 550 ierr = PetscMemcpy(indices,irow_indices,sizeof(PetscInt)*irow_n);CHKERRQ(ierr); 551 ierr = ISRestoreIndices(irow[i],&irow_indices);CHKERRQ(ierr); 552 ierr = ISGetIndices(icol[i],&icol_indices);CHKERRQ(ierr); 553 ierr = PetscMemcpy(indices+irow_n,icol_indices,sizeof(PetscInt)*icol_n);CHKERRQ(ierr); 554 ierr = ISRestoreIndices(icol[i],&icol_indices);CHKERRQ(ierr); 555 nindx = irow_n+icol_n; 556 ierr = PetscSortRemoveDupsInt(&nindx,indices);CHKERRQ(ierr); 557 /* renumber columns */ 558 for(j=0; j<irow_n; j++){ 559 for(k=sxadj[j]; k<sxadj[j+1]; k++){ 560 ierr = PetscFindInt(sadjncy[k],nindx,indices,&loc);CHKERRQ(ierr); 561 #if PETSC_USE_DEBUG 562 if(loc<0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"can not find col %d \n",sadjncy[k]); 563 #endif 564 sadjncy[k] = loc; 565 } 566 } 567 if(scall==MAT_INITIAL_MATRIX){ 568 ierr = MatCreateMPIAdj(scomm_row,irow_n,icol_n,sxadj,sadjncy,svalues,submat[i]);CHKERRQ(ierr); 569 }else{ 570 Mat sadj = *(submat[i]); 571 ierr = PetscObjectGetComm((PetscObject)sadj,&scomm_mat);CHKERRQ(ierr); 572 ierr = MPI_Comm_compare(scomm_row,scomm_mat,&issame);CHKERRQ(ierr); 573 if(issame != MPI_IDENT) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"submatrix must have the same comm as the col index set\n"); 574 Mat_MPIAdj *sa = (Mat_MPIAdj*)((sadj)->data); 575 ierr = PetscMemcpy(sa->i,sxadj,sizeof(PetscInt)*(irow_n+1));CHKERRQ(ierr); 576 ierr = PetscMemcpy(sa->j,sadjncy,sizeof(PetscInt)*sxadj[irow_n]);CHKERRQ(ierr); 577 if(svalues){ierr = PetscMemcpy(sa->values,svalues,sizeof(PetscInt)*sxadj[irow_n]);CHKERRQ(ierr);} 578 ierr = PetscFree(sxadj);CHKERRQ(ierr); 579 ierr = PetscFree(sadjncy);CHKERRQ(ierr); 580 if(svalues) {ierr = PetscFree(svalues);CHKERRQ(ierr);} 581 } 582 } 583 ierr = PetscFree(indices);CHKERRQ(ierr); 584 PetscFunctionReturn(0); 585 } 586 587 588 #undef __FUNCT__ 589 #define __FUNCT__ "MatGetSubMatrix_MPIAdj_data" 590 /* 591 * The interface should be easy to use for both MatGetSubMatrix (parallel sub-matrix) and MatGetSubMatrices (sequential sub-matrices) 592 * */ 593 static PetscErrorCode MatGetSubMatrix_MPIAdj_data(Mat adj,IS irows, IS icols, PetscInt **sadj_xadj,PetscInt **sadj_adjncy,PetscInt **sadj_values) 594 { 595 PetscInt nlrows_is,icols_n,i,j,nroots,nleaves,owner,rlocalindex,*ncols_send,*ncols_recv; 596 PetscInt nlrows_mat,*adjncy_recv,Ncols_recv,Ncols_send,*xadj_recv,*values_recv; 597 PetscInt *ncols_recv_offsets,loc,rnclos,*sadjncy,*sxadj,*svalues,isvalue; 598 const PetscInt *irows_indices,*icols_indices,*xadj, *adjncy; 599 Mat_MPIAdj *a = (Mat_MPIAdj*)adj->data; 600 PetscLayout rmap; 601 MPI_Comm comm; 602 PetscSF sf; 603 PetscSFNode *iremote; 604 PetscBool done; 605 PetscErrorCode ierr; 606 607 PetscFunctionBegin; 608 /* communicator */ 609 ierr = PetscObjectGetComm((PetscObject)adj,&comm);CHKERRQ(ierr); 610 /* Layouts */ 611 ierr = MatGetLayouts(adj,&rmap,PETSC_NULL);CHKERRQ(ierr); 612 /* get rows information */ 613 ierr = ISGetLocalSize(irows,&nlrows_is);CHKERRQ(ierr); 614 ierr = ISGetIndices(irows,&irows_indices);CHKERRQ(ierr); 615 ierr = PetscCalloc1(nlrows_is,&iremote);CHKERRQ(ierr); 616 /* construct sf graph*/ 617 nleaves = nlrows_is; 618 for(i=0; i<nlrows_is; i++){ 619 ierr = PetscLayoutFindOwnerIndex(rmap,irows_indices[i],&owner,&rlocalindex);CHKERRQ(ierr); 620 iremote[i].rank = owner; 621 iremote[i].index = rlocalindex; 622 } 623 ierr = MatGetRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlrows_mat,&xadj,&adjncy,&done);CHKERRQ(ierr); 624 ierr = PetscCalloc4(nlrows_mat,&ncols_send,nlrows_is,&xadj_recv,nlrows_is+1,&ncols_recv_offsets,nlrows_is,&ncols_recv);CHKERRQ(ierr); 625 nroots = nlrows_mat; 626 for(i=0; i<nlrows_mat; i++){ 627 ncols_send[i] = xadj[i+1]-xadj[i]; 628 } 629 ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 630 ierr = PetscSFSetGraph(sf,nroots,nleaves,PETSC_NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 631 ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr); 632 ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 633 ierr = PetscSFBcastBegin(sf,MPIU_INT,ncols_send,ncols_recv);CHKERRQ(ierr); 634 ierr = PetscSFBcastEnd(sf,MPIU_INT,ncols_send,ncols_recv);CHKERRQ(ierr); 635 ierr = PetscSFBcastBegin(sf,MPIU_INT,xadj,xadj_recv);CHKERRQ(ierr); 636 ierr = PetscSFBcastEnd(sf,MPIU_INT,xadj,xadj_recv);CHKERRQ(ierr); 637 ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 638 Ncols_recv =0; 639 for(i=0; i<nlrows_is; i++){ 640 Ncols_recv += ncols_recv[i]; 641 ncols_recv_offsets[i+1] = ncols_recv[i]+ncols_recv_offsets[i]; 642 } 643 Ncols_send = 0; 644 for(i=0; i<nlrows_mat; i++){ 645 Ncols_send += ncols_send[i]; 646 } 647 ierr = PetscCalloc1(Ncols_recv,&iremote);CHKERRQ(ierr); 648 ierr = PetscCalloc1(Ncols_recv,&adjncy_recv);CHKERRQ(ierr); 649 nleaves = Ncols_recv; 650 Ncols_recv = 0; 651 for(i=0; i<nlrows_is; i++){ 652 ierr = PetscLayoutFindOwner(rmap,irows_indices[i],&owner);CHKERRQ(ierr); 653 for(j=0; j<ncols_recv[i]; j++){ 654 iremote[Ncols_recv].rank = owner; 655 iremote[Ncols_recv++].index = xadj_recv[i]+j; 656 } 657 } 658 ierr = ISRestoreIndices(irows,&irows_indices);CHKERRQ(ierr); 659 /*if we need to deal with edge weights ???*/ 660 if(a->values){isvalue=1;}else{isvalue=0;} 661 /*involve a global communication */ 662 /*ierr = MPI_Allreduce(&isvalue,&isvalue,1,MPIU_INT,MPI_SUM,comm);CHKERRQ(ierr);*/ 663 if(isvalue){ierr = PetscCalloc1(Ncols_recv,&values_recv);CHKERRQ(ierr);} 664 nroots = Ncols_send; 665 ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 666 ierr = PetscSFSetGraph(sf,nroots,nleaves,PETSC_NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 667 ierr = PetscSFSetType(sf,PETSCSFBASIC);CHKERRQ(ierr); 668 ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 669 ierr = PetscSFBcastBegin(sf,MPIU_INT,adjncy,adjncy_recv);CHKERRQ(ierr); 670 ierr = PetscSFBcastEnd(sf,MPIU_INT,adjncy,adjncy_recv);CHKERRQ(ierr); 671 if(isvalue){ 672 ierr = PetscSFBcastBegin(sf,MPIU_INT,a->values,values_recv);CHKERRQ(ierr); 673 ierr = PetscSFBcastEnd(sf,MPIU_INT,a->values,values_recv);CHKERRQ(ierr); 674 } 675 ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 676 ierr = MatRestoreRowIJ(adj,0,PETSC_FALSE,PETSC_FALSE,&nlrows_mat,&xadj,&adjncy,&done);CHKERRQ(ierr); 677 ierr = ISGetLocalSize(icols,&icols_n);CHKERRQ(ierr); 678 ierr = ISGetIndices(icols,&icols_indices);CHKERRQ(ierr); 679 rnclos = 0; 680 for(i=0; i<nlrows_is; i++){ 681 for(j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++){ 682 ierr = PetscFindInt(adjncy_recv[j], icols_n, icols_indices, &loc);CHKERRQ(ierr); 683 if(loc<0){ 684 adjncy_recv[j] = -1; 685 if(isvalue) values_recv[j] = -1; 686 ncols_recv[i]--; 687 }else{ 688 rnclos++; 689 } 690 } 691 } 692 ierr = ISRestoreIndices(icols,&icols_indices);CHKERRQ(ierr); 693 ierr = PetscCalloc1(rnclos,&sadjncy);CHKERRQ(ierr); 694 if(isvalue) {ierr = PetscCalloc1(rnclos,&svalues);CHKERRQ(ierr);} 695 ierr = PetscCalloc1(nlrows_is+1,&sxadj);CHKERRQ(ierr); 696 rnclos = 0; 697 for(i=0; i<nlrows_is; i++){ 698 for(j=ncols_recv_offsets[i]; j<ncols_recv_offsets[i+1]; j++){ 699 if(adjncy_recv[j]<0) continue; 700 sadjncy[rnclos] = adjncy_recv[j]; 701 if(isvalue) svalues[rnclos] = values_recv[j]; 702 rnclos++; 703 } 704 } 705 for(i=0; i<nlrows_is; i++){ 706 sxadj[i+1] = sxadj[i]+ncols_recv[i]; 707 } 708 if(sadj_xadj) { *sadj_xadj = sxadj;}else { ierr = PetscFree(sxadj);CHKERRQ(ierr);} 709 if(sadj_adjncy){ *sadj_adjncy = sadjncy;}else{ ierr = PetscFree(sadjncy);CHKERRQ(ierr);} 710 if(sadj_values){ 711 if(isvalue) *sadj_values = svalues; else *sadj_values=0; 712 }else{ 713 if(isvalue) {ierr = PetscFree(svalues);CHKERRQ(ierr);} 714 } 715 ierr = PetscFree4(ncols_send,xadj_recv,ncols_recv_offsets,ncols_recv);CHKERRQ(ierr); 716 ierr = PetscFree(adjncy_recv);CHKERRQ(ierr); 717 if(isvalue) {ierr = PetscFree(values_recv);CHKERRQ(ierr);} 718 PetscFunctionReturn(0); 719 } 720 721 722 #undef __FUNCT__ 723 #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat_MPIAdj" 724 static PetscErrorCode MatMPIAdjCreateNonemptySubcommMat_MPIAdj(Mat A,Mat *B) 725 { 726 Mat_MPIAdj *a = (Mat_MPIAdj*)A->data; 727 PetscErrorCode ierr; 728 const PetscInt *ranges; 729 MPI_Comm acomm,bcomm; 730 MPI_Group agroup,bgroup; 731 PetscMPIInt i,rank,size,nranks,*ranks; 732 733 PetscFunctionBegin; 734 *B = NULL; 735 ierr = PetscObjectGetComm((PetscObject)A,&acomm);CHKERRQ(ierr); 736 ierr = MPI_Comm_size(acomm,&size);CHKERRQ(ierr); 737 ierr = MPI_Comm_size(acomm,&rank);CHKERRQ(ierr); 738 ierr = MatGetOwnershipRanges(A,&ranges);CHKERRQ(ierr); 739 for (i=0,nranks=0; i<size; i++) { 740 if (ranges[i+1] - ranges[i] > 0) nranks++; 741 } 742 if (nranks == size) { /* All ranks have a positive number of rows, so we do not need to create a subcomm; */ 743 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 744 *B = A; 745 PetscFunctionReturn(0); 746 } 747 748 ierr = PetscMalloc1(nranks,&ranks);CHKERRQ(ierr); 749 for (i=0,nranks=0; i<size; i++) { 750 if (ranges[i+1] - ranges[i] > 0) ranks[nranks++] = i; 751 } 752 ierr = MPI_Comm_group(acomm,&agroup);CHKERRQ(ierr); 753 ierr = MPI_Group_incl(agroup,nranks,ranks,&bgroup);CHKERRQ(ierr); 754 ierr = PetscFree(ranks);CHKERRQ(ierr); 755 ierr = MPI_Comm_create(acomm,bgroup,&bcomm);CHKERRQ(ierr); 756 ierr = MPI_Group_free(&agroup);CHKERRQ(ierr); 757 ierr = MPI_Group_free(&bgroup);CHKERRQ(ierr); 758 if (bcomm != MPI_COMM_NULL) { 759 PetscInt m,N; 760 Mat_MPIAdj *b; 761 ierr = MatGetLocalSize(A,&m,NULL);CHKERRQ(ierr); 762 ierr = MatGetSize(A,NULL,&N);CHKERRQ(ierr); 763 ierr = MatCreateMPIAdj(bcomm,m,N,a->i,a->j,a->values,B);CHKERRQ(ierr); 764 b = (Mat_MPIAdj*)(*B)->data; 765 b->freeaij = PETSC_FALSE; 766 ierr = MPI_Comm_free(&bcomm);CHKERRQ(ierr); 767 } 768 PetscFunctionReturn(0); 769 } 770 771 #undef __FUNCT__ 772 #define __FUNCT__ "MatMPIAdjCreateNonemptySubcommMat" 773 /*@ 774 MatMPIAdjCreateNonemptySubcommMat - create the same MPIAdj matrix on a subcommunicator containing only processes owning a positive number of rows 775 776 Collective 777 778 Input Arguments: 779 . A - original MPIAdj matrix 780 781 Output Arguments: 782 . B - matrix on subcommunicator, NULL on ranks that owned zero rows of A 783 784 Level: developer 785 786 Note: 787 This function is mostly useful for internal use by mesh partitioning packages that require that every process owns at least one row. 788 789 The matrix B should be destroyed with MatDestroy(). The arrays are not copied, so B should be destroyed before A is destroyed. 790 791 .seealso: MatCreateMPIAdj() 792 @*/ 793 PetscErrorCode MatMPIAdjCreateNonemptySubcommMat(Mat A,Mat *B) 794 { 795 PetscErrorCode ierr; 796 797 PetscFunctionBegin; 798 PetscValidHeaderSpecific(A,MAT_CLASSID,1); 799 ierr = PetscUseMethod(A,"MatMPIAdjCreateNonemptySubcommMat_C",(Mat,Mat*),(A,B));CHKERRQ(ierr); 800 PetscFunctionReturn(0); 801 } 802 803 /*MC 804 MATMPIADJ - MATMPIADJ = "mpiadj" - A matrix type to be used for distributed adjacency matrices, 805 intended for use constructing orderings and partitionings. 806 807 Level: beginner 808 809 .seealso: MatCreateMPIAdj 810 M*/ 811 812 #undef __FUNCT__ 813 #define __FUNCT__ "MatCreate_MPIAdj" 814 PETSC_EXTERN PetscErrorCode MatCreate_MPIAdj(Mat B) 815 { 816 Mat_MPIAdj *b; 817 PetscErrorCode ierr; 818 819 PetscFunctionBegin; 820 ierr = PetscNewLog(B,&b);CHKERRQ(ierr); 821 B->data = (void*)b; 822 ierr = PetscMemcpy(B->ops,&MatOps_Values,sizeof(struct _MatOps));CHKERRQ(ierr); 823 B->assembled = PETSC_FALSE; 824 825 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjSetPreallocation_C",MatMPIAdjSetPreallocation_MPIAdj);CHKERRQ(ierr); 826 ierr = PetscObjectComposeFunction((PetscObject)B,"MatMPIAdjCreateNonemptySubcommMat_C",MatMPIAdjCreateNonemptySubcommMat_MPIAdj);CHKERRQ(ierr); 827 ierr = PetscObjectChangeTypeName((PetscObject)B,MATMPIADJ);CHKERRQ(ierr); 828 PetscFunctionReturn(0); 829 } 830 831 #undef __FUNCT__ 832 #define __FUNCT__ "MatMPIAdjSetPreallocation" 833 /*@C 834 MatMPIAdjSetPreallocation - Sets the array used for storing the matrix elements 835 836 Logically Collective on MPI_Comm 837 838 Input Parameters: 839 + A - the matrix 840 . i - the indices into j for the start of each row 841 . j - the column indices for each row (sorted for each row). 842 The indices in i and j start with zero (NOT with one). 843 - values - [optional] edge weights 844 845 Level: intermediate 846 847 .seealso: MatCreate(), MatCreateMPIAdj(), MatSetValues() 848 @*/ 849 PetscErrorCode MatMPIAdjSetPreallocation(Mat B,PetscInt *i,PetscInt *j,PetscInt *values) 850 { 851 PetscErrorCode ierr; 852 853 PetscFunctionBegin; 854 ierr = PetscTryMethod(B,"MatMPIAdjSetPreallocation_C",(Mat,PetscInt*,PetscInt*,PetscInt*),(B,i,j,values));CHKERRQ(ierr); 855 PetscFunctionReturn(0); 856 } 857 858 #undef __FUNCT__ 859 #define __FUNCT__ "MatCreateMPIAdj" 860 /*@C 861 MatCreateMPIAdj - Creates a sparse matrix representing an adjacency list. 862 The matrix does not have numerical values associated with it, but is 863 intended for ordering (to reduce bandwidth etc) and partitioning. 864 865 Collective on MPI_Comm 866 867 Input Parameters: 868 + comm - MPI communicator 869 . m - number of local rows 870 . N - number of global columns 871 . i - the indices into j for the start of each row 872 . j - the column indices for each row (sorted for each row). 873 The indices in i and j start with zero (NOT with one). 874 - values -[optional] edge weights 875 876 Output Parameter: 877 . A - the matrix 878 879 Level: intermediate 880 881 Notes: This matrix object does not support most matrix operations, include 882 MatSetValues(). 883 You must NOT free the ii, values and jj arrays yourself. PETSc will free them 884 when the matrix is destroyed; you must allocate them with PetscMalloc(). If you 885 call from Fortran you need not create the arrays with PetscMalloc(). 886 Should not include the matrix diagonals. 887 888 If you already have a matrix, you can create its adjacency matrix by a call 889 to MatConvert, specifying a type of MATMPIADJ. 890 891 Possible values for MatSetOption() - MAT_STRUCTURALLY_SYMMETRIC 892 893 .seealso: MatCreate(), MatConvert(), MatGetOrdering() 894 @*/ 895 PetscErrorCode MatCreateMPIAdj(MPI_Comm comm,PetscInt m,PetscInt N,PetscInt *i,PetscInt *j,PetscInt *values,Mat *A) 896 { 897 PetscErrorCode ierr; 898 899 PetscFunctionBegin; 900 ierr = MatCreate(comm,A);CHKERRQ(ierr); 901 ierr = MatSetSizes(*A,m,PETSC_DETERMINE,PETSC_DETERMINE,N);CHKERRQ(ierr); 902 ierr = MatSetType(*A,MATMPIADJ);CHKERRQ(ierr); 903 ierr = MatMPIAdjSetPreallocation(*A,i,j,values);CHKERRQ(ierr); 904 PetscFunctionReturn(0); 905 } 906 907 908 909 910 911 912 913