1 /* 2 Defines matrix-matrix product routines for pairs of SeqAIJ 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 9 static int logkey_matmatmult = 0; 10 static int logkey_matmatmult_symbolic = 0; 11 static int logkey_matmatmult_numeric = 0; 12 13 #undef __FUNCT__ 14 #define __FUNCT__ "MatMatMult" 15 /*@ 16 MatMatMult - Performs Matrix-Matrix Multiplication C=A*B. 17 18 Collective on Mat 19 20 Input Parameters: 21 + A - the left matrix 22 - B - the right matrix 23 24 Output Parameters: 25 . C - the product matrix 26 27 Notes: 28 C will be created and must be destroyed by the user with MatDestroy(). 29 30 This routine is currently only implemented for pairs of SeqAIJ matrices and classes 31 which inherit from SeqAIJ. C will be of type MATSEQAIJ. 32 33 Level: intermediate 34 35 .seealso: MatMatMultSymbolic(),MatMatMultNumeric() 36 @*/ 37 int MatMatMult(Mat A,Mat B, Mat *C) { 38 /* Perhaps this "interface" routine should be moved into the interface directory.*/ 39 /* To facilitate implementations with varying types, QueryFunction is used.*/ 40 /* It is assumed that implementations will be composed as "MatMatMult_<type of A><type of B>". */ 41 int ierr; 42 char funct[80]; 43 int (*mult)(Mat,Mat,Mat*); 44 45 PetscFunctionBegin; 46 PetscValidHeaderSpecific(A,MAT_COOKIE,1); 47 PetscValidType(A,1); 48 MatPreallocated(A); 49 if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 50 if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 51 52 PetscValidHeaderSpecific(B,MAT_COOKIE,2); 53 PetscValidType(B,2); 54 MatPreallocated(B); 55 if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 56 if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 57 58 PetscValidPointer(C,3); 59 60 if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N); 61 62 /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */ 63 /* When other implementations exist, attack the multiple dispatch problem. */ 64 ierr = PetscStrcpy(funct,"MatMatMult_seqaijseqaij");CHKERRQ(ierr); 65 ierr = PetscObjectQueryFunction((PetscObject)B,funct,(PetscVoidFunction)&mult);CHKERRQ(ierr); 66 if (!mult) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name); 67 ierr = PetscObjectQueryFunction((PetscObject)A,funct,(PetscVoidFunction)&mult);CHKERRQ(ierr); 68 if (!mult) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name); 69 70 ierr = (*mult)(A,B,C);CHKERRQ(ierr); 71 PetscFunctionReturn(0); 72 } 73 74 #undef __FUNCT__ 75 #define __FUNCT__ "MatMatMult_SeqAIJ_SeqAIJ" 76 int MatMatMult_SeqAIJ_SeqAIJ(Mat A,Mat B, Mat *C) { 77 int ierr; 78 char symfunct[80],numfunct[80]; 79 int (*symbolic)(Mat,Mat,Mat*),(*numeric)(Mat,Mat,Mat); 80 81 PetscFunctionBegin; 82 83 /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */ 84 /* When other implementations exist, attack the multiple dispatch problem. */ 85 ierr = PetscStrcpy(symfunct,"MatMatMultSymbolic_seqaijseqaij");CHKERRQ(ierr); 86 ierr = PetscObjectQueryFunction((PetscObject)B,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr); 87 if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name); 88 ierr = PetscObjectQueryFunction((PetscObject)A,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr); 89 if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name); 90 91 ierr = PetscStrcpy(numfunct,"MatMatMultNumeric_seqaijseqaij");CHKERRQ(ierr); 92 ierr = PetscObjectQueryFunction((PetscObject)A,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr); 93 if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name); 94 ierr = PetscObjectQueryFunction((PetscObject)B,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr); 95 if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name); 96 97 ierr = PetscLogEventBegin(logkey_matmatmult,A,B,0,0);CHKERRQ(ierr); 98 ierr = (*symbolic)(A,B,C);CHKERRQ(ierr); 99 ierr = (*numeric)(A,B,*C);CHKERRQ(ierr); 100 ierr = PetscLogEventEnd(logkey_matmatmult,A,B,0,0);CHKERRQ(ierr); 101 PetscFunctionReturn(0); 102 } 103 104 #undef __FUNCT__ 105 #define __FUNCT__ "MatMatMultSymbolic" 106 /*@ 107 MatMatMultSymbolic - Performs construction, preallocation, and computes the ij structure 108 of the matrix-matrix product C=A*B. Call this routine before calling MatMatMultNumeric(). 109 110 Collective on Mat 111 112 Input Parameters: 113 + A - the left matrix 114 - B - the right matrix 115 116 Output Parameters: 117 . C - the matrix containing the ij structure of product matrix 118 119 Notes: 120 C will be created as a MATSEQAIJ matrix and must be destroyed by the user with MatDestroy(). 121 122 This routine is currently only implemented for SeqAIJ matrices and classes which inherit from SeqAIJ. 123 124 Level: intermediate 125 126 .seealso: MatMatMult(),MatMatMultNumeric() 127 @*/ 128 int MatMatMultSymbolic(Mat A,Mat B,Mat *C) { 129 /* Perhaps this "interface" routine should be moved into the interface directory.*/ 130 /* To facilitate implementations with varying types, QueryFunction is used.*/ 131 /* It is assumed that implementations will be composed as "MatMatMultSymbolic_<type of A><type of B>". */ 132 int ierr; 133 char symfunct[80]; 134 int (*symbolic)(Mat,Mat,Mat *); 135 136 PetscFunctionBegin; 137 PetscValidHeaderSpecific(A,MAT_COOKIE,1); 138 PetscValidType(A,1); 139 MatPreallocated(A); 140 if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 141 if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 142 143 PetscValidHeaderSpecific(B,MAT_COOKIE,2); 144 PetscValidType(B,2); 145 MatPreallocated(B); 146 if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 147 if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 148 PetscValidPointer(C,3); 149 150 151 if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N); 152 153 /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */ 154 /* When other implementations exist, attack the multiple dispatch problem. */ 155 ierr = PetscStrcpy(symfunct,"MatMatMultSymbolic_seqaijseqaij");CHKERRQ(ierr); 156 ierr = PetscObjectQueryFunction((PetscObject)B,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr); 157 if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name); 158 ierr = PetscObjectQueryFunction((PetscObject)A,symfunct,(PetscVoidFunction)&symbolic);CHKERRQ(ierr); 159 if (!symbolic) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name); 160 ierr = (*symbolic)(A,B,C);CHKERRQ(ierr); 161 162 PetscFunctionReturn(0); 163 } 164 165 #undef __FUNCT__ 166 #define __FUNCT__ "MatMatMult_Symbolic_SeqAIJ_SeqAIJ" 167 int MatMatMult_Symbolic_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat *C) 168 { 169 int ierr; 170 FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 171 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ*)B->data,*c; 172 int *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj; 173 int *ci,*cj,*lnk,idx0,idx,bcol; 174 int am=A->M,bn=B->N,bm=B->M; 175 int i,j,k,anzi,brow,bnzj,cnzi; 176 MatScalar *ca; 177 178 PetscFunctionBegin; 179 /* Start timers */ 180 ierr = PetscLogEventBegin(logkey_matmatmult_symbolic,A,B,0,0);CHKERRQ(ierr); 181 182 /* Set up */ 183 /* Allocate ci array, arrays for fill computation and */ 184 /* free space for accumulating nonzero column info */ 185 ierr = PetscMalloc(((am+1)+1)*sizeof(int),&ci);CHKERRQ(ierr); 186 ci[0] = 0; 187 188 ierr = PetscMalloc((bn+1)*sizeof(int),&lnk);CHKERRQ(ierr); 189 for (i=0; i<bn; i++) lnk[i] = -1; 190 191 /* Initial FreeSpace size is nnz(B)=4*bi[bm] */ 192 ierr = GetMoreSpace(4*bi[bm],&free_space);CHKERRQ(ierr); 193 current_space = free_space; 194 195 /* Determine symbolic info for each row of the product: */ 196 for (i=0;i<am;i++) { 197 anzi = ai[i+1] - ai[i]; 198 cnzi = 0; 199 lnk[bn] = bn; 200 for (j=0;j<anzi;j++) { 201 brow = *aj++; 202 bnzj = bi[brow+1] - bi[brow]; 203 bjj = bj + bi[brow]; 204 idx = bn; 205 for (k=0;k<bnzj;k++) { 206 bcol = bjj[k]; 207 if (lnk[bcol] == -1) { /* new col */ 208 if (k>0) idx = bjj[k-1]; 209 do { 210 idx0 = idx; 211 idx = lnk[idx0]; 212 } while (bcol > idx); 213 lnk[idx0] = bcol; 214 lnk[bcol] = idx; 215 cnzi++; 216 } 217 } 218 } 219 220 /* If free space is not available, make more free space */ 221 /* Double the amount of total space in the list */ 222 if (current_space->local_remaining<cnzi) { 223 printf("...%d -th row, double space ...\n",i); 224 ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 225 } 226 227 /* Copy data into free space, and zero out denserow and lnk */ 228 idx = bn; 229 for (j=0; j<cnzi; j++){ 230 idx0 = idx; 231 idx = lnk[idx0]; 232 *current_space->array++ = idx; 233 lnk[idx0] = -1; 234 } 235 lnk[idx] = -1; 236 237 current_space->local_used += cnzi; 238 current_space->local_remaining -= cnzi; 239 240 ci[i+1] = ci[i] + cnzi; 241 } 242 243 /* Column indices are in the list of free space */ 244 /* Allocate space for cj, initialize cj, and */ 245 /* destroy list of free space and other temporary array(s) */ 246 ierr = PetscMalloc((ci[am]+1)*sizeof(int),&cj);CHKERRQ(ierr); 247 ierr = MakeSpaceContiguous(&free_space,cj);CHKERRQ(ierr); 248 ierr = PetscFree(lnk);CHKERRQ(ierr); 249 250 /* Allocate space for ca */ 251 ierr = PetscMalloc((ci[am]+1)*sizeof(MatScalar),&ca);CHKERRQ(ierr); 252 ierr = PetscMemzero(ca,(ci[am]+1)*sizeof(MatScalar));CHKERRQ(ierr); 253 254 /* put together the new matrix */ 255 ierr = MatCreateSeqAIJWithArrays(A->comm,am,bn,ci,cj,ca,C);CHKERRQ(ierr); 256 257 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 258 /* These are PETSc arrays, so change flags so arrays can be deleted by PETSc */ 259 c = (Mat_SeqAIJ *)((*C)->data); 260 c->freedata = PETSC_TRUE; 261 c->nonew = 0; 262 263 ierr = PetscLogEventEnd(logkey_matmatmult_symbolic,A,B,0,0);CHKERRQ(ierr); 264 PetscFunctionReturn(0); 265 } 266 267 #undef __FUNCT__ 268 #define __FUNCT__ "MatMatMultNumeric" 269 /*@ 270 MatMatMultNumeric - Performs the numeric matrix-matrix product. 271 Call this routine after first calling MatMatMultSymbolic(). 272 273 Collective on Mat 274 275 Input Parameters: 276 + A - the left matrix 277 - B - the right matrix 278 279 Output Parameters: 280 . C - the product matrix, whose ij structure was defined from MatMatMultSymbolic(). 281 282 Notes: 283 C must have been created with MatMatMultSymbolic. 284 285 This routine is currently only implemented for SeqAIJ type matrices. 286 287 Level: intermediate 288 289 .seealso: MatMatMult(),MatMatMultSymbolic() 290 @*/ 291 int MatMatMultNumeric(Mat A,Mat B,Mat C){ 292 /* Perhaps this "interface" routine should be moved into the interface directory.*/ 293 /* To facilitate implementations with varying types, QueryFunction is used.*/ 294 /* It is assumed that implementations will be composed as "MatMatMultNumeric_<type of A><type of B>". */ 295 int ierr; 296 char numfunct[80]; 297 int (*numeric)(Mat,Mat,Mat); 298 299 PetscFunctionBegin; 300 301 PetscValidHeaderSpecific(A,MAT_COOKIE,1); 302 PetscValidType(A,1); 303 MatPreallocated(A); 304 if (!A->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 305 if (A->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 306 307 PetscValidHeaderSpecific(B,MAT_COOKIE,2); 308 PetscValidType(B,2); 309 MatPreallocated(B); 310 if (!B->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 311 if (B->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 312 313 PetscValidHeaderSpecific(C,MAT_COOKIE,3); 314 PetscValidType(C,3); 315 MatPreallocated(C); 316 if (!C->assembled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for unassembled matrix"); 317 if (C->factor) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Not for factored matrix"); 318 319 if (B->N!=C->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->N,C->N); 320 if (B->M!=A->N) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",B->M,A->N); 321 if (A->M!=C->M) SETERRQ2(PETSC_ERR_ARG_SIZ,"Matrix dimensions are incompatible, %d != %d",A->M,C->M); 322 323 /* Currently only _seqaijseqaij is implemented, so just query for it in A and B. */ 324 /* When other implementations exist, attack the multiple dispatch problem. */ 325 ierr = PetscStrcpy(numfunct,"MatMatMultNumeric_seqaijseqaij");CHKERRQ(ierr); 326 ierr = PetscObjectQueryFunction((PetscObject)A,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr); 327 if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for A of type %s",A->type_name); 328 ierr = PetscObjectQueryFunction((PetscObject)B,numfunct,(PetscVoidFunction)&numeric);CHKERRQ(ierr); 329 if (!numeric) SETERRQ1(PETSC_ERR_SUP,"C=A*B not implemented for B of type %s",B->type_name); 330 331 ierr = (*numeric)(A,B,C);CHKERRQ(ierr); 332 333 PetscFunctionReturn(0); 334 } 335 336 #undef __FUNCT__ 337 #define __FUNCT__ "MatMatMult_Numeric_SeqAIJ_SeqAIJ" 338 int MatMatMult_Numeric_SeqAIJ_SeqAIJ(Mat A,Mat B,Mat C) 339 { 340 int ierr,flops=0; 341 Mat_SeqAIJ *a = (Mat_SeqAIJ *)A->data; 342 Mat_SeqAIJ *b = (Mat_SeqAIJ *)B->data; 343 Mat_SeqAIJ *c = (Mat_SeqAIJ *)C->data; 344 int *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bjj,*ci=c->i,*cj=c->j; 345 int am=A->M,cn=C->N; 346 int i,j,k,anzi,bnzi,cnzi,brow; 347 MatScalar *aa=a->a,*ba=b->a,*baj,*ca=c->a,*temp; 348 349 PetscFunctionBegin; 350 351 /* Start timers */ 352 ierr = PetscLogEventBegin(logkey_matmatmult_numeric,A,B,C,0);CHKERRQ(ierr); 353 354 /* Allocate temp accumulation space to avoid searching for nonzero columns in C */ 355 ierr = PetscMalloc((cn+1)*sizeof(MatScalar),&temp);CHKERRQ(ierr); 356 ierr = PetscMemzero(temp,cn*sizeof(MatScalar));CHKERRQ(ierr); 357 /* Traverse A row-wise. */ 358 /* Build the ith row in C by summing over nonzero columns in A, */ 359 /* the rows of B corresponding to nonzeros of A. */ 360 for (i=0;i<am;i++) { 361 anzi = ai[i+1] - ai[i]; 362 for (j=0;j<anzi;j++) { 363 brow = *aj++; 364 bnzi = bi[brow+1] - bi[brow]; 365 bjj = bj + bi[brow]; 366 baj = ba + bi[brow]; 367 for (k=0;k<bnzi;k++) { 368 temp[bjj[k]] += (*aa)*baj[k]; 369 } 370 flops += 2*bnzi; 371 aa++; 372 } 373 /* Store row back into C, and re-zero temp */ 374 cnzi = ci[i+1] - ci[i]; 375 for (j=0;j<cnzi;j++) { 376 ca[j] = temp[cj[j]]; 377 temp[cj[j]] = 0.0; 378 } 379 ca += cnzi; 380 cj += cnzi; 381 } 382 ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 383 ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 384 385 /* Free temp */ 386 ierr = PetscFree(temp);CHKERRQ(ierr); 387 ierr = PetscLogFlops(flops);CHKERRQ(ierr); 388 ierr = PetscLogEventEnd(logkey_matmatmult_numeric,A,B,C,0);CHKERRQ(ierr); 389 PetscFunctionReturn(0); 390 } 391 392 #undef __FUNCT__ 393 #define __FUNCT__ "RegisterMatMatMultRoutines_Private" 394 int RegisterMatMatMultRoutines_Private(Mat A) { 395 int ierr; 396 397 PetscFunctionBegin; 398 if (!logkey_matmatmult) { 399 ierr = PetscLogEventRegister(&logkey_matmatmult,"MatMatMult",MAT_COOKIE);CHKERRQ(ierr); 400 } 401 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMatMult_seqaijseqaij", 402 "MatMatMult_SeqAIJ_SeqAIJ", 403 MatMatMult_SeqAIJ_SeqAIJ);CHKERRQ(ierr); 404 if (!logkey_matmatmult_symbolic) { 405 ierr = PetscLogEventRegister(&logkey_matmatmult_symbolic,"MatMatMult_Symbolic",MAT_COOKIE);CHKERRQ(ierr); 406 } 407 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMatMultSymbolic_seqaijseqaij", 408 "MatMatMult_Symbolic_SeqAIJ_SeqAIJ", 409 MatMatMult_Symbolic_SeqAIJ_SeqAIJ);CHKERRQ(ierr); 410 if (!logkey_matmatmult_numeric) { 411 ierr = PetscLogEventRegister(&logkey_matmatmult_numeric,"MatMatMult_Numeric",MAT_COOKIE);CHKERRQ(ierr); 412 } 413 ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatMatMultNumeric_seqaijseqaij", 414 "MatMatMult_Numeric_SeqAIJ_SeqAIJ", 415 MatMatMult_Numeric_SeqAIJ_SeqAIJ);CHKERRQ(ierr); 416 PetscFunctionReturn(0); 417 } 418