1 /* 2 Defines matrix-matrix product routines for pairs of AIJ matrices 3 C = A * B 4 */ 5 6 #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/ 7 #include "src/mat/utils/freespace.h" 8 #include "src/mat/impls/aij/mpi/mpiaij.h" 9 10 /* 11 Initialize a linked list 12 Input Parameters: 13 lnk_init - the initial index value indicating the entry in the list is not set yet 14 nlnk - max length of the list 15 lnk - linked list(an integer array) that is allocated 16 output Parameters: 17 lnk - the linked list with all values set as lnk_int 18 */ 19 #define LNKLISTINITIALIZE(lnk_init,nlnk,lnk){\ 20 int _i;\ 21 for (_i=0; _i<nlnk; _i++) lnk[_i] = lnk_init;\ 22 } 23 24 /* 25 Add a index set into a sorted linked list 26 Input Parameters: 27 nidx - number of input indices 28 indices - interger array 29 lnk_head - the header of the list 30 lnk_init - the initial index value indicating the entry in the list is not set yet 31 lnk - linked list(an integer array) that is created 32 output Parameters: 33 nlnk - number of newly added indices 34 lnk - the sorted(increasing order) linked list containing new and non-redundate entries from indices 35 */ 36 #define LNKLISTADD(nidx,indices,lnk_head,lnk_init,nlnk,lnk){\ 37 int _k,_entry,_lidx=lnk_head,_idx;\ 38 nlnk = 0;\ 39 _k=nidx;\ 40 while (_k){/* assume indices are almost in increasing order, starting from its end saves computation */\ 41 _entry = indices[--_k];\ 42 if (lnk[_entry] == lnk_init) { /* new col */\ 43 do {\ 44 _idx = _lidx;\ 45 _lidx = lnk[_idx];\ 46 } while (_entry > _lidx);\ 47 lnk[_idx] = _entry;\ 48 lnk[_entry] = _lidx;\ 49 nlnk++;\ 50 }\ 51 }\ 52 } 53 /* 54 Copy data on the list into an array, then initialize the list 55 Input Parameters: 56 lnk_head - the header of the list 57 lnk_init - the initial index value indicating the entry in the list is not set yet 58 nlnk - number of data on the list to be copied 59 lnk - linked list 60 output Parameters: 61 indices - array that contains the copied data 62 */ 63 #define LNKLISTCLEAR(lnk_head,lnk_init,nlnk,lnk,indices){\ 64 int _j,_idx=lnk_head,_idx0;\ 65 for (_j=0; _j<nlnk; _j++){\ 66 _idx0 = _idx; _idx = lnk[_idx0];\ 67 *(indices+_j) = _idx;\ 68 lnk[_idx0] = lnk_init;\ 69 }\ 70 lnk[_idx] = lnk_init;\ 71 } 72 73 typedef struct { /* used by MatMatMult_MPIAIJ_MPIAIJ for reusing symbolic mat product */ 74 IS isrowa,isrowb,iscolb; 75 Mat *aseq,*bseq,C_seq; 76 } Mat_MatMatMultMPI; 77 78 static int logkey_matmatmult_symbolic = 0; 79 static int logkey_matmatmult_numeric = 0; 80 81 #undef __FUNCT__ 82 #define __FUNCT__ "MatMatMult" 83 /*@ 84 MatMatMult - Performs Matrix-Matrix Multiplication C=A*B. 85 86 Collective on Mat 87 88 Input Parameters: 89 + A - the left matrix 90 . B - the right matrix 91 . scall - either MAT_INITIAL_MATRIX or MAT_REUSE_MATRIX 92 - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)) 93 94 Output Parameters: 95 . C - the product matrix 96 97 Notes: 98 C will be created and must be destroyed by the user with MatDestroy(). 99 100 This routine is currently only implemented for pairs of SeqAIJ matrices and classes 101 which inherit from SeqAIJ. C will be of type MATSEQAIJ. 102 103 Level: intermediate 104 105 .seealso: MatMatMultSymbolic(),MatMatMultNumeric() 106 @*/ 107 int MatMatMult(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) 108 { 109 int ierr; 110 111 PetscFunctionBegin; 112 PetscValidHeaderSpecific(A,MAT_COOKIE,1); 113 PetscValidType(A,1); 114 MatPreallocated(A); 115 if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 116 if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 117 PetscValidHeaderSpecific(B,MAT_COOKIE,2); 118 PetscValidType(B,2); 119 MatPreallocated(B); 120 if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 121 if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 122 PetscValidPointer(C,3); 123 if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N); 124 125 if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill); 126 127 ierr = PetscLogEventBegin(MAT_MatMult,A,B,0,0);CHKERRQ(ierr); 128 ierr = (*A->ops->matmult)(A,B,scall,fill,C);CHKERRQ(ierr); 129 ierr = PetscLogEventEnd(MAT_MatMult,A,B,0,0);CHKERRQ(ierr); 130 131 PetscFunctionReturn(0); 132 } 133 134 #undef __FUNCT__ 135 #define __FUNCT__ "MatMatMult_MPIAIJ_MPIAIJ" 136 int MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C) 137 { 138 int ierr; 139 140 PetscFunctionBegin; 141 if (scall == MAT_INITIAL_MATRIX){ 142 ierr = MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);CHKERRQ(ierr); 143 } 144 ierr = MatMatMultNumeric_MPIAIJ_MPIAIJ(A,B,*C);CHKERRQ(ierr); 145 PetscFunctionReturn(0); 146 } 147 148 #undef __FUNCT__ 149 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ" 150 int MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C) { 151 int ierr; 152 153 PetscFunctionBegin; 154 if (scall == MAT_INITIAL_MATRIX){ 155 ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(A,B,fill,C);CHKERRQ(ierr); 156 } 157 ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ(A,B,*C);CHKERRQ(ierr); 158 PetscFunctionReturn(0); 159 } 160 161 #undef __FUNCT__ 162 #define __FUNCT__ "MatMatMultSymbolic" 163 /*@ 164 MatMatMultSymbolic - Performs construction, preallocation, and computes the ij structure 165 of the matrix-matrix product C=A*B. Call this routine before calling MatMatMultNumeric(). 166 167 Collective on Mat 168 169 Input Parameters: 170 + A - the left matrix 171 . B - the right matrix 172 - fill - expected fill as ratio of nnz(C)/(nnz(A) + nnz(B)) 173 174 Output Parameters: 175 . C - the matrix containing the ij structure of product matrix 176 177 Notes: 178 C will be created as a MATSEQAIJ matrix and must be destroyed by the user with MatDestroy(). 179 180 This routine is currently only implemented for SeqAIJ matrices and classes which inherit from SeqAIJ. 181 182 Level: intermediate 183 184 .seealso: MatMatMult(),MatMatMultNumeric() 185 @*/ 186 int MatMatMultSymbolic(Mat A,Mat B,PetscReal fill,Mat *C) { 187 /* Perhaps this "interface" routine should be moved into the interface directory.*/ 188 /* To facilitate implementations with varying types, QueryFunction is used.*/ 189 /* It is assumed that implementations will be composed as "MatMatMultSymbolic_<type of A><type of B>". */ 190 int ierr; 191 char symfunct[80]; 192 int (*symbolic)(Mat,Mat,PetscReal,Mat *); 193 194 PetscFunctionBegin; 195 PetscValidHeaderSpecific(A,MAT_COOKIE,1); 196 PetscValidType(A,1); 197 MatPreallocated(A); 198 if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 199 if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 200 201 PetscValidHeaderSpecific(B,MAT_COOKIE,2); 202 PetscValidType(B,2); 203 MatPreallocated(B); 204 if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 205 if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 206 PetscValidPointer(C,3); 207 208 if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N); 209 if (fill <=0.0) SETERRQ1(PETSC_ERR_ARG_SIZ,"fill=%g must be > 0.0",fill); 210 211 /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */ 212 /* When other implementations exist, attack the multiple dispatch problem. */ 213 ierr = PetscStrcpy(symfunct,"MatMatMultSymbolic_seqaijseqaij");CHKERRQ(ierr); 214 ierr = PetscObjectQueryFunction((PetscObject)B,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr); 215 if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name); 216 ierr = PetscObjectQueryFunction((PetscObject)A,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr); 217 if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name); 218 ierr = (*symbolic)(A,B,fill,C);CHKERRQ(ierr); 219 220 PetscFunctionReturn(0); 221 } 222 223 EXTERN int MatDestroy_MPIAIJ(Mat); 224 #undef __FUNCT__ 225 #define __FUNCT__ "MatDestroy_MPIAIJ_MatMatMult" 226 int MatDestroy_MPIAIJ_MatMatMult(Mat A) 227 { 228 int ierr; 229 Mat_MatMatMultMPI *mult = ( Mat_MatMatMultMPI*)A->spptr; 230 231 PetscFunctionBegin; 232 /* printf(" ...MatDestroy_MPIAIJ_MatMatMult is called\n"); */ 233 ierr = PetscFree(mult);CHKERRQ(ierr); 234 ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr); 235 236 PetscFunctionReturn(0); 237 } 238 239 #undef __FUNCT__ 240 #define __FUNCT__ "MatMatMultSymbolic_MPIAIJ_MPIAIJ" 241 int MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 242 { 243 Mat *aseq,*bseq,A_seq=PETSC_NULL,B_seq=PETSC_NULL,C_seq; 244 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 245 int ierr,*idx,i,start,end,ncols,imark,nzA,nzB,*cmap; 246 IS isrow,iscol; 247 Mat_MatMatMultMPI *mult; 248 249 PetscFunctionBegin; 250 /* create a seq matrix B_seq = submatrix of B by taking rows of B that equal to nonzero col of A */ 251 start = a->cstart; 252 cmap = a->garray; 253 nzA = a->A->n; 254 nzB = a->B->n; 255 ierr = PetscMalloc((nzA+nzB)*sizeof(int), &idx);CHKERRQ(ierr); 256 ncols = 0; 257 for (i=0; i<nzB; i++) { 258 if (cmap[i] < start) idx[ncols++] = cmap[i]; 259 else break; 260 } 261 imark = i; 262 for (i=0; i<nzA; i++) idx[ncols++] = start + i; 263 for (i=imark; i<nzB; i++) idx[ncols++] = cmap[i]; 264 ierr = ISCreateGeneral(PETSC_COMM_SELF,ncols,idx,&isrow);CHKERRQ(ierr); /* isrow of B = iscol of A! */ 265 ierr = PetscFree(idx);CHKERRQ(ierr); 266 ierr = ISCreateStride(PETSC_COMM_SELF,B->N,0,1,&iscol);CHKERRQ(ierr); 267 ierr = MatGetSubMatrices(B,1,&isrow,&iscol,MAT_INITIAL_MATRIX,&bseq);CHKERRQ(ierr); 268 ierr = ISDestroy(iscol);CHKERRQ(ierr); 269 B_seq = bseq[0]; 270 271 /* create a seq matrix A_seq = submatrix of A by taking all local rows of A */ 272 start = a->rstart; end = a->rend; 273 ierr = ISCreateStride(PETSC_COMM_SELF,end-start,start,1,&iscol);CHKERRQ(ierr); /* isrow of A = iscol */ 274 ierr = MatGetSubMatrices(A,1,&iscol,&isrow,MAT_INITIAL_MATRIX,&aseq);CHKERRQ(ierr); 275 ierr = ISDestroy(isrow);CHKERRQ(ierr); 276 ierr = ISDestroy(iscol);CHKERRQ(ierr); 277 A_seq = aseq[0]; 278 279 /* compute C_seq = A_seq * B_seq */ 280 ierr = MatMatMult_SeqAIJ_SeqAIJ(A_seq,B_seq,MAT_INITIAL_MATRIX,fill,&C_seq);CHKERRQ(ierr); 281 /* 282 int rank; 283 ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr); 284 ierr = PetscPrintf(PETSC_COMM_SELF," [%d] A_seq: %d, %d; B_seq: %d, %d; C_seq: %d, %d;\n",rank,A_seq->m,A_seq->n,B_seq->m,B_seq->n,C_seq->m,C_seq->n); 285 */ 286 ierr = MatDestroyMatrices(1,&aseq);CHKERRQ(ierr); 287 ierr = MatDestroyMatrices(1,&bseq);CHKERRQ(ierr); 288 289 /* create mpi matrix C by concatinating C_seq */ 290 ierr = MatMerge(A->comm,C_seq,C);CHKERRQ(ierr); 291 292 /* attach the supporting struct to C for reuse of symbolic C */ 293 ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr); 294 (*C)->spptr = (void*)mult; 295 (*C)->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 296 297 PetscFunctionReturn(0); 298 } 299 300 #undef __FUNCT__ 301 #define __FUNCT__ "MatMatMultSymbolic_SeqAIJ_SeqAIJ" 302 int MatMatMultSymbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,PetscReal fill,Mat *C) 303 { 304 int ierr; 305 FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 306 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 307 int *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj; 308 int *ci,*cj,*lnk; 309 int am=A->M,bn=B->N,bm=B->M; 310 int i,j,anzi,brow,bnzj,cnzi,nlnk,lnk_init=-1,nspacedouble=0; 311 MatScalar *ca; 312 313 PetscFunctionBegin; 314 /* Start timers */ 315 ierr = PetscLogEventBegin(logkey_matmatmult_symbolic,A,B,0,0);CHKERRQ(ierr); 316 /* Set up */ 317 /* Allocate ci array, arrays for fill computation and */ 318 /* free space for accumulating nonzero column info */ 319 ierr = PetscMalloc(((am+1)+1)*sizeof(int),&ci);CHKERRQ(ierr); 320 ci[0] = 0; 321 322 ierr = PetscMalloc((bn+1)*sizeof(int),&lnk);CHKERRQ(ierr); 323 LNKLISTINITIALIZE(lnk_init,bn,lnk); 324 325 /* Initial FreeSpace size is fill*(nnz(A)+nnz(B)) */ 326 ierr = GetMoreSpace(int(fill*(ai[am]+bi[bm])),&free_space);CHKERRQ(ierr); 327 current_space = free_space; 328 329 /* Determine symbolic info for each row of the product: */ 330 for (i=0;i<am;i++) { 331 anzi = ai[i+1] - ai[i]; 332 cnzi = 0; 333 lnk[bn] = bn; 334 j = anzi; 335 aj = a->j + ai[i]; 336 while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 337 j--; 338 brow = *(aj + j); 339 bnzj = bi[brow+1] - bi[brow]; 340 bjj = bj + bi[brow]; 341 /* add non-zero cols of B into the sorted linked list lnk */ 342 LNKLISTADD(bnzj,bjj,bn,lnk_init,nlnk,lnk); 343 cnzi += nlnk; 344 } 345 346 /* If free space is not available, make more free space */ 347 /* Double the amount of total space in the list */ 348 if (current_space->local_remaining<cnzi) { 349 /* printf("MatMatMultSymbolic_SeqAIJ_SeqAIJ()...%d -th row, double space ...\n",i);*/ 350 ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 351 nspacedouble++; 352 } 353 354 /* Copy data into free space, then initialize lnk */ 355 LNKLISTCLEAR(bn,lnk_init,cnzi,lnk,current_space->array); 356 current_space->array += cnzi; 357 358 current_space->local_used += cnzi; 359 current_space->local_remaining -= cnzi; 360 361 ci[i+1] = ci[i] + cnzi; 362 } 363 364 /* Column indices are in the list of free space */ 365 /* Allocate space for cj, initialize cj, and */ 366 /* destroy list of free space and other temporary array(s) */ 367 ierr = PetscMalloc((ci[am]+1)*sizeof(int),&cj);CHKERRQ(ierr); 368 ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 369 ierr = PetscFree(lnk);CHKERRQ(ierr); 370 371 /* Allocate space for ca */ 372 ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 373 ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 374 375 /* put together the new symbolic matrix */ 376 ierr = MatCreateSeqAIJWithArrays(A->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr); 377 378 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 379 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 380 c = (Mat_SeqAIJ *)((*C)->data); 381 c->freedata = PETSC_TRUE; 382 c->nonew = 0; 383 384 ierr = PetscLogEventEnd(logkey_matmatmult_symbolic,A,B,0,0);CHKERRQ(ierr); 385 PetscLogInfo((PetscObject)(*C),"Number of calls to GetMoreSpace(): %d\n",nspacedouble); 386 PetscFunctionReturn(0); 387 } 388 389 #undef __FUNCT__ 390 #define __FUNCT__ "MatMatMultNumeric" 391 /*@ 392 MatMatMultNumeric - Performs the numeric matrix-matrix product. 393 Call this routine after first calling MatMatMultSymbolic(). 394 395 Collective on Mat 396 397 Input Parameters: 398 + A - the left matrix 399 - B - the right matrix 400 401 Output Parameters: 402 . C - the product matrix, whose ij structure was defined from MatMatMultSymbolic(). 403 404 Notes: 405 C must have been created with MatMatMultSymbolic. 406 407 This routine is currently only implemented for SeqAIJ type matrices. 408 409 Level: intermediate 410 411 .seealso: MatMatMult(),MatMatMultSymbolic() 412 @*/ 413 int MatMatMultNumeric(Mat A,Mat B,Mat C){ 414 /* Perhaps this "interface" routine should be moved into the interface directory.*/ 415 /* To facilitate implementations with varying types, QueryFunction is used.*/ 416 /* It is assumed that implementations will be composed as "MatMatMultNumeric_<type of A><type of B>". */ 417 int ierr; 418 char numfunct[80]; 419 int (*numeric)(Mat,Mat,Mat); 420 421 PetscFunctionBegin; 422 423 PetscValidHeaderSpecific(A,MAT_COOKIE,1); 424 PetscValidType(A,1); 425 MatPreallocated(A); 426 if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 427 if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 428 429 PetscValidHeaderSpecific(B,MAT_COOKIE,2); 430 PetscValidType(B,2); 431 MatPreallocated(B); 432 if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 433 if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 434 435 PetscValidHeaderSpecific(C,MAT_COOKIE,3); 436 PetscValidType(C,3); 437 MatPreallocated(C); 438 if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 439 if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 440 441 if (B->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->N,C->N); 442 if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N); 443 if (A->M!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->M,C->M); 444 445 /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */ 446 /* When other implementations exist, attack the multiple dispatch problem. */ 447 ierr = PetscStrcpy(numfunct,"MatMatMultNumeric_seqaijseqaij");CHKERRQ(ierr); 448 ierr = PetscObjectQueryFunction((PetscObject)A,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr); 449 if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name); 450 ierr = PetscObjectQueryFunction((PetscObject)B,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr); 451 if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name); 452 453 ierr = (*numeric)(A,B,C);CHKERRQ(ierr); 454 455 PetscFunctionReturn(0); 456 } 457 458 #undef __FUNCT__ 459 #define __FUNCT__ "MatMatMultNumeric_MPIAIJ_MPIAIJ" 460 int MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat B,Mat C) 461 { 462 PetscFunctionBegin; 463 PetscFunctionReturn(0); 464 } 465 466 #undef __FUNCT__ 467 #define __FUNCT__ "MatMatMultNumeric_SeqAIJ_SeqAIJ" 468 int MatMatMultNumeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 469 { 470 int ierr,flops=0; 471 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 472 Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 473 Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 474 int *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 475 int am=A->M,cn=C->N; 476 int i,j,k,anzi,bnzi,cnzi,brow; 477 MatScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,*temp; 478 479 PetscFunctionBegin; 480 481 /* Start timers */ 482 ierr = PetscLogEventBegin(logkey_matmatmult_numeric,A,B,C,0);CHKERRQ(ierr); 483 484 /* Allocate temp accumulation space to avoid searching for nonzero columns in C */ 485 ierr = PetscMalloc((cn+1)*sizeof(MatScalar),&temp);CHKERRQ(ierr); 486 ierr = PetscMemzero(temp,cn*sizeof(MatScalar));CHKERRQ(ierr); 487 /* Traverse A row-wise. */ 488 /* Build the ith row in C by summing over nonzero columns in A, */ 489 /* the rows of B corresponding to nonzeros of A. */ 490 for (i=0;i<am;i++) { 491 anzi = ai[i+1] - ai[i]; 492 for (j=0;j<anzi;j++) { 493 brow = *aj++; 494 bnzi = bi[brow+1] - bi[brow]; 495 bjj = bj + bi[brow]; 496 baj = ba + bi[brow]; 497 for (k=0;k<bnzi;k++) { 498 temp[bjj[k]] += (*aa)*baj[k]; 499 } 500 flops += 2*bnzi; 501 aa++; 502 } 503 /* Store row back into C, and re-zero temp */ 504 cnzi = ci[i+1] - ci[i]; 505 for (j=0;j<cnzi;j++) { 506 ca[j] = temp[cj[j]]; 507 temp[cj[j]] = 0.0; 508 } 509 ca += cnzi; 510 cj += cnzi; 511 } 512 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 513 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 514 515 /* Free temp */ 516 ierr = PetscFree(temp);CHKERRQ(ierr); 517 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 518 ierr = PetscLogEventEnd(logkey_matmatmult_numeric,A,B,C,0);CHKERRQ(ierr); 519 PetscFunctionReturn(0); 520 } 521 522 #undef __FUNCT__ 523 #define __FUNCT__ "RegisterMatMatMultRoutines_Private" 524 int RegisterMatMatMultRoutines_Private(Mat A) { 525 int ierr; 526 527 PetscFunctionBegin; 528 if (!logkey_matmatmult_symbolic) { 529 ierr = PetscLogEventRegister(&logkey_matmatmult_symbolic,"MatMatMultSymbolic",MAT_COOKIE);CHKERRQ(ierr); 530 } 531 /* 532 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMatMultSymbolic_seqaijseqaij", 533 "MatMatMultSymbolic_SeqAIJ_SeqAIJ", 534 MatMatMultSymbolic_SeqAIJ_SeqAIJ);CHKERRQ(ierr);*/ 535 if (!logkey_matmatmult_numeric) { 536 ierr = PetscLogEventRegister(&logkey_matmatmult_numeric,"MatMatMultNumeric",MAT_COOKIE);CHKERRQ(ierr); 537 } 538 /* 539 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMatMultNumeric_seqaijseqaij", 540 "MatMatMultNumeric_SeqAIJ_SeqAIJ", 541 MatMatMultNumeric_SeqAIJ_SeqAIJ);CHKERRQ(ierr);*/ 542 PetscFunctionReturn(0); 543 } 544