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