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