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