1 2 /* 3 Defines matrix-matrix product routines for pairs of SeqAIJ matrices 4 C = A * B 5 */ 6 7 #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 8 #include <../src/mat/utils/freespace.h> 9 #include <../src/mat/utils/petscheap.h> 10 #include <petscbt.h> 11 #include <../src/mat/impls/dense/seq/dense.h> 12 #include <petsctime.h> 13 14 static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat*); 15 16 #undef __FUNCT__ 17 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ" 18 PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 19 { 20 PetscErrorCode ierr; 21 const char *algTypes[6] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed"}; 22 PetscInt alg=0; /* set default algorithm */ 23 24 PetscFunctionBegin; 25 if (scall == MAT_INITIAL_MATRIX) { 26 ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 27 ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,6,algTypes[0],&alg,NULL);CHKERRQ(ierr); 28 ierr = PetscOptionsEnd();CHKERRQ(ierr); 29 ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 30 switch (alg) { 31 case 1: 32 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr); 33 break; 34 case 2: 35 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr); 36 break; 37 case 3: 38 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr); 39 break; 40 case 4: 41 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);CHKERRQ(ierr); 42 break; 43 case 5: 44 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(A,B,fill,C);CHKERRQ(ierr); 45 break; 46 default: 47 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 48 break; 49 } 50 ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 51 } 52 53 ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 54 ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 55 ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 56 PetscFunctionReturn(0); 57 } 58 59 #undef __FUNCT__ 60 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed" 61 static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat *C) 62 { 63 PetscErrorCode ierr; 64 Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 65 PetscInt *ai=a->i,*bi=b->i,*ci,*cj; 66 PetscInt am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 67 PetscReal afill; 68 PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0; 69 PetscBT lnkbt; 70 PetscFreeSpaceList free_space=NULL,current_space=NULL; 71 72 PetscFunctionBegin; 73 /* Get ci and cj */ 74 /*---------------*/ 75 /* Allocate ci array, arrays for fill computation and */ 76 /* free space for accumulating nonzero column info */ 77 ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 78 ci[0] = 0; 79 80 /* create and initialize a linked list */ 81 nlnk_max = a->rmax*b->rmax; 82 if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; 83 ierr = PetscLLCondensedCreate(nlnk_max,bn,&lnk,&lnkbt);CHKERRQ(ierr); 84 85 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 86 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 87 88 current_space = free_space; 89 90 /* Determine ci and cj */ 91 for (i=0; i<am; i++) { 92 anzi = ai[i+1] - ai[i]; 93 aj = a->j + ai[i]; 94 for (j=0; j<anzi; j++) { 95 brow = aj[j]; 96 bnzj = bi[brow+1] - bi[brow]; 97 bj = b->j + bi[brow]; 98 /* add non-zero cols of B into the sorted linked list lnk */ 99 ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr); 100 } 101 cnzi = lnk[0]; 102 103 /* If free space is not available, make more free space */ 104 /* Double the amount of total space in the list */ 105 if (current_space->local_remaining<cnzi) { 106 ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 107 ndouble++; 108 } 109 110 /* Copy data into free space, then initialize lnk */ 111 ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 112 113 current_space->array += cnzi; 114 current_space->local_used += cnzi; 115 current_space->local_remaining -= cnzi; 116 117 ci[i+1] = ci[i] + cnzi; 118 } 119 120 /* Column indices are in the list of free space */ 121 /* Allocate space for cj, initialize cj, and */ 122 /* destroy list of free space and other temporary array(s) */ 123 ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 124 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 125 ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr); 126 127 /* put together the new symbolic matrix */ 128 ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr); 129 130 (*C)->rmap->bs = A->rmap->bs; 131 (*C)->cmap->bs = B->cmap->bs; 132 133 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 134 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 135 c = (Mat_SeqAIJ*)((*C)->data); 136 c->free_a = PETSC_FALSE; 137 c->free_ij = PETSC_TRUE; 138 c->nonew = 0; 139 (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */ 140 141 /* set MatInfo */ 142 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 143 if (afill < 1.0) afill = 1.0; 144 c->maxnz = ci[am]; 145 c->nz = ci[am]; 146 (*C)->info.mallocs = ndouble; 147 (*C)->info.fill_ratio_given = fill; 148 (*C)->info.fill_ratio_needed = afill; 149 150 #if defined(PETSC_USE_INFO) 151 if (ci[am]) { 152 ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 153 ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 154 } else { 155 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 156 } 157 #endif 158 PetscFunctionReturn(0); 159 } 160 161 #undef __FUNCT__ 162 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ" 163 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 164 { 165 PetscErrorCode ierr; 166 PetscLogDouble flops=0.0; 167 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 168 Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 169 Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 170 PetscInt *ai =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 171 PetscInt am =A->rmap->n,cm=C->rmap->n; 172 PetscInt i,j,k,anzi,bnzi,cnzi,brow; 173 PetscScalar *aa=a->a,*ba=b->a,*baj,*ca,valtmp; 174 PetscScalar *ab_dense; 175 176 PetscFunctionBegin; 177 if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */ 178 ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 179 c->a = ca; 180 c->free_a = PETSC_TRUE; 181 182 ierr = PetscMalloc(B->cmap->N*sizeof(PetscScalar),&ab_dense);CHKERRQ(ierr); 183 ierr = PetscMemzero(ab_dense,B->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr); 184 185 c->matmult_abdense = ab_dense; 186 } else { 187 ca = c->a; 188 ab_dense = c->matmult_abdense; 189 } 190 191 /* clean old values in C */ 192 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 193 /* Traverse A row-wise. */ 194 /* Build the ith row in C by summing over nonzero columns in A, */ 195 /* the rows of B corresponding to nonzeros of A. */ 196 for (i=0; i<am; i++) { 197 anzi = ai[i+1] - ai[i]; 198 for (j=0; j<anzi; j++) { 199 brow = aj[j]; 200 bnzi = bi[brow+1] - bi[brow]; 201 bjj = bj + bi[brow]; 202 baj = ba + bi[brow]; 203 /* perform dense axpy */ 204 valtmp = aa[j]; 205 for (k=0; k<bnzi; k++) { 206 ab_dense[bjj[k]] += valtmp*baj[k]; 207 } 208 flops += 2*bnzi; 209 } 210 aj += anzi; aa += anzi; 211 212 cnzi = ci[i+1] - ci[i]; 213 for (k=0; k<cnzi; k++) { 214 ca[k] += ab_dense[cj[k]]; 215 ab_dense[cj[k]] = 0.0; /* zero ab_dense */ 216 } 217 flops += cnzi; 218 cj += cnzi; ca += cnzi; 219 } 220 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 221 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 222 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 223 PetscFunctionReturn(0); 224 } 225 226 #undef __FUNCT__ 227 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable" 228 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C) 229 { 230 PetscErrorCode ierr; 231 PetscLogDouble flops=0.0; 232 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 233 Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 234 Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 235 PetscInt *ai = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 236 PetscInt am = A->rmap->N,cm=C->rmap->N; 237 PetscInt i,j,k,anzi,bnzi,cnzi,brow; 238 PetscScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp; 239 PetscInt nextb; 240 241 PetscFunctionBegin; 242 /* clean old values in C */ 243 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 244 /* Traverse A row-wise. */ 245 /* Build the ith row in C by summing over nonzero columns in A, */ 246 /* the rows of B corresponding to nonzeros of A. */ 247 for (i=0; i<am; i++) { 248 anzi = ai[i+1] - ai[i]; 249 cnzi = ci[i+1] - ci[i]; 250 for (j=0; j<anzi; j++) { 251 brow = aj[j]; 252 bnzi = bi[brow+1] - bi[brow]; 253 bjj = bj + bi[brow]; 254 baj = ba + bi[brow]; 255 /* perform sparse axpy */ 256 valtmp = aa[j]; 257 nextb = 0; 258 for (k=0; nextb<bnzi; k++) { 259 if (cj[k] == bjj[nextb]) { /* ccol == bcol */ 260 ca[k] += valtmp*baj[nextb++]; 261 } 262 } 263 flops += 2*bnzi; 264 } 265 aj += anzi; aa += anzi; 266 cj += cnzi; ca += cnzi; 267 } 268 269 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 270 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 271 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 272 PetscFunctionReturn(0); 273 } 274 275 #undef __FUNCT__ 276 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast" 277 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat *C) 278 { 279 PetscErrorCode ierr; 280 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 281 PetscInt *ai = a->i,*bi=b->i,*ci,*cj; 282 PetscInt am = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 283 MatScalar *ca; 284 PetscReal afill; 285 PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0; 286 PetscFreeSpaceList free_space=NULL,current_space=NULL; 287 288 PetscFunctionBegin; 289 /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */ 290 /*-----------------------------------------------------------------------------------------*/ 291 /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 292 ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 293 ci[0] = 0; 294 295 /* create and initialize a linked list */ 296 nlnk_max = a->rmax*b->rmax; 297 if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */ 298 ierr = PetscLLCondensedCreate_fast(nlnk_max,&lnk);CHKERRQ(ierr); 299 300 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 301 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 302 current_space = free_space; 303 304 /* Determine ci and cj */ 305 for (i=0; i<am; i++) { 306 anzi = ai[i+1] - ai[i]; 307 aj = a->j + ai[i]; 308 for (j=0; j<anzi; j++) { 309 brow = aj[j]; 310 bnzj = bi[brow+1] - bi[brow]; 311 bj = b->j + bi[brow]; 312 /* add non-zero cols of B into the sorted linked list lnk */ 313 ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr); 314 } 315 cnzi = lnk[1]; 316 317 /* If free space is not available, make more free space */ 318 /* Double the amount of total space in the list */ 319 if (current_space->local_remaining<cnzi) { 320 ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 321 ndouble++; 322 } 323 324 /* Copy data into free space, then initialize lnk */ 325 ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr); 326 327 current_space->array += cnzi; 328 current_space->local_used += cnzi; 329 current_space->local_remaining -= cnzi; 330 331 ci[i+1] = ci[i] + cnzi; 332 } 333 334 /* Column indices are in the list of free space */ 335 /* Allocate space for cj, initialize cj, and */ 336 /* destroy list of free space and other temporary array(s) */ 337 ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 338 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 339 ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr); 340 341 /* Allocate space for ca */ 342 ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 343 ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 344 345 /* put together the new symbolic matrix */ 346 ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr); 347 348 (*C)->rmap->bs = A->rmap->bs; 349 (*C)->cmap->bs = B->cmap->bs; 350 351 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 352 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 353 c = (Mat_SeqAIJ*)((*C)->data); 354 c->free_a = PETSC_TRUE; 355 c->free_ij = PETSC_TRUE; 356 c->nonew = 0; 357 358 (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */ 359 360 /* set MatInfo */ 361 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 362 if (afill < 1.0) afill = 1.0; 363 c->maxnz = ci[am]; 364 c->nz = ci[am]; 365 (*C)->info.mallocs = ndouble; 366 (*C)->info.fill_ratio_given = fill; 367 (*C)->info.fill_ratio_needed = afill; 368 369 #if defined(PETSC_USE_INFO) 370 if (ci[am]) { 371 ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 372 ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 373 } else { 374 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 375 } 376 #endif 377 PetscFunctionReturn(0); 378 } 379 380 381 #undef __FUNCT__ 382 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable" 383 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C) 384 { 385 PetscErrorCode ierr; 386 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 387 PetscInt *ai = a->i,*bi=b->i,*ci,*cj; 388 PetscInt am = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 389 MatScalar *ca; 390 PetscReal afill; 391 PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0; 392 PetscFreeSpaceList free_space=NULL,current_space=NULL; 393 394 PetscFunctionBegin; 395 /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */ 396 /*---------------------------------------------------------------------------------------------*/ 397 /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 398 ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 399 ci[0] = 0; 400 401 /* create and initialize a linked list */ 402 nlnk_max = a->rmax*b->rmax; 403 if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; /* in case rmax is not defined for A or B */ 404 ierr = PetscLLCondensedCreate_Scalable(nlnk_max,&lnk);CHKERRQ(ierr); 405 406 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 407 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 408 current_space = free_space; 409 410 /* Determine ci and cj */ 411 for (i=0; i<am; i++) { 412 anzi = ai[i+1] - ai[i]; 413 aj = a->j + ai[i]; 414 for (j=0; j<anzi; j++) { 415 brow = aj[j]; 416 bnzj = bi[brow+1] - bi[brow]; 417 bj = b->j + bi[brow]; 418 /* add non-zero cols of B into the sorted linked list lnk */ 419 ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr); 420 } 421 cnzi = lnk[0]; 422 423 /* If free space is not available, make more free space */ 424 /* Double the amount of total space in the list */ 425 if (current_space->local_remaining<cnzi) { 426 ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 427 ndouble++; 428 } 429 430 /* Copy data into free space, then initialize lnk */ 431 ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr); 432 433 current_space->array += cnzi; 434 current_space->local_used += cnzi; 435 current_space->local_remaining -= cnzi; 436 437 ci[i+1] = ci[i] + cnzi; 438 } 439 440 /* Column indices are in the list of free space */ 441 /* Allocate space for cj, initialize cj, and */ 442 /* destroy list of free space and other temporary array(s) */ 443 ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 444 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 445 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 446 447 /* Allocate space for ca */ 448 /*-----------------------*/ 449 ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 450 ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 451 452 /* put together the new symbolic matrix */ 453 ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,C);CHKERRQ(ierr); 454 455 (*C)->rmap->bs = A->rmap->bs; 456 (*C)->cmap->bs = B->cmap->bs; 457 458 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 459 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 460 c = (Mat_SeqAIJ*)((*C)->data); 461 c->free_a = PETSC_TRUE; 462 c->free_ij = PETSC_TRUE; 463 c->nonew = 0; 464 465 (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */ 466 467 /* set MatInfo */ 468 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 469 if (afill < 1.0) afill = 1.0; 470 c->maxnz = ci[am]; 471 c->nz = ci[am]; 472 (*C)->info.mallocs = ndouble; 473 (*C)->info.fill_ratio_given = fill; 474 (*C)->info.fill_ratio_needed = afill; 475 476 #if defined(PETSC_USE_INFO) 477 if (ci[am]) { 478 ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 479 ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 480 } else { 481 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 482 } 483 #endif 484 PetscFunctionReturn(0); 485 } 486 487 #undef __FUNCT__ 488 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap" 489 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat *C) 490 { 491 PetscErrorCode ierr; 492 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 493 const PetscInt *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j; 494 PetscInt *ci,*cj,*bb; 495 PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 496 PetscReal afill; 497 PetscInt i,j,col,ndouble = 0; 498 PetscFreeSpaceList free_space=NULL,current_space=NULL; 499 PetscHeap h; 500 501 PetscFunctionBegin; 502 /* Get ci and cj - by merging sorted rows using a heap */ 503 /*---------------------------------------------------------------------------------------------*/ 504 /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 505 ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 506 ci[0] = 0; 507 508 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 509 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 510 current_space = free_space; 511 512 ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr); 513 ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr); 514 515 /* Determine ci and cj */ 516 for (i=0; i<am; i++) { 517 const PetscInt anzi = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */ 518 const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 519 ci[i+1] = ci[i]; 520 /* Populate the min heap */ 521 for (j=0; j<anzi; j++) { 522 bb[j] = bi[acol[j]]; /* bb points at the start of the row */ 523 if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */ 524 ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr); 525 } 526 } 527 /* Pick off the min element, adding it to free space */ 528 ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 529 while (j >= 0) { 530 if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */ 531 ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),¤t_space);CHKERRQ(ierr); 532 ndouble++; 533 } 534 *(current_space->array++) = col; 535 current_space->local_used++; 536 current_space->local_remaining--; 537 ci[i+1]++; 538 539 /* stash if anything else remains in this row of B */ 540 if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);} 541 while (1) { /* pop and stash any other rows of B that also had an entry in this column */ 542 PetscInt j2,col2; 543 ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr); 544 if (col2 != col) break; 545 ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr); 546 if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);} 547 } 548 /* Put any stashed elements back into the min heap */ 549 ierr = PetscHeapUnstash(h);CHKERRQ(ierr); 550 ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 551 } 552 } 553 ierr = PetscFree(bb);CHKERRQ(ierr); 554 ierr = PetscHeapDestroy(&h);CHKERRQ(ierr); 555 556 /* Column indices are in the list of free space */ 557 /* Allocate space for cj, initialize cj, and */ 558 /* destroy list of free space and other temporary array(s) */ 559 ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr); 560 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 561 562 /* put together the new symbolic matrix */ 563 ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr); 564 565 (*C)->rmap->bs = A->rmap->bs; 566 (*C)->cmap->bs = B->cmap->bs; 567 568 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 569 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 570 c = (Mat_SeqAIJ*)((*C)->data); 571 c->free_a = PETSC_TRUE; 572 c->free_ij = PETSC_TRUE; 573 c->nonew = 0; 574 575 (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; 576 577 /* set MatInfo */ 578 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 579 if (afill < 1.0) afill = 1.0; 580 c->maxnz = ci[am]; 581 c->nz = ci[am]; 582 (*C)->info.mallocs = ndouble; 583 (*C)->info.fill_ratio_given = fill; 584 (*C)->info.fill_ratio_needed = afill; 585 586 #if defined(PETSC_USE_INFO) 587 if (ci[am]) { 588 ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 589 ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 590 } else { 591 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 592 } 593 #endif 594 PetscFunctionReturn(0); 595 } 596 597 #undef __FUNCT__ 598 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap" 599 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C) 600 { 601 PetscErrorCode ierr; 602 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 603 const PetscInt *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j; 604 PetscInt *ci,*cj,*bb; 605 PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 606 PetscReal afill; 607 PetscInt i,j,col,ndouble = 0; 608 PetscFreeSpaceList free_space=NULL,current_space=NULL; 609 PetscHeap h; 610 PetscBT bt; 611 612 PetscFunctionBegin; 613 /* Get ci and cj - using a heap for the sorted rows, but use BT so that each index is only added once */ 614 /*---------------------------------------------------------------------------------------------*/ 615 /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 616 ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 617 ci[0] = 0; 618 619 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 620 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 621 622 current_space = free_space; 623 624 ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr); 625 ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr); 626 ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr); 627 628 /* Determine ci and cj */ 629 for (i=0; i<am; i++) { 630 const PetscInt anzi = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */ 631 const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 632 const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */ 633 ci[i+1] = ci[i]; 634 /* Populate the min heap */ 635 for (j=0; j<anzi; j++) { 636 PetscInt brow = acol[j]; 637 for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) { 638 PetscInt bcol = bj[bb[j]]; 639 if (!PetscBTLookupSet(bt,bcol)) { /* new entry */ 640 ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr); 641 bb[j]++; 642 break; 643 } 644 } 645 } 646 /* Pick off the min element, adding it to free space */ 647 ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 648 while (j >= 0) { 649 if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */ 650 fptr = NULL; /* need PetscBTMemzero */ 651 ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),¤t_space);CHKERRQ(ierr); 652 ndouble++; 653 } 654 *(current_space->array++) = col; 655 current_space->local_used++; 656 current_space->local_remaining--; 657 ci[i+1]++; 658 659 /* stash if anything else remains in this row of B */ 660 for (; bb[j] < bi[acol[j]+1]; bb[j]++) { 661 PetscInt bcol = bj[bb[j]]; 662 if (!PetscBTLookupSet(bt,bcol)) { /* new entry */ 663 ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr); 664 bb[j]++; 665 break; 666 } 667 } 668 ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 669 } 670 if (fptr) { /* Clear the bits for this row */ 671 for (; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);} 672 } else { /* We reallocated so we don't remember (easily) how to clear only the bits we changed */ 673 ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr); 674 } 675 } 676 ierr = PetscFree(bb);CHKERRQ(ierr); 677 ierr = PetscHeapDestroy(&h);CHKERRQ(ierr); 678 ierr = PetscBTDestroy(&bt);CHKERRQ(ierr); 679 680 /* Column indices are in the list of free space */ 681 /* Allocate space for cj, initialize cj, and */ 682 /* destroy list of free space and other temporary array(s) */ 683 ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr); 684 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 685 686 /* put together the new symbolic matrix */ 687 ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr); 688 689 (*C)->rmap->bs = A->rmap->bs; 690 (*C)->cmap->bs = B->cmap->bs; 691 692 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 693 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 694 c = (Mat_SeqAIJ*)((*C)->data); 695 c->free_a = PETSC_TRUE; 696 c->free_ij = PETSC_TRUE; 697 c->nonew = 0; 698 699 (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; 700 701 /* set MatInfo */ 702 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 703 if (afill < 1.0) afill = 1.0; 704 c->maxnz = ci[am]; 705 c->nz = ci[am]; 706 (*C)->info.mallocs = ndouble; 707 (*C)->info.fill_ratio_given = fill; 708 (*C)->info.fill_ratio_needed = afill; 709 710 #if defined(PETSC_USE_INFO) 711 if (ci[am]) { 712 ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 713 ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 714 } else { 715 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 716 } 717 #endif 718 PetscFunctionReturn(0); 719 } 720 721 #undef __FUNCT__ 722 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ" 723 /* concatenate unique entries and then sort */ 724 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 725 { 726 PetscErrorCode ierr; 727 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 728 const PetscInt *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j; 729 PetscInt *ci,*cj; 730 PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 731 PetscReal afill; 732 PetscInt i,j,ndouble = 0; 733 PetscSegBuffer seg,segrow; 734 char *seen; 735 736 PetscFunctionBegin; 737 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 738 ci[0] = 0; 739 740 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 741 ierr = PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);CHKERRQ(ierr); 742 ierr = PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);CHKERRQ(ierr); 743 ierr = PetscMalloc(bn*sizeof(char),&seen);CHKERRQ(ierr); 744 ierr = PetscMemzero(seen,bn*sizeof(char));CHKERRQ(ierr); 745 746 /* Determine ci and cj */ 747 for (i=0; i<am; i++) { 748 const PetscInt anzi = ai[i+1] - ai[i]; /* number of nonzeros in this row of A, this is the number of rows of B that we merge */ 749 const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 750 PetscInt packlen = 0,*PETSC_RESTRICT crow; 751 /* Pack segrow */ 752 for (j=0; j<anzi; j++) { 753 PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k; 754 for (k=bjstart; k<bjend; k++) { 755 PetscInt bcol = bj[k]; 756 if (!seen[bcol]) { /* new entry */ 757 PetscInt *PETSC_RESTRICT slot; 758 ierr = PetscSegBufferGetInts(segrow,1,&slot);CHKERRQ(ierr); 759 *slot = bcol; 760 seen[bcol] = 1; 761 packlen++; 762 } 763 } 764 } 765 ierr = PetscSegBufferGetInts(seg,packlen,&crow);CHKERRQ(ierr); 766 ierr = PetscSegBufferExtractTo(segrow,crow);CHKERRQ(ierr); 767 ierr = PetscSortInt(packlen,crow);CHKERRQ(ierr); 768 ci[i+1] = ci[i] + packlen; 769 for (j=0; j<packlen; j++) seen[crow[j]] = 0; 770 } 771 ierr = PetscSegBufferDestroy(&segrow);CHKERRQ(ierr); 772 ierr = PetscFree(seen);CHKERRQ(ierr); 773 774 /* Column indices are in the segmented buffer */ 775 ierr = PetscSegBufferExtractAlloc(seg,&cj);CHKERRQ(ierr); 776 ierr = PetscSegBufferDestroy(&seg);CHKERRQ(ierr); 777 778 /* put together the new symbolic matrix */ 779 ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr); 780 781 (*C)->rmap->bs = A->rmap->bs; 782 (*C)->cmap->bs = B->cmap->bs; 783 784 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 785 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 786 c = (Mat_SeqAIJ*)((*C)->data); 787 c->free_a = PETSC_TRUE; 788 c->free_ij = PETSC_TRUE; 789 c->nonew = 0; 790 791 (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; 792 793 /* set MatInfo */ 794 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 795 if (afill < 1.0) afill = 1.0; 796 c->maxnz = ci[am]; 797 c->nz = ci[am]; 798 (*C)->info.mallocs = ndouble; 799 (*C)->info.fill_ratio_given = fill; 800 (*C)->info.fill_ratio_needed = afill; 801 802 #if defined(PETSC_USE_INFO) 803 if (ci[am]) { 804 ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 805 ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 806 } else { 807 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 808 } 809 #endif 810 PetscFunctionReturn(0); 811 } 812 813 /* This routine is not used. Should be removed! */ 814 #undef __FUNCT__ 815 #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ" 816 PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 817 { 818 PetscErrorCode ierr; 819 820 PetscFunctionBegin; 821 if (scall == MAT_INITIAL_MATRIX) { 822 ierr = PetscLogEventBegin(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr); 823 ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 824 ierr = PetscLogEventEnd(MAT_MatTransposeMultSymbolic,A,B,0,0);CHKERRQ(ierr); 825 } 826 ierr = PetscLogEventBegin(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr); 827 ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 828 ierr = PetscLogEventEnd(MAT_MatTransposeMultNumeric,A,B,0,0);CHKERRQ(ierr); 829 PetscFunctionReturn(0); 830 } 831 832 #undef __FUNCT__ 833 #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans" 834 PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A) 835 { 836 PetscErrorCode ierr; 837 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 838 Mat_MatMatTransMult *abt=a->abt; 839 840 PetscFunctionBegin; 841 ierr = (abt->destroy)(A);CHKERRQ(ierr); 842 ierr = MatTransposeColoringDestroy(&abt->matcoloring);CHKERRQ(ierr); 843 ierr = MatDestroy(&abt->Bt_den);CHKERRQ(ierr); 844 ierr = MatDestroy(&abt->ABt_den);CHKERRQ(ierr); 845 ierr = PetscFree(abt);CHKERRQ(ierr); 846 PetscFunctionReturn(0); 847 } 848 849 #undef __FUNCT__ 850 #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ" 851 PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 852 { 853 PetscErrorCode ierr; 854 Mat Bt; 855 PetscInt *bti,*btj; 856 Mat_MatMatTransMult *abt; 857 Mat_SeqAIJ *c; 858 859 PetscFunctionBegin; 860 /* create symbolic Bt */ 861 ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 862 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);CHKERRQ(ierr); 863 864 Bt->rmap->bs = A->cmap->bs; 865 Bt->cmap->bs = B->cmap->bs; 866 867 /* get symbolic C=A*Bt */ 868 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr); 869 870 /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */ 871 ierr = PetscNew(Mat_MatMatTransMult,&abt);CHKERRQ(ierr); 872 c = (Mat_SeqAIJ*)(*C)->data; 873 c->abt = abt; 874 875 abt->usecoloring = PETSC_FALSE; 876 abt->destroy = (*C)->ops->destroy; 877 (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans; 878 879 ierr = PetscOptionsGetBool(NULL,"-matmattransmult_color",&abt->usecoloring,NULL);CHKERRQ(ierr); 880 if (abt->usecoloring) { 881 /* Create MatTransposeColoring from symbolic C=A*B^T */ 882 MatTransposeColoring matcoloring; 883 ISColoring iscoloring; 884 Mat Bt_dense,C_dense; 885 Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data; 886 /* inode causes memory problem, don't know why */ 887 if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'"); 888 889 ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 890 ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 891 892 abt->matcoloring = matcoloring; 893 894 ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 895 896 /* Create Bt_dense and C_dense = A*Bt_dense */ 897 ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr); 898 ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 899 ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr); 900 ierr = MatSeqDenseSetPreallocation(Bt_dense,NULL);CHKERRQ(ierr); 901 902 Bt_dense->assembled = PETSC_TRUE; 903 abt->Bt_den = Bt_dense; 904 905 ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr); 906 ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 907 ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr); 908 ierr = MatSeqDenseSetPreallocation(C_dense,NULL);CHKERRQ(ierr); 909 910 Bt_dense->assembled = PETSC_TRUE; 911 abt->ABt_den = C_dense; 912 913 #if defined(PETSC_USE_INFO) 914 { 915 Mat_SeqAIJ *c = (Mat_SeqAIJ*)(*C)->data; 916 ierr = PetscInfo7(*C,"Use coloring of C=A*B^T; B^T: %D %D, Bt_dense: %D,%D; Cnz %D / (cm*ncolors %D) = %g\n",B->cmap->n,B->rmap->n,Bt_dense->rmap->n,Bt_dense->cmap->n,c->nz,A->rmap->n*matcoloring->ncolors,(PetscReal)(c->nz)/(A->rmap->n*matcoloring->ncolors));CHKERRQ(ierr); 917 } 918 #endif 919 } 920 /* clean up */ 921 ierr = MatDestroy(&Bt);CHKERRQ(ierr); 922 ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 923 PetscFunctionReturn(0); 924 } 925 926 #undef __FUNCT__ 927 #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ" 928 PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 929 { 930 PetscErrorCode ierr; 931 Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 932 PetscInt *ai =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow; 933 PetscInt cm =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol; 934 PetscLogDouble flops=0.0; 935 MatScalar *aa =a->a,*aval,*ba=b->a,*bval,*ca,*cval; 936 Mat_MatMatTransMult *abt = c->abt; 937 938 PetscFunctionBegin; 939 /* clear old values in C */ 940 if (!c->a) { 941 ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 942 c->a = ca; 943 c->free_a = PETSC_TRUE; 944 } else { 945 ca = c->a; 946 } 947 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 948 949 if (abt->usecoloring) { 950 MatTransposeColoring matcoloring = abt->matcoloring; 951 Mat Bt_dense,C_dense = abt->ABt_den; 952 953 /* Get Bt_dense by Apply MatTransposeColoring to B */ 954 Bt_dense = abt->Bt_den; 955 ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr); 956 957 /* C_dense = A*Bt_dense */ 958 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr); 959 960 /* Recover C from C_dense */ 961 ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr); 962 PetscFunctionReturn(0); 963 } 964 965 for (i=0; i<cm; i++) { 966 anzi = ai[i+1] - ai[i]; 967 acol = aj + ai[i]; 968 aval = aa + ai[i]; 969 cnzi = ci[i+1] - ci[i]; 970 ccol = cj + ci[i]; 971 cval = ca + ci[i]; 972 for (j=0; j<cnzi; j++) { 973 brow = ccol[j]; 974 bnzj = bi[brow+1] - bi[brow]; 975 bcol = bj + bi[brow]; 976 bval = ba + bi[brow]; 977 978 /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 979 nexta = 0; nextb = 0; 980 while (nexta<anzi && nextb<bnzj) { 981 while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++; 982 if (nexta == anzi) break; 983 while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++; 984 if (nextb == bnzj) break; 985 if (acol[nexta] == bcol[nextb]) { 986 cval[j] += aval[nexta]*bval[nextb]; 987 nexta++; nextb++; 988 flops += 2; 989 } 990 } 991 } 992 } 993 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 994 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 995 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 996 PetscFunctionReturn(0); 997 } 998 999 #undef __FUNCT__ 1000 #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ" 1001 PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1002 { 1003 PetscErrorCode ierr; 1004 1005 PetscFunctionBegin; 1006 if (scall == MAT_INITIAL_MATRIX) { 1007 ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1008 ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 1009 ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1010 } 1011 ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1012 ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 1013 ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1014 PetscFunctionReturn(0); 1015 } 1016 1017 #undef __FUNCT__ 1018 #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ" 1019 PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 1020 { 1021 PetscErrorCode ierr; 1022 Mat At; 1023 PetscInt *ati,*atj; 1024 1025 PetscFunctionBegin; 1026 /* create symbolic At */ 1027 ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 1028 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);CHKERRQ(ierr); 1029 1030 At->rmap->bs = A->cmap->bs; 1031 At->cmap->bs = B->cmap->bs; 1032 1033 /* get symbolic C=At*B */ 1034 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 1035 1036 /* clean up */ 1037 ierr = MatDestroy(&At);CHKERRQ(ierr); 1038 ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 1039 PetscFunctionReturn(0); 1040 } 1041 1042 #undef __FUNCT__ 1043 #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ" 1044 PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 1045 { 1046 PetscErrorCode ierr; 1047 Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 1048 PetscInt am =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 1049 PetscInt cm =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 1050 PetscLogDouble flops=0.0; 1051 MatScalar *aa =a->a,*ba,*ca,*caj; 1052 1053 PetscFunctionBegin; 1054 if (!c->a) { 1055 ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 1056 1057 c->a = ca; 1058 c->free_a = PETSC_TRUE; 1059 } else { 1060 ca = c->a; 1061 } 1062 /* clear old values in C */ 1063 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 1064 1065 /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 1066 for (i=0; i<am; i++) { 1067 bj = b->j + bi[i]; 1068 ba = b->a + bi[i]; 1069 bnzi = bi[i+1] - bi[i]; 1070 anzi = ai[i+1] - ai[i]; 1071 for (j=0; j<anzi; j++) { 1072 nextb = 0; 1073 crow = *aj++; 1074 cjj = cj + ci[crow]; 1075 caj = ca + ci[crow]; 1076 /* perform sparse axpy operation. Note cjj includes bj. */ 1077 for (k=0; nextb<bnzi; k++) { 1078 if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 1079 caj[k] += (*aa)*(*(ba+nextb)); 1080 nextb++; 1081 } 1082 } 1083 flops += 2*bnzi; 1084 aa++; 1085 } 1086 } 1087 1088 /* Assemble the final matrix and clean up */ 1089 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1090 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1091 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 1092 PetscFunctionReturn(0); 1093 } 1094 1095 #undef __FUNCT__ 1096 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense" 1097 PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1098 { 1099 PetscErrorCode ierr; 1100 1101 PetscFunctionBegin; 1102 if (scall == MAT_INITIAL_MATRIX) { 1103 ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1104 ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1105 ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1106 } 1107 ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1108 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 1109 ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1110 PetscFunctionReturn(0); 1111 } 1112 1113 #undef __FUNCT__ 1114 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense" 1115 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1116 { 1117 PetscErrorCode ierr; 1118 1119 PetscFunctionBegin; 1120 ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr); 1121 1122 (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense; 1123 PetscFunctionReturn(0); 1124 } 1125 1126 #undef __FUNCT__ 1127 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense" 1128 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 1129 { 1130 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 1131 PetscErrorCode ierr; 1132 PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 1133 MatScalar *aa; 1134 PetscInt cm = C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n; 1135 PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam; 1136 1137 PetscFunctionBegin; 1138 if (!cm || !cn) PetscFunctionReturn(0); 1139 if (bm != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in A %D not equal rows in B %D\n",A->cmap->n,bm); 1140 if (A->rmap->n != C->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number rows in C %D not equal rows in A %D\n",C->rmap->n,A->rmap->n); 1141 if (B->cmap->n != C->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Number columns in B %D not equal columns in C %D\n",B->cmap->n,C->cmap->n); 1142 ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 1143 ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr); 1144 b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 1145 for (col=0; col<cn-4; col += 4) { /* over columns of C */ 1146 colam = col*am; 1147 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1148 r1 = r2 = r3 = r4 = 0.0; 1149 n = a->i[i+1] - a->i[i]; 1150 aj = a->j + a->i[i]; 1151 aa = a->a + a->i[i]; 1152 for (j=0; j<n; j++) { 1153 r1 += (*aa)*b1[*aj]; 1154 r2 += (*aa)*b2[*aj]; 1155 r3 += (*aa)*b3[*aj]; 1156 r4 += (*aa++)*b4[*aj++]; 1157 } 1158 c[colam + i] = r1; 1159 c[colam + am + i] = r2; 1160 c[colam + am2 + i] = r3; 1161 c[colam + am3 + i] = r4; 1162 } 1163 b1 += bm4; 1164 b2 += bm4; 1165 b3 += bm4; 1166 b4 += bm4; 1167 } 1168 for (; col<cn; col++) { /* over extra columns of C */ 1169 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1170 r1 = 0.0; 1171 n = a->i[i+1] - a->i[i]; 1172 aj = a->j + a->i[i]; 1173 aa = a->a + a->i[i]; 1174 1175 for (j=0; j<n; j++) { 1176 r1 += (*aa++)*b1[*aj++]; 1177 } 1178 c[col*am + i] = r1; 1179 } 1180 b1 += bm; 1181 } 1182 ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 1183 ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 1184 ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr); 1185 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1186 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1187 PetscFunctionReturn(0); 1188 } 1189 1190 /* 1191 Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 1192 */ 1193 #undef __FUNCT__ 1194 #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense" 1195 PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 1196 { 1197 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 1198 PetscErrorCode ierr; 1199 PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 1200 MatScalar *aa; 1201 PetscInt cm = C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n,*ii,arm; 1202 PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 1203 1204 PetscFunctionBegin; 1205 if (!cm || !cn) PetscFunctionReturn(0); 1206 ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 1207 ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr); 1208 b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 1209 1210 if (a->compressedrow.use) { /* use compressed row format */ 1211 for (col=0; col<cn-4; col += 4) { /* over columns of C */ 1212 colam = col*am; 1213 arm = a->compressedrow.nrows; 1214 ii = a->compressedrow.i; 1215 ridx = a->compressedrow.rindex; 1216 for (i=0; i<arm; i++) { /* over rows of C in those columns */ 1217 r1 = r2 = r3 = r4 = 0.0; 1218 n = ii[i+1] - ii[i]; 1219 aj = a->j + ii[i]; 1220 aa = a->a + ii[i]; 1221 for (j=0; j<n; j++) { 1222 r1 += (*aa)*b1[*aj]; 1223 r2 += (*aa)*b2[*aj]; 1224 r3 += (*aa)*b3[*aj]; 1225 r4 += (*aa++)*b4[*aj++]; 1226 } 1227 c[colam + ridx[i]] += r1; 1228 c[colam + am + ridx[i]] += r2; 1229 c[colam + am2 + ridx[i]] += r3; 1230 c[colam + am3 + ridx[i]] += r4; 1231 } 1232 b1 += bm4; 1233 b2 += bm4; 1234 b3 += bm4; 1235 b4 += bm4; 1236 } 1237 for (; col<cn; col++) { /* over extra columns of C */ 1238 colam = col*am; 1239 arm = a->compressedrow.nrows; 1240 ii = a->compressedrow.i; 1241 ridx = a->compressedrow.rindex; 1242 for (i=0; i<arm; i++) { /* over rows of C in those columns */ 1243 r1 = 0.0; 1244 n = ii[i+1] - ii[i]; 1245 aj = a->j + ii[i]; 1246 aa = a->a + ii[i]; 1247 1248 for (j=0; j<n; j++) { 1249 r1 += (*aa++)*b1[*aj++]; 1250 } 1251 c[colam + ridx[i]] += r1; 1252 } 1253 b1 += bm; 1254 } 1255 } else { 1256 for (col=0; col<cn-4; col += 4) { /* over columns of C */ 1257 colam = col*am; 1258 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1259 r1 = r2 = r3 = r4 = 0.0; 1260 n = a->i[i+1] - a->i[i]; 1261 aj = a->j + a->i[i]; 1262 aa = a->a + a->i[i]; 1263 for (j=0; j<n; j++) { 1264 r1 += (*aa)*b1[*aj]; 1265 r2 += (*aa)*b2[*aj]; 1266 r3 += (*aa)*b3[*aj]; 1267 r4 += (*aa++)*b4[*aj++]; 1268 } 1269 c[colam + i] += r1; 1270 c[colam + am + i] += r2; 1271 c[colam + am2 + i] += r3; 1272 c[colam + am3 + i] += r4; 1273 } 1274 b1 += bm4; 1275 b2 += bm4; 1276 b3 += bm4; 1277 b4 += bm4; 1278 } 1279 for (; col<cn; col++) { /* over extra columns of C */ 1280 colam = col*am; 1281 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1282 r1 = 0.0; 1283 n = a->i[i+1] - a->i[i]; 1284 aj = a->j + a->i[i]; 1285 aa = a->a + a->i[i]; 1286 1287 for (j=0; j<n; j++) { 1288 r1 += (*aa++)*b1[*aj++]; 1289 } 1290 c[colam + i] += r1; 1291 } 1292 b1 += bm; 1293 } 1294 } 1295 ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr); 1296 ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 1297 ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr); 1298 PetscFunctionReturn(0); 1299 } 1300 1301 #undef __FUNCT__ 1302 #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ" 1303 PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense) 1304 { 1305 PetscErrorCode ierr; 1306 Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 1307 Mat_SeqDense *btdense = (Mat_SeqDense*)Btdense->data; 1308 PetscInt *bi = b->i,*bj=b->j; 1309 PetscInt m = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns; 1310 MatScalar *btval,*btval_den,*ba=b->a; 1311 PetscInt *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors; 1312 1313 PetscFunctionBegin; 1314 btval_den=btdense->v; 1315 ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr); 1316 for (k=0; k<ncolors; k++) { 1317 ncolumns = coloring->ncolumns[k]; 1318 for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */ 1319 col = *(columns + colorforcol[k] + l); 1320 btcol = bj + bi[col]; 1321 btval = ba + bi[col]; 1322 anz = bi[col+1] - bi[col]; 1323 for (j=0; j<anz; j++) { 1324 brow = btcol[j]; 1325 btval_den[brow] = btval[j]; 1326 } 1327 } 1328 btval_den += m; 1329 } 1330 PetscFunctionReturn(0); 1331 } 1332 1333 #undef __FUNCT__ 1334 #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ" 1335 PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp) 1336 { 1337 PetscErrorCode ierr; 1338 Mat_SeqAIJ *csp = (Mat_SeqAIJ*)Csp->data; 1339 PetscScalar *ca_den,*ca_den_ptr,*ca=csp->a; 1340 PetscInt k,l,m=Cden->rmap->n,ncolors=matcoloring->ncolors; 1341 PetscInt brows=matcoloring->brows,*den2sp=matcoloring->den2sp; 1342 PetscInt nrows,*row,*idx; 1343 PetscInt *rows=matcoloring->rows,*colorforrow=matcoloring->colorforrow; 1344 1345 PetscFunctionBegin; 1346 ierr = MatDenseGetArray(Cden,&ca_den);CHKERRQ(ierr); 1347 1348 if (brows > 0) { 1349 PetscInt *lstart,row_end,row_start; 1350 lstart = matcoloring->lstart; 1351 ierr = PetscMemzero(lstart,ncolors*sizeof(PetscInt));CHKERRQ(ierr); 1352 1353 row_end = brows; 1354 if (row_end > m) row_end = m; 1355 for (row_start=0; row_start<m; row_start+=brows) { /* loop over row blocks of Csp */ 1356 ca_den_ptr = ca_den; 1357 for (k=0; k<ncolors; k++) { /* loop over colors (columns of Cden) */ 1358 nrows = matcoloring->nrows[k]; 1359 row = rows + colorforrow[k]; 1360 idx = den2sp + colorforrow[k]; 1361 for (l=lstart[k]; l<nrows; l++) { 1362 if (row[l] >= row_end) { 1363 lstart[k] = l; 1364 break; 1365 } else { 1366 ca[idx[l]] = ca_den_ptr[row[l]]; 1367 } 1368 } 1369 ca_den_ptr += m; 1370 } 1371 row_end += brows; 1372 if (row_end > m) row_end = m; 1373 } 1374 } else { /* non-blocked impl: loop over columns of Csp - slow if Csp is large */ 1375 ca_den_ptr = ca_den; 1376 for (k=0; k<ncolors; k++) { 1377 nrows = matcoloring->nrows[k]; 1378 row = rows + colorforrow[k]; 1379 idx = den2sp + colorforrow[k]; 1380 for (l=0; l<nrows; l++) { 1381 ca[idx[l]] = ca_den_ptr[row[l]]; 1382 } 1383 ca_den_ptr += m; 1384 } 1385 } 1386 1387 ierr = MatDenseRestoreArray(Cden,&ca_den);CHKERRQ(ierr); 1388 #if defined(PETSC_USE_INFO) 1389 if (matcoloring->brows > 0) { 1390 ierr = PetscInfo1(Csp,"Loop over %D row blocks for den2sp\n",brows);CHKERRQ(ierr); 1391 } else { 1392 ierr = PetscInfo(Csp,"Loop over colors/columns of Cden, inefficient for large sparse matrix product \n");CHKERRQ(ierr); 1393 } 1394 #endif 1395 PetscFunctionReturn(0); 1396 } 1397 1398 /* 1399 MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from 1400 MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output 1401 spidx[], index of a->a, to be used in MatTransposeColoringCreate_SeqAIJ(). 1402 */ 1403 #undef __FUNCT__ 1404 #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color" 1405 PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *nn,const PetscInt *ia[],const PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1406 { 1407 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1408 PetscErrorCode ierr; 1409 PetscInt i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n; 1410 PetscInt nz = a->i[m],row,*jj,mr,col; 1411 PetscInt *cspidx; 1412 1413 PetscFunctionBegin; 1414 *nn = n; 1415 if (!ia) PetscFunctionReturn(0); 1416 if (symmetric) { 1417 SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric"); 1418 ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,(PetscInt**)ia,(PetscInt**)ja);CHKERRQ(ierr); 1419 } else { 1420 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr); 1421 ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 1422 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr); 1423 ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr); 1424 ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr); 1425 jj = a->j; 1426 for (i=0; i<nz; i++) { 1427 collengths[jj[i]]++; 1428 } 1429 cia[0] = oshift; 1430 for (i=0; i<n; i++) { 1431 cia[i+1] = cia[i] + collengths[i]; 1432 } 1433 ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 1434 jj = a->j; 1435 for (row=0; row<m; row++) { 1436 mr = a->i[row+1] - a->i[row]; 1437 for (i=0; i<mr; i++) { 1438 col = *jj++; 1439 1440 cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */ 1441 cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 1442 } 1443 } 1444 ierr = PetscFree(collengths);CHKERRQ(ierr); 1445 *ia = cia; *ja = cja; 1446 *spidx = cspidx; 1447 } 1448 PetscFunctionReturn(0); 1449 } 1450 1451 #undef __FUNCT__ 1452 #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color" 1453 PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1454 { 1455 PetscErrorCode ierr; 1456 1457 PetscFunctionBegin; 1458 if (!ia) PetscFunctionReturn(0); 1459 1460 ierr = PetscFree(*ia);CHKERRQ(ierr); 1461 ierr = PetscFree(*ja);CHKERRQ(ierr); 1462 ierr = PetscFree(*spidx);CHKERRQ(ierr); 1463 PetscFunctionReturn(0); 1464 } 1465 1466 #undef __FUNCT__ 1467 #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ" 1468 PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c) 1469 { 1470 PetscErrorCode ierr; 1471 PetscInt i,n,nrows,Nbs,j,k,m,ncols,col,cm; 1472 const PetscInt *is,*ci,*cj,*row_idx; 1473 PetscInt nis = iscoloring->n,*rowhit,bs = 1; 1474 IS *isa; 1475 Mat_SeqAIJ *csp = (Mat_SeqAIJ*)mat->data; 1476 PetscInt *colorforrow,*rows,*rows_i,*idxhit,*spidx,*den2sp,*den2sp_i; 1477 PetscInt *colorforcol,*columns,*columns_i,brows; 1478 PetscBool flg; 1479 1480 PetscFunctionBegin; 1481 ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 1482 1483 /* bs >1 is not being tested yet! */ 1484 Nbs = mat->cmap->N/bs; 1485 c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 1486 c->N = Nbs; 1487 c->m = c->M; 1488 c->rstart = 0; 1489 c->brows = 100; 1490 1491 c->ncolors = nis; 1492 ierr = PetscMalloc3(nis,PetscInt,&c->ncolumns,nis,PetscInt,&c->nrows,nis+1,PetscInt,&colorforrow);CHKERRQ(ierr); 1493 ierr = PetscMalloc((csp->nz+1)*sizeof(PetscInt),&rows);CHKERRQ(ierr); 1494 ierr = PetscMalloc((csp->nz+1)*sizeof(PetscInt),&den2sp);CHKERRQ(ierr); 1495 1496 brows = c->brows; 1497 ierr = PetscOptionsGetInt(NULL,"-matden2sp_brows",&brows,&flg);CHKERRQ(ierr); 1498 if (flg) c->brows = brows; 1499 if (brows > 0) { 1500 ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&c->lstart);CHKERRQ(ierr); 1501 } 1502 1503 colorforrow[0] = 0; 1504 rows_i = rows; 1505 den2sp_i = den2sp; 1506 1507 ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr); 1508 ierr = PetscMalloc((Nbs+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr); 1509 1510 colorforcol[0] = 0; 1511 columns_i = columns; 1512 1513 /* get column-wise storage of mat */ 1514 ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 1515 1516 cm = c->m; 1517 ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 1518 ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr); 1519 for (i=0; i<nis; i++) { /* loop over color */ 1520 ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 1521 ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 1522 1523 c->ncolumns[i] = n; 1524 if (n) { 1525 ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr); 1526 } 1527 colorforcol[i+1] = colorforcol[i] + n; 1528 columns_i += n; 1529 1530 /* fast, crude version requires O(N*N) work */ 1531 ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr); 1532 1533 /* loop over columns*/ 1534 for (j=0; j<n; j++) { 1535 col = is[j]; 1536 row_idx = cj + ci[col]; 1537 m = ci[col+1] - ci[col]; 1538 /* loop over columns marking them in rowhit */ 1539 for (k=0; k<m; k++) { 1540 idxhit[*row_idx] = spidx[ci[col] + k]; 1541 rowhit[*row_idx++] = col + 1; 1542 } 1543 } 1544 /* count the number of hits */ 1545 nrows = 0; 1546 for (j=0; j<cm; j++) { 1547 if (rowhit[j]) nrows++; 1548 } 1549 c->nrows[i] = nrows; 1550 colorforrow[i+1] = colorforrow[i] + nrows; 1551 1552 nrows = 0; 1553 for (j=0; j<cm; j++) { 1554 if (rowhit[j]) { 1555 rows_i[nrows] = j; 1556 den2sp_i[nrows] = idxhit[j]; 1557 nrows++; 1558 } 1559 } 1560 den2sp_i += nrows; 1561 1562 ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 1563 rows_i += nrows; 1564 } 1565 ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 1566 ierr = PetscFree(rowhit);CHKERRQ(ierr); 1567 ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 1568 if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]); 1569 1570 c->colorforrow = colorforrow; 1571 c->rows = rows; 1572 c->den2sp = den2sp; 1573 c->colorforcol = colorforcol; 1574 c->columns = columns; 1575 1576 ierr = PetscFree(idxhit);CHKERRQ(ierr); 1577 PetscFunctionReturn(0); 1578 } 1579