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