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