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