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