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 <petscbt.h> 10 #include <petsc/private/isimpl.h> 11 #include <../src/mat/impls/dense/seq/dense.h> 12 13 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 14 { 15 PetscErrorCode ierr; 16 17 PetscFunctionBegin; 18 if (C->ops->matmultnumeric) { 19 if (C->ops->matmultnumeric == MatMatMultNumeric_SeqAIJ_SeqAIJ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Recursive call"); 20 ierr = (*C->ops->matmultnumeric)(A,B,C);CHKERRQ(ierr); 21 } else { 22 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(A,B,C);CHKERRQ(ierr); 23 } 24 PetscFunctionReturn(0); 25 } 26 27 /* Modified from MatCreateSeqAIJWithArrays() */ 28 PETSC_INTERN PetscErrorCode MatSetSeqAIJWithArrays_private(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt i[],PetscInt j[],PetscScalar a[],MatType mtype,Mat mat) 29 { 30 PetscErrorCode ierr; 31 PetscInt ii; 32 Mat_SeqAIJ *aij; 33 PetscBool isseqaij; 34 35 PetscFunctionBegin; 36 if (m > 0 && i[0]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"i (row indices) must start with 0"); 37 ierr = MatSetSizes(mat,m,n,m,n);CHKERRQ(ierr); 38 39 if (!mtype) { 40 ierr = PetscObjectBaseTypeCompare((PetscObject)mat,MATSEQAIJ,&isseqaij);CHKERRQ(ierr); 41 if (!isseqaij) { ierr = MatSetType(mat,MATSEQAIJ);CHKERRQ(ierr); } 42 } else { 43 ierr = MatSetType(mat,mtype);CHKERRQ(ierr); 44 } 45 ierr = MatSeqAIJSetPreallocation_SeqAIJ(mat,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 46 aij = (Mat_SeqAIJ*)(mat)->data; 47 ierr = PetscMalloc1(m,&aij->imax);CHKERRQ(ierr); 48 ierr = PetscMalloc1(m,&aij->ilen);CHKERRQ(ierr); 49 50 aij->i = i; 51 aij->j = j; 52 aij->a = a; 53 aij->singlemalloc = PETSC_FALSE; 54 aij->nonew = -1; /*this indicates that inserting a new value in the matrix that generates a new nonzero is an error*/ 55 aij->free_a = PETSC_FALSE; 56 aij->free_ij = PETSC_FALSE; 57 58 for (ii=0; ii<m; ii++) { 59 aij->ilen[ii] = aij->imax[ii] = i[ii+1] - i[ii]; 60 } 61 62 PetscFunctionReturn(0); 63 } 64 65 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat C) 66 { 67 PetscErrorCode ierr; 68 Mat_Product *product = C->product; 69 MatProductAlgorithm alg; 70 PetscBool flg; 71 72 PetscFunctionBegin; 73 if (product) { 74 alg = product->alg; 75 } else { 76 alg = "sorted"; 77 } 78 /* sorted */ 79 ierr = PetscStrcmp(alg,"sorted",&flg);CHKERRQ(ierr); 80 if (flg) { 81 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(A,B,fill,C);CHKERRQ(ierr); 82 PetscFunctionReturn(0); 83 } 84 85 /* scalable */ 86 ierr = PetscStrcmp(alg,"scalable",&flg);CHKERRQ(ierr); 87 if (flg) { 88 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr); 89 PetscFunctionReturn(0); 90 } 91 92 /* scalable_fast */ 93 ierr = PetscStrcmp(alg,"scalable_fast",&flg);CHKERRQ(ierr); 94 if (flg) { 95 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(A,B,fill,C);CHKERRQ(ierr); 96 PetscFunctionReturn(0); 97 } 98 99 /* heap */ 100 ierr = PetscStrcmp(alg,"heap",&flg);CHKERRQ(ierr); 101 if (flg) { 102 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(A,B,fill,C);CHKERRQ(ierr); 103 PetscFunctionReturn(0); 104 } 105 106 /* btheap */ 107 ierr = PetscStrcmp(alg,"btheap",&flg);CHKERRQ(ierr); 108 if (flg) { 109 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(A,B,fill,C);CHKERRQ(ierr); 110 PetscFunctionReturn(0); 111 } 112 113 /* llcondensed */ 114 ierr = PetscStrcmp(alg,"llcondensed",&flg);CHKERRQ(ierr); 115 if (flg) { 116 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(A,B,fill,C);CHKERRQ(ierr); 117 PetscFunctionReturn(0); 118 } 119 120 /* rowmerge */ 121 ierr = PetscStrcmp(alg,"rowmerge",&flg);CHKERRQ(ierr); 122 if (flg) { 123 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(A,B,fill,C);CHKERRQ(ierr); 124 PetscFunctionReturn(0); 125 } 126 127 #if defined(PETSC_HAVE_HYPRE) 128 ierr = PetscStrcmp(alg,"hypre",&flg);CHKERRQ(ierr); 129 if (flg) { 130 ierr = MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);CHKERRQ(ierr); 131 PetscFunctionReturn(0); 132 } 133 #endif 134 135 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat Product Algorithm is not supported"); 136 PetscFunctionReturn(0); 137 } 138 139 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_LLCondensed(Mat A,Mat B,PetscReal fill,Mat C) 140 { 141 PetscErrorCode ierr; 142 Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 143 PetscInt *ai=a->i,*bi=b->i,*ci,*cj; 144 PetscInt am =A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 145 PetscReal afill; 146 PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax; 147 PetscTable ta; 148 PetscBT lnkbt; 149 PetscFreeSpaceList free_space=NULL,current_space=NULL; 150 151 PetscFunctionBegin; 152 /* Get ci and cj */ 153 /*---------------*/ 154 /* Allocate ci array, arrays for fill computation and */ 155 /* free space for accumulating nonzero column info */ 156 ierr = PetscMalloc1(am+2,&ci);CHKERRQ(ierr); 157 ci[0] = 0; 158 159 /* create and initialize a linked list */ 160 ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr); 161 MatRowMergeMax_SeqAIJ(b,bm,ta); 162 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 163 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 164 165 ierr = PetscLLCondensedCreate(Crmax,bn,&lnk,&lnkbt);CHKERRQ(ierr); 166 167 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 168 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr); 169 170 current_space = free_space; 171 172 /* Determine ci and cj */ 173 for (i=0; i<am; i++) { 174 anzi = ai[i+1] - ai[i]; 175 aj = a->j + ai[i]; 176 for (j=0; j<anzi; j++) { 177 brow = aj[j]; 178 bnzj = bi[brow+1] - bi[brow]; 179 bj = b->j + bi[brow]; 180 /* add non-zero cols of B into the sorted linked list lnk */ 181 ierr = PetscLLCondensedAddSorted(bnzj,bj,lnk,lnkbt);CHKERRQ(ierr); 182 } 183 cnzi = lnk[0]; 184 185 /* If free space is not available, make more free space */ 186 /* Double the amount of total space in the list */ 187 if (current_space->local_remaining<cnzi) { 188 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 189 ndouble++; 190 } 191 192 /* Copy data into free space, then initialize lnk */ 193 ierr = PetscLLCondensedClean(bn,cnzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 194 195 current_space->array += cnzi; 196 current_space->local_used += cnzi; 197 current_space->local_remaining -= cnzi; 198 199 ci[i+1] = ci[i] + cnzi; 200 } 201 202 /* Column indices are in the list of free space */ 203 /* Allocate space for cj, initialize cj, and */ 204 /* destroy list of free space and other temporary array(s) */ 205 ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr); 206 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 207 ierr = PetscLLCondensedDestroy(lnk,lnkbt);CHKERRQ(ierr); 208 209 /* put together the new symbolic matrix */ 210 ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);CHKERRQ(ierr); 211 ierr = MatSetBlockSizesFromMats(C,A,B);CHKERRQ(ierr); 212 213 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 214 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 215 c = (Mat_SeqAIJ*)(C->data); 216 c->free_a = PETSC_FALSE; 217 c->free_ij = PETSC_TRUE; 218 c->nonew = 0; 219 220 /* fast, needs non-scalable O(bn) array 'abdense' */ 221 C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted; 222 223 /* set MatInfo */ 224 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 225 if (afill < 1.0) afill = 1.0; 226 c->maxnz = ci[am]; 227 c->nz = ci[am]; 228 C->info.mallocs = ndouble; 229 C->info.fill_ratio_given = fill; 230 C->info.fill_ratio_needed = afill; 231 232 #if defined(PETSC_USE_INFO) 233 if (ci[am]) { 234 ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr); 235 ierr = PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 236 } else { 237 ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr); 238 } 239 #endif 240 PetscFunctionReturn(0); 241 } 242 243 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted(Mat A,Mat B,Mat C) 244 { 245 PetscErrorCode ierr; 246 PetscLogDouble flops=0.0; 247 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 248 Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 249 Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 250 PetscInt *ai =a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 251 PetscInt am =A->rmap->n,cm=C->rmap->n; 252 PetscInt i,j,k,anzi,bnzi,cnzi,brow; 253 PetscScalar *aa=a->a,*ba=b->a,*baj,*ca,valtmp; 254 PetscScalar *ab_dense; 255 PetscContainer cab_dense; 256 257 PetscFunctionBegin; 258 if (!c->a) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */ 259 ierr = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr); 260 c->a = ca; 261 c->free_a = PETSC_TRUE; 262 } else ca = c->a; 263 264 /* TODO this should be done in the symbolic phase */ 265 /* However, this function is so heavily used (sometimes in an hidden way through multnumeric function pointers 266 that is hard to eradicate) */ 267 ierr = PetscObjectQuery((PetscObject)C,"__PETSc__ab_dense",(PetscObject*)&cab_dense);CHKERRQ(ierr); 268 if (!cab_dense) { 269 ierr = PetscMalloc1(B->cmap->N,&ab_dense);CHKERRQ(ierr); 270 ierr = PetscContainerCreate(PETSC_COMM_SELF,&cab_dense);CHKERRQ(ierr); 271 ierr = PetscContainerSetPointer(cab_dense,ab_dense);CHKERRQ(ierr); 272 ierr = PetscContainerSetUserDestroy(cab_dense,PetscContainerUserDestroyDefault);CHKERRQ(ierr); 273 ierr = PetscObjectCompose((PetscObject)C,"__PETSc__ab_dense",(PetscObject)cab_dense);CHKERRQ(ierr); 274 ierr = PetscObjectDereference((PetscObject)cab_dense);CHKERRQ(ierr); 275 } 276 ierr = PetscContainerGetPointer(cab_dense,(void**)&ab_dense);CHKERRQ(ierr); 277 ierr = PetscArrayzero(ab_dense,B->cmap->N);CHKERRQ(ierr); 278 279 /* clean old values in C */ 280 ierr = PetscArrayzero(ca,ci[cm]);CHKERRQ(ierr); 281 /* Traverse A row-wise. */ 282 /* Build the ith row in C by summing over nonzero columns in A, */ 283 /* the rows of B corresponding to nonzeros of A. */ 284 for (i=0; i<am; i++) { 285 anzi = ai[i+1] - ai[i]; 286 for (j=0; j<anzi; j++) { 287 brow = aj[j]; 288 bnzi = bi[brow+1] - bi[brow]; 289 bjj = bj + bi[brow]; 290 baj = ba + bi[brow]; 291 /* perform dense axpy */ 292 valtmp = aa[j]; 293 for (k=0; k<bnzi; k++) { 294 ab_dense[bjj[k]] += valtmp*baj[k]; 295 } 296 flops += 2*bnzi; 297 } 298 aj += anzi; aa += anzi; 299 300 cnzi = ci[i+1] - ci[i]; 301 for (k=0; k<cnzi; k++) { 302 ca[k] += ab_dense[cj[k]]; 303 ab_dense[cj[k]] = 0.0; /* zero ab_dense */ 304 } 305 flops += cnzi; 306 cj += cnzi; ca += cnzi; 307 } 308 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 309 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 310 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 311 PetscFunctionReturn(0); 312 } 313 314 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C) 315 { 316 PetscErrorCode ierr; 317 PetscLogDouble flops=0.0; 318 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 319 Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 320 Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 321 PetscInt *ai = a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 322 PetscInt am = A->rmap->N,cm=C->rmap->N; 323 PetscInt i,j,k,anzi,bnzi,cnzi,brow; 324 PetscScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp; 325 PetscInt nextb; 326 327 PetscFunctionBegin; 328 if (!ca) { /* first call of MatMatMultNumeric_SeqAIJ_SeqAIJ, allocate ca and matmult_abdense */ 329 ierr = PetscMalloc1(ci[cm]+1,&ca);CHKERRQ(ierr); 330 c->a = ca; 331 c->free_a = PETSC_TRUE; 332 } 333 334 /* clean old values in C */ 335 ierr = PetscArrayzero(ca,ci[cm]);CHKERRQ(ierr); 336 /* Traverse A row-wise. */ 337 /* Build the ith row in C by summing over nonzero columns in A, */ 338 /* the rows of B corresponding to nonzeros of A. */ 339 for (i=0; i<am; i++) { 340 anzi = ai[i+1] - ai[i]; 341 cnzi = ci[i+1] - ci[i]; 342 for (j=0; j<anzi; j++) { 343 brow = aj[j]; 344 bnzi = bi[brow+1] - bi[brow]; 345 bjj = bj + bi[brow]; 346 baj = ba + bi[brow]; 347 /* perform sparse axpy */ 348 valtmp = aa[j]; 349 nextb = 0; 350 for (k=0; nextb<bnzi; k++) { 351 if (cj[k] == bjj[nextb]) { /* ccol == bcol */ 352 ca[k] += valtmp*baj[nextb++]; 353 } 354 } 355 flops += 2*bnzi; 356 } 357 aj += anzi; aa += anzi; 358 cj += cnzi; ca += cnzi; 359 } 360 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 361 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 362 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 363 PetscFunctionReturn(0); 364 } 365 366 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable_fast(Mat A,Mat B,PetscReal fill,Mat C) 367 { 368 PetscErrorCode ierr; 369 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 370 PetscInt *ai = a->i,*bi=b->i,*ci,*cj; 371 PetscInt am = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 372 MatScalar *ca; 373 PetscReal afill; 374 PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax; 375 PetscTable ta; 376 PetscFreeSpaceList free_space=NULL,current_space=NULL; 377 378 PetscFunctionBegin; 379 /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_fast() */ 380 /*-----------------------------------------------------------------------------------------*/ 381 /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 382 ierr = PetscMalloc1(am+2,&ci);CHKERRQ(ierr); 383 ci[0] = 0; 384 385 /* create and initialize a linked list */ 386 ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr); 387 MatRowMergeMax_SeqAIJ(b,bm,ta); 388 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 389 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 390 391 ierr = PetscLLCondensedCreate_fast(Crmax,&lnk);CHKERRQ(ierr); 392 393 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 394 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr); 395 current_space = free_space; 396 397 /* Determine ci and cj */ 398 for (i=0; i<am; i++) { 399 anzi = ai[i+1] - ai[i]; 400 aj = a->j + ai[i]; 401 for (j=0; j<anzi; j++) { 402 brow = aj[j]; 403 bnzj = bi[brow+1] - bi[brow]; 404 bj = b->j + bi[brow]; 405 /* add non-zero cols of B into the sorted linked list lnk */ 406 ierr = PetscLLCondensedAddSorted_fast(bnzj,bj,lnk);CHKERRQ(ierr); 407 } 408 cnzi = lnk[1]; 409 410 /* If free space is not available, make more free space */ 411 /* Double the amount of total space in the list */ 412 if (current_space->local_remaining<cnzi) { 413 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 414 ndouble++; 415 } 416 417 /* Copy data into free space, then initialize lnk */ 418 ierr = PetscLLCondensedClean_fast(cnzi,current_space->array,lnk);CHKERRQ(ierr); 419 420 current_space->array += cnzi; 421 current_space->local_used += cnzi; 422 current_space->local_remaining -= cnzi; 423 424 ci[i+1] = ci[i] + cnzi; 425 } 426 427 /* Column indices are in the list of free space */ 428 /* Allocate space for cj, initialize cj, and */ 429 /* destroy list of free space and other temporary array(s) */ 430 ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr); 431 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 432 ierr = PetscLLCondensedDestroy_fast(lnk);CHKERRQ(ierr); 433 434 /* Allocate space for ca */ 435 ierr = PetscCalloc1(ci[am]+1,&ca);CHKERRQ(ierr); 436 437 /* put together the new symbolic matrix */ 438 ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,((PetscObject)A)->type_name,C);CHKERRQ(ierr); 439 ierr = MatSetBlockSizesFromMats(C,A,B);CHKERRQ(ierr); 440 441 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 442 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 443 c = (Mat_SeqAIJ*)(C->data); 444 c->free_a = PETSC_TRUE; 445 c->free_ij = PETSC_TRUE; 446 c->nonew = 0; 447 448 /* slower, less memory */ 449 C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; 450 451 /* set MatInfo */ 452 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 453 if (afill < 1.0) afill = 1.0; 454 c->maxnz = ci[am]; 455 c->nz = ci[am]; 456 C->info.mallocs = ndouble; 457 C->info.fill_ratio_given = fill; 458 C->info.fill_ratio_needed = afill; 459 460 #if defined(PETSC_USE_INFO) 461 if (ci[am]) { 462 ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr); 463 ierr = PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 464 } else { 465 ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr); 466 } 467 #endif 468 PetscFunctionReturn(0); 469 } 470 471 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat C) 472 { 473 PetscErrorCode ierr; 474 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 475 PetscInt *ai = a->i,*bi=b->i,*ci,*cj; 476 PetscInt am = A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 477 MatScalar *ca; 478 PetscReal afill; 479 PetscInt i,j,anzi,brow,bnzj,cnzi,*bj,*aj,*lnk,ndouble=0,Crmax; 480 PetscTable ta; 481 PetscFreeSpaceList free_space=NULL,current_space=NULL; 482 483 PetscFunctionBegin; 484 /* Get ci and cj - same as MatMatMultSymbolic_SeqAIJ_SeqAIJ except using PetscLLxxx_Scalalbe() */ 485 /*---------------------------------------------------------------------------------------------*/ 486 /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 487 ierr = PetscMalloc1(am+2,&ci);CHKERRQ(ierr); 488 ci[0] = 0; 489 490 /* create and initialize a linked list */ 491 ierr = PetscTableCreate(bn,bn,&ta);CHKERRQ(ierr); 492 MatRowMergeMax_SeqAIJ(b,bm,ta); 493 ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 494 ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 495 ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 496 497 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 498 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr); 499 current_space = free_space; 500 501 /* Determine ci and cj */ 502 for (i=0; i<am; i++) { 503 anzi = ai[i+1] - ai[i]; 504 aj = a->j + ai[i]; 505 for (j=0; j<anzi; j++) { 506 brow = aj[j]; 507 bnzj = bi[brow+1] - bi[brow]; 508 bj = b->j + bi[brow]; 509 /* add non-zero cols of B into the sorted linked list lnk */ 510 ierr = PetscLLCondensedAddSorted_Scalable(bnzj,bj,lnk);CHKERRQ(ierr); 511 } 512 cnzi = lnk[0]; 513 514 /* If free space is not available, make more free space */ 515 /* Double the amount of total space in the list */ 516 if (current_space->local_remaining<cnzi) { 517 ierr = PetscFreeSpaceGet(PetscIntSumTruncate(cnzi,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 518 ndouble++; 519 } 520 521 /* Copy data into free space, then initialize lnk */ 522 ierr = PetscLLCondensedClean_Scalable(cnzi,current_space->array,lnk);CHKERRQ(ierr); 523 524 current_space->array += cnzi; 525 current_space->local_used += cnzi; 526 current_space->local_remaining -= cnzi; 527 528 ci[i+1] = ci[i] + cnzi; 529 } 530 531 /* Column indices are in the list of free space */ 532 /* Allocate space for cj, initialize cj, and */ 533 /* destroy list of free space and other temporary array(s) */ 534 ierr = PetscMalloc1(ci[am]+1,&cj);CHKERRQ(ierr); 535 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 536 ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 537 538 /* Allocate space for ca */ 539 /*-----------------------*/ 540 ierr = PetscCalloc1(ci[am]+1,&ca);CHKERRQ(ierr); 541 542 /* put together the new symbolic matrix */ 543 ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,ca,((PetscObject)A)->type_name,C);CHKERRQ(ierr); 544 ierr = MatSetBlockSizesFromMats(C,A,B);CHKERRQ(ierr); 545 546 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 547 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 548 c = (Mat_SeqAIJ*)(C->data); 549 c->free_a = PETSC_TRUE; 550 c->free_ij = PETSC_TRUE; 551 c->nonew = 0; 552 553 /* slower, less memory */ 554 C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; 555 556 /* set MatInfo */ 557 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 558 if (afill < 1.0) afill = 1.0; 559 c->maxnz = ci[am]; 560 c->nz = ci[am]; 561 C->info.mallocs = ndouble; 562 C->info.fill_ratio_given = fill; 563 C->info.fill_ratio_needed = afill; 564 565 #if defined(PETSC_USE_INFO) 566 if (ci[am]) { 567 ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr); 568 ierr = PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 569 } else { 570 ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr); 571 } 572 #endif 573 PetscFunctionReturn(0); 574 } 575 576 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Heap(Mat A,Mat B,PetscReal fill,Mat C) 577 { 578 PetscErrorCode ierr; 579 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 580 const PetscInt *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j; 581 PetscInt *ci,*cj,*bb; 582 PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 583 PetscReal afill; 584 PetscInt i,j,col,ndouble = 0; 585 PetscFreeSpaceList free_space=NULL,current_space=NULL; 586 PetscHeap h; 587 588 PetscFunctionBegin; 589 /* Get ci and cj - by merging sorted rows using a heap */ 590 /*---------------------------------------------------------------------------------------------*/ 591 /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 592 ierr = PetscMalloc1(am+2,&ci);CHKERRQ(ierr); 593 ci[0] = 0; 594 595 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 596 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr); 597 current_space = free_space; 598 599 ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr); 600 ierr = PetscMalloc1(a->rmax,&bb);CHKERRQ(ierr); 601 602 /* Determine ci and cj */ 603 for (i=0; i<am; i++) { 604 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 */ 605 const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 606 ci[i+1] = ci[i]; 607 /* Populate the min heap */ 608 for (j=0; j<anzi; j++) { 609 bb[j] = bi[acol[j]]; /* bb points at the start of the row */ 610 if (bb[j] < bi[acol[j]+1]) { /* Add if row is nonempty */ 611 ierr = PetscHeapAdd(h,j,bj[bb[j]++]);CHKERRQ(ierr); 612 } 613 } 614 /* Pick off the min element, adding it to free space */ 615 ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 616 while (j >= 0) { 617 if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */ 618 ierr = PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),¤t_space);CHKERRQ(ierr); 619 ndouble++; 620 } 621 *(current_space->array++) = col; 622 current_space->local_used++; 623 current_space->local_remaining--; 624 ci[i+1]++; 625 626 /* stash if anything else remains in this row of B */ 627 if (bb[j] < bi[acol[j]+1]) {ierr = PetscHeapStash(h,j,bj[bb[j]++]);CHKERRQ(ierr);} 628 while (1) { /* pop and stash any other rows of B that also had an entry in this column */ 629 PetscInt j2,col2; 630 ierr = PetscHeapPeek(h,&j2,&col2);CHKERRQ(ierr); 631 if (col2 != col) break; 632 ierr = PetscHeapPop(h,&j2,&col2);CHKERRQ(ierr); 633 if (bb[j2] < bi[acol[j2]+1]) {ierr = PetscHeapStash(h,j2,bj[bb[j2]++]);CHKERRQ(ierr);} 634 } 635 /* Put any stashed elements back into the min heap */ 636 ierr = PetscHeapUnstash(h);CHKERRQ(ierr); 637 ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 638 } 639 } 640 ierr = PetscFree(bb);CHKERRQ(ierr); 641 ierr = PetscHeapDestroy(&h);CHKERRQ(ierr); 642 643 /* Column indices are in the list of free space */ 644 /* Allocate space for cj, initialize cj, and */ 645 /* destroy list of free space and other temporary array(s) */ 646 ierr = PetscMalloc1(ci[am],&cj);CHKERRQ(ierr); 647 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 648 649 /* put together the new symbolic matrix */ 650 ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);CHKERRQ(ierr); 651 ierr = MatSetBlockSizesFromMats(C,A,B);CHKERRQ(ierr); 652 653 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 654 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 655 c = (Mat_SeqAIJ*)(C->data); 656 c->free_a = PETSC_TRUE; 657 c->free_ij = PETSC_TRUE; 658 c->nonew = 0; 659 660 C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted; 661 662 /* set MatInfo */ 663 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 664 if (afill < 1.0) afill = 1.0; 665 c->maxnz = ci[am]; 666 c->nz = ci[am]; 667 C->info.mallocs = ndouble; 668 C->info.fill_ratio_given = fill; 669 C->info.fill_ratio_needed = afill; 670 671 #if defined(PETSC_USE_INFO) 672 if (ci[am]) { 673 ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr); 674 ierr = PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 675 } else { 676 ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr); 677 } 678 #endif 679 PetscFunctionReturn(0); 680 } 681 682 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_BTHeap(Mat A,Mat B,PetscReal fill,Mat C) 683 { 684 PetscErrorCode ierr; 685 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 686 const PetscInt *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j; 687 PetscInt *ci,*cj,*bb; 688 PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 689 PetscReal afill; 690 PetscInt i,j,col,ndouble = 0; 691 PetscFreeSpaceList free_space=NULL,current_space=NULL; 692 PetscHeap h; 693 PetscBT bt; 694 695 PetscFunctionBegin; 696 /* Get ci and cj - using a heap for the sorted rows, but use BT so that each index is only added once */ 697 /*---------------------------------------------------------------------------------------------*/ 698 /* Allocate arrays for fill computation and free space for accumulating nonzero column */ 699 ierr = PetscMalloc1(am+2,&ci);CHKERRQ(ierr); 700 ci[0] = 0; 701 702 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 703 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ai[am],bi[bm])),&free_space);CHKERRQ(ierr); 704 705 current_space = free_space; 706 707 ierr = PetscHeapCreate(a->rmax,&h);CHKERRQ(ierr); 708 ierr = PetscMalloc1(a->rmax,&bb);CHKERRQ(ierr); 709 ierr = PetscBTCreate(bn,&bt);CHKERRQ(ierr); 710 711 /* Determine ci and cj */ 712 for (i=0; i<am; i++) { 713 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 */ 714 const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 715 const PetscInt *fptr = current_space->array; /* Save beginning of the row so we can clear the BT later */ 716 ci[i+1] = ci[i]; 717 /* Populate the min heap */ 718 for (j=0; j<anzi; j++) { 719 PetscInt brow = acol[j]; 720 for (bb[j] = bi[brow]; bb[j] < bi[brow+1]; bb[j]++) { 721 PetscInt bcol = bj[bb[j]]; 722 if (!PetscBTLookupSet(bt,bcol)) { /* new entry */ 723 ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr); 724 bb[j]++; 725 break; 726 } 727 } 728 } 729 /* Pick off the min element, adding it to free space */ 730 ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 731 while (j >= 0) { 732 if (current_space->local_remaining < 1) { /* double the size, but don't exceed 16 MiB */ 733 fptr = NULL; /* need PetscBTMemzero */ 734 ierr = PetscFreeSpaceGet(PetscMin(PetscIntMultTruncate(2,current_space->total_array_size),16 << 20),¤t_space);CHKERRQ(ierr); 735 ndouble++; 736 } 737 *(current_space->array++) = col; 738 current_space->local_used++; 739 current_space->local_remaining--; 740 ci[i+1]++; 741 742 /* stash if anything else remains in this row of B */ 743 for (; bb[j] < bi[acol[j]+1]; bb[j]++) { 744 PetscInt bcol = bj[bb[j]]; 745 if (!PetscBTLookupSet(bt,bcol)) { /* new entry */ 746 ierr = PetscHeapAdd(h,j,bcol);CHKERRQ(ierr); 747 bb[j]++; 748 break; 749 } 750 } 751 ierr = PetscHeapPop(h,&j,&col);CHKERRQ(ierr); 752 } 753 if (fptr) { /* Clear the bits for this row */ 754 for (; fptr<current_space->array; fptr++) {ierr = PetscBTClear(bt,*fptr);CHKERRQ(ierr);} 755 } else { /* We reallocated so we don't remember (easily) how to clear only the bits we changed */ 756 ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr); 757 } 758 } 759 ierr = PetscFree(bb);CHKERRQ(ierr); 760 ierr = PetscHeapDestroy(&h);CHKERRQ(ierr); 761 ierr = PetscBTDestroy(&bt);CHKERRQ(ierr); 762 763 /* Column indices are in the list of free space */ 764 /* Allocate space for cj, initialize cj, and */ 765 /* destroy list of free space and other temporary array(s) */ 766 ierr = PetscMalloc1(ci[am],&cj);CHKERRQ(ierr); 767 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 768 769 /* put together the new symbolic matrix */ 770 ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);CHKERRQ(ierr); 771 ierr = MatSetBlockSizesFromMats(C,A,B);CHKERRQ(ierr); 772 773 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 774 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 775 c = (Mat_SeqAIJ*)(C->data); 776 c->free_a = PETSC_TRUE; 777 c->free_ij = PETSC_TRUE; 778 c->nonew = 0; 779 780 C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted; 781 782 /* set MatInfo */ 783 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 784 if (afill < 1.0) afill = 1.0; 785 c->maxnz = ci[am]; 786 c->nz = ci[am]; 787 C->info.mallocs = ndouble; 788 C->info.fill_ratio_given = fill; 789 C->info.fill_ratio_needed = afill; 790 791 #if defined(PETSC_USE_INFO) 792 if (ci[am]) { 793 ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr); 794 ierr = PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 795 } else { 796 ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr); 797 } 798 #endif 799 PetscFunctionReturn(0); 800 } 801 802 803 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_RowMerge(Mat A,Mat B,PetscReal fill,Mat C) 804 { 805 PetscErrorCode ierr; 806 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 807 const PetscInt *ai=a->i,*bi=b->i,*aj=a->j,*bj=b->j,*inputi,*inputj,*inputcol,*inputcol_L1; 808 PetscInt *ci,*cj,*outputj,worki_L1[9],worki_L2[9]; 809 PetscInt c_maxmem,a_maxrownnz=0,a_rownnz; 810 const PetscInt workcol[8]={0,1,2,3,4,5,6,7}; 811 const PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 812 const PetscInt *brow_ptr[8],*brow_end[8]; 813 PetscInt window[8]; 814 PetscInt window_min,old_window_min,ci_nnz,outputi_nnz=0,L1_nrows,L2_nrows; 815 PetscInt i,k,ndouble=0,L1_rowsleft,rowsleft; 816 PetscReal afill; 817 PetscInt *workj_L1,*workj_L2,*workj_L3; 818 PetscInt L1_nnz,L2_nnz; 819 820 /* Step 1: Get upper bound on memory required for allocation. 821 Because of the way virtual memory works, 822 only the memory pages that are actually needed will be physically allocated. */ 823 PetscFunctionBegin; 824 ierr = PetscMalloc1(am+1,&ci);CHKERRQ(ierr); 825 for (i=0; i<am; i++) { 826 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 */ 827 const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 828 a_rownnz = 0; 829 for (k=0; k<anzi; ++k) { 830 a_rownnz += bi[acol[k]+1] - bi[acol[k]]; 831 if (a_rownnz > bn) { 832 a_rownnz = bn; 833 break; 834 } 835 } 836 a_maxrownnz = PetscMax(a_maxrownnz, a_rownnz); 837 } 838 /* temporary work areas for merging rows */ 839 ierr = PetscMalloc1(a_maxrownnz*8,&workj_L1);CHKERRQ(ierr); 840 ierr = PetscMalloc1(a_maxrownnz*8,&workj_L2);CHKERRQ(ierr); 841 ierr = PetscMalloc1(a_maxrownnz,&workj_L3);CHKERRQ(ierr); 842 843 /* This should be enough for almost all matrices. If not, memory is reallocated later. */ 844 c_maxmem = 8*(ai[am]+bi[bm]); 845 /* Step 2: Populate pattern for C */ 846 ierr = PetscMalloc1(c_maxmem,&cj);CHKERRQ(ierr); 847 848 ci_nnz = 0; 849 ci[0] = 0; 850 worki_L1[0] = 0; 851 worki_L2[0] = 0; 852 for (i=0; i<am; i++) { 853 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 */ 854 const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 855 rowsleft = anzi; 856 inputcol_L1 = acol; 857 L2_nnz = 0; 858 L2_nrows = 1; /* Number of rows to be merged on Level 3. output of L3 already exists -> initial value 1 */ 859 worki_L2[1] = 0; 860 outputi_nnz = 0; 861 862 /* If the number of indices in C so far + the max number of columns in the next row > c_maxmem -> allocate more memory */ 863 while (ci_nnz+a_maxrownnz > c_maxmem) { 864 c_maxmem *= 2; 865 ndouble++; 866 ierr = PetscRealloc(sizeof(PetscInt)*c_maxmem,&cj);CHKERRQ(ierr); 867 } 868 869 while (rowsleft) { 870 L1_rowsleft = PetscMin(64, rowsleft); /* In the inner loop max 64 rows of B can be merged */ 871 L1_nrows = 0; 872 L1_nnz = 0; 873 inputcol = inputcol_L1; 874 inputi = bi; 875 inputj = bj; 876 877 /* The following macro is used to specialize for small rows in A. 878 This helps with compiler unrolling, improving performance substantially. 879 Input: inputj inputi inputcol bn 880 Output: outputj outputi_nnz */ 881 #define MatMatMultSymbolic_RowMergeMacro(ANNZ) \ 882 window_min = bn; \ 883 outputi_nnz = 0; \ 884 for (k=0; k<ANNZ; ++k) { \ 885 brow_ptr[k] = inputj + inputi[inputcol[k]]; \ 886 brow_end[k] = inputj + inputi[inputcol[k]+1]; \ 887 window[k] = (brow_ptr[k] != brow_end[k]) ? *brow_ptr[k] : bn; \ 888 window_min = PetscMin(window[k], window_min); \ 889 } \ 890 while (window_min < bn) { \ 891 outputj[outputi_nnz++] = window_min; \ 892 /* advance front and compute new minimum */ \ 893 old_window_min = window_min; \ 894 window_min = bn; \ 895 for (k=0; k<ANNZ; ++k) { \ 896 if (window[k] == old_window_min) { \ 897 brow_ptr[k]++; \ 898 window[k] = (brow_ptr[k] != brow_end[k]) ? *brow_ptr[k] : bn; \ 899 } \ 900 window_min = PetscMin(window[k], window_min); \ 901 } \ 902 } 903 904 /************** L E V E L 1 ***************/ 905 /* Merge up to 8 rows of B to L1 work array*/ 906 while (L1_rowsleft) { 907 outputi_nnz = 0; 908 if (anzi > 8) outputj = workj_L1 + L1_nnz; /* Level 1 rowmerge*/ 909 else outputj = cj + ci_nnz; /* Merge directly to C */ 910 911 switch (L1_rowsleft) { 912 case 1: brow_ptr[0] = inputj + inputi[inputcol[0]]; 913 brow_end[0] = inputj + inputi[inputcol[0]+1]; 914 for (; brow_ptr[0] != brow_end[0]; ++brow_ptr[0]) outputj[outputi_nnz++] = *brow_ptr[0]; /* copy row in b over */ 915 inputcol += L1_rowsleft; 916 rowsleft -= L1_rowsleft; 917 L1_rowsleft = 0; 918 break; 919 case 2: MatMatMultSymbolic_RowMergeMacro(2); 920 inputcol += L1_rowsleft; 921 rowsleft -= L1_rowsleft; 922 L1_rowsleft = 0; 923 break; 924 case 3: MatMatMultSymbolic_RowMergeMacro(3); 925 inputcol += L1_rowsleft; 926 rowsleft -= L1_rowsleft; 927 L1_rowsleft = 0; 928 break; 929 case 4: MatMatMultSymbolic_RowMergeMacro(4); 930 inputcol += L1_rowsleft; 931 rowsleft -= L1_rowsleft; 932 L1_rowsleft = 0; 933 break; 934 case 5: MatMatMultSymbolic_RowMergeMacro(5); 935 inputcol += L1_rowsleft; 936 rowsleft -= L1_rowsleft; 937 L1_rowsleft = 0; 938 break; 939 case 6: MatMatMultSymbolic_RowMergeMacro(6); 940 inputcol += L1_rowsleft; 941 rowsleft -= L1_rowsleft; 942 L1_rowsleft = 0; 943 break; 944 case 7: MatMatMultSymbolic_RowMergeMacro(7); 945 inputcol += L1_rowsleft; 946 rowsleft -= L1_rowsleft; 947 L1_rowsleft = 0; 948 break; 949 default: MatMatMultSymbolic_RowMergeMacro(8); 950 inputcol += 8; 951 rowsleft -= 8; 952 L1_rowsleft -= 8; 953 break; 954 } 955 inputcol_L1 = inputcol; 956 L1_nnz += outputi_nnz; 957 worki_L1[++L1_nrows] = L1_nnz; 958 } 959 960 /********************** L E V E L 2 ************************/ 961 /* Merge from L1 work array to either C or to L2 work array */ 962 if (anzi > 8) { 963 inputi = worki_L1; 964 inputj = workj_L1; 965 inputcol = workcol; 966 outputi_nnz = 0; 967 968 if (anzi <= 64) outputj = cj + ci_nnz; /* Merge from L1 work array to C */ 969 else outputj = workj_L2 + L2_nnz; /* Merge from L1 work array to L2 work array */ 970 971 switch (L1_nrows) { 972 case 1: brow_ptr[0] = inputj + inputi[inputcol[0]]; 973 brow_end[0] = inputj + inputi[inputcol[0]+1]; 974 for (; brow_ptr[0] != brow_end[0]; ++brow_ptr[0]) outputj[outputi_nnz++] = *brow_ptr[0]; /* copy row in b over */ 975 break; 976 case 2: MatMatMultSymbolic_RowMergeMacro(2); break; 977 case 3: MatMatMultSymbolic_RowMergeMacro(3); break; 978 case 4: MatMatMultSymbolic_RowMergeMacro(4); break; 979 case 5: MatMatMultSymbolic_RowMergeMacro(5); break; 980 case 6: MatMatMultSymbolic_RowMergeMacro(6); break; 981 case 7: MatMatMultSymbolic_RowMergeMacro(7); break; 982 case 8: MatMatMultSymbolic_RowMergeMacro(8); break; 983 default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatMult logic error: Not merging 1-8 rows from L1 work array!"); 984 } 985 L2_nnz += outputi_nnz; 986 worki_L2[++L2_nrows] = L2_nnz; 987 988 /************************ L E V E L 3 **********************/ 989 /* Merge from L2 work array to either C or to L2 work array */ 990 if (anzi > 64 && (L2_nrows == 8 || rowsleft == 0)) { 991 inputi = worki_L2; 992 inputj = workj_L2; 993 inputcol = workcol; 994 outputi_nnz = 0; 995 if (rowsleft) outputj = workj_L3; 996 else outputj = cj + ci_nnz; 997 switch (L2_nrows) { 998 case 1: brow_ptr[0] = inputj + inputi[inputcol[0]]; 999 brow_end[0] = inputj + inputi[inputcol[0]+1]; 1000 for (; brow_ptr[0] != brow_end[0]; ++brow_ptr[0]) outputj[outputi_nnz++] = *brow_ptr[0]; /* copy row in b over */ 1001 break; 1002 case 2: MatMatMultSymbolic_RowMergeMacro(2); break; 1003 case 3: MatMatMultSymbolic_RowMergeMacro(3); break; 1004 case 4: MatMatMultSymbolic_RowMergeMacro(4); break; 1005 case 5: MatMatMultSymbolic_RowMergeMacro(5); break; 1006 case 6: MatMatMultSymbolic_RowMergeMacro(6); break; 1007 case 7: MatMatMultSymbolic_RowMergeMacro(7); break; 1008 case 8: MatMatMultSymbolic_RowMergeMacro(8); break; 1009 default: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatMatMult logic error: Not merging 1-8 rows from L2 work array!"); 1010 } 1011 L2_nrows = 1; 1012 L2_nnz = outputi_nnz; 1013 worki_L2[1] = outputi_nnz; 1014 /* Copy to workj_L2 */ 1015 if (rowsleft) { 1016 for (k=0; k<outputi_nnz; ++k) workj_L2[k] = outputj[k]; 1017 } 1018 } 1019 } 1020 } /* while (rowsleft) */ 1021 #undef MatMatMultSymbolic_RowMergeMacro 1022 1023 /* terminate current row */ 1024 ci_nnz += outputi_nnz; 1025 ci[i+1] = ci_nnz; 1026 } 1027 1028 /* Step 3: Create the new symbolic matrix */ 1029 ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);CHKERRQ(ierr); 1030 ierr = MatSetBlockSizesFromMats(C,A,B);CHKERRQ(ierr); 1031 1032 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 1033 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 1034 c = (Mat_SeqAIJ*)(C->data); 1035 c->free_a = PETSC_TRUE; 1036 c->free_ij = PETSC_TRUE; 1037 c->nonew = 0; 1038 1039 C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted; 1040 1041 /* set MatInfo */ 1042 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 1043 if (afill < 1.0) afill = 1.0; 1044 c->maxnz = ci[am]; 1045 c->nz = ci[am]; 1046 C->info.mallocs = ndouble; 1047 C->info.fill_ratio_given = fill; 1048 C->info.fill_ratio_needed = afill; 1049 1050 #if defined(PETSC_USE_INFO) 1051 if (ci[am]) { 1052 ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr); 1053 ierr = PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 1054 } else { 1055 ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr); 1056 } 1057 #endif 1058 1059 /* Step 4: Free temporary work areas */ 1060 ierr = PetscFree(workj_L1);CHKERRQ(ierr); 1061 ierr = PetscFree(workj_L2);CHKERRQ(ierr); 1062 ierr = PetscFree(workj_L3);CHKERRQ(ierr); 1063 PetscFunctionReturn(0); 1064 } 1065 1066 /* concatenate unique entries and then sort */ 1067 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Sorted(Mat A,Mat B,PetscReal fill,Mat C) 1068 { 1069 PetscErrorCode ierr; 1070 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 1071 const PetscInt *ai = a->i,*bi=b->i,*aj=a->j,*bj=b->j; 1072 PetscInt *ci,*cj; 1073 PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N; 1074 PetscReal afill; 1075 PetscInt i,j,ndouble = 0; 1076 PetscSegBuffer seg,segrow; 1077 char *seen; 1078 1079 PetscFunctionBegin; 1080 ierr = PetscMalloc1(am+1,&ci);CHKERRQ(ierr); 1081 ci[0] = 0; 1082 1083 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 1084 ierr = PetscSegBufferCreate(sizeof(PetscInt),(PetscInt)(fill*(ai[am]+bi[bm])),&seg);CHKERRQ(ierr); 1085 ierr = PetscSegBufferCreate(sizeof(PetscInt),100,&segrow);CHKERRQ(ierr); 1086 ierr = PetscCalloc1(bn,&seen);CHKERRQ(ierr); 1087 1088 /* Determine ci and cj */ 1089 for (i=0; i<am; i++) { 1090 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 */ 1091 const PetscInt *acol = aj + ai[i]; /* column indices of nonzero entries in this row */ 1092 PetscInt packlen = 0,*PETSC_RESTRICT crow; 1093 /* Pack segrow */ 1094 for (j=0; j<anzi; j++) { 1095 PetscInt brow = acol[j],bjstart = bi[brow],bjend = bi[brow+1],k; 1096 for (k=bjstart; k<bjend; k++) { 1097 PetscInt bcol = bj[k]; 1098 if (!seen[bcol]) { /* new entry */ 1099 PetscInt *PETSC_RESTRICT slot; 1100 ierr = PetscSegBufferGetInts(segrow,1,&slot);CHKERRQ(ierr); 1101 *slot = bcol; 1102 seen[bcol] = 1; 1103 packlen++; 1104 } 1105 } 1106 } 1107 ierr = PetscSegBufferGetInts(seg,packlen,&crow);CHKERRQ(ierr); 1108 ierr = PetscSegBufferExtractTo(segrow,crow);CHKERRQ(ierr); 1109 ierr = PetscSortInt(packlen,crow);CHKERRQ(ierr); 1110 ci[i+1] = ci[i] + packlen; 1111 for (j=0; j<packlen; j++) seen[crow[j]] = 0; 1112 } 1113 ierr = PetscSegBufferDestroy(&segrow);CHKERRQ(ierr); 1114 ierr = PetscFree(seen);CHKERRQ(ierr); 1115 1116 /* Column indices are in the segmented buffer */ 1117 ierr = PetscSegBufferExtractAlloc(seg,&cj);CHKERRQ(ierr); 1118 ierr = PetscSegBufferDestroy(&seg);CHKERRQ(ierr); 1119 1120 /* put together the new symbolic matrix */ 1121 ierr = MatSetSeqAIJWithArrays_private(PetscObjectComm((PetscObject)A),am,bn,ci,cj,NULL,((PetscObject)A)->type_name,C);CHKERRQ(ierr); 1122 ierr = MatSetBlockSizesFromMats(C,A,B);CHKERRQ(ierr); 1123 1124 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 1125 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 1126 c = (Mat_SeqAIJ*)(C->data); 1127 c->free_a = PETSC_TRUE; 1128 c->free_ij = PETSC_TRUE; 1129 c->nonew = 0; 1130 1131 C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Sorted; 1132 1133 /* set MatInfo */ 1134 afill = (PetscReal)ci[am]/PetscMax(ai[am]+bi[bm],1) + 1.e-5; 1135 if (afill < 1.0) afill = 1.0; 1136 c->maxnz = ci[am]; 1137 c->nz = ci[am]; 1138 C->info.mallocs = ndouble; 1139 C->info.fill_ratio_given = fill; 1140 C->info.fill_ratio_needed = afill; 1141 1142 #if defined(PETSC_USE_INFO) 1143 if (ci[am]) { 1144 ierr = PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",ndouble,(double)fill,(double)afill);CHKERRQ(ierr); 1145 ierr = PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);CHKERRQ(ierr); 1146 } else { 1147 ierr = PetscInfo(C,"Empty matrix product\n");CHKERRQ(ierr); 1148 } 1149 #endif 1150 PetscFunctionReturn(0); 1151 } 1152 1153 PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(void *data) 1154 { 1155 PetscErrorCode ierr; 1156 Mat_MatMatTransMult *abt=(Mat_MatMatTransMult *)data; 1157 1158 PetscFunctionBegin; 1159 ierr = MatTransposeColoringDestroy(&abt->matcoloring);CHKERRQ(ierr); 1160 ierr = MatDestroy(&abt->Bt_den);CHKERRQ(ierr); 1161 ierr = MatDestroy(&abt->ABt_den);CHKERRQ(ierr); 1162 ierr = PetscFree(abt);CHKERRQ(ierr); 1163 PetscFunctionReturn(0); 1164 } 1165 1166 PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat C) 1167 { 1168 PetscErrorCode ierr; 1169 Mat Bt; 1170 PetscInt *bti,*btj; 1171 Mat_MatMatTransMult *abt; 1172 Mat_Product *product = C->product; 1173 char *alg; 1174 1175 PetscFunctionBegin; 1176 if (!product) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct"); 1177 if (product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty"); 1178 1179 /* create symbolic Bt */ 1180 ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 1181 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,NULL,&Bt);CHKERRQ(ierr); 1182 ierr = MatSetBlockSizes(Bt,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));CHKERRQ(ierr); 1183 ierr = MatSetType(Bt,((PetscObject)A)->type_name);CHKERRQ(ierr); 1184 1185 /* get symbolic C=A*Bt */ 1186 ierr = PetscStrallocpy(product->alg,&alg);CHKERRQ(ierr); 1187 ierr = MatProductSetAlgorithm(C,"sorted");CHKERRQ(ierr); /* set algorithm for C = A*Bt */ 1188 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr); 1189 ierr = MatProductSetAlgorithm(C,alg);CHKERRQ(ierr); /* resume original algorithm for ABt product */ 1190 ierr = PetscFree(alg);CHKERRQ(ierr); 1191 1192 /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */ 1193 ierr = PetscNew(&abt);CHKERRQ(ierr); 1194 1195 product->data = abt; 1196 product->destroy = MatDestroy_SeqAIJ_MatMatMultTrans; 1197 1198 C->ops->mattransposemultnumeric = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ; 1199 1200 abt->usecoloring = PETSC_FALSE; 1201 ierr = PetscStrcmp(product->alg,"color",&abt->usecoloring);CHKERRQ(ierr); 1202 if (abt->usecoloring) { 1203 /* Create MatTransposeColoring from symbolic C=A*B^T */ 1204 MatTransposeColoring matcoloring; 1205 MatColoring coloring; 1206 ISColoring iscoloring; 1207 Mat Bt_dense,C_dense; 1208 1209 /* inode causes memory problem */ 1210 ierr = MatSetOption(C,MAT_USE_INODES,PETSC_FALSE);CHKERRQ(ierr); 1211 1212 ierr = MatColoringCreate(C,&coloring);CHKERRQ(ierr); 1213 ierr = MatColoringSetDistance(coloring,2);CHKERRQ(ierr); 1214 ierr = MatColoringSetType(coloring,MATCOLORINGSL);CHKERRQ(ierr); 1215 ierr = MatColoringSetFromOptions(coloring);CHKERRQ(ierr); 1216 ierr = MatColoringApply(coloring,&iscoloring);CHKERRQ(ierr); 1217 ierr = MatColoringDestroy(&coloring);CHKERRQ(ierr); 1218 ierr = MatTransposeColoringCreate(C,iscoloring,&matcoloring);CHKERRQ(ierr); 1219 1220 abt->matcoloring = matcoloring; 1221 1222 ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 1223 1224 /* Create Bt_dense and C_dense = A*Bt_dense */ 1225 ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr); 1226 ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 1227 ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr); 1228 ierr = MatSeqDenseSetPreallocation(Bt_dense,NULL);CHKERRQ(ierr); 1229 1230 Bt_dense->assembled = PETSC_TRUE; 1231 abt->Bt_den = Bt_dense; 1232 1233 ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr); 1234 ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 1235 ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr); 1236 ierr = MatSeqDenseSetPreallocation(C_dense,NULL);CHKERRQ(ierr); 1237 1238 Bt_dense->assembled = PETSC_TRUE; 1239 abt->ABt_den = C_dense; 1240 1241 #if defined(PETSC_USE_INFO) 1242 { 1243 Mat_SeqAIJ *c = (Mat_SeqAIJ*)C->data; 1244 ierr = PetscInfo7(C,"Use coloring of C=A*B^T; B^T: %D %D, Bt_dense: %D,%D; Cnz %D / (cm*ncolors %D) = %g\n",B->cmap->n,B->rmap->n,Bt_dense->rmap->n,Bt_dense->cmap->n,c->nz,A->rmap->n*matcoloring->ncolors,(PetscReal)(c->nz)/(A->rmap->n*matcoloring->ncolors));CHKERRQ(ierr); 1245 } 1246 #endif 1247 } 1248 /* clean up */ 1249 ierr = MatDestroy(&Bt);CHKERRQ(ierr); 1250 ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 1251 PetscFunctionReturn(0); 1252 } 1253 1254 PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 1255 { 1256 PetscErrorCode ierr; 1257 Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 1258 PetscInt *ai =a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow; 1259 PetscInt cm =C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol; 1260 PetscLogDouble flops=0.0; 1261 MatScalar *aa =a->a,*aval,*ba=b->a,*bval,*ca,*cval; 1262 Mat_MatMatTransMult *abt; 1263 Mat_Product *product = C->product; 1264 1265 PetscFunctionBegin; 1266 if (!product) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct"); 1267 abt = (Mat_MatMatTransMult *)product->data; 1268 if (!abt) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct"); 1269 /* clear old values in C */ 1270 if (!c->a) { 1271 ierr = PetscCalloc1(ci[cm]+1,&ca);CHKERRQ(ierr); 1272 c->a = ca; 1273 c->free_a = PETSC_TRUE; 1274 } else { 1275 ca = c->a; 1276 ierr = PetscArrayzero(ca,ci[cm]+1);CHKERRQ(ierr); 1277 } 1278 1279 if (abt->usecoloring) { 1280 MatTransposeColoring matcoloring = abt->matcoloring; 1281 Mat Bt_dense,C_dense = abt->ABt_den; 1282 1283 /* Get Bt_dense by Apply MatTransposeColoring to B */ 1284 Bt_dense = abt->Bt_den; 1285 ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr); 1286 1287 /* C_dense = A*Bt_dense */ 1288 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr); 1289 1290 /* Recover C from C_dense */ 1291 ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr); 1292 PetscFunctionReturn(0); 1293 } 1294 1295 for (i=0; i<cm; i++) { 1296 anzi = ai[i+1] - ai[i]; 1297 acol = aj + ai[i]; 1298 aval = aa + ai[i]; 1299 cnzi = ci[i+1] - ci[i]; 1300 ccol = cj + ci[i]; 1301 cval = ca + ci[i]; 1302 for (j=0; j<cnzi; j++) { 1303 brow = ccol[j]; 1304 bnzj = bi[brow+1] - bi[brow]; 1305 bcol = bj + bi[brow]; 1306 bval = ba + bi[brow]; 1307 1308 /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 1309 nexta = 0; nextb = 0; 1310 while (nexta<anzi && nextb<bnzj) { 1311 while (nexta < anzi && acol[nexta] < bcol[nextb]) nexta++; 1312 if (nexta == anzi) break; 1313 while (nextb < bnzj && acol[nexta] > bcol[nextb]) nextb++; 1314 if (nextb == bnzj) break; 1315 if (acol[nexta] == bcol[nextb]) { 1316 cval[j] += aval[nexta]*bval[nextb]; 1317 nexta++; nextb++; 1318 flops += 2; 1319 } 1320 } 1321 } 1322 } 1323 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1324 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1325 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 1326 PetscFunctionReturn(0); 1327 } 1328 1329 PetscErrorCode MatDestroy_SeqAIJ_MatTransMatMult(void *data) 1330 { 1331 PetscErrorCode ierr; 1332 Mat_MatTransMatMult *atb = (Mat_MatTransMatMult*)data; 1333 1334 PetscFunctionBegin; 1335 ierr = MatDestroy(&atb->At);CHKERRQ(ierr); 1336 if (atb->destroy) { 1337 ierr = (*atb->destroy)(atb->data);CHKERRQ(ierr); 1338 } 1339 ierr = PetscFree(atb);CHKERRQ(ierr); 1340 PetscFunctionReturn(0); 1341 } 1342 1343 PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat C) 1344 { 1345 PetscErrorCode ierr; 1346 Mat At; 1347 PetscInt *ati,*atj; 1348 Mat_Product *product = C->product; 1349 MatProductAlgorithm alg; 1350 PetscBool flg; 1351 1352 PetscFunctionBegin; 1353 if (product) { 1354 alg = product->alg; 1355 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"!product, not supported yet"); 1356 1357 /* outerproduct */ 1358 ierr = PetscStrcmp(alg,"outerproduct",&flg);CHKERRQ(ierr); 1359 if (flg) { 1360 /* create symbolic At */ 1361 ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 1362 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,NULL,&At);CHKERRQ(ierr); 1363 ierr = MatSetBlockSizes(At,PetscAbs(A->cmap->bs),PetscAbs(B->cmap->bs));CHKERRQ(ierr); 1364 ierr = MatSetType(At,((PetscObject)A)->type_name);CHKERRQ(ierr); 1365 1366 /* get symbolic C=At*B */ 1367 ierr = MatProductSetAlgorithm(C,"sorted");CHKERRQ(ierr); 1368 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 1369 1370 /* clean up */ 1371 ierr = MatDestroy(&At);CHKERRQ(ierr); 1372 ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 1373 1374 C->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ; /* outerproduct */ 1375 ierr = MatProductSetAlgorithm(C,"outerproduct");CHKERRQ(ierr); 1376 PetscFunctionReturn(0); 1377 } 1378 1379 /* matmatmult */ 1380 ierr = PetscStrcmp(alg,"at*b",&flg);CHKERRQ(ierr); 1381 if (flg) { 1382 Mat_MatTransMatMult *atb; 1383 1384 if (product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty"); 1385 ierr = PetscNew(&atb);CHKERRQ(ierr); 1386 ierr = MatTranspose_SeqAIJ(A,MAT_INITIAL_MATRIX,&At);CHKERRQ(ierr); 1387 ierr = MatProductSetAlgorithm(C,"sorted");CHKERRQ(ierr); 1388 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 1389 ierr = MatProductSetAlgorithm(C,"at*b");CHKERRQ(ierr); 1390 product->data = atb; 1391 product->destroy = MatDestroy_SeqAIJ_MatTransMatMult; 1392 atb->At = At; 1393 atb->updateAt = PETSC_FALSE; /* because At is computed here */ 1394 1395 C->ops->mattransposemultnumeric = NULL; /* see MatProductNumeric_AtB_SeqAIJ_SeqAIJ */ 1396 PetscFunctionReturn(0); 1397 } 1398 1399 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Mat Product Algorithm is not supported"); 1400 PetscFunctionReturn(0); 1401 } 1402 1403 PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 1404 { 1405 PetscErrorCode ierr; 1406 Mat_SeqAIJ *a =(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 1407 PetscInt am =A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 1408 PetscInt cm =C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 1409 PetscLogDouble flops=0.0; 1410 MatScalar *aa =a->a,*ba,*ca,*caj; 1411 1412 PetscFunctionBegin; 1413 if (!c->a) { 1414 ierr = PetscCalloc1(ci[cm]+1,&ca);CHKERRQ(ierr); 1415 1416 c->a = ca; 1417 c->free_a = PETSC_TRUE; 1418 } else { 1419 ca = c->a; 1420 ierr = PetscArrayzero(ca,ci[cm]);CHKERRQ(ierr); 1421 } 1422 1423 /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 1424 for (i=0; i<am; i++) { 1425 bj = b->j + bi[i]; 1426 ba = b->a + bi[i]; 1427 bnzi = bi[i+1] - bi[i]; 1428 anzi = ai[i+1] - ai[i]; 1429 for (j=0; j<anzi; j++) { 1430 nextb = 0; 1431 crow = *aj++; 1432 cjj = cj + ci[crow]; 1433 caj = ca + ci[crow]; 1434 /* perform sparse axpy operation. Note cjj includes bj. */ 1435 for (k=0; nextb<bnzi; k++) { 1436 if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 1437 caj[k] += (*aa)*(*(ba+nextb)); 1438 nextb++; 1439 } 1440 } 1441 flops += 2*bnzi; 1442 aa++; 1443 } 1444 } 1445 1446 /* Assemble the final matrix and clean up */ 1447 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1448 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1449 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 1450 PetscFunctionReturn(0); 1451 } 1452 1453 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat C) 1454 { 1455 PetscErrorCode ierr; 1456 1457 PetscFunctionBegin; 1458 ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr); 1459 C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqDense; 1460 PetscFunctionReturn(0); 1461 } 1462 1463 PETSC_INTERN PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C,const PetscBool add) 1464 { 1465 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 1466 Mat_SeqDense *bd=(Mat_SeqDense*)B->data; 1467 Mat_SeqDense *cd=(Mat_SeqDense*)C->data; 1468 PetscErrorCode ierr; 1469 PetscScalar *c,r1,r2,r3,r4,*c1,*c2,*c3,*c4; 1470 const PetscScalar *aa,*b,*b1,*b2,*b3,*b4,*av; 1471 const PetscInt *aj; 1472 PetscInt cm=C->rmap->n,cn=B->cmap->n,bm=bd->lda,am=A->rmap->n; 1473 PetscInt clda=cd->lda; 1474 PetscInt am4=4*clda,bm4=4*bm,col,i,j,n; 1475 1476 PetscFunctionBegin; 1477 if (!cm || !cn) PetscFunctionReturn(0); 1478 ierr = MatSeqAIJGetArrayRead(A,&av);CHKERRQ(ierr); 1479 if (add) { 1480 ierr = MatDenseGetArray(C,&c);CHKERRQ(ierr); 1481 } else { 1482 ierr = MatDenseGetArrayWrite(C,&c);CHKERRQ(ierr); 1483 } 1484 ierr = MatDenseGetArrayRead(B,&b);CHKERRQ(ierr); 1485 b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 1486 c1 = c; c2 = c1 + clda; c3 = c2 + clda; c4 = c3 + clda; 1487 for (col=0; col<(cn/4)*4; col += 4) { /* over columns of C */ 1488 for (i=0; i<am; i++) { /* over rows of A in those columns */ 1489 r1 = r2 = r3 = r4 = 0.0; 1490 n = a->i[i+1] - a->i[i]; 1491 aj = a->j + a->i[i]; 1492 aa = av + a->i[i]; 1493 for (j=0; j<n; j++) { 1494 const PetscScalar aatmp = aa[j]; 1495 const PetscInt ajtmp = aj[j]; 1496 r1 += aatmp*b1[ajtmp]; 1497 r2 += aatmp*b2[ajtmp]; 1498 r3 += aatmp*b3[ajtmp]; 1499 r4 += aatmp*b4[ajtmp]; 1500 } 1501 if (add) { 1502 c1[i] += r1; 1503 c2[i] += r2; 1504 c3[i] += r3; 1505 c4[i] += r4; 1506 } else { 1507 c1[i] = r1; 1508 c2[i] = r2; 1509 c3[i] = r3; 1510 c4[i] = r4; 1511 } 1512 } 1513 b1 += bm4; b2 += bm4; b3 += bm4; b4 += bm4; 1514 c1 += am4; c2 += am4; c3 += am4; c4 += am4; 1515 } 1516 /* process remaining columns */ 1517 if (col != cn) { 1518 PetscInt rc = cn-col; 1519 1520 if (rc == 1) { 1521 for (i=0; i<am; i++) { 1522 r1 = 0.0; 1523 n = a->i[i+1] - a->i[i]; 1524 aj = a->j + a->i[i]; 1525 aa = av + a->i[i]; 1526 for (j=0; j<n; j++) r1 += aa[j]*b1[aj[j]]; 1527 if (add) c1[i] += r1; 1528 else c1[i] = r1; 1529 } 1530 } else if (rc == 2) { 1531 for (i=0; i<am; i++) { 1532 r1 = r2 = 0.0; 1533 n = a->i[i+1] - a->i[i]; 1534 aj = a->j + a->i[i]; 1535 aa = av + a->i[i]; 1536 for (j=0; j<n; j++) { 1537 const PetscScalar aatmp = aa[j]; 1538 const PetscInt ajtmp = aj[j]; 1539 r1 += aatmp*b1[ajtmp]; 1540 r2 += aatmp*b2[ajtmp]; 1541 } 1542 if (add) { 1543 c1[i] += r1; 1544 c2[i] += r2; 1545 } else { 1546 c1[i] = r1; 1547 c2[i] = r2; 1548 } 1549 } 1550 } else { 1551 for (i=0; i<am; i++) { 1552 r1 = r2 = r3 = 0.0; 1553 n = a->i[i+1] - a->i[i]; 1554 aj = a->j + a->i[i]; 1555 aa = av + a->i[i]; 1556 for (j=0; j<n; j++) { 1557 const PetscScalar aatmp = aa[j]; 1558 const PetscInt ajtmp = aj[j]; 1559 r1 += aatmp*b1[ajtmp]; 1560 r2 += aatmp*b2[ajtmp]; 1561 r3 += aatmp*b3[ajtmp]; 1562 } 1563 if (add) { 1564 c1[i] += r1; 1565 c2[i] += r2; 1566 c3[i] += r3; 1567 } else { 1568 c1[i] = r1; 1569 c2[i] = r2; 1570 c3[i] = r3; 1571 } 1572 } 1573 } 1574 } 1575 ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 1576 if (add) { 1577 ierr = MatDenseRestoreArray(C,&c);CHKERRQ(ierr); 1578 } else { 1579 ierr = MatDenseRestoreArrayWrite(C,&c);CHKERRQ(ierr); 1580 } 1581 ierr = MatDenseRestoreArrayRead(B,&b);CHKERRQ(ierr); 1582 ierr = MatSeqAIJRestoreArrayRead(A,&av);CHKERRQ(ierr); 1583 PetscFunctionReturn(0); 1584 } 1585 1586 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 1587 { 1588 PetscErrorCode ierr; 1589 1590 PetscFunctionBegin; 1591 if (B->rmap->n != 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,B->rmap->n); 1592 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); 1593 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); 1594 1595 ierr = MatMatMultNumericAdd_SeqAIJ_SeqDense(A,B,C,PETSC_FALSE);CHKERRQ(ierr); 1596 PetscFunctionReturn(0); 1597 } 1598 1599 /* ------------------------------------------------------- */ 1600 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense_AB(Mat C) 1601 { 1602 PetscFunctionBegin; 1603 C->ops->matmultsymbolic = MatMatMultSymbolic_SeqAIJ_SeqDense; 1604 C->ops->productsymbolic = MatProductSymbolic_AB; 1605 PetscFunctionReturn(0); 1606 } 1607 1608 PETSC_INTERN PetscErrorCode MatTMatTMultSymbolic_SeqAIJ_SeqDense(Mat,Mat,PetscReal,Mat); 1609 1610 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense_AtB(Mat C) 1611 { 1612 PetscFunctionBegin; 1613 C->ops->transposematmultsymbolic = MatTMatTMultSymbolic_SeqAIJ_SeqDense; 1614 C->ops->productsymbolic = MatProductSymbolic_AtB; 1615 PetscFunctionReturn(0); 1616 } 1617 1618 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense_ABt(Mat C) 1619 { 1620 PetscFunctionBegin; 1621 C->ops->mattransposemultsymbolic = MatTMatTMultSymbolic_SeqAIJ_SeqDense; 1622 C->ops->productsymbolic = MatProductSymbolic_ABt; 1623 PetscFunctionReturn(0); 1624 } 1625 1626 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqAIJ_SeqDense(Mat C) 1627 { 1628 PetscErrorCode ierr; 1629 Mat_Product *product = C->product; 1630 1631 PetscFunctionBegin; 1632 switch (product->type) { 1633 case MATPRODUCT_AB: 1634 ierr = MatProductSetFromOptions_SeqAIJ_SeqDense_AB(C);CHKERRQ(ierr); 1635 break; 1636 case MATPRODUCT_AtB: 1637 ierr = MatProductSetFromOptions_SeqAIJ_SeqDense_AtB(C);CHKERRQ(ierr); 1638 break; 1639 case MATPRODUCT_ABt: 1640 ierr = MatProductSetFromOptions_SeqAIJ_SeqDense_ABt(C);CHKERRQ(ierr); 1641 break; 1642 default: 1643 break; 1644 } 1645 PetscFunctionReturn(0); 1646 } 1647 /* ------------------------------------------------------- */ 1648 static PetscErrorCode MatProductSetFromOptions_SeqXBAIJ_SeqDense_AB(Mat C) 1649 { 1650 PetscErrorCode ierr; 1651 Mat_Product *product = C->product; 1652 Mat A = product->A; 1653 PetscBool baij; 1654 1655 PetscFunctionBegin; 1656 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQBAIJ,&baij);CHKERRQ(ierr); 1657 if (!baij) { /* A is seqsbaij */ 1658 PetscBool sbaij; 1659 ierr = PetscObjectTypeCompare((PetscObject)A,MATSEQSBAIJ,&sbaij);CHKERRQ(ierr); 1660 if (!sbaij) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"Mat must be either seqbaij or seqsbaij format"); 1661 1662 C->ops->matmultsymbolic = MatMatMultSymbolic_SeqSBAIJ_SeqDense; 1663 } else { /* A is seqbaij */ 1664 C->ops->matmultsymbolic = MatMatMultSymbolic_SeqBAIJ_SeqDense; 1665 } 1666 1667 C->ops->productsymbolic = MatProductSymbolic_AB; 1668 PetscFunctionReturn(0); 1669 } 1670 1671 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqXBAIJ_SeqDense(Mat C) 1672 { 1673 PetscErrorCode ierr; 1674 Mat_Product *product = C->product; 1675 1676 PetscFunctionBegin; 1677 MatCheckProduct(C,1); 1678 if (!product->A) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing A"); 1679 if (product->type == MATPRODUCT_AB || (product->type == MATPRODUCT_AtB && product->A->symmetric)) { 1680 ierr = MatProductSetFromOptions_SeqXBAIJ_SeqDense_AB(C);CHKERRQ(ierr); 1681 } 1682 PetscFunctionReturn(0); 1683 } 1684 1685 /* ------------------------------------------------------- */ 1686 static PetscErrorCode MatProductSetFromOptions_SeqDense_SeqAIJ_AB(Mat C) 1687 { 1688 PetscFunctionBegin; 1689 C->ops->matmultsymbolic = MatMatMultSymbolic_SeqDense_SeqAIJ; 1690 C->ops->productsymbolic = MatProductSymbolic_AB; 1691 PetscFunctionReturn(0); 1692 } 1693 1694 PETSC_INTERN PetscErrorCode MatProductSetFromOptions_SeqDense_SeqAIJ(Mat C) 1695 { 1696 PetscErrorCode ierr; 1697 Mat_Product *product = C->product; 1698 1699 PetscFunctionBegin; 1700 if (product->type == MATPRODUCT_AB) { 1701 ierr = MatProductSetFromOptions_SeqDense_SeqAIJ_AB(C);CHKERRQ(ierr); 1702 } 1703 PetscFunctionReturn(0); 1704 } 1705 /* ------------------------------------------------------- */ 1706 1707 PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense) 1708 { 1709 PetscErrorCode ierr; 1710 Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 1711 Mat_SeqDense *btdense = (Mat_SeqDense*)Btdense->data; 1712 PetscInt *bi = b->i,*bj=b->j; 1713 PetscInt m = Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns; 1714 MatScalar *btval,*btval_den,*ba=b->a; 1715 PetscInt *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors; 1716 1717 PetscFunctionBegin; 1718 btval_den=btdense->v; 1719 ierr = PetscArrayzero(btval_den,m*n);CHKERRQ(ierr); 1720 for (k=0; k<ncolors; k++) { 1721 ncolumns = coloring->ncolumns[k]; 1722 for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */ 1723 col = *(columns + colorforcol[k] + l); 1724 btcol = bj + bi[col]; 1725 btval = ba + bi[col]; 1726 anz = bi[col+1] - bi[col]; 1727 for (j=0; j<anz; j++) { 1728 brow = btcol[j]; 1729 btval_den[brow] = btval[j]; 1730 } 1731 } 1732 btval_den += m; 1733 } 1734 PetscFunctionReturn(0); 1735 } 1736 1737 PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp) 1738 { 1739 PetscErrorCode ierr; 1740 Mat_SeqAIJ *csp = (Mat_SeqAIJ*)Csp->data; 1741 const PetscScalar *ca_den,*ca_den_ptr; 1742 PetscScalar *ca=csp->a; 1743 PetscInt k,l,m=Cden->rmap->n,ncolors=matcoloring->ncolors; 1744 PetscInt brows=matcoloring->brows,*den2sp=matcoloring->den2sp; 1745 PetscInt nrows,*row,*idx; 1746 PetscInt *rows=matcoloring->rows,*colorforrow=matcoloring->colorforrow; 1747 1748 PetscFunctionBegin; 1749 ierr = MatDenseGetArrayRead(Cden,&ca_den);CHKERRQ(ierr); 1750 1751 if (brows > 0) { 1752 PetscInt *lstart,row_end,row_start; 1753 lstart = matcoloring->lstart; 1754 ierr = PetscArrayzero(lstart,ncolors);CHKERRQ(ierr); 1755 1756 row_end = brows; 1757 if (row_end > m) row_end = m; 1758 for (row_start=0; row_start<m; row_start+=brows) { /* loop over row blocks of Csp */ 1759 ca_den_ptr = ca_den; 1760 for (k=0; k<ncolors; k++) { /* loop over colors (columns of Cden) */ 1761 nrows = matcoloring->nrows[k]; 1762 row = rows + colorforrow[k]; 1763 idx = den2sp + colorforrow[k]; 1764 for (l=lstart[k]; l<nrows; l++) { 1765 if (row[l] >= row_end) { 1766 lstart[k] = l; 1767 break; 1768 } else { 1769 ca[idx[l]] = ca_den_ptr[row[l]]; 1770 } 1771 } 1772 ca_den_ptr += m; 1773 } 1774 row_end += brows; 1775 if (row_end > m) row_end = m; 1776 } 1777 } else { /* non-blocked impl: loop over columns of Csp - slow if Csp is large */ 1778 ca_den_ptr = ca_den; 1779 for (k=0; k<ncolors; k++) { 1780 nrows = matcoloring->nrows[k]; 1781 row = rows + colorforrow[k]; 1782 idx = den2sp + colorforrow[k]; 1783 for (l=0; l<nrows; l++) { 1784 ca[idx[l]] = ca_den_ptr[row[l]]; 1785 } 1786 ca_den_ptr += m; 1787 } 1788 } 1789 1790 ierr = MatDenseRestoreArrayRead(Cden,&ca_den);CHKERRQ(ierr); 1791 #if defined(PETSC_USE_INFO) 1792 if (matcoloring->brows > 0) { 1793 ierr = PetscInfo1(Csp,"Loop over %D row blocks for den2sp\n",brows);CHKERRQ(ierr); 1794 } else { 1795 ierr = PetscInfo(Csp,"Loop over colors/columns of Cden, inefficient for large sparse matrix product \n");CHKERRQ(ierr); 1796 } 1797 #endif 1798 PetscFunctionReturn(0); 1799 } 1800 1801 PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c) 1802 { 1803 PetscErrorCode ierr; 1804 PetscInt i,n,nrows,Nbs,j,k,m,ncols,col,cm; 1805 const PetscInt *is,*ci,*cj,*row_idx; 1806 PetscInt nis = iscoloring->n,*rowhit,bs = 1; 1807 IS *isa; 1808 Mat_SeqAIJ *csp = (Mat_SeqAIJ*)mat->data; 1809 PetscInt *colorforrow,*rows,*rows_i,*idxhit,*spidx,*den2sp,*den2sp_i; 1810 PetscInt *colorforcol,*columns,*columns_i,brows; 1811 PetscBool flg; 1812 1813 PetscFunctionBegin; 1814 ierr = ISColoringGetIS(iscoloring,PETSC_USE_POINTER,PETSC_IGNORE,&isa);CHKERRQ(ierr); 1815 1816 /* bs >1 is not being tested yet! */ 1817 Nbs = mat->cmap->N/bs; 1818 c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 1819 c->N = Nbs; 1820 c->m = c->M; 1821 c->rstart = 0; 1822 c->brows = 100; 1823 1824 c->ncolors = nis; 1825 ierr = PetscMalloc3(nis,&c->ncolumns,nis,&c->nrows,nis+1,&colorforrow);CHKERRQ(ierr); 1826 ierr = PetscMalloc1(csp->nz+1,&rows);CHKERRQ(ierr); 1827 ierr = PetscMalloc1(csp->nz+1,&den2sp);CHKERRQ(ierr); 1828 1829 brows = c->brows; 1830 ierr = PetscOptionsGetInt(NULL,NULL,"-matden2sp_brows",&brows,&flg);CHKERRQ(ierr); 1831 if (flg) c->brows = brows; 1832 if (brows > 0) { 1833 ierr = PetscMalloc1(nis+1,&c->lstart);CHKERRQ(ierr); 1834 } 1835 1836 colorforrow[0] = 0; 1837 rows_i = rows; 1838 den2sp_i = den2sp; 1839 1840 ierr = PetscMalloc1(nis+1,&colorforcol);CHKERRQ(ierr); 1841 ierr = PetscMalloc1(Nbs+1,&columns);CHKERRQ(ierr); 1842 1843 colorforcol[0] = 0; 1844 columns_i = columns; 1845 1846 /* get column-wise storage of mat */ 1847 ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 1848 1849 cm = c->m; 1850 ierr = PetscMalloc1(cm+1,&rowhit);CHKERRQ(ierr); 1851 ierr = PetscMalloc1(cm+1,&idxhit);CHKERRQ(ierr); 1852 for (i=0; i<nis; i++) { /* loop over color */ 1853 ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 1854 ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 1855 1856 c->ncolumns[i] = n; 1857 if (n) { 1858 ierr = PetscArraycpy(columns_i,is,n);CHKERRQ(ierr); 1859 } 1860 colorforcol[i+1] = colorforcol[i] + n; 1861 columns_i += n; 1862 1863 /* fast, crude version requires O(N*N) work */ 1864 ierr = PetscArrayzero(rowhit,cm);CHKERRQ(ierr); 1865 1866 for (j=0; j<n; j++) { /* loop over columns*/ 1867 col = is[j]; 1868 row_idx = cj + ci[col]; 1869 m = ci[col+1] - ci[col]; 1870 for (k=0; k<m; k++) { /* loop over columns marking them in rowhit */ 1871 idxhit[*row_idx] = spidx[ci[col] + k]; 1872 rowhit[*row_idx++] = col + 1; 1873 } 1874 } 1875 /* count the number of hits */ 1876 nrows = 0; 1877 for (j=0; j<cm; j++) { 1878 if (rowhit[j]) nrows++; 1879 } 1880 c->nrows[i] = nrows; 1881 colorforrow[i+1] = colorforrow[i] + nrows; 1882 1883 nrows = 0; 1884 for (j=0; j<cm; j++) { /* loop over rows */ 1885 if (rowhit[j]) { 1886 rows_i[nrows] = j; 1887 den2sp_i[nrows] = idxhit[j]; 1888 nrows++; 1889 } 1890 } 1891 den2sp_i += nrows; 1892 1893 ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 1894 rows_i += nrows; 1895 } 1896 ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,NULL);CHKERRQ(ierr); 1897 ierr = PetscFree(rowhit);CHKERRQ(ierr); 1898 ierr = ISColoringRestoreIS(iscoloring,PETSC_USE_POINTER,&isa);CHKERRQ(ierr); 1899 if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]); 1900 1901 c->colorforrow = colorforrow; 1902 c->rows = rows; 1903 c->den2sp = den2sp; 1904 c->colorforcol = colorforcol; 1905 c->columns = columns; 1906 1907 ierr = PetscFree(idxhit);CHKERRQ(ierr); 1908 PetscFunctionReturn(0); 1909 } 1910 1911 /* --------------------------------------------------------------- */ 1912 static PetscErrorCode MatProductNumeric_AtB_SeqAIJ_SeqAIJ(Mat C) 1913 { 1914 PetscErrorCode ierr; 1915 Mat_Product *product = C->product; 1916 Mat A=product->A,B=product->B; 1917 1918 PetscFunctionBegin; 1919 if (C->ops->mattransposemultnumeric) { 1920 /* Alg: "outerproduct" */ 1921 ierr = (*C->ops->mattransposemultnumeric)(A,B,C);CHKERRQ(ierr); 1922 } else { 1923 /* Alg: "matmatmult" -- C = At*B */ 1924 Mat_MatTransMatMult *atb = (Mat_MatTransMatMult *)product->data; 1925 Mat At; 1926 1927 if (!atb) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Missing product struct"); 1928 At = atb->At; 1929 if (atb->updateAt) { /* At is computed in MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ() */ 1930 ierr = MatTranspose_SeqAIJ(A,MAT_REUSE_MATRIX,&At);CHKERRQ(ierr); 1931 } 1932 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(At,B,C);CHKERRQ(ierr); 1933 atb->updateAt = PETSC_TRUE; 1934 } 1935 PetscFunctionReturn(0); 1936 } 1937 1938 static PetscErrorCode MatProductSymbolic_AtB_SeqAIJ_SeqAIJ(Mat C) 1939 { 1940 PetscErrorCode ierr; 1941 Mat_Product *product = C->product; 1942 Mat A=product->A,B=product->B; 1943 PetscReal fill=product->fill; 1944 1945 PetscFunctionBegin; 1946 ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 1947 1948 C->ops->productnumeric = MatProductNumeric_AtB_SeqAIJ_SeqAIJ; 1949 PetscFunctionReturn(0); 1950 } 1951 1952 /* --------------------------------------------------------------- */ 1953 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_AB(Mat C) 1954 { 1955 PetscErrorCode ierr; 1956 Mat_Product *product = C->product; 1957 PetscInt alg = 0; /* default algorithm */ 1958 PetscBool flg = PETSC_FALSE; 1959 #if !defined(PETSC_HAVE_HYPRE) 1960 const char *algTypes[7] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","rowmerge"}; 1961 PetscInt nalg = 7; 1962 #else 1963 const char *algTypes[8] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","rowmerge","hypre"}; 1964 PetscInt nalg = 8; 1965 #endif 1966 1967 PetscFunctionBegin; 1968 /* Set default algorithm */ 1969 ierr = PetscStrcmp(C->product->alg,"default",&flg);CHKERRQ(ierr); 1970 if (flg) { 1971 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 1972 } 1973 1974 /* Get runtime option */ 1975 if (product->api_user) { 1976 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMult","Mat");CHKERRQ(ierr); 1977 ierr = PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[0],&alg,&flg);CHKERRQ(ierr); 1978 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1979 } else { 1980 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AB","Mat");CHKERRQ(ierr); 1981 ierr = PetscOptionsEList("-matproduct_ab_via","Algorithmic approach","MatProduct_AB",algTypes,nalg,algTypes[0],&alg,&flg);CHKERRQ(ierr); 1982 ierr = PetscOptionsEnd();CHKERRQ(ierr); 1983 } 1984 if (flg) { 1985 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 1986 } 1987 1988 C->ops->productsymbolic = MatProductSymbolic_AB; 1989 C->ops->matmultsymbolic = MatMatMultSymbolic_SeqAIJ_SeqAIJ; 1990 PetscFunctionReturn(0); 1991 } 1992 1993 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_AtB(Mat C) 1994 { 1995 PetscErrorCode ierr; 1996 Mat_Product *product = C->product; 1997 PetscInt alg = 0; /* default algorithm */ 1998 PetscBool flg = PETSC_FALSE; 1999 const char *algTypes[2] = {"at*b","outerproduct"}; 2000 PetscInt nalg = 2; 2001 2002 PetscFunctionBegin; 2003 /* Set default algorithm */ 2004 ierr = PetscStrcmp(product->alg,"default",&flg);CHKERRQ(ierr); 2005 if (flg) { 2006 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2007 } 2008 2009 /* Get runtime option */ 2010 if (product->api_user) { 2011 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatTransposeMatMult","Mat");CHKERRQ(ierr); 2012 ierr = PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 2013 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2014 } else { 2015 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AtB","Mat");CHKERRQ(ierr); 2016 ierr = PetscOptionsEList("-matproduct_atb_via","Algorithmic approach","MatProduct_AtB",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 2017 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2018 } 2019 if (flg) { 2020 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2021 } 2022 2023 C->ops->productsymbolic = MatProductSymbolic_AtB_SeqAIJ_SeqAIJ; 2024 PetscFunctionReturn(0); 2025 } 2026 2027 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_ABt(Mat C) 2028 { 2029 PetscErrorCode ierr; 2030 Mat_Product *product = C->product; 2031 PetscInt alg = 0; /* default algorithm */ 2032 PetscBool flg = PETSC_FALSE; 2033 const char *algTypes[2] = {"default","color"}; 2034 PetscInt nalg = 2; 2035 2036 PetscFunctionBegin; 2037 /* Set default algorithm */ 2038 ierr = PetscStrcmp(C->product->alg,"default",&flg);CHKERRQ(ierr); 2039 if (!flg) { 2040 alg = 1; 2041 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2042 } 2043 2044 /* Get runtime option */ 2045 if (product->api_user) { 2046 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatTransposeMult","Mat");CHKERRQ(ierr); 2047 ierr = PetscOptionsEList("-matmattransmult_via","Algorithmic approach","MatMatTransposeMult",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 2048 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2049 } else { 2050 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_ABt","Mat");CHKERRQ(ierr); 2051 ierr = PetscOptionsEList("-matproduct_abt_via","Algorithmic approach","MatProduct_ABt",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 2052 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2053 } 2054 if (flg) { 2055 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2056 } 2057 2058 C->ops->mattransposemultsymbolic = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ; 2059 C->ops->productsymbolic = MatProductSymbolic_ABt; 2060 PetscFunctionReturn(0); 2061 } 2062 2063 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_PtAP(Mat C) 2064 { 2065 PetscErrorCode ierr; 2066 Mat_Product *product = C->product; 2067 PetscBool flg = PETSC_FALSE; 2068 PetscInt alg = 0; /* default algorithm -- alg=1 should be default!!! */ 2069 #if !defined(PETSC_HAVE_HYPRE) 2070 const char *algTypes[2] = {"scalable","rap"}; 2071 PetscInt nalg = 2; 2072 #else 2073 const char *algTypes[3] = {"scalable","rap","hypre"}; 2074 PetscInt nalg = 3; 2075 #endif 2076 2077 PetscFunctionBegin; 2078 /* Set default algorithm */ 2079 ierr = PetscStrcmp(product->alg,"default",&flg);CHKERRQ(ierr); 2080 if (flg) { 2081 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2082 } 2083 2084 /* Get runtime option */ 2085 if (product->api_user) { 2086 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatPtAP","Mat");CHKERRQ(ierr); 2087 ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[0],&alg,&flg);CHKERRQ(ierr); 2088 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2089 } else { 2090 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_PtAP","Mat");CHKERRQ(ierr); 2091 ierr = PetscOptionsEList("-matproduct_ptap_via","Algorithmic approach","MatProduct_PtAP",algTypes,nalg,algTypes[0],&alg,&flg);CHKERRQ(ierr); 2092 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2093 } 2094 if (flg) { 2095 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2096 } 2097 2098 C->ops->productsymbolic = MatProductSymbolic_PtAP_SeqAIJ_SeqAIJ; 2099 PetscFunctionReturn(0); 2100 } 2101 2102 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_RARt(Mat C) 2103 { 2104 PetscErrorCode ierr; 2105 Mat_Product *product = C->product; 2106 PetscBool flg = PETSC_FALSE; 2107 PetscInt alg = 0; /* default algorithm */ 2108 const char *algTypes[3] = {"r*a*rt","r*art","coloring_rart"}; 2109 PetscInt nalg = 3; 2110 2111 PetscFunctionBegin; 2112 /* Set default algorithm */ 2113 ierr = PetscStrcmp(product->alg,"default",&flg);CHKERRQ(ierr); 2114 if (flg) { 2115 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2116 } 2117 2118 /* Get runtime option */ 2119 if (product->api_user) { 2120 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatRARt","Mat");CHKERRQ(ierr); 2121 ierr = PetscOptionsEList("-matrart_via","Algorithmic approach","MatRARt",algTypes,nalg,algTypes[0],&alg,&flg);CHKERRQ(ierr); 2122 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2123 } else { 2124 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_RARt","Mat");CHKERRQ(ierr); 2125 ierr = PetscOptionsEList("-matproduct_rart_via","Algorithmic approach","MatProduct_RARt",algTypes,nalg,algTypes[0],&alg,&flg);CHKERRQ(ierr); 2126 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2127 } 2128 if (flg) { 2129 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2130 } 2131 2132 C->ops->productsymbolic = MatProductSymbolic_RARt_SeqAIJ_SeqAIJ; 2133 PetscFunctionReturn(0); 2134 } 2135 2136 /* ABC = A*B*C = A*(B*C); ABC's algorithm must be chosen from AB's algorithm */ 2137 static PetscErrorCode MatProductSetFromOptions_SeqAIJ_ABC(Mat C) 2138 { 2139 PetscErrorCode ierr; 2140 Mat_Product *product = C->product; 2141 PetscInt alg = 0; /* default algorithm */ 2142 PetscBool flg = PETSC_FALSE; 2143 const char *algTypes[7] = {"sorted","scalable","scalable_fast","heap","btheap","llcondensed","rowmerge"}; 2144 PetscInt nalg = 7; 2145 2146 PetscFunctionBegin; 2147 /* Set default algorithm */ 2148 ierr = PetscStrcmp(product->alg,"default",&flg);CHKERRQ(ierr); 2149 if (flg) { 2150 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2151 } 2152 2153 /* Get runtime option */ 2154 if (product->api_user) { 2155 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMatMult","Mat");CHKERRQ(ierr); 2156 ierr = PetscOptionsEList("-matmatmatmult_via","Algorithmic approach","MatMatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 2157 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2158 } else { 2159 ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_ABC","Mat");CHKERRQ(ierr); 2160 ierr = PetscOptionsEList("-matproduct_abc_via","Algorithmic approach","MatProduct_ABC",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 2161 ierr = PetscOptionsEnd();CHKERRQ(ierr); 2162 } 2163 if (flg) { 2164 ierr = MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);CHKERRQ(ierr); 2165 } 2166 2167 C->ops->matmatmultsymbolic = MatMatMatMultSymbolic_SeqAIJ_SeqAIJ_SeqAIJ; 2168 C->ops->productsymbolic = MatProductSymbolic_ABC; 2169 PetscFunctionReturn(0); 2170 } 2171 2172 PetscErrorCode MatProductSetFromOptions_SeqAIJ(Mat C) 2173 { 2174 PetscErrorCode ierr; 2175 Mat_Product *product = C->product; 2176 2177 PetscFunctionBegin; 2178 switch (product->type) { 2179 case MATPRODUCT_AB: 2180 ierr = MatProductSetFromOptions_SeqAIJ_AB(C);CHKERRQ(ierr); 2181 break; 2182 case MATPRODUCT_AtB: 2183 ierr = MatProductSetFromOptions_SeqAIJ_AtB(C);CHKERRQ(ierr); 2184 break; 2185 case MATPRODUCT_ABt: 2186 ierr = MatProductSetFromOptions_SeqAIJ_ABt(C);CHKERRQ(ierr); 2187 break; 2188 case MATPRODUCT_PtAP: 2189 ierr = MatProductSetFromOptions_SeqAIJ_PtAP(C);CHKERRQ(ierr); 2190 break; 2191 case MATPRODUCT_RARt: 2192 ierr = MatProductSetFromOptions_SeqAIJ_RARt(C);CHKERRQ(ierr); 2193 break; 2194 case MATPRODUCT_ABC: 2195 ierr = MatProductSetFromOptions_SeqAIJ_ABC(C);CHKERRQ(ierr); 2196 break; 2197 default: 2198 break; 2199 } 2200 PetscFunctionReturn(0); 2201 } 2202