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