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