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