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 PetscReal afill; 481 PetscInt i,j,col,ndouble = 0; 482 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 483 PetscHeap h; 484 485 PetscFunctionBegin; 486 /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */ 487 /*---------------------------------------------------------------------------------------------*/ 488 /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 489 ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 490 ci[0] = 0; 491 492 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 493 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 494 current_space = free_space; 495 496 ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr); 497 ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr); 498 499 /* Determine ci and cj */ 500 for (i=0; i<am; i++) { 501 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 */ 502 const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 503 ci[i+1] = ci[i]; 504 /* Populate the min heap */ 505 for (j=0; j<anzi; j++) { 506 bb[j] = bi[acol[j]]; /* bb points at the start of the row */ 507 if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */ 508 ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr); 509 } 510 } 511 /* Pick off the min element, adding it to free space */ 512 ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 513 while (j >= 0) { 514 if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */ 515 ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),¤t_space);CHKERRQ(ierr); 516 ndouble++; 517 } 518 *(current_space->array++) = col; 519 current_space->local_used++; 520 current_space->local_remaining--; 521 ci[i+1]++; 522 523 /* stash if anything else remains in this row of B */ 524 if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);} 525 while (1) { /* pop and stash any other rows of B that also had an entry in this column */ 526 PetscInt j2,col2; 527 ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr); 528 if (col2 != col) break; 529 ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr); 530 if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);} 531 } 532 /* Put any stashed elements back into the min heap */ 533 ierr = PetscHeapUnstash(h);CHKERRQ(ierr); 534 ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 535 } 536 } 537 ierr = PetscFree(bb);CHKERRQ(ierr); 538 ierr = PetscHeapDestroy(&h);CHKERRQ(ierr); 539 540 /* Column indices are in the list of free space */ 541 /* Allocate space for cj, initialize cj, and */ 542 /* destroy list of free space and other temporary array(s) */ 543 ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr); 544 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 545 546 /* put together the new symbolic matrix */ 547 ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,PETSC_NULL,C);CHKERRQ(ierr); 548 (*C)->rmap->bs = A->rmap->bs; 549 (*C)->cmap->bs = B->cmap->bs; 550 551 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 552 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 553 c = (Mat_SeqAIJ *)((*C)->data); 554 c->free_a = PETSC_TRUE; 555 c->free_ij = PETSC_TRUE; 556 c->nonew = 0; 557 (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; 558 559 /* set MatInfo */ 560 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 561 if (afill < 1.0) afill = 1.0; 562 c->maxnz = ci[am]; 563 c->nz = ci[am]; 564 (*C)->info.mallocs = ndouble; 565 (*C)->info.fill_ratio_given = fill; 566 (*C)->info.fill_ratio_needed = afill; 567 568 #if defined(PETSC_USE_INFO) 569 if (ci[am]) { 570 ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 571 ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 572 } else { 573 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 574 } 575 #endif 576 PetscFunctionReturn(0); 577 } 578 579 #undef __FUNCT__ 580 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap" 581 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat *C) 582 { 583 PetscErrorCode ierr; 584 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 585 const PetscInt *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j; 586 PetscInt *ci,*cj,*bb; 587 PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 588 PetscReal afill; 589 PetscInt i,j,col,ndouble = 0; 590 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 591 PetscHeap h; 592 PetscBT bt; 593 594 PetscFunctionBegin; 595 /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */ 596 /*---------------------------------------------------------------------------------------------*/ 597 /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 598 ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 599 ci[0] = 0; 600 601 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 602 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 603 current_space = free_space; 604 605 ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr); 606 ierr = PetscMalloc(a->rmax*sizeof(PetscInt),&bb);CHKERRQ(ierr); 607 ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr); 608 609 /* Determine ci and cj */ 610 for (i=0; i<am; i++) { 611 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 */ 612 const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 613 const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */ 614 ci[i+1] = ci[i]; 615 /* Populate the min heap */ 616 for (j=0; j<anzi; j++) { 617 PetscInt brow = acol[j]; 618 for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) { 619 PetscInt bcol = bj[bb[j]]; 620 if (!PetscBTLookupSet(bt,bcol)) { /* new entry */ 621 ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr); 622 bb[j]++; 623 break; 624 } 625 } 626 } 627 /* Pick off the min element, adding it to free space */ 628 ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 629 while (j >= 0) { 630 if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */ 631 fptr = PETSC_NULL; /* need PetscBTMemzero */ 632 ierr = PetscFreeSpaceGet(PetscMin(2*current_space->total_array_size,16 << 20),¤t_space);CHKERRQ(ierr); 633 ndouble++; 634 } 635 *(current_space->array++) = col; 636 current_space->local_used++; 637 current_space->local_remaining--; 638 ci[i+1]++; 639 640 /* stash if anything else remains in this row of B */ 641 for ( ; bb[j] < bi[acol[j]+1]; bb[j]++) { 642 PetscInt bcol = bj[bb[j]]; 643 if (!PetscBTLookupSet(bt,bcol)) { /* new entry */ 644 ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr); 645 bb[j]++; 646 break; 647 } 648 } 649 ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 650 } 651 if (fptr) { /* Clear the bits for this row */ 652 for ( ; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);} 653 } else { /* We reallocated so we don't remember (easily) how to clear only the bits we changed */ 654 ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr); 655 } 656 } 657 ierr = PetscFree(bb);CHKERRQ(ierr); 658 ierr = PetscHeapDestroy(&h);CHKERRQ(ierr); 659 ierr = PetscBTDestroy(&bt);CHKERRQ(ierr); 660 661 /* Column indices are in the list of free space */ 662 /* Allocate space for cj, initialize cj, and */ 663 /* destroy list of free space and other temporary array(s) */ 664 ierr = PetscMalloc(ci[am]*sizeof(PetscInt),&cj);CHKERRQ(ierr); 665 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 666 667 /* put together the new symbolic matrix */ 668 ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,PETSC_NULL,C);CHKERRQ(ierr); 669 (*C)->rmap->bs = A->rmap->bs; 670 (*C)->cmap->bs = B->cmap->bs; 671 672 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 673 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 674 c = (Mat_SeqAIJ *)((*C)->data); 675 c->free_a = PETSC_TRUE; 676 c->free_ij = PETSC_TRUE; 677 c->nonew = 0; 678 (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; 679 680 /* set MatInfo */ 681 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 682 if (afill < 1.0) afill = 1.0; 683 c->maxnz = ci[am]; 684 c->nz = ci[am]; 685 (*C)->info.mallocs = ndouble; 686 (*C)->info.fill_ratio_given = fill; 687 (*C)->info.fill_ratio_needed = afill; 688 689 #if defined(PETSC_USE_INFO) 690 if (ci[am]) { 691 ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",ndouble,fill,afill);CHKERRQ(ierr); 692 ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 693 } else { 694 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 695 } 696 #endif 697 PetscFunctionReturn(0); 698 } 699 700 /* This routine is not used. Should be removed! */ 701 #undef __FUNCT__ 702 #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ" 703 PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 704 { 705 PetscErrorCode ierr; 706 707 PetscFunctionBegin; 708 if (scall == MAT_INITIAL_MATRIX){ 709 ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 710 } 711 ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 712 PetscFunctionReturn(0); 713 } 714 715 #undef __FUNCT__ 716 #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult" 717 PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr) 718 { 719 PetscErrorCode ierr; 720 Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr; 721 722 PetscFunctionBegin; 723 ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr); 724 ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr); 725 ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr); 726 ierr = PetscFree(multtrans);CHKERRQ(ierr); 727 PetscFunctionReturn(0); 728 } 729 730 #undef __FUNCT__ 731 #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans" 732 PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A) 733 { 734 PetscErrorCode ierr; 735 PetscContainer container; 736 Mat_MatMatTransMult *multtrans=PETSC_NULL; 737 738 PetscFunctionBegin; 739 ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr); 740 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 741 ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr); 742 A->ops->destroy = multtrans->destroy; 743 if (A->ops->destroy) { 744 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 745 } 746 ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr); 747 PetscFunctionReturn(0); 748 } 749 750 #undef __FUNCT__ 751 #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ" 752 PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 753 { 754 PetscErrorCode ierr; 755 Mat Bt; 756 PetscInt *bti,*btj; 757 Mat_MatMatTransMult *multtrans; 758 PetscContainer container; 759 PetscLogDouble t0,tf,etime2=0.0; 760 761 PetscFunctionBegin; 762 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 763 /* create symbolic Bt */ 764 ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 765 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr); 766 Bt->rmap->bs = A->cmap->bs; 767 Bt->cmap->bs = B->cmap->bs; 768 769 /* get symbolic C=A*Bt */ 770 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr); 771 772 /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */ 773 ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr); 774 775 /* attach the supporting struct to C */ 776 ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 777 ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr); 778 ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr); 779 ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr); 780 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 781 782 multtrans->usecoloring = PETSC_FALSE; 783 multtrans->destroy = (*C)->ops->destroy; 784 (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans; 785 786 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 787 etime2 += tf - t0; 788 789 ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr); 790 if (multtrans->usecoloring){ 791 /* Create MatTransposeColoring from symbolic C=A*B^T */ 792 MatTransposeColoring matcoloring; 793 ISColoring iscoloring; 794 Mat Bt_dense,C_dense; 795 PetscLogDouble etime0=0.0,etime01=0.0,etime1=0.0; 796 797 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 798 ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 799 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 800 etime0 += tf - t0; 801 802 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 803 ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 804 multtrans->matcoloring = matcoloring; 805 ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 806 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 807 etime01 += tf - t0; 808 809 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 810 /* Create Bt_dense and C_dense = A*Bt_dense */ 811 ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr); 812 ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 813 ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr); 814 ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr); 815 Bt_dense->assembled = PETSC_TRUE; 816 multtrans->Bt_den = Bt_dense; 817 818 ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr); 819 ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 820 ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr); 821 ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr); 822 Bt_dense->assembled = PETSC_TRUE; 823 multtrans->ABt_den = C_dense; 824 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 825 etime1 += tf - t0; 826 827 #if defined(PETSC_USE_INFO) 828 { 829 Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data; 830 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)); 831 ierr = PetscInfo5(*C,"Sym = GetColor %g + ColorCreate %g + MatDenseCreate %g + non-colorSym %g = %g\n",etime0,etime01,etime1,etime2,etime0+etime01+etime1+etime2); 832 } 833 #endif 834 } 835 /* clean up */ 836 ierr = MatDestroy(&Bt);CHKERRQ(ierr); 837 ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 838 839 840 841 #if defined(INEFFICIENT_ALGORITHM) 842 /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */ 843 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 844 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 845 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol; 846 PetscInt am=A->rmap->N,bm=B->rmap->N; 847 PetscInt i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1]; 848 MatScalar *ca; 849 PetscBT lnkbt; 850 PetscReal afill; 851 852 /* Allocate row pointer array ci */ 853 ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 854 ci[0] = 0; 855 856 /* Create and initialize a linked list for C columns */ 857 nlnk = bm+1; 858 ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 859 860 /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */ 861 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 862 current_space = free_space; 863 864 /* Determine symbolic info for each row of the product A*B^T: */ 865 for (i=0; i<am; i++) { 866 anzi = ai[i+1] - ai[i]; 867 cnzi = 0; 868 acol = aj + ai[i]; 869 for (j=0; j<bm; j++){ 870 bnzj = bi[j+1] - bi[j]; 871 bcol= bj + bi[j]; 872 /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 873 ka = 0; kb = 0; 874 while (ka < anzi && kb < bnzj){ 875 while (acol[ka] < bcol[kb] && ka < anzi) ka++; 876 if (ka == anzi) break; 877 while (acol[ka] > bcol[kb] && kb < bnzj) kb++; 878 if (kb == bnzj) break; 879 if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */ 880 index[0] = j; 881 ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 882 cnzi++; 883 break; 884 } 885 } 886 } 887 888 /* If free space is not available, make more free space */ 889 /* Double the amount of total space in the list */ 890 if (current_space->local_remaining<cnzi) { 891 ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 892 nspacedouble++; 893 } 894 895 /* Copy data into free space, then initialize lnk */ 896 ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 897 current_space->array += cnzi; 898 current_space->local_used += cnzi; 899 current_space->local_remaining -= cnzi; 900 901 ci[i+1] = ci[i] + cnzi; 902 } 903 904 905 /* Column indices are in the list of free space. 906 Allocate array cj, copy column indices to cj, and destroy list of free space */ 907 ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 908 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 909 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 910 911 /* Allocate space for ca */ 912 ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 913 ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 914 915 /* put together the new symbolic matrix */ 916 ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr); 917 (*C)->rmap->bs = A->cmap->bs; 918 (*C)->cmap->bs = B->cmap->bs; 919 920 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 921 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 922 c = (Mat_SeqAIJ *)((*C)->data); 923 c->free_a = PETSC_TRUE; 924 c->free_ij = PETSC_TRUE; 925 c->nonew = 0; 926 927 /* set MatInfo */ 928 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 929 if (afill < 1.0) afill = 1.0; 930 c->maxnz = ci[am]; 931 c->nz = ci[am]; 932 (*C)->info.mallocs = nspacedouble; 933 (*C)->info.fill_ratio_given = fill; 934 (*C)->info.fill_ratio_needed = afill; 935 936 #if defined(PETSC_USE_INFO) 937 if (ci[am]) { 938 ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 939 ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 940 } else { 941 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 942 } 943 #endif 944 #endif 945 PetscFunctionReturn(0); 946 } 947 948 /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */ 949 #undef __FUNCT__ 950 #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ" 951 PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 952 { 953 PetscErrorCode ierr; 954 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 955 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow; 956 PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol; 957 PetscLogDouble flops=0.0; 958 MatScalar *aa=a->a,*aval,*ba=b->a,*bval,*ca,*cval; 959 Mat_MatMatTransMult *multtrans; 960 PetscContainer container; 961 #if defined(USE_ARRAY) 962 MatScalar *spdot; 963 #endif 964 965 PetscFunctionBegin; 966 ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr); 967 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 968 ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr); 969 if (multtrans->usecoloring){ 970 MatTransposeColoring matcoloring = multtrans->matcoloring; 971 Mat Bt_dense; 972 PetscInt m,n; 973 PetscLogDouble t0,tf,etime0=0.0,etime1=0.0,etime2=0.0; 974 Mat C_dense = multtrans->ABt_den; 975 976 Bt_dense = multtrans->Bt_den; 977 ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr); 978 979 /* Get Bt_dense by Apply MatTransposeColoring to B */ 980 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 981 ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr); 982 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 983 etime0 += tf - t0; 984 985 /* C_dense = A*Bt_dense */ 986 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 987 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr); 988 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 989 etime2 += tf - t0; 990 991 /* Recover C from C_dense */ 992 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 993 ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr); 994 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 995 etime1 += tf - t0; 996 #if defined(PETSC_USE_INFO) 997 ierr = PetscInfo4(C,"Num = ColoringApply: %g %g + Mult_sp_dense %g = %g\n",etime0,etime1,etime2,etime0+etime1+etime2); 998 #endif 999 PetscFunctionReturn(0); 1000 } 1001 1002 #if defined(USE_ARRAY) 1003 /* allocate an array for implementing sparse inner-product */ 1004 ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr); 1005 ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr); 1006 #endif 1007 1008 /* clear old values in C */ 1009 if (!c->a){ 1010 ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 1011 c->a = ca; 1012 c->free_a = PETSC_TRUE; 1013 } else { 1014 ca = c->a; 1015 } 1016 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 1017 1018 for (i=0; i<cm; i++) { 1019 anzi = ai[i+1] - ai[i]; 1020 acol = aj + ai[i]; 1021 aval = aa + ai[i]; 1022 cnzi = ci[i+1] - ci[i]; 1023 ccol = cj + ci[i]; 1024 cval = ca + ci[i]; 1025 for (j=0; j<cnzi; j++){ 1026 brow = ccol[j]; 1027 bnzj = bi[brow+1] - bi[brow]; 1028 bcol = bj + bi[brow]; 1029 bval = ba + bi[brow]; 1030 1031 /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 1032 #if defined(USE_ARRAY) 1033 /* put ba to spdot array */ 1034 for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb]; 1035 /* c(i,j)=A[i,:]*B[j,:]^T */ 1036 for (nexta=0; nexta<anzi; nexta++){ 1037 cval[j] += spdot[acol[nexta]]*aval[nexta]; 1038 } 1039 /* zero spdot array */ 1040 for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0; 1041 #else 1042 nexta = 0; nextb = 0; 1043 while (nexta<anzi && nextb<bnzj){ 1044 while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++; 1045 if (nexta == anzi) break; 1046 while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++; 1047 if (nextb == bnzj) break; 1048 if (acol[nexta] == bcol[nextb]){ 1049 cval[j] += aval[nexta]*bval[nextb]; 1050 nexta++; nextb++; 1051 flops += 2; 1052 } 1053 } 1054 #endif 1055 } 1056 } 1057 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1058 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1059 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 1060 #if defined(USE_ARRAY) 1061 ierr = PetscFree(spdot); 1062 #endif 1063 PetscFunctionReturn(0); 1064 } 1065 1066 #undef __FUNCT__ 1067 #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ" 1068 PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) { 1069 PetscErrorCode ierr; 1070 1071 PetscFunctionBegin; 1072 if (scall == MAT_INITIAL_MATRIX){ 1073 ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 1074 } 1075 ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 1076 PetscFunctionReturn(0); 1077 } 1078 1079 #undef __FUNCT__ 1080 #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ" 1081 PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 1082 { 1083 PetscErrorCode ierr; 1084 Mat At; 1085 PetscInt *ati,*atj; 1086 1087 PetscFunctionBegin; 1088 /* create symbolic At */ 1089 ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 1090 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr); 1091 At->rmap->bs = A->cmap->bs; 1092 At->cmap->bs = B->cmap->bs; 1093 1094 /* get symbolic C=At*B */ 1095 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 1096 1097 /* clean up */ 1098 ierr = MatDestroy(&At);CHKERRQ(ierr); 1099 ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 1100 PetscFunctionReturn(0); 1101 } 1102 1103 #undef __FUNCT__ 1104 #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ" 1105 PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 1106 { 1107 PetscErrorCode ierr; 1108 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 1109 PetscInt am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 1110 PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 1111 PetscLogDouble flops=0.0; 1112 MatScalar *aa=a->a,*ba,*ca,*caj; 1113 1114 PetscFunctionBegin; 1115 if (!c->a){ 1116 ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 1117 c->a = ca; 1118 c->free_a = PETSC_TRUE; 1119 } else { 1120 ca = c->a; 1121 } 1122 /* clear old values in C */ 1123 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 1124 1125 /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 1126 for (i=0;i<am;i++) { 1127 bj = b->j + bi[i]; 1128 ba = b->a + bi[i]; 1129 bnzi = bi[i+1] - bi[i]; 1130 anzi = ai[i+1] - ai[i]; 1131 for (j=0; j<anzi; j++) { 1132 nextb = 0; 1133 crow = *aj++; 1134 cjj = cj + ci[crow]; 1135 caj = ca + ci[crow]; 1136 /* perform sparse axpy operation. Note cjj includes bj. */ 1137 for (k=0; nextb<bnzi; k++) { 1138 if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 1139 caj[k] += (*aa)*(*(ba+nextb)); 1140 nextb++; 1141 } 1142 } 1143 flops += 2*bnzi; 1144 aa++; 1145 } 1146 } 1147 1148 /* Assemble the final matrix and clean up */ 1149 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1150 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1151 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 1152 PetscFunctionReturn(0); 1153 } 1154 1155 EXTERN_C_BEGIN 1156 #undef __FUNCT__ 1157 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense" 1158 PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1159 { 1160 PetscErrorCode ierr; 1161 1162 PetscFunctionBegin; 1163 if (scall == MAT_INITIAL_MATRIX){ 1164 ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1165 } 1166 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 1167 PetscFunctionReturn(0); 1168 } 1169 EXTERN_C_END 1170 1171 #undef __FUNCT__ 1172 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense" 1173 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1174 { 1175 PetscErrorCode ierr; 1176 1177 PetscFunctionBegin; 1178 ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr); 1179 (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqDense; 1180 PetscFunctionReturn(0); 1181 } 1182 1183 #undef __FUNCT__ 1184 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense" 1185 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 1186 { 1187 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 1188 PetscErrorCode ierr; 1189 PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 1190 MatScalar *aa; 1191 PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n; 1192 PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam; 1193 1194 PetscFunctionBegin; 1195 if (!cm || !cn) PetscFunctionReturn(0); 1196 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); 1197 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); 1198 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); 1199 ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 1200 ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr); 1201 b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 1202 for (col=0; col<cn-4; col += 4){ /* over columns of C */ 1203 colam = col*am; 1204 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1205 r1 = r2 = r3 = r4 = 0.0; 1206 n = a->i[i+1] - a->i[i]; 1207 aj = a->j + a->i[i]; 1208 aa = a->a + a->i[i]; 1209 for (j=0; j<n; j++) { 1210 r1 += (*aa)*b1[*aj]; 1211 r2 += (*aa)*b2[*aj]; 1212 r3 += (*aa)*b3[*aj]; 1213 r4 += (*aa++)*b4[*aj++]; 1214 } 1215 c[colam + i] = r1; 1216 c[colam + am + i] = r2; 1217 c[colam + am2 + i] = r3; 1218 c[colam + am3 + i] = r4; 1219 } 1220 b1 += bm4; 1221 b2 += bm4; 1222 b3 += bm4; 1223 b4 += bm4; 1224 } 1225 for (;col<cn; col++){ /* over extra columns of C */ 1226 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1227 r1 = 0.0; 1228 n = a->i[i+1] - a->i[i]; 1229 aj = a->j + a->i[i]; 1230 aa = a->a + a->i[i]; 1231 1232 for (j=0; j<n; j++) { 1233 r1 += (*aa++)*b1[*aj++]; 1234 } 1235 c[col*am + i] = r1; 1236 } 1237 b1 += bm; 1238 } 1239 ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 1240 ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 1241 ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr); 1242 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1243 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1244 PetscFunctionReturn(0); 1245 } 1246 1247 /* 1248 Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 1249 */ 1250 #undef __FUNCT__ 1251 #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense" 1252 PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 1253 { 1254 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 1255 PetscErrorCode ierr; 1256 PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 1257 MatScalar *aa; 1258 PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n,*ii,arm; 1259 PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 1260 1261 PetscFunctionBegin; 1262 if (!cm || !cn) PetscFunctionReturn(0); 1263 ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 1264 ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr); 1265 b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 1266 1267 if (a->compressedrow.use){ /* use compressed row format */ 1268 for (col=0; col<cn-4; col += 4){ /* over columns of C */ 1269 colam = col*am; 1270 arm = a->compressedrow.nrows; 1271 ii = a->compressedrow.i; 1272 ridx = a->compressedrow.rindex; 1273 for (i=0; i<arm; i++) { /* over rows of C in those columns */ 1274 r1 = r2 = r3 = r4 = 0.0; 1275 n = ii[i+1] - ii[i]; 1276 aj = a->j + ii[i]; 1277 aa = a->a + ii[i]; 1278 for (j=0; j<n; j++) { 1279 r1 += (*aa)*b1[*aj]; 1280 r2 += (*aa)*b2[*aj]; 1281 r3 += (*aa)*b3[*aj]; 1282 r4 += (*aa++)*b4[*aj++]; 1283 } 1284 c[colam + ridx[i]] += r1; 1285 c[colam + am + ridx[i]] += r2; 1286 c[colam + am2 + ridx[i]] += r3; 1287 c[colam + am3 + ridx[i]] += r4; 1288 } 1289 b1 += bm4; 1290 b2 += bm4; 1291 b3 += bm4; 1292 b4 += bm4; 1293 } 1294 for (;col<cn; col++){ /* over extra columns of C */ 1295 colam = col*am; 1296 arm = a->compressedrow.nrows; 1297 ii = a->compressedrow.i; 1298 ridx = a->compressedrow.rindex; 1299 for (i=0; i<arm; i++) { /* over rows of C in those columns */ 1300 r1 = 0.0; 1301 n = ii[i+1] - ii[i]; 1302 aj = a->j + ii[i]; 1303 aa = a->a + ii[i]; 1304 1305 for (j=0; j<n; j++) { 1306 r1 += (*aa++)*b1[*aj++]; 1307 } 1308 c[col*am + ridx[i]] += r1; 1309 } 1310 b1 += bm; 1311 } 1312 } else { 1313 for (col=0; col<cn-4; col += 4){ /* over columns of C */ 1314 colam = col*am; 1315 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1316 r1 = r2 = r3 = r4 = 0.0; 1317 n = a->i[i+1] - a->i[i]; 1318 aj = a->j + a->i[i]; 1319 aa = a->a + a->i[i]; 1320 for (j=0; j<n; j++) { 1321 r1 += (*aa)*b1[*aj]; 1322 r2 += (*aa)*b2[*aj]; 1323 r3 += (*aa)*b3[*aj]; 1324 r4 += (*aa++)*b4[*aj++]; 1325 } 1326 c[colam + i] += r1; 1327 c[colam + am + i] += r2; 1328 c[colam + am2 + i] += r3; 1329 c[colam + am3 + i] += r4; 1330 } 1331 b1 += bm4; 1332 b2 += bm4; 1333 b3 += bm4; 1334 b4 += bm4; 1335 } 1336 for (;col<cn; col++){ /* over extra columns of C */ 1337 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1338 r1 = 0.0; 1339 n = a->i[i+1] - a->i[i]; 1340 aj = a->j + a->i[i]; 1341 aa = a->a + a->i[i]; 1342 1343 for (j=0; j<n; j++) { 1344 r1 += (*aa++)*b1[*aj++]; 1345 } 1346 c[col*am + i] += r1; 1347 } 1348 b1 += bm; 1349 } 1350 } 1351 ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr); 1352 ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 1353 ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr); 1354 PetscFunctionReturn(0); 1355 } 1356 1357 #undef __FUNCT__ 1358 #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ" 1359 PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense) 1360 { 1361 PetscErrorCode ierr; 1362 Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 1363 Mat_SeqDense *btdense = (Mat_SeqDense*)Btdense->data; 1364 PetscInt *bi=b->i,*bj=b->j; 1365 PetscInt m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns; 1366 MatScalar *btval,*btval_den,*ba=b->a; 1367 PetscInt *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors; 1368 1369 PetscFunctionBegin; 1370 btval_den=btdense->v; 1371 ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr); 1372 for (k=0; k<ncolors; k++) { 1373 ncolumns = coloring->ncolumns[k]; 1374 for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */ 1375 col = *(columns + colorforcol[k] + l); 1376 btcol = bj + bi[col]; 1377 btval = ba + bi[col]; 1378 anz = bi[col+1] - bi[col]; 1379 for (j=0; j<anz; j++){ 1380 brow = btcol[j]; 1381 btval_den[brow] = btval[j]; 1382 } 1383 } 1384 btval_den += m; 1385 } 1386 PetscFunctionReturn(0); 1387 } 1388 1389 #undef __FUNCT__ 1390 #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ" 1391 PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp) 1392 { 1393 PetscErrorCode ierr; 1394 Mat_SeqAIJ *csp = (Mat_SeqAIJ*)Csp->data; 1395 PetscInt k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows; 1396 PetscScalar *ca_den,*cp_den,*ca=csp->a; 1397 PetscInt *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow; 1398 1399 PetscFunctionBegin; 1400 ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr); 1401 ierr = MatDenseGetArray(Cden,&ca_den);CHKERRQ(ierr); 1402 cp_den = ca_den; 1403 for (k=0; k<ncolors; k++) { 1404 nrows = matcoloring->nrows[k]; 1405 row = rows + colorforrow[k]; 1406 idx = spidx + colorforrow[k]; 1407 for (l=0; l<nrows; l++){ 1408 ca[idx[l]] = cp_den[row[l]]; 1409 } 1410 cp_den += m; 1411 } 1412 ierr = MatDenseRestoreArray(Cden,&ca_den);CHKERRQ(ierr); 1413 PetscFunctionReturn(0); 1414 } 1415 1416 /* 1417 MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from 1418 MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output 1419 spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ(). 1420 */ 1421 #undef __FUNCT__ 1422 #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color" 1423 PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1424 { 1425 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1426 PetscErrorCode ierr; 1427 PetscInt i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n; 1428 PetscInt nz = a->i[m],row,*jj,mr,col; 1429 PetscInt *cspidx; 1430 1431 PetscFunctionBegin; 1432 *nn = n; 1433 if (!ia) PetscFunctionReturn(0); 1434 if (symmetric) { 1435 SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric"); 1436 ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 1437 } else { 1438 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr); 1439 ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 1440 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr); 1441 ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr); 1442 ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr); 1443 jj = a->j; 1444 for (i=0; i<nz; i++) { 1445 collengths[jj[i]]++; 1446 } 1447 cia[0] = oshift; 1448 for (i=0; i<n; i++) { 1449 cia[i+1] = cia[i] + collengths[i]; 1450 } 1451 ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 1452 jj = a->j; 1453 for (row=0; row<m; row++) { 1454 mr = a->i[row+1] - a->i[row]; 1455 for (i=0; i<mr; i++) { 1456 col = *jj++; 1457 cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */ 1458 cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 1459 } 1460 } 1461 ierr = PetscFree(collengths);CHKERRQ(ierr); 1462 *ia = cia; *ja = cja; 1463 *spidx = cspidx; 1464 } 1465 PetscFunctionReturn(0); 1466 } 1467 1468 #undef __FUNCT__ 1469 #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color" 1470 PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1471 { 1472 PetscErrorCode ierr; 1473 1474 PetscFunctionBegin; 1475 if (!ia) PetscFunctionReturn(0); 1476 1477 ierr = PetscFree(*ia);CHKERRQ(ierr); 1478 ierr = PetscFree(*ja);CHKERRQ(ierr); 1479 ierr = PetscFree(*spidx);CHKERRQ(ierr); 1480 PetscFunctionReturn(0); 1481 } 1482 1483 #undef __FUNCT__ 1484 #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ" 1485 PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c) 1486 { 1487 PetscErrorCode ierr; 1488 PetscInt i,n,nrows,N,j,k,m,*row_idx,*ci,*cj,ncols,col,cm; 1489 const PetscInt *is; 1490 PetscInt nis = iscoloring->n,*rowhit,bs = 1; 1491 IS *isa; 1492 PetscBool flg1,flg2; 1493 Mat_SeqAIJ *csp = (Mat_SeqAIJ*)mat->data; 1494 PetscInt *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx; 1495 PetscInt *colorforcol,*columns,*columns_i; 1496 1497 PetscFunctionBegin; 1498 ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 1499 1500 /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */ 1501 ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr); 1502 ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr); 1503 if (flg1 || flg2) { 1504 ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 1505 } 1506 1507 N = mat->cmap->N/bs; 1508 c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 1509 c->N = mat->cmap->N/bs; 1510 c->m = mat->rmap->N/bs; 1511 c->rstart = 0; 1512 1513 c->ncolors = nis; 1514 ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 1515 ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 1516 ierr = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr); 1517 ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr); 1518 colorforrow[0] = 0; 1519 rows_i = rows; 1520 columnsforspidx_i = columnsforspidx; 1521 1522 ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr); 1523 ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr); 1524 colorforcol[0] = 0; 1525 columns_i = columns; 1526 1527 ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,PETSC_NULL);CHKERRQ(ierr); /* column-wise storage of mat */ 1528 1529 cm = c->m; 1530 ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 1531 ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr); 1532 for (i=0; i<nis; i++) { 1533 ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 1534 ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 1535 c->ncolumns[i] = n; 1536 if (n) { 1537 ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr); 1538 } 1539 colorforcol[i+1] = colorforcol[i] + n; 1540 columns_i += n; 1541 1542 /* fast, crude version requires O(N*N) work */ 1543 ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr); 1544 1545 /* loop over columns*/ 1546 for (j=0; j<n; j++) { 1547 col = is[j]; 1548 row_idx = cj + ci[col]; 1549 m = ci[col+1] - ci[col]; 1550 /* loop over columns marking them in rowhit */ 1551 for (k=0; k<m; k++) { 1552 idxhit[*row_idx] = spidx[ci[col] + k]; 1553 rowhit[*row_idx++] = col + 1; 1554 } 1555 } 1556 /* count the number of hits */ 1557 nrows = 0; 1558 for (j=0; j<cm; j++) { 1559 if (rowhit[j]) nrows++; 1560 } 1561 c->nrows[i] = nrows; 1562 colorforrow[i+1] = colorforrow[i] + nrows; 1563 1564 nrows = 0; 1565 for (j=0; j<cm; j++) { 1566 if (rowhit[j]) { 1567 rows_i[nrows] = j; 1568 columnsforspidx_i[nrows] = idxhit[j]; 1569 nrows++; 1570 } 1571 } 1572 ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 1573 rows_i += nrows; columnsforspidx_i += nrows; 1574 } 1575 ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,PETSC_NULL);CHKERRQ(ierr); 1576 ierr = PetscFree(rowhit);CHKERRQ(ierr); 1577 ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 1578 #if defined(PETSC_USE_DEBUG) 1579 if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]); 1580 #endif 1581 1582 c->colorforrow = colorforrow; 1583 c->rows = rows; 1584 c->columnsforspidx = columnsforspidx; 1585 c->colorforcol = colorforcol; 1586 c->columns = columns; 1587 ierr = PetscFree(idxhit);CHKERRQ(ierr); 1588 PetscFunctionReturn(0); 1589 } 1590