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 <../src/mat/impls/dense/seq/dense.h> /*I "petscmat.h" I*/ 11 /* 12 #define DEBUG_MATMATMULT 13 */ 14 EXTERN_C_BEGIN 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; 21 22 PetscFunctionBegin; 23 if (scall == MAT_INITIAL_MATRIX){ 24 ierr = PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 25 ierr = PetscOptionsGetBool(PETSC_NULL,"-matmatmult_scalable",&scalable,PETSC_NULL);CHKERRQ(ierr); 26 if (scalable){ 27 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(A,B,fill,C);CHKERRQ(ierr); 28 } else { 29 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 30 } 31 ierr = PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);CHKERRQ(ierr); 32 } 33 ierr = PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 34 ierr = (*(*C)->ops->matmultnumeric)(A,B,*C);CHKERRQ(ierr); 35 ierr = PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);CHKERRQ(ierr); 36 PetscFunctionReturn(0); 37 } 38 EXTERN_C_END 39 40 /* 41 MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ - Get symbolic structure of C=A*B 42 Input Parameter: 43 . am, Ai, Aj - number of rows and structure of A 44 . bm, bn, Bi, Bj - number of rows, columns, and structure of B 45 . fill - filll ratio See MatMatMult() 46 47 Output Parameter: 48 . Ci, Cj - structure of C = A*B 49 . nspacedouble - number of extra mallocs 50 */ 51 #undef __FUNCT__ 52 #define __FUNCT__ "MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ" 53 PetscErrorCode MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(PetscInt am,PetscInt *Ai,PetscInt *Aj,PetscInt bm,PetscInt bn,PetscInt *Bi,PetscInt *Bj,PetscReal fill,PetscInt *Ci[],PetscInt *Cj[],PetscInt *nspacedouble) 54 { 55 PetscErrorCode ierr; 56 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 57 PetscInt *ai=Ai,*aj=Aj,*bi=Bi,*bj=Bj,*bjj,*ci,*cj; 58 PetscInt i,j,anzi,brow,bnzj,cnzi,nlnk,*lnk,ndouble=0; 59 PetscBT lnkbt; 60 61 PetscFunctionBegin; 62 /* Allocate ci array, arrays for fill computation and */ 63 /* free space for accumulating nonzero column info */ 64 ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 65 ci[0] = 0; 66 67 /* create and initialize a linked list */ 68 nlnk = bn+1; 69 ierr = PetscLLCreate(bn,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 70 71 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 72 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 73 current_space = free_space; 74 75 /* Determine symbolic info for each row of the product: */ 76 for (i=0; i<am; i++) { 77 anzi = ai[i+1] - ai[i]; 78 cnzi = 0; 79 aj = Aj + ai[i]; 80 for (j=0; j<anzi; j++){ 81 brow = aj[j]; 82 bnzj = bi[brow+1] - bi[brow]; 83 bjj = bj + bi[brow]; 84 /* add non-zero cols of B into the sorted linked list lnk */ 85 ierr = PetscLLAddSorted(bnzj,bjj,bn,nlnk,lnk,lnkbt);CHKERRQ(ierr); 86 cnzi += nlnk; 87 } 88 89 /* If free space is not available, make more free space */ 90 /* Double the amount of total space in the list */ 91 if (current_space->local_remaining<cnzi) { 92 ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 93 ndouble++; 94 } 95 96 /* Copy data into free space, then initialize lnk */ 97 ierr = PetscLLClean(bn,bn,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 98 current_space->array += cnzi; 99 current_space->local_used += cnzi; 100 current_space->local_remaining -= cnzi; 101 ci[i+1] = ci[i] + cnzi; 102 } 103 104 /* Column indices are in the list of free space */ 105 /* Allocate space for cj, initialize cj, and */ 106 /* destroy list of free space and other temporary array(s) */ 107 ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 108 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 109 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 110 111 *Ci = ci; 112 *Cj = cj; 113 *nspacedouble = ndouble; 114 PetscFunctionReturn(0); 115 } 116 117 /* 118 MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable - same as MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ() 119 but replaces O(bn) array 'lnkbt' with a scalable array 'abj' of size Crmax. 120 */ 121 #undef __FUNCT__ 122 #define __FUNCT__ "MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable" 123 PetscErrorCode MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable(PetscInt am,PetscInt *Ai,PetscInt *Aj,PetscInt bm,PetscInt bn,PetscInt *Bi,PetscInt *Bj,PetscReal fill,PetscInt *Ci[],PetscInt *Cj[],PetscInt *nspacedouble) 124 { 125 PetscErrorCode ierr; 126 PetscInt *aj=Aj,*bi=Bi,*bj,*ci,*cj,rmax=0,*abj,*cj_tmp,nextabj; 127 PetscInt i,j,anzi,brow,bnzj,cnzi,k; 128 PetscBT bt; 129 130 PetscFunctionBegin; 131 /* Allocate ci array and bt */ 132 ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 133 ci[0] = 0; 134 135 /* Get ci and count rmax of C */ 136 ierr = PetscBTCreate(bn,bt);CHKERRQ(ierr); 137 ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr); 138 for (i=0; i<am; i++) { 139 anzi = Ai[i+1] - Ai[i]; 140 cnzi = 0; 141 aj = Aj + Ai[i]; 142 for (j=0; j<anzi; j++){ 143 brow = aj[j]; 144 bnzj = bi[brow+1] - bi[brow]; 145 bj = Bj + bi[brow]; 146 for (k=0; k<bnzj; k++){ 147 if (!PetscBTLookupSet(bt,bj[k])) cnzi++; /* new entry */ 148 } 149 } 150 ierr = PetscBTMemzero(bn,bt);CHKERRQ(ierr); 151 ci[i+1] = ci[i] + cnzi; 152 if (rmax < cnzi) rmax = cnzi; 153 } 154 155 /* Allocate space for cj */ 156 ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 157 158 /* Allocate a temp array for storing column indices of A*B */ 159 ierr = PetscMalloc((rmax+1)*sizeof(PetscInt),&abj);CHKERRQ(ierr); 160 161 /* Determine cj */ 162 for (i=0; i<am; i++) { 163 anzi = Ai[i+1] - Ai[i]; 164 cnzi = 0; 165 nextabj = 0; 166 aj = Aj + Ai[i]; 167 for (j=0; j<anzi; j++){ 168 brow = aj[j]; 169 bnzj = bi[brow+1] - bi[brow]; 170 bj = Bj + bi[brow]; 171 for (k=0; k<bnzj; k++){ 172 if (!PetscBTLookupSet(bt,bj[k])){ /* new entry */ 173 abj[nextabj] = bj[k]; nextabj++; 174 } 175 } 176 } 177 178 /* Sort abj, then copy it to cj */ 179 cnzi = ci[i+1] - ci[i]; 180 ierr = PetscSortInt(cnzi,abj);CHKERRQ(ierr); 181 cj_tmp = cj + ci[i]; 182 for (k=0; k< cnzi; k++){ 183 cj_tmp[k] = abj[k]; 184 ierr = PetscBTClear(bt,abj[k]);CHKERRQ(ierr); 185 } 186 } 187 188 ierr = PetscBTDestroy(bt);CHKERRQ(ierr); 189 ierr = PetscFree(abj);CHKERRQ(ierr); 190 191 *Ci = ci; 192 *Cj = cj; 193 *nspacedouble = 0; 194 PetscFunctionReturn(0); 195 } 196 197 #undef __FUNCT__ 198 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ" 199 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 200 { 201 PetscErrorCode ierr; 202 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 203 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj; 204 PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N,nspacedouble; 205 MatScalar *ca; 206 PetscReal afill; 207 208 PetscFunctionBegin; 209 /* Get ci and cj */ 210 ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ(am,ai,aj,bm,bn,bi,bj,fill,&ci,&cj,&nspacedouble);CHKERRQ(ierr); 211 #if defined(DEBUG_MATMATMULT) 212 ierr = PetscPrintf(PETSC_COMM_SELF,"MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ() is done \n");CHKERRQ(ierr); 213 #endif 214 215 /* Allocate space for ca */ 216 ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 217 ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 218 219 /* put together the new symbolic matrix */ 220 ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr); 221 222 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 223 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 224 c = (Mat_SeqAIJ *)((*C)->data); 225 c->free_a = PETSC_TRUE; 226 c->free_ij = PETSC_TRUE; 227 c->nonew = 0; 228 (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqAIJ; 229 230 ierr = PetscMalloc(bn*sizeof(PetscScalar),&c->matmult_abdense);CHKERRQ(ierr); 231 ierr = PetscMemzero(c->matmult_abdense,bn*sizeof(PetscScalar));CHKERRQ(ierr); 232 (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; /* fast, takes additional dense_axpy*bn*sizeof(PetscScalar) space */ 233 234 /* set MatInfo */ 235 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 236 if (afill < 1.0) afill = 1.0; 237 c->maxnz = ci[am]; 238 c->nz = ci[am]; 239 (*C)->info.mallocs = nspacedouble; 240 (*C)->info.fill_ratio_given = fill; 241 (*C)->info.fill_ratio_needed = afill; 242 243 #if defined(PETSC_USE_INFO) 244 if (ci[am]) { 245 ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 246 ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 247 } else { 248 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 249 } 250 #endif 251 PetscFunctionReturn(0); 252 } 253 254 #undef __FUNCT__ 255 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ" 256 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 257 { 258 PetscErrorCode ierr; 259 PetscLogDouble flops=0.0; 260 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 261 Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 262 Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 263 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 264 PetscInt am=A->rmap->n,cm=C->rmap->n; 265 PetscInt i,j,k,anzi,bnzi,cnzi,brow; 266 PetscScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp; 267 PetscScalar *ab_dense=c->matmult_abdense; 268 269 PetscFunctionBegin; 270 /* clean old values in C */ 271 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 272 /* Traverse A row-wise. */ 273 /* Build the ith row in C by summing over nonzero columns in A, */ 274 /* the rows of B corresponding to nonzeros of A. */ 275 for (i=0; i<am; i++) { 276 anzi = ai[i+1] - ai[i]; 277 for (j=0; j<anzi; j++) { 278 brow = aj[j]; 279 bnzi = bi[brow+1] - bi[brow]; 280 bjj = bj + bi[brow]; 281 baj = ba + bi[brow]; 282 /* perform dense axpy */ 283 valtmp = aa[j]; 284 for (k=0; k<bnzi; k++) { 285 ab_dense[bjj[k]] += valtmp*baj[k]; 286 } 287 flops += 2*bnzi; 288 } 289 aj += anzi; aa += anzi; 290 291 cnzi = ci[i+1] - ci[i]; 292 for (k=0; k<cnzi; k++) { 293 ca[k] += ab_dense[cj[k]]; 294 ab_dense[cj[k]] = 0.0; /* zero ab_dense */ 295 } 296 flops += cnzi; 297 cj += cnzi; ca += cnzi; 298 } 299 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 300 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 301 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 302 C->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ; 303 PetscFunctionReturn(0); 304 } 305 306 #undef __FUNCT__ 307 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable" 308 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,Mat C) 309 { 310 PetscErrorCode ierr; 311 PetscLogDouble flops=0.0; 312 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 313 Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 314 Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 315 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 316 PetscInt am=A->rmap->N,cm=C->rmap->N; 317 PetscInt i,j,k,anzi,bnzi,cnzi,brow; 318 PetscScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,valtmp; 319 PetscInt nextb; 320 321 PetscFunctionBegin; 322 #if defined(DEBUG_MATMATMULT) 323 //ierr = PetscPrintf(PETSC_COMM_SELF,"MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable...\n");CHKERRQ(ierr); 324 #endif 325 /* clean old values in C */ 326 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 327 /* Traverse A row-wise. */ 328 /* Build the ith row in C by summing over nonzero columns in A, */ 329 /* the rows of B corresponding to nonzeros of A. */ 330 for (i=0;i<am;i++) { 331 anzi = ai[i+1] - ai[i]; 332 cnzi = ci[i+1] - ci[i]; 333 for (j=0;j<anzi;j++) { 334 brow = aj[j]; 335 bnzi = bi[brow+1] - bi[brow]; 336 bjj = bj + bi[brow]; 337 baj = ba + bi[brow]; 338 /* perform sparse axpy */ 339 valtmp = aa[j]; 340 nextb = 0; 341 for (k=0; nextb<bnzi; k++) { 342 if (cj[k] == bjj[nextb]){ /* ccol == bcol */ 343 ca[k] += valtmp*baj[nextb++]; 344 } 345 } 346 flops += 2*bnzi; 347 } 348 aj += anzi; aa += anzi; 349 cj += cnzi; ca += cnzi; 350 } 351 352 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 353 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 354 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 355 PetscFunctionReturn(0); 356 } 357 358 #undef __FUNCT__ 359 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable" 360 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable(Mat A,Mat B,PetscReal fill,Mat *C) 361 { 362 PetscErrorCode ierr; 363 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 364 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj; 365 PetscInt am=A->rmap->N,bn=B->cmap->N,bm=B->rmap->N,nspacedouble; 366 MatScalar *ca; 367 PetscReal afill; 368 369 PetscFunctionBegin; 370 #if defined(DEBUG_MATMATMULT) 371 ierr = PetscPrintf(PETSC_COMM_SELF,"MatMatMultSymbolic_SeqAIJ_SeqAIJ_Scalable \n");CHKERRQ(ierr); 372 #endif 373 /* Get ci and cj */ 374 ierr = MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable(am,ai,aj,bm,bn,bi,bj,fill,&ci,&cj,&nspacedouble);CHKERRQ(ierr); 375 #if defined(DEBUG_MATMATMULT) 376 ierr = PetscPrintf(PETSC_COMM_SELF,"MatGetSymbolicMatMatMult_SeqAIJ_SeqAIJ_Scalable() is done \n");CHKERRQ(ierr); 377 #endif 378 379 /* Allocate space for ca */ 380 ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 381 ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 382 383 /* put together the new symbolic matrix */ 384 ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr); 385 386 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 387 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 388 c = (Mat_SeqAIJ *)((*C)->data); 389 c->free_a = PETSC_TRUE; 390 c->free_ij = PETSC_TRUE; 391 c->nonew = 0; 392 (*C)->ops->matmultnumeric = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable; /* slower, less memory */ 393 394 /* set MatInfo */ 395 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 396 if (afill < 1.0) afill = 1.0; 397 c->maxnz = ci[am]; 398 c->nz = ci[am]; 399 (*C)->info.mallocs = nspacedouble; 400 (*C)->info.fill_ratio_given = fill; 401 (*C)->info.fill_ratio_needed = afill; 402 403 #if defined(PETSC_USE_INFO) 404 if (ci[am]) { 405 ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 406 ierr = PetscInfo1((*C),"Use MatMatMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 407 } else { 408 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 409 } 410 #endif 411 PetscFunctionReturn(0); 412 } 413 414 415 /* This routine is not used. Should be removed! */ 416 #undef __FUNCT__ 417 #define __FUNCT__ "MatMatTransposeMult_SeqAIJ_SeqAIJ" 418 PetscErrorCode MatMatTransposeMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 419 { 420 PetscErrorCode ierr; 421 422 PetscFunctionBegin; 423 if (scall == MAT_INITIAL_MATRIX){ 424 ierr = MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 425 } 426 ierr = MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 427 PetscFunctionReturn(0); 428 } 429 430 #undef __FUNCT__ 431 #define __FUNCT__ "PetscContainerDestroy_Mat_MatMatTransMult" 432 PetscErrorCode PetscContainerDestroy_Mat_MatMatTransMult(void *ptr) 433 { 434 PetscErrorCode ierr; 435 Mat_MatMatTransMult *multtrans=(Mat_MatMatTransMult*)ptr; 436 437 PetscFunctionBegin; 438 ierr = MatTransposeColoringDestroy(&multtrans->matcoloring);CHKERRQ(ierr); 439 ierr = MatDestroy(&multtrans->Bt_den);CHKERRQ(ierr); 440 ierr = MatDestroy(&multtrans->ABt_den);CHKERRQ(ierr); 441 ierr = PetscFree(multtrans);CHKERRQ(ierr); 442 PetscFunctionReturn(0); 443 } 444 445 #undef __FUNCT__ 446 #define __FUNCT__ "MatDestroy_SeqAIJ_MatMatMultTrans" 447 PetscErrorCode MatDestroy_SeqAIJ_MatMatMultTrans(Mat A) 448 { 449 PetscErrorCode ierr; 450 PetscContainer container; 451 Mat_MatMatTransMult *multtrans=PETSC_NULL; 452 453 PetscFunctionBegin; 454 ierr = PetscObjectQuery((PetscObject)A,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr); 455 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 456 ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr); 457 A->ops->destroy = multtrans->destroy; 458 if (A->ops->destroy) { 459 ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 460 } 461 ierr = PetscObjectCompose((PetscObject)A,"Mat_MatMatTransMult",0);CHKERRQ(ierr); 462 PetscFunctionReturn(0); 463 } 464 465 #undef __FUNCT__ 466 #define __FUNCT__ "MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ" 467 PetscErrorCode MatMatTransposeMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 468 { 469 PetscErrorCode ierr; 470 Mat Bt; 471 PetscInt *bti,*btj; 472 Mat_MatMatTransMult *multtrans; 473 PetscContainer container; 474 PetscLogDouble t0,tf,etime2=0.0; 475 476 PetscFunctionBegin; 477 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 478 /* create symbolic Bt */ 479 ierr = MatGetSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 480 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,B->cmap->n,B->rmap->n,bti,btj,PETSC_NULL,&Bt);CHKERRQ(ierr); 481 482 /* get symbolic C=A*Bt */ 483 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,Bt,fill,C);CHKERRQ(ierr); 484 485 /* create a supporting struct for reuse intermidiate dense matrices with matcoloring */ 486 ierr = PetscNew(Mat_MatMatTransMult,&multtrans);CHKERRQ(ierr); 487 488 /* attach the supporting struct to C */ 489 ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 490 ierr = PetscContainerSetPointer(container,multtrans);CHKERRQ(ierr); 491 ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatTransMult);CHKERRQ(ierr); 492 ierr = PetscObjectCompose((PetscObject)(*C),"Mat_MatMatTransMult",(PetscObject)container);CHKERRQ(ierr); 493 ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 494 495 multtrans->usecoloring = PETSC_FALSE; 496 multtrans->destroy = (*C)->ops->destroy; 497 (*C)->ops->destroy = MatDestroy_SeqAIJ_MatMatMultTrans; 498 499 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 500 etime2 += tf - t0; 501 502 ierr = PetscOptionsGetBool(PETSC_NULL,"-matmattransmult_color",&multtrans->usecoloring,PETSC_NULL);CHKERRQ(ierr); 503 if (multtrans->usecoloring){ 504 /* Create MatTransposeColoring from symbolic C=A*B^T */ 505 MatTransposeColoring matcoloring; 506 ISColoring iscoloring; 507 Mat Bt_dense,C_dense; 508 PetscLogDouble etime0=0.0,etime01=0.0,etime1=0.0; 509 510 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 511 ierr = MatGetColoring(*C,MATCOLORINGLF,&iscoloring);CHKERRQ(ierr); 512 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 513 etime0 += tf - t0; 514 515 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 516 ierr = MatTransposeColoringCreate(*C,iscoloring,&matcoloring);CHKERRQ(ierr); 517 multtrans->matcoloring = matcoloring; 518 ierr = ISColoringDestroy(&iscoloring);CHKERRQ(ierr); 519 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 520 etime01 += tf - t0; 521 522 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 523 /* Create Bt_dense and C_dense = A*Bt_dense */ 524 ierr = MatCreate(PETSC_COMM_SELF,&Bt_dense);CHKERRQ(ierr); 525 ierr = MatSetSizes(Bt_dense,A->cmap->n,matcoloring->ncolors,A->cmap->n,matcoloring->ncolors);CHKERRQ(ierr); 526 ierr = MatSetType(Bt_dense,MATSEQDENSE);CHKERRQ(ierr); 527 ierr = MatSeqDenseSetPreallocation(Bt_dense,PETSC_NULL);CHKERRQ(ierr); 528 Bt_dense->assembled = PETSC_TRUE; 529 multtrans->Bt_den = Bt_dense; 530 531 ierr = MatCreate(PETSC_COMM_SELF,&C_dense);CHKERRQ(ierr); 532 ierr = MatSetSizes(C_dense,A->rmap->n,matcoloring->ncolors,A->rmap->n,matcoloring->ncolors);CHKERRQ(ierr); 533 ierr = MatSetType(C_dense,MATSEQDENSE);CHKERRQ(ierr); 534 ierr = MatSeqDenseSetPreallocation(C_dense,PETSC_NULL);CHKERRQ(ierr); 535 Bt_dense->assembled = PETSC_TRUE; 536 multtrans->ABt_den = C_dense; 537 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 538 etime1 += tf - t0; 539 540 #if defined(PETSC_USE_INFO) 541 { 542 Mat_SeqAIJ *c=(Mat_SeqAIJ*)(*C)->data; 543 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)); 544 ierr = PetscInfo5(*C,"Sym = GetColor %g + ColorCreate %g + MatDenseCreate %g + non-colorSym %g = %g\n",etime0,etime01,etime1,etime2,etime0+etime01+etime1+etime2); 545 } 546 #endif 547 } 548 /* clean up */ 549 ierr = MatDestroy(&Bt);CHKERRQ(ierr); 550 ierr = MatRestoreSymbolicTranspose_SeqAIJ(B,&bti,&btj);CHKERRQ(ierr); 551 552 553 554 #if defined(INEFFICIENT_ALGORITHM) 555 /* The algorithm below computes am*bm sparse inner-product - inefficient! It will be deleted later. */ 556 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 557 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 558 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ci,*cj,*acol,*bcol; 559 PetscInt am=A->rmap->N,bm=B->rmap->N; 560 PetscInt i,j,anzi,bnzj,cnzi,nlnk,*lnk,nspacedouble=0,ka,kb,index[1]; 561 MatScalar *ca; 562 PetscBT lnkbt; 563 PetscReal afill; 564 565 /* Allocate row pointer array ci */ 566 ierr = PetscMalloc(((am+1)+1)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 567 ci[0] = 0; 568 569 /* Create and initialize a linked list for C columns */ 570 nlnk = bm+1; 571 ierr = PetscLLCreate(bm,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 572 573 /* Initial FreeSpace with size fill*(nnz(A)+nnz(B)) */ 574 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 575 current_space = free_space; 576 577 /* Determine symbolic info for each row of the product A*B^T: */ 578 for (i=0; i<am; i++) { 579 anzi = ai[i+1] - ai[i]; 580 cnzi = 0; 581 acol = aj + ai[i]; 582 for (j=0; j<bm; j++){ 583 bnzj = bi[j+1] - bi[j]; 584 bcol= bj + bi[j]; 585 /* sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 586 ka = 0; kb = 0; 587 while (ka < anzi && kb < bnzj){ 588 while (acol[ka] < bcol[kb] && ka < anzi) ka++; 589 if (ka == anzi) break; 590 while (acol[ka] > bcol[kb] && kb < bnzj) kb++; 591 if (kb == bnzj) break; 592 if (acol[ka] == bcol[kb]){ /* add nonzero c(i,j) to lnk */ 593 index[0] = j; 594 ierr = PetscLLAdd(1,index,bm,nlnk,lnk,lnkbt);CHKERRQ(ierr); 595 cnzi++; 596 break; 597 } 598 } 599 } 600 601 /* If free space is not available, make more free space */ 602 /* Double the amount of total space in the list */ 603 if (current_space->local_remaining<cnzi) { 604 ierr = PetscFreeSpaceGet(cnzi+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 605 nspacedouble++; 606 } 607 608 /* Copy data into free space, then initialize lnk */ 609 ierr = PetscLLClean(bm,bm,cnzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 610 current_space->array += cnzi; 611 current_space->local_used += cnzi; 612 current_space->local_remaining -= cnzi; 613 614 ci[i+1] = ci[i] + cnzi; 615 } 616 617 618 /* Column indices are in the list of free space. 619 Allocate array cj, copy column indices to cj, and destroy list of free space */ 620 ierr = PetscMalloc((ci[am]+1)*sizeof(PetscInt),&cj);CHKERRQ(ierr); 621 ierr = PetscFreeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 622 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 623 624 /* Allocate space for ca */ 625 ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 626 ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 627 628 /* put together the new symbolic matrix */ 629 ierr = MatCreateSeqAIJWithArrays(((PetscObject)A)->comm,am,bm,ci,cj,ca,C);CHKERRQ(ierr); 630 631 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 632 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 633 c = (Mat_SeqAIJ *)((*C)->data); 634 c->free_a = PETSC_TRUE; 635 c->free_ij = PETSC_TRUE; 636 c->nonew = 0; 637 638 /* set MatInfo */ 639 afill = (PetscReal)ci[am]/(ai[am]+bi[bm]) + 1.e-5; 640 if (afill < 1.0) afill = 1.0; 641 c->maxnz = ci[am]; 642 c->nz = ci[am]; 643 (*C)->info.mallocs = nspacedouble; 644 (*C)->info.fill_ratio_given = fill; 645 (*C)->info.fill_ratio_needed = afill; 646 647 #if defined(PETSC_USE_INFO) 648 if (ci[am]) { 649 ierr = PetscInfo3((*C),"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 650 ierr = PetscInfo1((*C),"Use MatMatTransposeMult(A,B,MatReuse,%G,&C) for best performance.;\n",afill);CHKERRQ(ierr); 651 } else { 652 ierr = PetscInfo((*C),"Empty matrix product\n");CHKERRQ(ierr); 653 } 654 #endif 655 #endif 656 PetscFunctionReturn(0); 657 } 658 659 /* #define USE_ARRAY - for sparse dot product. Slower than !USE_ARRAY */ 660 #undef __FUNCT__ 661 #define __FUNCT__ "MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ" 662 PetscErrorCode MatMatTransposeMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 663 { 664 PetscErrorCode ierr; 665 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 666 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,anzi,bnzj,nexta,nextb,*acol,*bcol,brow; 667 PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,i,j,cnzi,*ccol; 668 PetscLogDouble flops=0.0; 669 MatScalar *aa=a->a,*aval,*ba=b->a,*bval,*ca=c->a,*cval; 670 Mat_MatMatTransMult *multtrans; 671 PetscContainer container; 672 #if defined(USE_ARRAY) 673 MatScalar *spdot; 674 #endif 675 676 PetscFunctionBegin; 677 ierr = PetscObjectQuery((PetscObject)C,"Mat_MatMatTransMult",(PetscObject *)&container);CHKERRQ(ierr); 678 if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 679 ierr = PetscContainerGetPointer(container,(void **)&multtrans);CHKERRQ(ierr); 680 if (multtrans->usecoloring){ 681 MatTransposeColoring matcoloring = multtrans->matcoloring; 682 Mat Bt_dense; 683 PetscInt m,n; 684 PetscLogDouble t0,tf,etime0=0.0,etime1=0.0,etime2=0.0; 685 Mat C_dense = multtrans->ABt_den; 686 687 Bt_dense = multtrans->Bt_den; 688 ierr = MatGetLocalSize(Bt_dense,&m,&n);CHKERRQ(ierr); 689 690 /* Get Bt_dense by Apply MatTransposeColoring to B */ 691 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 692 ierr = MatTransColoringApplySpToDen(matcoloring,B,Bt_dense);CHKERRQ(ierr); 693 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 694 etime0 += tf - t0; 695 696 /* C_dense = A*Bt_dense */ 697 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 698 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,Bt_dense,C_dense);CHKERRQ(ierr); 699 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 700 etime2 += tf - t0; 701 702 /* Recover C from C_dense */ 703 ierr = PetscGetTime(&t0);CHKERRQ(ierr); 704 ierr = MatTransColoringApplyDenToSp(matcoloring,C_dense,C);CHKERRQ(ierr); 705 ierr = PetscGetTime(&tf);CHKERRQ(ierr); 706 etime1 += tf - t0; 707 #if defined(PETSC_USE_INFO) 708 ierr = PetscInfo4(C,"Num = ColoringApply: %g %g + Mult_sp_dense %g = %g\n",etime0,etime1,etime2,etime0+etime1+etime2); 709 #endif 710 PetscFunctionReturn(0); 711 } 712 713 #if defined(USE_ARRAY) 714 /* allocate an array for implementing sparse inner-product */ 715 ierr = PetscMalloc((A->cmap->n+1)*sizeof(MatScalar),&spdot);CHKERRQ(ierr); 716 ierr = PetscMemzero(spdot,(A->cmap->n+1)*sizeof(MatScalar));CHKERRQ(ierr); 717 #endif 718 719 /* clear old values in C */ 720 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 721 722 for (i=0; i<cm; i++) { 723 anzi = ai[i+1] - ai[i]; 724 acol = aj + ai[i]; 725 aval = aa + ai[i]; 726 cnzi = ci[i+1] - ci[i]; 727 ccol = cj + ci[i]; 728 cval = ca + ci[i]; 729 for (j=0; j<cnzi; j++){ 730 brow = ccol[j]; 731 bnzj = bi[brow+1] - bi[brow]; 732 bcol = bj + bi[brow]; 733 bval = ba + bi[brow]; 734 735 /* perform sparse inner-product c(i,j)=A[i,:]*B[j,:]^T */ 736 #if defined(USE_ARRAY) 737 /* put ba to spdot array */ 738 for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = bval[nextb]; 739 /* c(i,j)=A[i,:]*B[j,:]^T */ 740 for (nexta=0; nexta<anzi; nexta++){ 741 cval[j] += spdot[acol[nexta]]*aval[nexta]; 742 } 743 /* zero spdot array */ 744 for (nextb=0; nextb<bnzj; nextb++) spdot[bcol[nextb]] = 0.0; 745 #else 746 nexta = 0; nextb = 0; 747 while (nexta<anzi && nextb<bnzj){ 748 while (acol[nexta] < bcol[nextb] && nexta < anzi) nexta++; 749 if (nexta == anzi) break; 750 while (acol[nexta] > bcol[nextb] && nextb < bnzj) nextb++; 751 if (nextb == bnzj) break; 752 if (acol[nexta] == bcol[nextb]){ 753 cval[j] += aval[nexta]*bval[nextb]; 754 nexta++; nextb++; 755 flops += 2; 756 } 757 } 758 #endif 759 } 760 } 761 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 762 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 763 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 764 #if defined(USE_ARRAY) 765 ierr = PetscFree(spdot); 766 #endif 767 PetscFunctionReturn(0); 768 } 769 770 #undef __FUNCT__ 771 #define __FUNCT__ "MatTransposeMatMult_SeqAIJ_SeqAIJ" 772 PetscErrorCode MatTransposeMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) { 773 PetscErrorCode ierr; 774 775 PetscFunctionBegin; 776 if (scall == MAT_INITIAL_MATRIX){ 777 ierr = MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 778 } 779 ierr = MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 780 PetscFunctionReturn(0); 781 } 782 783 #undef __FUNCT__ 784 #define __FUNCT__ "MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ" 785 PetscErrorCode MatTransposeMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 786 { 787 PetscErrorCode ierr; 788 Mat At; 789 PetscInt *ati,*atj; 790 791 PetscFunctionBegin; 792 /* create symbolic At */ 793 ierr = MatGetSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 794 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,A->cmap->n,A->rmap->n,ati,atj,PETSC_NULL,&At);CHKERRQ(ierr); 795 796 /* get symbolic C=At*B */ 797 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(At,B,fill,C);CHKERRQ(ierr); 798 799 /* clean up */ 800 ierr = MatDestroy(&At);CHKERRQ(ierr); 801 ierr = MatRestoreSymbolicTranspose_SeqAIJ(A,&ati,&atj);CHKERRQ(ierr); 802 PetscFunctionReturn(0); 803 } 804 805 #undef __FUNCT__ 806 #define __FUNCT__ "MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ" 807 PetscErrorCode MatTransposeMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 808 { 809 PetscErrorCode ierr; 810 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c=(Mat_SeqAIJ*)C->data; 811 PetscInt am=A->rmap->n,anzi,*ai=a->i,*aj=a->j,*bi=b->i,*bj,bnzi,nextb; 812 PetscInt cm=C->rmap->n,*ci=c->i,*cj=c->j,crow,*cjj,i,j,k; 813 PetscLogDouble flops=0.0; 814 MatScalar *aa=a->a,*ba,*ca=c->a,*caj; 815 816 PetscFunctionBegin; 817 /* clear old values in C */ 818 ierr = PetscMemzero(ca,ci[cm]*sizeof(MatScalar));CHKERRQ(ierr); 819 820 /* compute A^T*B using outer product (A^T)[:,i]*B[i,:] */ 821 for (i=0;i<am;i++) { 822 bj = b->j + bi[i]; 823 ba = b->a + bi[i]; 824 bnzi = bi[i+1] - bi[i]; 825 anzi = ai[i+1] - ai[i]; 826 for (j=0; j<anzi; j++) { 827 nextb = 0; 828 crow = *aj++; 829 cjj = cj + ci[crow]; 830 caj = ca + ci[crow]; 831 /* perform sparse axpy operation. Note cjj includes bj. */ 832 for (k=0; nextb<bnzi; k++) { 833 if (cjj[k] == *(bj+nextb)) { /* ccol == bcol */ 834 caj[k] += (*aa)*(*(ba+nextb)); 835 nextb++; 836 } 837 } 838 flops += 2*bnzi; 839 aa++; 840 } 841 } 842 843 /* Assemble the final matrix and clean up */ 844 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 845 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 846 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 847 PetscFunctionReturn(0); 848 } 849 850 EXTERN_C_BEGIN 851 #undef __FUNCT__ 852 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqDense" 853 PetscErrorCode MatMatMult_SeqAIJ_SeqDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 854 { 855 PetscErrorCode ierr; 856 857 PetscFunctionBegin; 858 if (scall == MAT_INITIAL_MATRIX){ 859 ierr = MatMatMultSymbolic_SeqAIJ_SeqDense(A,B,fill,C);CHKERRQ(ierr); 860 } 861 ierr = MatMatMultNumeric_SeqAIJ_SeqDense(A,B,*C);CHKERRQ(ierr); 862 PetscFunctionReturn(0); 863 } 864 EXTERN_C_END 865 866 #undef __FUNCT__ 867 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqDense" 868 PetscErrorCode MatMatMultSymbolic_SeqAIJ_SeqDense(Mat A,Mat B,PetscReal fill,Mat *C) 869 { 870 PetscErrorCode ierr; 871 872 PetscFunctionBegin; 873 ierr = MatMatMultSymbolic_SeqDense_SeqDense(A,B,0.0,C);CHKERRQ(ierr); 874 (*C)->ops->matmult = MatMatMult_SeqAIJ_SeqDense; 875 PetscFunctionReturn(0); 876 } 877 878 #undef __FUNCT__ 879 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqDense" 880 PetscErrorCode MatMatMultNumeric_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 881 { 882 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 883 PetscErrorCode ierr; 884 PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 885 MatScalar *aa; 886 PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n; 887 PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam; 888 889 PetscFunctionBegin; 890 if (!cm || !cn) PetscFunctionReturn(0); 891 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); 892 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); 893 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); 894 ierr = MatGetArray(B,&b);CHKERRQ(ierr); 895 ierr = MatGetArray(C,&c);CHKERRQ(ierr); 896 b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 897 for (col=0; col<cn-4; col += 4){ /* over columns of C */ 898 colam = col*am; 899 for (i=0; i<am; i++) { /* over rows of C in those columns */ 900 r1 = r2 = r3 = r4 = 0.0; 901 n = a->i[i+1] - a->i[i]; 902 aj = a->j + a->i[i]; 903 aa = a->a + a->i[i]; 904 for (j=0; j<n; j++) { 905 r1 += (*aa)*b1[*aj]; 906 r2 += (*aa)*b2[*aj]; 907 r3 += (*aa)*b3[*aj]; 908 r4 += (*aa++)*b4[*aj++]; 909 } 910 c[colam + i] = r1; 911 c[colam + am + i] = r2; 912 c[colam + am2 + i] = r3; 913 c[colam + am3 + i] = r4; 914 } 915 b1 += bm4; 916 b2 += bm4; 917 b3 += bm4; 918 b4 += bm4; 919 } 920 for (;col<cn; col++){ /* over extra columns of C */ 921 for (i=0; i<am; i++) { /* over rows of C in those columns */ 922 r1 = 0.0; 923 n = a->i[i+1] - a->i[i]; 924 aj = a->j + a->i[i]; 925 aa = a->a + a->i[i]; 926 927 for (j=0; j<n; j++) { 928 r1 += (*aa++)*b1[*aj++]; 929 } 930 c[col*am + i] = r1; 931 } 932 b1 += bm; 933 } 934 ierr = PetscLogFlops(cn*(2.0*a->nz));CHKERRQ(ierr); 935 ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 936 ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 937 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 938 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 939 PetscFunctionReturn(0); 940 } 941 942 /* 943 Note very similar to MatMult_SeqAIJ(), should generate both codes from same base 944 */ 945 #undef __FUNCT__ 946 #define __FUNCT__ "MatMatMultNumericAdd_SeqAIJ_SeqDense" 947 PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat A,Mat B,Mat C) 948 { 949 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 950 PetscErrorCode ierr; 951 PetscScalar *b,*c,r1,r2,r3,r4,*b1,*b2,*b3,*b4; 952 MatScalar *aa; 953 PetscInt cm=C->rmap->n, cn=B->cmap->n, bm=B->rmap->n, col, i,j,n,*aj, am = A->rmap->n,*ii,arm; 954 PetscInt am2 = 2*am, am3 = 3*am, bm4 = 4*bm,colam,*ridx; 955 956 PetscFunctionBegin; 957 if (!cm || !cn) PetscFunctionReturn(0); 958 ierr = MatGetArray(B,&b);CHKERRQ(ierr); 959 ierr = MatGetArray(C,&c);CHKERRQ(ierr); 960 b1 = b; b2 = b1 + bm; b3 = b2 + bm; b4 = b3 + bm; 961 962 if (a->compressedrow.use){ /* use compressed row format */ 963 for (col=0; col<cn-4; col += 4){ /* over columns of C */ 964 colam = col*am; 965 arm = a->compressedrow.nrows; 966 ii = a->compressedrow.i; 967 ridx = a->compressedrow.rindex; 968 for (i=0; i<arm; i++) { /* over rows of C in those columns */ 969 r1 = r2 = r3 = r4 = 0.0; 970 n = ii[i+1] - ii[i]; 971 aj = a->j + ii[i]; 972 aa = a->a + ii[i]; 973 for (j=0; j<n; j++) { 974 r1 += (*aa)*b1[*aj]; 975 r2 += (*aa)*b2[*aj]; 976 r3 += (*aa)*b3[*aj]; 977 r4 += (*aa++)*b4[*aj++]; 978 } 979 c[colam + ridx[i]] += r1; 980 c[colam + am + ridx[i]] += r2; 981 c[colam + am2 + ridx[i]] += r3; 982 c[colam + am3 + ridx[i]] += r4; 983 } 984 b1 += bm4; 985 b2 += bm4; 986 b3 += bm4; 987 b4 += bm4; 988 } 989 for (;col<cn; col++){ /* over extra columns of C */ 990 colam = col*am; 991 arm = a->compressedrow.nrows; 992 ii = a->compressedrow.i; 993 ridx = a->compressedrow.rindex; 994 for (i=0; i<arm; i++) { /* over rows of C in those columns */ 995 r1 = 0.0; 996 n = ii[i+1] - ii[i]; 997 aj = a->j + ii[i]; 998 aa = a->a + ii[i]; 999 1000 for (j=0; j<n; j++) { 1001 r1 += (*aa++)*b1[*aj++]; 1002 } 1003 c[col*am + ridx[i]] += r1; 1004 } 1005 b1 += bm; 1006 } 1007 } else { 1008 for (col=0; col<cn-4; col += 4){ /* over columns of C */ 1009 colam = col*am; 1010 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1011 r1 = r2 = r3 = r4 = 0.0; 1012 n = a->i[i+1] - a->i[i]; 1013 aj = a->j + a->i[i]; 1014 aa = a->a + a->i[i]; 1015 for (j=0; j<n; j++) { 1016 r1 += (*aa)*b1[*aj]; 1017 r2 += (*aa)*b2[*aj]; 1018 r3 += (*aa)*b3[*aj]; 1019 r4 += (*aa++)*b4[*aj++]; 1020 } 1021 c[colam + i] += r1; 1022 c[colam + am + i] += r2; 1023 c[colam + am2 + i] += r3; 1024 c[colam + am3 + i] += r4; 1025 } 1026 b1 += bm4; 1027 b2 += bm4; 1028 b3 += bm4; 1029 b4 += bm4; 1030 } 1031 for (;col<cn; col++){ /* over extra columns of C */ 1032 for (i=0; i<am; i++) { /* over rows of C in those columns */ 1033 r1 = 0.0; 1034 n = a->i[i+1] - a->i[i]; 1035 aj = a->j + a->i[i]; 1036 aa = a->a + a->i[i]; 1037 1038 for (j=0; j<n; j++) { 1039 r1 += (*aa++)*b1[*aj++]; 1040 } 1041 c[col*am + i] += r1; 1042 } 1043 b1 += bm; 1044 } 1045 } 1046 ierr = PetscLogFlops(cn*2.0*a->nz);CHKERRQ(ierr); 1047 ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 1048 ierr = MatRestoreArray(C,&c);CHKERRQ(ierr); 1049 PetscFunctionReturn(0); 1050 } 1051 1052 #undef __FUNCT__ 1053 #define __FUNCT__ "MatTransColoringApplySpToDen_SeqAIJ" 1054 PetscErrorCode MatTransColoringApplySpToDen_SeqAIJ(MatTransposeColoring coloring,Mat B,Mat Btdense) 1055 { 1056 PetscErrorCode ierr; 1057 Mat_SeqAIJ *b = (Mat_SeqAIJ*)B->data; 1058 Mat_SeqDense *btdense = (Mat_SeqDense*)Btdense->data; 1059 PetscInt *bi=b->i,*bj=b->j; 1060 PetscInt m=Btdense->rmap->n,n=Btdense->cmap->n,j,k,l,col,anz,*btcol,brow,ncolumns; 1061 MatScalar *btval,*btval_den,*ba=b->a; 1062 PetscInt *columns=coloring->columns,*colorforcol=coloring->colorforcol,ncolors=coloring->ncolors; 1063 1064 PetscFunctionBegin; 1065 btval_den=btdense->v; 1066 ierr = PetscMemzero(btval_den,(m*n)*sizeof(MatScalar));CHKERRQ(ierr); 1067 for (k=0; k<ncolors; k++) { 1068 ncolumns = coloring->ncolumns[k]; 1069 for (l=0; l<ncolumns; l++) { /* insert a row of B to a column of Btdense */ 1070 col = *(columns + colorforcol[k] + l); 1071 btcol = bj + bi[col]; 1072 btval = ba + bi[col]; 1073 anz = bi[col+1] - bi[col]; 1074 for (j=0; j<anz; j++){ 1075 brow = btcol[j]; 1076 btval_den[brow] = btval[j]; 1077 } 1078 } 1079 btval_den += m; 1080 } 1081 PetscFunctionReturn(0); 1082 } 1083 1084 #undef __FUNCT__ 1085 #define __FUNCT__ "MatTransColoringApplyDenToSp_SeqAIJ" 1086 PetscErrorCode MatTransColoringApplyDenToSp_SeqAIJ(MatTransposeColoring matcoloring,Mat Cden,Mat Csp) 1087 { 1088 PetscErrorCode ierr; 1089 Mat_SeqAIJ *csp = (Mat_SeqAIJ*)Csp->data; 1090 PetscInt k,l,*row,*idx,m,ncolors=matcoloring->ncolors,nrows; 1091 PetscScalar *ca_den,*cp_den,*ca=csp->a; 1092 PetscInt *rows=matcoloring->rows,*spidx=matcoloring->columnsforspidx,*colorforrow=matcoloring->colorforrow; 1093 1094 PetscFunctionBegin; 1095 ierr = MatGetLocalSize(Csp,&m,PETSC_NULL);CHKERRQ(ierr); 1096 ierr = MatGetArray(Cden,&ca_den);CHKERRQ(ierr); 1097 cp_den = ca_den; 1098 for (k=0; k<ncolors; k++) { 1099 nrows = matcoloring->nrows[k]; 1100 row = rows + colorforrow[k]; 1101 idx = spidx + colorforrow[k]; 1102 for (l=0; l<nrows; l++){ 1103 ca[idx[l]] = cp_den[row[l]]; 1104 } 1105 cp_den += m; 1106 } 1107 ierr = MatRestoreArray(Cden,&ca_den);CHKERRQ(ierr); 1108 PetscFunctionReturn(0); 1109 } 1110 1111 /* 1112 MatGetColumnIJ_SeqAIJ_Color() and MatRestoreColumnIJ_SeqAIJ_Color() are customized from 1113 MatGetColumnIJ_SeqAIJ() and MatRestoreColumnIJ_SeqAIJ() by adding an output 1114 spidx[], index of a->j, to be used for setting 'columnsforspidx' in MatTransposeColoringCreate_SeqAIJ(). 1115 */ 1116 #undef __FUNCT__ 1117 #define __FUNCT__ "MatGetColumnIJ_SeqAIJ_Color" 1118 PetscErrorCode MatGetColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *nn,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1119 { 1120 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1121 PetscErrorCode ierr; 1122 PetscInt i,*collengths,*cia,*cja,n = A->cmap->n,m = A->rmap->n; 1123 PetscInt nz = a->i[m],row,*jj,mr,col; 1124 PetscInt *cspidx; 1125 1126 PetscFunctionBegin; 1127 *nn = n; 1128 if (!ia) PetscFunctionReturn(0); 1129 if (symmetric) { 1130 SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"MatGetColumnIJ_SeqAIJ_Color() not supported for the case symmetric"); 1131 ierr = MatToSymmetricIJ_SeqAIJ(A->rmap->n,a->i,a->j,0,oshift,ia,ja);CHKERRQ(ierr); 1132 } else { 1133 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&collengths);CHKERRQ(ierr); 1134 ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 1135 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&cia);CHKERRQ(ierr); 1136 ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cja);CHKERRQ(ierr); 1137 ierr = PetscMalloc((nz+1)*sizeof(PetscInt),&cspidx);CHKERRQ(ierr); 1138 jj = a->j; 1139 for (i=0; i<nz; i++) { 1140 collengths[jj[i]]++; 1141 } 1142 cia[0] = oshift; 1143 for (i=0; i<n; i++) { 1144 cia[i+1] = cia[i] + collengths[i]; 1145 } 1146 ierr = PetscMemzero(collengths,n*sizeof(PetscInt));CHKERRQ(ierr); 1147 jj = a->j; 1148 for (row=0; row<m; row++) { 1149 mr = a->i[row+1] - a->i[row]; 1150 for (i=0; i<mr; i++) { 1151 col = *jj++; 1152 cspidx[cia[col] + collengths[col] - oshift] = a->i[row] + i; /* index of a->j */ 1153 cja[cia[col] + collengths[col]++ - oshift] = row + oshift; 1154 } 1155 } 1156 ierr = PetscFree(collengths);CHKERRQ(ierr); 1157 *ia = cia; *ja = cja; 1158 *spidx = cspidx; 1159 } 1160 PetscFunctionReturn(0); 1161 } 1162 1163 #undef __FUNCT__ 1164 #define __FUNCT__ "MatRestoreColumnIJ_SeqAIJ_Color" 1165 PetscErrorCode MatRestoreColumnIJ_SeqAIJ_Color(Mat A,PetscInt oshift,PetscBool symmetric,PetscBool inodecompressed,PetscInt *n,PetscInt *ia[],PetscInt *ja[],PetscInt *spidx[],PetscBool *done) 1166 { 1167 PetscErrorCode ierr; 1168 1169 PetscFunctionBegin; 1170 if (!ia) PetscFunctionReturn(0); 1171 1172 ierr = PetscFree(*ia);CHKERRQ(ierr); 1173 ierr = PetscFree(*ja);CHKERRQ(ierr); 1174 ierr = PetscFree(*spidx);CHKERRQ(ierr); 1175 PetscFunctionReturn(0); 1176 } 1177 1178 #undef __FUNCT__ 1179 #define __FUNCT__ "MatTransposeColoringCreate_SeqAIJ" 1180 PetscErrorCode MatTransposeColoringCreate_SeqAIJ(Mat mat,ISColoring iscoloring,MatTransposeColoring c) 1181 { 1182 PetscErrorCode ierr; 1183 PetscInt i,n,nrows,N,j,k,m,*row_idx,*ci,*cj,ncols,col,cm; 1184 const PetscInt *is; 1185 PetscInt nis = iscoloring->n,*rowhit,bs = 1; 1186 IS *isa; 1187 PetscBool done; 1188 PetscBool flg1,flg2; 1189 Mat_SeqAIJ *csp = (Mat_SeqAIJ*)mat->data; 1190 PetscInt *colorforrow,*rows,*rows_i,*columnsforspidx,*columnsforspidx_i,*idxhit,*spidx; 1191 PetscInt *colorforcol,*columns,*columns_i; 1192 1193 PetscFunctionBegin; 1194 ierr = ISColoringGetIS(iscoloring,PETSC_IGNORE,&isa);CHKERRQ(ierr); 1195 1196 /* this is ugly way to get blocksize but cannot call MatGetBlockSize() because AIJ can have bs > 1 */ 1197 ierr = PetscTypeCompare((PetscObject)mat,MATSEQBAIJ,&flg1);CHKERRQ(ierr); 1198 ierr = PetscTypeCompare((PetscObject)mat,MATMPIBAIJ,&flg2);CHKERRQ(ierr); 1199 if (flg1 || flg2) { 1200 ierr = MatGetBlockSize(mat,&bs);CHKERRQ(ierr); 1201 } 1202 1203 N = mat->cmap->N/bs; 1204 c->M = mat->rmap->N/bs; /* set total rows, columns and local rows */ 1205 c->N = mat->cmap->N/bs; 1206 c->m = mat->rmap->N/bs; 1207 c->rstart = 0; 1208 1209 c->ncolors = nis; 1210 ierr = PetscMalloc(nis*sizeof(PetscInt),&c->ncolumns);CHKERRQ(ierr); 1211 ierr = PetscMalloc(nis*sizeof(PetscInt),&c->nrows);CHKERRQ(ierr); 1212 ierr = PetscMalloc2(csp->nz+1,PetscInt,&rows,csp->nz+1,PetscInt,&columnsforspidx);CHKERRQ(ierr); 1213 ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforrow);CHKERRQ(ierr); 1214 colorforrow[0] = 0; 1215 rows_i = rows; 1216 columnsforspidx_i = columnsforspidx; 1217 1218 ierr = PetscMalloc((nis+1)*sizeof(PetscInt),&colorforcol);CHKERRQ(ierr); 1219 ierr = PetscMalloc((N+1)*sizeof(PetscInt),&columns);CHKERRQ(ierr); 1220 colorforcol[0] = 0; 1221 columns_i = columns; 1222 1223 ierr = MatGetColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); /* column-wise storage of mat */ 1224 if (!done) SETERRQ1(((PetscObject)mat)->comm,PETSC_ERR_SUP,"MatGetColumnIJ() not supported for matrix type %s",((PetscObject)mat)->type_name); 1225 1226 cm = c->m; 1227 ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&rowhit);CHKERRQ(ierr); 1228 ierr = PetscMalloc((cm+1)*sizeof(PetscInt),&idxhit);CHKERRQ(ierr); 1229 for (i=0; i<nis; i++) { 1230 ierr = ISGetLocalSize(isa[i],&n);CHKERRQ(ierr); 1231 ierr = ISGetIndices(isa[i],&is);CHKERRQ(ierr); 1232 c->ncolumns[i] = n; 1233 if (n) { 1234 ierr = PetscMemcpy(columns_i,is,n*sizeof(PetscInt));CHKERRQ(ierr); 1235 } 1236 colorforcol[i+1] = colorforcol[i] + n; 1237 columns_i += n; 1238 1239 /* fast, crude version requires O(N*N) work */ 1240 ierr = PetscMemzero(rowhit,cm*sizeof(PetscInt));CHKERRQ(ierr); 1241 1242 /* loop over columns*/ 1243 for (j=0; j<n; j++) { 1244 col = is[j]; 1245 row_idx = cj + ci[col]; 1246 m = ci[col+1] - ci[col]; 1247 /* loop over columns marking them in rowhit */ 1248 for (k=0; k<m; k++) { 1249 idxhit[*row_idx] = spidx[ci[col] + k]; 1250 rowhit[*row_idx++] = col + 1; 1251 } 1252 } 1253 /* count the number of hits */ 1254 nrows = 0; 1255 for (j=0; j<cm; j++) { 1256 if (rowhit[j]) nrows++; 1257 } 1258 c->nrows[i] = nrows; 1259 colorforrow[i+1] = colorforrow[i] + nrows; 1260 1261 nrows = 0; 1262 for (j=0; j<cm; j++) { 1263 if (rowhit[j]) { 1264 rows_i[nrows] = j; 1265 columnsforspidx_i[nrows] = idxhit[j]; 1266 nrows++; 1267 } 1268 } 1269 ierr = ISRestoreIndices(isa[i],&is);CHKERRQ(ierr); 1270 rows_i += nrows; columnsforspidx_i += nrows; 1271 } 1272 ierr = MatRestoreColumnIJ_SeqAIJ_Color(mat,0,PETSC_FALSE,PETSC_FALSE,&ncols,&ci,&cj,&spidx,&done);CHKERRQ(ierr); 1273 ierr = PetscFree(rowhit);CHKERRQ(ierr); 1274 ierr = ISColoringRestoreIS(iscoloring,&isa);CHKERRQ(ierr); 1275 #if defined(PETSC_USE_DEBUG) 1276 if (csp->nz != colorforrow[nis]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"csp->nz %d != colorforrow[nis] %d",csp->nz,colorforrow[nis]); 1277 #endif 1278 1279 c->colorforrow = colorforrow; 1280 c->rows = rows; 1281 c->columnsforspidx = columnsforspidx; 1282 c->colorforcol = colorforcol; 1283 c->columns = columns; 1284 ierr = PetscFree(idxhit);CHKERRQ(ierr); 1285 PetscFunctionReturn(0); 1286 } 1287