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