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