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