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