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