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