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