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