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