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