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