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 13 static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat,Mat,PetscReal,Mat*); 14 15 #undef __FUNCT__ 16 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ" 17 PetscErrorCode MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 18 { 19 PetscErrorCode ierr; 20 const char *algTypes[6] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed"}; 21 PetscInt alg=0; /* set default algorithm */ 22 23 PetscFunctionBegin; 24 if (scall == MAT_INITIAL_MATRIX) { 25 ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 26 ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,6,algTypes[0],&alg,NULL);CHKERRQ(ierr); 27 ierr = PetscOptionsEnd();CHKERRQ(ierr); 28 ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 29 switch (alg) { 30 case 1: 31 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr); 32 break; 33 case 2: 34 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr); 35 break; 36 case 3: 37 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr); 38 break; 39 case 4: 40 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);CHKERRQ(ierr); 41 break; 42 case 5: 43 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(A,B,fill,C);CHKERRQ(ierr); 44 break; 45 default: 46 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 47 break; 48 } 49 ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 50 } 51 52 ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 53 ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 54 ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 55 PetscFunctionReturn(0); 56 } 57 58 #undef __FUNCT__ 59 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed" 60 static PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat *C) 61 { 62 PetscErrorCode ierr; 63 Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 64 PetscInt *ai=a->i,*bi=b->i,*ci,*cj; 65 PetscInt am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 66 PetscReal afill; 67 PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,nlnk_max,*lnk,ndouble=0; 68 PetscBT lnkbt; 69 PetscFreeSpaceList free_space=NULL,current_space=NULL; 70 71 PetscFunctionBegin; 72 /* Get ci and cj */ 73 /*---------------*/ 74 /* Allocate ci array, arrays for fill computation and */ 75 /* free space for accumulating nonzero column info */ 76 ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 77 ci[0] = 0; 78 79 /* create and initialize a linked list */ 80 nlnk_max = a->rmax*b->rmax; 81 if (!nlnk_max || nlnk_max > bn) nlnk_max = bn; 82 ierr = PetscLLCondensedCreate(nlnk_max,bn,&lnk,&lnkbt);CHKERRQ(ierr); 83 84 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 85 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 86 87 current_space = free_space; 88 89 /* Determine ci and cj */ 90 for (i=0; i<am; i++) { 91 anzi = ai[i+1] - ai[i]; 92 aj = a->j + ai[i]; 93 for (j=0; j<anzi; j++) { 94 brow = aj[j]; 95 bnzj = bi[brow+1] - bi[brow]; 96 bj = b->j + bi[brow]; 97 /* add non-zero cols of B into the sorted linked list lnk */ 98 ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr); 99 } 100 cnzi = lnk[0]; 101 102 /* If free space is not available, make more free space */ 103 /* Double the amount of total space in the list */ 104 if (current_space->local_remaining<cnzi) { 105 ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 106 ndouble++; 107 } 108 109 /* Copy data into free space, then initialize lnk */ 110 ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 111 112 current_space->array += cnzi; 113 current_space->local_used += cnzi; 114 current_space->local_remaining -= cnzi; 115 116 ci[i+1] = ci[i] + cnzi; 117 } 118 119 /* Column indices are in the list of free space */ 120 /* Allocate space for cj, initialize cj, and */ 121 /* destroy list of free space and other temporary array(s) */ 122 ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 123 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 124 ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr); 125 126 /* put together the new symbolic matrix */ 127 ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,C);CHKERRQ(ierr); 128 129 (*C)->rmap->bs = A->rmap->bs; 130 (*C)->cmap->bs = B->cmap->bs; 131 132 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 133 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 134 c = (Mat_SeqAIJ*)((*C)->data); 135 c->free_a = PETSC_FALSE; 136 c->free_ij = PETSC_TRUE; 137 c->nonew = 0; 138 (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, needs non-scalable O(bn) array 'abdense' */ 139 140 /* set MatInfo */ 141 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 142 if (afill < 1.0) afill = 1.0; 143 c->maxnz = ci[am]; 144 c->nz = ci[am]; 145 (*C)->info.mallocs = ndouble; 146 (*C)->info.fill_ratio_given = fill; 147 (*C)->info.fill_ratio_needed = afill; 148 149 #if defined(PETSC_USE_INFO) 150 if (ci[am]) { 151 ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 152 ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 153 } else { 154 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 155 } 156 #endif 157 PetscFunctionReturn(0); 158 } 159 160 #undef __FUNCT__ 161 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ" 162 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 163 { 164 PetscErrorCode ierr; 165 PetscLogDouble flops=0.0; 166 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 167 Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 168 Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 169 PetscInt *ai =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 170 PetscInt am =A->rmap->n,cm=C->rmap->n; 171 PetscInt i,j,k,anzi,bnzi,cnzi,brow; 172 PetscScalar *aa=a->a,*ba=b->a,*baj,*ca,valtmp; 173 PetscScalar *ab_dense; 174 175 PetscFunctionBegin; 176 if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */ 177 ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 178 c->a = ca; 179 c->free_a = PETSC_TRUE; 180 } else { 181 ca = c->a; 182 } 183 if (!c->matmult_abdense) { 184 ierr = PetscMalloc(B->cmap->N*sizeof(PetscScalar),&ab_dense);CHKERRQ(ierr); 185 ierr = PetscMemzero(ab_dense,B->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr); 186 c->matmult_abdense = ab_dense; 187 } else { 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 MatColoring coloring; 884 ISColoring iscoloring; 885 Mat Bt_dense,C_dense; 886 Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data; 887 /* inode causes memory problem, don't know why */ 888 if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'"); 889 890 ierr = MatColoringCreate(*C,&coloring);CHKERRQ(ierr); 891 ierr = MatColoringSetDistance(coloring,2);CHKERRQ(ierr); 892 ierr = MatColoringSetType(coloring,MATCOLORINGSL);CHKERRQ(ierr); 893 ierr = MatColoringSetFromOptions(coloring);CHKERRQ(ierr); 894 ierr = MatColoringApply(coloring,&iscoloring);CHKERRQ(ierr); 895 ierr = MatColoringDestroy(&coloring);CHKERRQ(ierr); 896 ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 897 898 abt->matcoloring = matcoloring; 899 900 ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 901 902 /* Create Bt_dense and C_dense = A*Bt_dense */ 903 ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr); 904 ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 905 ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr); 906 ierr = MatSeqDenseSetPreallocation(Bt_dense,NULL);CHKERRQ(ierr); 907 908 Bt_dense->assembled = PETSC_TRUE; 909 abt->Bt_den = Bt_dense; 910 911 ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr); 912 ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 913 ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr); 914 ierr = MatSeqDenseSetPreallocation(C_dense,NULL);CHKERRQ(ierr); 915 916 Bt_dense->assembled = PETSC_TRUE; 917 abt->ABt_den = C_dense; 918 919 #if defined(PETSC_USE_INFO) 920 { 921 Mat_SeqAIJ *c = (Mat_SeqAIJ*)(*C)->data; 922 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); 923 } 924 #endif 925 } 926 /* clean up */ 927 ierr = MatDestroy(&Bt);CHKERRQ(ierr); 928 ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 929 PetscFunctionReturn(0); 930 } 931 932 #undef __FUNCT__ 933 #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ" 934 PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 935 { 936 PetscErrorCode ierr; 937 Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 938 PetscInt *ai =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow; 939 PetscInt cm =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol; 940 PetscLogDouble flops=0.0; 941 MatScalar *aa =a->a,*aval,*ba=b->a,*bval,*ca,*cval; 942 Mat_MatMatTransMult *abt = c->abt; 943 944 PetscFunctionBegin; 945 /* clear old values in C */ 946 if (!c->a) { 947 ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 948 c->a = ca; 949 c->free_a = PETSC_TRUE; 950 } else { 951 ca = c->a; 952 } 953 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 954 955 if (abt->usecoloring) { 956 MatTransposeColoring matcoloring = abt->matcoloring; 957 Mat Bt_dense,C_dense = abt->ABt_den; 958 959 /* Get Bt_dense by Apply MatTransposeColoring to B */ 960 Bt_dense = abt->Bt_den; 961 ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr); 962 963 /* C_dense = A*Bt_dense */ 964 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr); 965 966 /* Recover C from C_dense */ 967 ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr); 968 PetscFunctionReturn(0); 969 } 970 971 for (i=0; i<cm; i++) { 972 anzi = ai[i+1] - ai[i]; 973 acol = aj + ai[i]; 974 aval = aa + ai[i]; 975 cnzi = ci[i+1] - ci[i]; 976 ccol = cj + ci[i]; 977 cval = ca + ci[i]; 978 for (j=0; j<cnzi; j++) { 979 brow = ccol[j]; 980 bnzj = bi[brow+1] - bi[brow]; 981 bcol = bj + bi[brow]; 982 bval = ba + bi[brow]; 983 984 /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 985 nexta = 0; nextb = 0; 986 while (nexta<anzi && nextb<bnzj) { 987 while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++; 988 if (nexta == anzi) break; 989 while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++; 990 if (nextb == bnzj) break; 991 if (acol[nexta] == bcol[nextb]) { 992 cval[j] += aval[nexta]*bval[nextb]; 993 nexta++; nextb++; 994 flops += 2; 995 } 996 } 997 } 998 } 999 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1000 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1001 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 1002 PetscFunctionReturn(0); 1003 } 1004 1005 #undef __FUNCT__ 1006 #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ" 1007 PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1008 { 1009 PetscErrorCode ierr; 1010 1011 PetscFunctionBegin; 1012 if (scall == MAT_INITIAL_MATRIX) { 1013 ierr = PetscLogEventBegin(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1014 ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 1015 ierr = PetscLogEventEnd(MAT_TransposeMatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1016 } 1017 ierr = PetscLogEventBegin(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1018 ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 1019 ierr = PetscLogEventEnd(MAT_TransposeMatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1020 PetscFunctionReturn(0); 1021 } 1022 1023 #undef __FUNCT__ 1024 #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ" 1025 PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 1026 { 1027 PetscErrorCode ierr; 1028 Mat At; 1029 PetscInt *ati,*atj; 1030 1031 PetscFunctionBegin; 1032 /* create symbolic At */ 1033 ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 1034 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);CHKERRQ(ierr); 1035 1036 At->rmap->bs = A->cmap->bs; 1037 At->cmap->bs = B->cmap->bs; 1038 1039 /* get symbolic C=At*B */ 1040 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 1041 1042 /* clean up */ 1043 ierr = MatDestroy(&At);CHKERRQ(ierr); 1044 ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 1045 PetscFunctionReturn(0); 1046 } 1047 1048 #undef __FUNCT__ 1049 #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ" 1050 PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 1051 { 1052 PetscErrorCode ierr; 1053 Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 1054 PetscInt am =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 1055 PetscInt cm =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 1056 PetscLogDouble flops=0.0; 1057 MatScalar *aa =a->a,*ba,*ca,*caj; 1058 1059 PetscFunctionBegin; 1060 if (!c->a) { 1061 ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 1062 1063 c->a = ca; 1064 c->free_a = PETSC_TRUE; 1065 } else { 1066 ca = c->a; 1067 } 1068 /* clear old values in C */ 1069 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 1070 1071 /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 1072 for (i=0; i<am; i++) { 1073 bj = b->j + bi[i]; 1074 ba = b->a + bi[i]; 1075 bnzi = bi[i+1] - bi[i]; 1076 anzi = ai[i+1] - ai[i]; 1077 for (j=0; j<anzi; j++) { 1078 nextb = 0; 1079 crow = *aj++; 1080 cjj = cj + ci[crow]; 1081 caj = ca + ci[crow]; 1082 /* perform sparse axpy operation. Note cjj includes bj. */ 1083 for (k=0; nextb<bnzi; k++) { 1084 if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 1085 caj[k] += (*aa)*(*(ba+nextb)); 1086 nextb++; 1087 } 1088 } 1089 flops += 2*bnzi; 1090 aa++; 1091 } 1092 } 1093 1094 /* Assemble the final matrix and clean up */ 1095 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1096 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1097 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 1098 PetscFunctionReturn(0); 1099 } 1100 1101 #undef __FUNCT__ 1102 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense" 1103 PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1104 { 1105 PetscErrorCode ierr; 1106 1107 PetscFunctionBegin; 1108 if (scall == MAT_INITIAL_MATRIX) { 1109 ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1110 ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1111 ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1112 } 1113 ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1114 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 1115 ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1116 PetscFunctionReturn(0); 1117 } 1118 1119 #undef __FUNCT__ 1120 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense" 1121 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1122 { 1123 PetscErrorCode ierr; 1124 1125 PetscFunctionBegin; 1126 ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr); 1127 1128 (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense; 1129 PetscFunctionReturn(0); 1130 } 1131 1132 #undef __FUNCT__ 1133 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense" 1134 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 1135 { 1136 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 1137 PetscErrorCode ierr; 1138 PetscScalar *c,*b,r1,r2,r3,r4,*c1,*c2,*c3,*c4,aatmp; 1139 const PetscScalar *aa,*b1,*b2,*b3,*b4; 1140 const PetscInt *aj; 1141 PetscInt cm=C->rmap->n,cn=B->cmap->n,bm=B->rmap->n,am=A->rmap->n; 1142 PetscInt am4=4*am,bm4=4*bm,col,i,j,n,ajtmp; 1143 1144 PetscFunctionBegin; 1145 if (!cm || !cn) PetscFunctionReturn(0); 1146 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); 1147 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); 1148 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); 1149 ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 1150 ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr); 1151 b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 1152 c1 = c; c2 = c1 + am; c3 = c2 + am; c4 = c3 + am; 1153 for (col=0; col<cn-4; col += 4) { /* over columns of C */ 1154 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1155 r1 = r2 = r3 = r4 = 0.0; 1156 n = a->i[i+1] - a->i[i]; 1157 aj = a->j + a->i[i]; 1158 aa = a->a + a->i[i]; 1159 for (j=0; j<n; j++) { 1160 aatmp = aa[j]; ajtmp = aj[j]; 1161 r1 += aatmp*b1[ajtmp]; 1162 r2 += aatmp*b2[ajtmp]; 1163 r3 += aatmp*b3[ajtmp]; 1164 r4 += aatmp*b4[ajtmp]; 1165 } 1166 c1[i] = r1; 1167 c2[i] = r2; 1168 c3[i] = r3; 1169 c4[i] = r4; 1170 } 1171 b1 += bm4; b2 += bm4; b3 += bm4; b4 += bm4; 1172 c1 += am4; c2 += am4; c3 += am4; c4 += am4; 1173 } 1174 for (; col<cn; col++) { /* over extra columns of C */ 1175 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1176 r1 = 0.0; 1177 n = a->i[i+1] - a->i[i]; 1178 aj = a->j + a->i[i]; 1179 aa = a->a + a->i[i]; 1180 for (j=0; j<n; j++) { 1181 r1 += aa[j]*b1[aj[j]]; 1182 } 1183 c1[i] = r1; 1184 } 1185 b1 += bm; 1186 c1 += am; 1187 } 1188 ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 1189 ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 1190 ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr); 1191 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1192 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1193 PetscFunctionReturn(0); 1194 } 1195 1196 /* 1197 Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 1198 */ 1199 #undef __FUNCT__ 1200 #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense" 1201 PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 1202 { 1203 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 1204 PetscErrorCode ierr; 1205 PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 1206 MatScalar *aa; 1207 PetscInt cm = C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n,*ii,arm; 1208 PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 1209 1210 PetscFunctionBegin; 1211 if (!cm || !cn) PetscFunctionReturn(0); 1212 ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 1213 ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr); 1214 b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 1215 1216 if (a->compressedrow.use) { /* use compressed row format */ 1217 for (col=0; col<cn-4; col += 4) { /* over columns of C */ 1218 colam = col*am; 1219 arm = a->compressedrow.nrows; 1220 ii = a->compressedrow.i; 1221 ridx = a->compressedrow.rindex; 1222 for (i=0; i<arm; i++) { /* over rows of C in those columns */ 1223 r1 = r2 = r3 = r4 = 0.0; 1224 n = ii[i+1] - ii[i]; 1225 aj = a->j + ii[i]; 1226 aa = a->a + ii[i]; 1227 for (j=0; j<n; j++) { 1228 r1 += (*aa)*b1[*aj]; 1229 r2 += (*aa)*b2[*aj]; 1230 r3 += (*aa)*b3[*aj]; 1231 r4 += (*aa++)*b4[*aj++]; 1232 } 1233 c[colam + ridx[i]] += r1; 1234 c[colam + am + ridx[i]] += r2; 1235 c[colam + am2 + ridx[i]] += r3; 1236 c[colam + am3 + ridx[i]] += r4; 1237 } 1238 b1 += bm4; 1239 b2 += bm4; 1240 b3 += bm4; 1241 b4 += bm4; 1242 } 1243 for (; col<cn; col++) { /* over extra columns of C */ 1244 colam = col*am; 1245 arm = a->compressedrow.nrows; 1246 ii = a->compressedrow.i; 1247 ridx = a->compressedrow.rindex; 1248 for (i=0; i<arm; i++) { /* over rows of C in those columns */ 1249 r1 = 0.0; 1250 n = ii[i+1] - ii[i]; 1251 aj = a->j + ii[i]; 1252 aa = a->a + ii[i]; 1253 1254 for (j=0; j<n; j++) { 1255 r1 += (*aa++)*b1[*aj++]; 1256 } 1257 c[colam + ridx[i]] += r1; 1258 } 1259 b1 += bm; 1260 } 1261 } else { 1262 for (col=0; col<cn-4; col += 4) { /* over columns of C */ 1263 colam = col*am; 1264 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1265 r1 = r2 = r3 = r4 = 0.0; 1266 n = a->i[i+1] - a->i[i]; 1267 aj = a->j + a->i[i]; 1268 aa = a->a + a->i[i]; 1269 for (j=0; j<n; j++) { 1270 r1 += (*aa)*b1[*aj]; 1271 r2 += (*aa)*b2[*aj]; 1272 r3 += (*aa)*b3[*aj]; 1273 r4 += (*aa++)*b4[*aj++]; 1274 } 1275 c[colam + i] += r1; 1276 c[colam + am + i] += r2; 1277 c[colam + am2 + i] += r3; 1278 c[colam + am3 + i] += r4; 1279 } 1280 b1 += bm4; 1281 b2 += bm4; 1282 b3 += bm4; 1283 b4 += bm4; 1284 } 1285 for (; col<cn; col++) { /* over extra columns of C */ 1286 colam = col*am; 1287 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1288 r1 = 0.0; 1289 n = a->i[i+1] - a->i[i]; 1290 aj = a->j + a->i[i]; 1291 aa = a->a + a->i[i]; 1292 1293 for (j=0; j<n; j++) { 1294 r1 += (*aa++)*b1[*aj++]; 1295 } 1296 c[colam + i] += r1; 1297 } 1298 b1 += bm; 1299 } 1300 } 1301 ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr); 1302 ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 1303 ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr); 1304 PetscFunctionReturn(0); 1305 } 1306 1307 #undef __FUNCT__ 1308 #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ" 1309 PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense) 1310 { 1311 PetscErrorCode ierr; 1312 Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 1313 Mat_SeqDense *btdense = (Mat_SeqDense*)Btdense->data; 1314 PetscInt *bi = b->i,*bj=b->j; 1315 PetscInt m = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns; 1316 MatScalar *btval,*btval_den,*ba=b->a; 1317 PetscInt *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors; 1318 1319 PetscFunctionBegin; 1320 btval_den=btdense->v; 1321 ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr); 1322 for (k=0; k<ncolors; k++) { 1323 ncolumns = coloring->ncolumns[k]; 1324 for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */ 1325 col = *(columns + colorforcol[k] + l); 1326 btcol = bj + bi[col]; 1327 btval = ba + bi[col]; 1328 anz = bi[col+1] - bi[col]; 1329 for (j=0; j<anz; j++) { 1330 brow = btcol[j]; 1331 btval_den[brow] = btval[j]; 1332 } 1333 } 1334 btval_den += m; 1335 } 1336 PetscFunctionReturn(0); 1337 } 1338 1339 #undef __FUNCT__ 1340 #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ" 1341 PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp) 1342 { 1343 PetscErrorCode ierr; 1344 Mat_SeqAIJ *csp = (Mat_SeqAIJ*)Csp->data; 1345 PetscScalar *ca_den,*ca_den_ptr,*ca=csp->a; 1346 PetscInt k,l,m=Cden->rmap->n,ncolors=matcoloring->ncolors; 1347 PetscInt brows=matcoloring->brows,*den2sp=matcoloring->den2sp; 1348 PetscInt nrows,*row,*idx; 1349 PetscInt *rows=matcoloring->rows,*colorforrow=matcoloring->colorforrow; 1350 1351 PetscFunctionBegin; 1352 ierr = MatDenseGetArray(Cden,&ca_den);CHKERRQ(ierr); 1353 1354 if (brows > 0) { 1355 PetscInt *lstart,row_end,row_start; 1356 lstart = matcoloring->lstart; 1357 ierr = PetscMemzero(lstart,ncolors*sizeof(PetscInt));CHKERRQ(ierr); 1358 1359 row_end = brows; 1360 if (row_end > m) row_end = m; 1361 for (row_start=0; row_start<m; row_start+=brows) { /* loop over row blocks of Csp */ 1362 ca_den_ptr = ca_den; 1363 for (k=0; k<ncolors; k++) { /* loop over colors (columns of Cden) */ 1364 nrows = matcoloring->nrows[k]; 1365 row = rows + colorforrow[k]; 1366 idx = den2sp + colorforrow[k]; 1367 for (l=lstart[k]; l<nrows; l++) { 1368 if (row[l] >= row_end) { 1369 lstart[k] = l; 1370 break; 1371 } else { 1372 ca[idx[l]] = ca_den_ptr[row[l]]; 1373 } 1374 } 1375 ca_den_ptr += m; 1376 } 1377 row_end += brows; 1378 if (row_end > m) row_end = m; 1379 } 1380 } else { /* non-blocked impl: loop over columns of Csp - slow if Csp is large */ 1381 ca_den_ptr = ca_den; 1382 for (k=0; k<ncolors; k++) { 1383 nrows = matcoloring->nrows[k]; 1384 row = rows + colorforrow[k]; 1385 idx = den2sp + colorforrow[k]; 1386 for (l=0; l<nrows; l++) { 1387 ca[idx[l]] = ca_den_ptr[row[l]]; 1388 } 1389 ca_den_ptr += m; 1390 } 1391 } 1392 1393 ierr = MatDenseRestoreArray(Cden,&ca_den);CHKERRQ(ierr); 1394 #if defined(PETSC_USE_INFO) 1395 if (matcoloring->brows > 0) { 1396 ierr = PetscInfo1(Csp,"Loop over %D row blocks for den2sp\n",brows);CHKERRQ(ierr); 1397 } else { 1398 ierr = PetscInfo(Csp,"Loop over colors/columns of Cden, inefficient for large sparse matrix product \n");CHKERRQ(ierr); 1399 } 1400 #endif 1401 PetscFunctionReturn(0); 1402 } 1403 1404 #undef __FUNCT__ 1405 #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ" 1406 PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c) 1407 { 1408 PetscErrorCode ierr; 1409 PetscInt i,n,nrows,Nbs,j,k,m,ncols,col,cm; 1410 const PetscInt *is,*ci,*cj,*row_idx; 1411 PetscInt nis = iscoloring->n,*rowhit,bs = 1; 1412 IS *isa; 1413 Mat_SeqAIJ *csp = (Mat_SeqAIJ*)mat->data; 1414 PetscInt *colorforrow,*rows,*rows_i,*idxhit,*spidx,*den2sp,*den2sp_i; 1415 PetscInt *colorforcol,*columns,*columns_i,brows; 1416 PetscBool flg; 1417 1418 PetscFunctionBegin; 1419 ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 1420 1421 /* bs >1 is not being tested yet! */ 1422 Nbs = mat->cmap->N/bs; 1423 c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 1424 c->N = Nbs; 1425 c->m = c->M; 1426 c->rstart = 0; 1427 c->brows = 100; 1428 1429 c->ncolors = nis; 1430 ierr = PetscMalloc3(nis,&c->ncolumns,nis,&c->nrows,nis+1,&colorforrow);CHKERRQ(ierr); 1431 ierr = PetscMalloc((csp->nz+1)*sizeof(PetscInt),&rows);CHKERRQ(ierr); 1432 ierr = PetscMalloc((csp->nz+1)*sizeof(PetscInt),&den2sp);CHKERRQ(ierr); 1433 1434 brows = c->brows; 1435 ierr = PetscOptionsGetInt(NULL,"-matden2sp_brows",&brows,&flg);CHKERRQ(ierr); 1436 if (flg) c->brows = brows; 1437 if (brows > 0) { 1438 ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&c->lstart);CHKERRQ(ierr); 1439 } 1440 1441 colorforrow[0] = 0; 1442 rows_i = rows; 1443 den2sp_i = den2sp; 1444 1445 ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr); 1446 ierr = PetscMalloc((Nbs+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr); 1447 1448 colorforcol[0] = 0; 1449 columns_i = columns; 1450 1451 /* get column-wise storage of mat */ 1452 ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 1453 1454 cm = c->m; 1455 ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 1456 ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr); 1457 for (i=0; i<nis; i++) { /* loop over color */ 1458 ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 1459 ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 1460 1461 c->ncolumns[i] = n; 1462 if (n) { 1463 ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr); 1464 } 1465 colorforcol[i+1] = colorforcol[i] + n; 1466 columns_i += n; 1467 1468 /* fast, crude version requires O(N*N) work */ 1469 ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr); 1470 1471 for (j=0; j<n; j++) { /* loop over columns*/ 1472 col = is[j]; 1473 row_idx = cj + ci[col]; 1474 m = ci[col+1] - ci[col]; 1475 for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */ 1476 idxhit[*row_idx] = spidx[ci[col] + k]; 1477 rowhit[*row_idx++] = col + 1; 1478 } 1479 } 1480 /* count the number of hits */ 1481 nrows = 0; 1482 for (j=0; j<cm; j++) { 1483 if (rowhit[j]) nrows++; 1484 } 1485 c->nrows[i] = nrows; 1486 colorforrow[i+1] = colorforrow[i] + nrows; 1487 1488 nrows = 0; 1489 for (j=0; j<cm; j++) { /* loop over rows */ 1490 if (rowhit[j]) { 1491 rows_i[nrows] = j; 1492 den2sp_i[nrows] = idxhit[j]; 1493 nrows++; 1494 } 1495 } 1496 den2sp_i += nrows; 1497 1498 ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 1499 rows_i += nrows; 1500 } 1501 ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 1502 ierr = PetscFree(rowhit);CHKERRQ(ierr); 1503 ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 1504 if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]); 1505 1506 c->colorforrow = colorforrow; 1507 c->rows = rows; 1508 c->den2sp = den2sp; 1509 c->colorforcol = colorforcol; 1510 c->columns = columns; 1511 1512 ierr = PetscFree(idxhit);CHKERRQ(ierr); 1513 PetscFunctionReturn(0); 1514 } 1515