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__ "PetscContainerDestroy_Mat_MatMatTransMult" 830 PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr) 831 { 832 PetscErrorCode ierr; 833 Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr; 834 835 PetscFunctionBegin; 836 ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr); 837 ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr); 838 ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr); 839 ierr = PetscFree(multtrans);CHKERRQ(ierr); 840 PetscFunctionReturn(0); 841 } 842 843 #undef __FUNCT__ 844 #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans" 845 PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A) 846 { 847 PetscErrorCode ierr; 848 PetscContainer container; 849 Mat_MatMatTransMult *multtrans=NULL; 850 851 PetscFunctionBegin; 852 ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject*)&container);CHKERRQ(ierr); 853 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 854 ierr = PetscContainerGetPointer(container,(void**)&multtrans);CHKERRQ(ierr); 855 856 A->ops->destroy = multtrans->destroy; 857 if (A->ops->destroy) { 858 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 859 } 860 ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr); 861 PetscFunctionReturn(0); 862 } 863 864 #undef __FUNCT__ 865 #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ" 866 PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 867 { 868 PetscErrorCode ierr; 869 Mat Bt; 870 PetscInt *bti,*btj; 871 Mat_MatMatTransMult *multtrans; 872 PetscContainer container; 873 874 PetscFunctionBegin; 875 /* create symbolic Bt */ 876 ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 877 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);CHKERRQ(ierr); 878 879 Bt->rmap->bs = A->cmap->bs; 880 Bt->cmap->bs = B->cmap->bs; 881 882 /* get symbolic C=A*Bt */ 883 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr); 884 885 /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */ 886 ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr); 887 888 /* attach the supporting struct to C */ 889 ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 890 ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr); 891 ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr); 892 ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr); 893 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 894 895 multtrans->usecoloring = PETSC_FALSE; 896 multtrans->destroy = (*C)->ops->destroy; 897 (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans; 898 899 ierr = PetscOptionsGetBool(NULL,"-matmattransmult_color",&multtrans->usecoloring,NULL);CHKERRQ(ierr); 900 if (multtrans->usecoloring) { 901 /* Create MatTransposeColoring from symbolic C=A*B^T */ 902 MatTransposeColoring matcoloring; 903 ISColoring iscoloring; 904 Mat Bt_dense,C_dense; 905 Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data; 906 /* inode causes memory problem, don't know why */ 907 if (c->inode.use) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MAT_USE_INODES is not supported. Use '-mat_no_inode'"); 908 909 ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 910 ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 911 912 multtrans->matcoloring = matcoloring; 913 914 ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 915 916 /* Create Bt_dense and C_dense = A*Bt_dense */ 917 ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr); 918 ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 919 ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr); 920 ierr = MatSeqDenseSetPreallocation(Bt_dense,NULL);CHKERRQ(ierr); 921 922 Bt_dense->assembled = PETSC_TRUE; 923 multtrans->Bt_den = Bt_dense; 924 925 ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr); 926 ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 927 ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr); 928 ierr = MatSeqDenseSetPreallocation(C_dense,NULL);CHKERRQ(ierr); 929 930 Bt_dense->assembled = PETSC_TRUE; 931 multtrans->ABt_den = C_dense; 932 933 #if defined(PETSC_USE_INFO) 934 { 935 Mat_SeqAIJ *c = (Mat_SeqAIJ*)(*C)->data; 936 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); 937 } 938 #endif 939 } 940 /* clean up */ 941 ierr = MatDestroy(&Bt);CHKERRQ(ierr); 942 ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 943 PetscFunctionReturn(0); 944 } 945 946 /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */ 947 #undef __FUNCT__ 948 #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ" 949 PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 950 { 951 PetscErrorCode ierr; 952 Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 953 PetscInt *ai =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow; 954 PetscInt cm =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol; 955 PetscLogDouble flops=0.0; 956 MatScalar *aa =a->a,*aval,*ba=b->a,*bval,*ca,*cval; 957 Mat_MatMatTransMult *multtrans; 958 PetscContainer container; 959 #if defined(USE_ARRAY) 960 MatScalar *spdot; 961 #endif 962 963 PetscFunctionBegin; 964 /* clear old values in C */ 965 if (!c->a) { 966 ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 967 c->a = ca; 968 c->free_a = PETSC_TRUE; 969 } else { 970 ca = c->a; 971 } 972 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 973 974 ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject*)&container);CHKERRQ(ierr); 975 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 976 ierr = PetscContainerGetPointer(container,(void**)&multtrans);CHKERRQ(ierr); 977 if (multtrans->usecoloring) { 978 MatTransposeColoring matcoloring = multtrans->matcoloring; 979 Mat Bt_dense; 980 PetscInt m,n; 981 Mat C_dense = multtrans->ABt_den; 982 PetscLogDouble t0,t1,t2,t3,Bt_den,C_den,C_sp; 983 984 /* Get Bt_dense by Apply MatTransposeColoring to B */ 985 ierr = PetscTime(&t0);CHKERRQ(ierr); 986 Bt_dense = multtrans->Bt_den; 987 ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr); 988 ierr = PetscTime(&t1);CHKERRQ(ierr); 989 Bt_den = t1 - t0; 990 991 /* C_dense = A*Bt_dense */ 992 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr); 993 ierr = PetscTime(&t2);CHKERRQ(ierr); 994 C_den = t2 - t1; 995 996 /* Recover C from C_dense */ 997 ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr); 998 ierr = PetscTime(&t3);CHKERRQ(ierr); 999 C_sp = t3 - t2; 1000 #if defined(PETSC_USE_INFO) 1001 ierr = PetscInfo4(C,"Coloring A*B^T = Bt_den %g + C_den %g + C_sp %g = %g;\n",Bt_den,C_den,C_sp,Bt_den+C_den+C_sp);CHKERRQ(ierr); 1002 ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr); 1003 ierr = PetscInfo4(C,"Bt_den: %d x %d, B: %d x %d\n",m,n,m,C->cmap->n);CHKERRQ(ierr); 1004 #endif 1005 PetscFunctionReturn(0); 1006 } 1007 1008 #if defined(USE_ARRAY) 1009 /* allocate an array for implementing sparse inner-product */ 1010 ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr); 1011 ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr); 1012 #endif 1013 1014 for (i=0; i<cm; i++) { 1015 anzi = ai[i+1] - ai[i]; 1016 acol = aj + ai[i]; 1017 aval = aa + ai[i]; 1018 cnzi = ci[i+1] - ci[i]; 1019 ccol = cj + ci[i]; 1020 cval = ca + ci[i]; 1021 for (j=0; j<cnzi; j++) { 1022 brow = ccol[j]; 1023 bnzj = bi[brow+1] - bi[brow]; 1024 bcol = bj + bi[brow]; 1025 bval = ba + bi[brow]; 1026 1027 /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 1028 #if defined(USE_ARRAY) 1029 /* put ba to spdot array */ 1030 for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb]; 1031 /* c(i,j)=A[i,:]*B[j,:]^T */ 1032 for (nexta=0; nexta<anzi; nexta++) { 1033 cval[j] += spdot[acol[nexta]]*aval[nexta]; 1034 } 1035 /* zero spdot array */ 1036 for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0; 1037 #else 1038 nexta = 0; nextb = 0; 1039 while (nexta<anzi && nextb<bnzj) { 1040 while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++; 1041 if (nexta == anzi) break; 1042 while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++; 1043 if (nextb == bnzj) break; 1044 if (acol[nexta] == bcol[nextb]) { 1045 cval[j] += aval[nexta]*bval[nextb]; 1046 nexta++; nextb++; 1047 flops += 2; 1048 } 1049 } 1050 #endif 1051 } 1052 } 1053 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1054 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1055 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 1056 #if defined(USE_ARRAY) 1057 ierr = PetscFree(spdot); 1058 #endif 1059 PetscFunctionReturn(0); 1060 } 1061 1062 #undef __FUNCT__ 1063 #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ" 1064 PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1065 { 1066 PetscErrorCode ierr; 1067 1068 PetscFunctionBegin; 1069 if (scall == MAT_INITIAL_MATRIX) { 1070 ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 1071 } 1072 ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 1073 PetscFunctionReturn(0); 1074 } 1075 1076 #undef __FUNCT__ 1077 #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ" 1078 PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 1079 { 1080 PetscErrorCode ierr; 1081 Mat At; 1082 PetscInt *ati,*atj; 1083 1084 PetscFunctionBegin; 1085 /* create symbolic At */ 1086 ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 1087 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);CHKERRQ(ierr); 1088 1089 At->rmap->bs = A->cmap->bs; 1090 At->cmap->bs = B->cmap->bs; 1091 1092 /* get symbolic C=At*B */ 1093 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 1094 1095 /* clean up */ 1096 ierr = MatDestroy(&At);CHKERRQ(ierr); 1097 ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 1098 PetscFunctionReturn(0); 1099 } 1100 1101 #undef __FUNCT__ 1102 #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ" 1103 PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 1104 { 1105 PetscErrorCode ierr; 1106 Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 1107 PetscInt am =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 1108 PetscInt cm =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 1109 PetscLogDouble flops=0.0; 1110 MatScalar *aa =a->a,*ba,*ca,*caj; 1111 1112 PetscFunctionBegin; 1113 if (!c->a) { 1114 ierr = PetscMalloc((ci[cm]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 1115 1116 c->a = ca; 1117 c->free_a = PETSC_TRUE; 1118 } else { 1119 ca = c->a; 1120 } 1121 /* clear old values in C */ 1122 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 1123 1124 /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 1125 for (i=0; i<am; i++) { 1126 bj = b->j + bi[i]; 1127 ba = b->a + bi[i]; 1128 bnzi = bi[i+1] - bi[i]; 1129 anzi = ai[i+1] - ai[i]; 1130 for (j=0; j<anzi; j++) { 1131 nextb = 0; 1132 crow = *aj++; 1133 cjj = cj + ci[crow]; 1134 caj = ca + ci[crow]; 1135 /* perform sparse axpy operation. Note cjj includes bj. */ 1136 for (k=0; nextb<bnzi; k++) { 1137 if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 1138 caj[k] += (*aa)*(*(ba+nextb)); 1139 nextb++; 1140 } 1141 } 1142 flops += 2*bnzi; 1143 aa++; 1144 } 1145 } 1146 1147 /* Assemble the final matrix and clean up */ 1148 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1149 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1150 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 1151 PetscFunctionReturn(0); 1152 } 1153 1154 #undef __FUNCT__ 1155 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense" 1156 PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 1157 { 1158 PetscErrorCode ierr; 1159 1160 PetscFunctionBegin; 1161 if (scall == MAT_INITIAL_MATRIX) { 1162 ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1163 ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 1164 ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 1165 } 1166 ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1167 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 1168 ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 1169 PetscFunctionReturn(0); 1170 } 1171 1172 #undef __FUNCT__ 1173 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense" 1174 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 1175 { 1176 PetscErrorCode ierr; 1177 1178 PetscFunctionBegin; 1179 ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr); 1180 1181 (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense; 1182 PetscFunctionReturn(0); 1183 } 1184 1185 #undef __FUNCT__ 1186 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense" 1187 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 1188 { 1189 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 1190 PetscErrorCode ierr; 1191 PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 1192 MatScalar *aa; 1193 PetscInt cm = C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n; 1194 PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam; 1195 1196 PetscFunctionBegin; 1197 if (!cm || !cn) PetscFunctionReturn(0); 1198 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); 1199 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); 1200 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); 1201 ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 1202 ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr); 1203 b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 1204 for (col=0; col<cn-4; col += 4) { /* over columns of C */ 1205 colam = col*am; 1206 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1207 r1 = r2 = r3 = r4 = 0.0; 1208 n = a->i[i+1] - a->i[i]; 1209 aj = a->j + a->i[i]; 1210 aa = a->a + a->i[i]; 1211 for (j=0; j<n; j++) { 1212 r1 += (*aa)*b1[*aj]; 1213 r2 += (*aa)*b2[*aj]; 1214 r3 += (*aa)*b3[*aj]; 1215 r4 += (*aa++)*b4[*aj++]; 1216 } 1217 c[colam + i] = r1; 1218 c[colam + am + i] = r2; 1219 c[colam + am2 + i] = r3; 1220 c[colam + am3 + i] = r4; 1221 } 1222 b1 += bm4; 1223 b2 += bm4; 1224 b3 += bm4; 1225 b4 += bm4; 1226 } 1227 for (; col<cn; col++) { /* over extra columns of C */ 1228 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1229 r1 = 0.0; 1230 n = a->i[i+1] - a->i[i]; 1231 aj = a->j + a->i[i]; 1232 aa = a->a + a->i[i]; 1233 1234 for (j=0; j<n; j++) { 1235 r1 += (*aa++)*b1[*aj++]; 1236 } 1237 c[col*am + i] = r1; 1238 } 1239 b1 += bm; 1240 } 1241 ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 1242 ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 1243 ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr); 1244 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1245 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1246 PetscFunctionReturn(0); 1247 } 1248 1249 /* 1250 Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 1251 */ 1252 #undef __FUNCT__ 1253 #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense" 1254 PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 1255 { 1256 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 1257 PetscErrorCode ierr; 1258 PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 1259 MatScalar *aa; 1260 PetscInt cm = C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n,*ii,arm; 1261 PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 1262 1263 PetscFunctionBegin; 1264 if (!cm || !cn) PetscFunctionReturn(0); 1265 ierr = MatDenseGetArray(B,&b);CHKERRQ(ierr); 1266 ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr); 1267 b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 1268 1269 if (a->compressedrow.use) { /* use compressed row format */ 1270 for (col=0; col<cn-4; col += 4) { /* over columns of C */ 1271 colam = col*am; 1272 arm = a->compressedrow.nrows; 1273 ii = a->compressedrow.i; 1274 ridx = a->compressedrow.rindex; 1275 for (i=0; i<arm; i++) { /* over rows of C in those columns */ 1276 r1 = r2 = r3 = r4 = 0.0; 1277 n = ii[i+1] - ii[i]; 1278 aj = a->j + ii[i]; 1279 aa = a->a + ii[i]; 1280 for (j=0; j<n; j++) { 1281 r1 += (*aa)*b1[*aj]; 1282 r2 += (*aa)*b2[*aj]; 1283 r3 += (*aa)*b3[*aj]; 1284 r4 += (*aa++)*b4[*aj++]; 1285 } 1286 c[colam + ridx[i]] += r1; 1287 c[colam + am + ridx[i]] += r2; 1288 c[colam + am2 + ridx[i]] += r3; 1289 c[colam + am3 + ridx[i]] += r4; 1290 } 1291 b1 += bm4; 1292 b2 += bm4; 1293 b3 += bm4; 1294 b4 += bm4; 1295 } 1296 for (; col<cn; col++) { /* over extra columns of C */ 1297 colam = col*am; 1298 arm = a->compressedrow.nrows; 1299 ii = a->compressedrow.i; 1300 ridx = a->compressedrow.rindex; 1301 for (i=0; i<arm; i++) { /* over rows of C in those columns */ 1302 r1 = 0.0; 1303 n = ii[i+1] - ii[i]; 1304 aj = a->j + ii[i]; 1305 aa = a->a + ii[i]; 1306 1307 for (j=0; j<n; j++) { 1308 r1 += (*aa++)*b1[*aj++]; 1309 } 1310 c[colam + ridx[i]] += r1; 1311 } 1312 b1 += bm; 1313 } 1314 } else { 1315 for (col=0; col<cn-4; col += 4) { /* over columns of C */ 1316 colam = col*am; 1317 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1318 r1 = r2 = r3 = r4 = 0.0; 1319 n = a->i[i+1] - a->i[i]; 1320 aj = a->j + a->i[i]; 1321 aa = a->a + a->i[i]; 1322 for (j=0; j<n; j++) { 1323 r1 += (*aa)*b1[*aj]; 1324 r2 += (*aa)*b2[*aj]; 1325 r3 += (*aa)*b3[*aj]; 1326 r4 += (*aa++)*b4[*aj++]; 1327 } 1328 c[colam + i] += r1; 1329 c[colam + am + i] += r2; 1330 c[colam + am2 + i] += r3; 1331 c[colam + am3 + i] += r4; 1332 } 1333 b1 += bm4; 1334 b2 += bm4; 1335 b3 += bm4; 1336 b4 += bm4; 1337 } 1338 for (; col<cn; col++) { /* over extra columns of C */ 1339 colam = col*am; 1340 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1341 r1 = 0.0; 1342 n = a->i[i+1] - a->i[i]; 1343 aj = a->j + a->i[i]; 1344 aa = a->a + a->i[i]; 1345 1346 for (j=0; j<n; j++) { 1347 r1 += (*aa++)*b1[*aj++]; 1348 } 1349 c[colam + i] += r1; 1350 } 1351 b1 += bm; 1352 } 1353 } 1354 ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr); 1355 ierr = MatDenseRestoreArray(B,&b);CHKERRQ(ierr); 1356 ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr); 1357 PetscFunctionReturn(0); 1358 } 1359 1360 #undef __FUNCT__ 1361 #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ" 1362 PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense) 1363 { 1364 PetscErrorCode ierr; 1365 Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 1366 Mat_SeqDense *btdense = (Mat_SeqDense*)Btdense->data; 1367 PetscInt *bi = b->i,*bj=b->j; 1368 PetscInt m = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns; 1369 MatScalar *btval,*btval_den,*ba=b->a; 1370 PetscInt *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors; 1371 1372 PetscFunctionBegin; 1373 btval_den=btdense->v; 1374 ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr); 1375 for (k=0; k<ncolors; k++) { 1376 ncolumns = coloring->ncolumns[k]; 1377 for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */ 1378 col = *(columns + colorforcol[k] + l); 1379 btcol = bj + bi[col]; 1380 btval = ba + bi[col]; 1381 anz = bi[col+1] - bi[col]; 1382 for (j=0; j<anz; j++) { 1383 brow = btcol[j]; 1384 btval_den[brow] = btval[j]; 1385 } 1386 } 1387 btval_den += m; 1388 } 1389 PetscFunctionReturn(0); 1390 } 1391 1392 #undef __FUNCT__ 1393 #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ" 1394 PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp) 1395 { 1396 PetscErrorCode ierr; 1397 Mat_SeqAIJ *csp = (Mat_SeqAIJ*)Csp->data; 1398 PetscInt k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows; 1399 PetscScalar *ca_den,*cp_den,*ca=csp->a; 1400 PetscInt *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow; 1401 1402 PetscFunctionBegin; 1403 ierr = MatGetLocalSize(Csp,&m,NULL);CHKERRQ(ierr); 1404 ierr = MatDenseGetArray(Cden,&ca_den);CHKERRQ(ierr); 1405 cp_den = ca_den; 1406 for (k=0; k<ncolors; k++) { 1407 nrows = matcoloring->nrows[k]; 1408 row = rows + colorforrow[k]; 1409 idx = spidx + colorforrow[k]; 1410 for (l=0; l<nrows; l++) { 1411 ca[idx[l]] = cp_den[row[l]]; 1412 } 1413 cp_den += m; 1414 } 1415 ierr = MatDenseRestoreArray(Cden,&ca_den);CHKERRQ(ierr); 1416 PetscFunctionReturn(0); 1417 } 1418 1419 /* 1420 MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from 1421 MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output 1422 spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ(). 1423 */ 1424 #undef __FUNCT__ 1425 #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color" 1426 PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *nn,const PetscInt *ia[],const PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1427 { 1428 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1429 PetscErrorCode ierr; 1430 PetscInt i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n; 1431 PetscInt nz = a->i[m],row,*jj,mr,col; 1432 PetscInt *cspidx; 1433 1434 PetscFunctionBegin; 1435 *nn = n; 1436 if (!ia) PetscFunctionReturn(0); 1437 if (symmetric) { 1438 SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric"); 1439 ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,(PetscInt**)ia,(PetscInt**)ja);CHKERRQ(ierr); 1440 } else { 1441 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr); 1442 ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 1443 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr); 1444 ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr); 1445 ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr); 1446 jj = a->j; 1447 for (i=0; i<nz; i++) { 1448 collengths[jj[i]]++; 1449 } 1450 cia[0] = oshift; 1451 for (i=0; i<n; i++) { 1452 cia[i+1] = cia[i] + collengths[i]; 1453 } 1454 ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 1455 jj = a->j; 1456 for (row=0; row<m; row++) { 1457 mr = a->i[row+1] - a->i[row]; 1458 for (i=0; i<mr; i++) { 1459 col = *jj++; 1460 1461 cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */ 1462 cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 1463 } 1464 } 1465 ierr = PetscFree(collengths);CHKERRQ(ierr); 1466 *ia = cia; *ja = cja; 1467 *spidx = cspidx; 1468 } 1469 PetscFunctionReturn(0); 1470 } 1471 1472 #undef __FUNCT__ 1473 #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color" 1474 PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,const PetscInt *ia[],const PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1475 { 1476 PetscErrorCode ierr; 1477 1478 PetscFunctionBegin; 1479 if (!ia) PetscFunctionReturn(0); 1480 1481 ierr = PetscFree(*ia);CHKERRQ(ierr); 1482 ierr = PetscFree(*ja);CHKERRQ(ierr); 1483 ierr = PetscFree(*spidx);CHKERRQ(ierr); 1484 PetscFunctionReturn(0); 1485 } 1486 1487 #undef __FUNCT__ 1488 #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ" 1489 PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c) 1490 { 1491 PetscErrorCode ierr; 1492 PetscInt i,n,nrows,N,j,k,m,ncols,col,cm; 1493 const PetscInt *is,*ci,*cj,*row_idx; 1494 PetscInt nis = iscoloring->n,*rowhit,bs = 1; 1495 IS *isa; 1496 PetscBool flg1,flg2; 1497 Mat_SeqAIJ *csp = (Mat_SeqAIJ*)mat->data; 1498 PetscInt *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx; 1499 PetscInt *colorforcol,*columns,*columns_i; 1500 1501 PetscFunctionBegin; 1502 ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 1503 1504 /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */ 1505 ierr = PetscObjectTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr); 1506 ierr = PetscObjectTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr); 1507 if (flg1 || flg2) { 1508 ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 1509 } 1510 1511 N = mat->cmap->N/bs; 1512 c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 1513 c->N = mat->cmap->N/bs; 1514 c->m = mat->rmap->N/bs; 1515 c->rstart = 0; 1516 1517 c->ncolors = nis; 1518 ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 1519 ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 1520 ierr = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr); 1521 ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr); 1522 1523 colorforrow[0] = 0; 1524 rows_i = rows; 1525 columnsforspidx_i = columnsforspidx; 1526 1527 ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr); 1528 ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr); 1529 1530 colorforcol[0] = 0; 1531 columns_i = columns; 1532 1533 ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); /* column-wise storage of mat */ 1534 1535 cm = c->m; 1536 ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 1537 ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr); 1538 for (i=0; i<nis; i++) { 1539 ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 1540 ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 1541 1542 c->ncolumns[i] = n; 1543 if (n) { 1544 ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr); 1545 } 1546 colorforcol[i+1] = colorforcol[i] + n; 1547 columns_i += n; 1548 1549 /* fast, crude version requires O(N*N) work */ 1550 ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr); 1551 1552 /* loop over columns*/ 1553 for (j=0; j<n; j++) { 1554 col = is[j]; 1555 row_idx = cj + ci[col]; 1556 m = ci[col+1] - ci[col]; 1557 /* loop over columns marking them in rowhit */ 1558 for (k=0; k<m; k++) { 1559 idxhit[*row_idx] = spidx[ci[col] + k]; 1560 rowhit[*row_idx++] = col + 1; 1561 } 1562 } 1563 /* count the number of hits */ 1564 nrows = 0; 1565 for (j=0; j<cm; j++) { 1566 if (rowhit[j]) nrows++; 1567 } 1568 c->nrows[i] = nrows; 1569 colorforrow[i+1] = colorforrow[i] + nrows; 1570 1571 nrows = 0; 1572 for (j=0; j<cm; j++) { 1573 if (rowhit[j]) { 1574 rows_i[nrows] = j; 1575 columnsforspidx_i[nrows] = idxhit[j]; 1576 nrows++; 1577 } 1578 } 1579 ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 1580 rows_i += nrows; columnsforspidx_i += nrows; 1581 } 1582 ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 1583 ierr = PetscFree(rowhit);CHKERRQ(ierr); 1584 ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 1585 #if defined(PETSC_USE_DEBUG) 1586 if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]); 1587 #endif 1588 1589 c->colorforrow = colorforrow; 1590 c->rows = rows; 1591 c->columnsforspidx = columnsforspidx; 1592 c->colorforcol = colorforcol; 1593 c->columns = columns; 1594 1595 ierr = PetscFree(idxhit);CHKERRQ(ierr); 1596 PetscFunctionReturn(0); 1597 } 1598