1 /* 2 Factorization code for BAIJ format. 3 */ 4 #include "src/mat/impls/baij/seq/baij.h" 5 #include "src/inline/ilu.h" 6 7 /* ------------------------------------------------------------*/ 8 /* 9 Version for when blocks are 2 by 2 10 */ 11 #undef __FUNCT__ 12 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2" 13 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat A,MatFactorInfo *info,Mat *B) 14 { 15 Mat C = *B; 16 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 17 IS isrow = b->row,isicol = b->icol; 18 PetscErrorCode ierr; 19 PetscInt *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j; 20 PetscInt *ajtmpold,*ajtmp,nz,row; 21 PetscInt *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj; 22 MatScalar *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4; 23 MatScalar p1,p2,p3,p4; 24 MatScalar *ba = b->a,*aa = a->a; 25 26 PetscFunctionBegin; 27 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 28 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 29 ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 30 31 for (i=0; i<n; i++) { 32 nz = bi[i+1] - bi[i]; 33 ajtmp = bj + bi[i]; 34 for (j=0; j<nz; j++) { 35 x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0; 36 } 37 /* load in initial (unfactored row) */ 38 idx = r[i]; 39 nz = ai[idx+1] - ai[idx]; 40 ajtmpold = aj + ai[idx]; 41 v = aa + 4*ai[idx]; 42 for (j=0; j<nz; j++) { 43 x = rtmp+4*ic[ajtmpold[j]]; 44 x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 45 v += 4; 46 } 47 row = *ajtmp++; 48 while (row < i) { 49 pc = rtmp + 4*row; 50 p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 51 if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) { 52 pv = ba + 4*diag_offset[row]; 53 pj = bj + diag_offset[row] + 1; 54 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 55 pc[0] = m1 = p1*x1 + p3*x2; 56 pc[1] = m2 = p2*x1 + p4*x2; 57 pc[2] = m3 = p1*x3 + p3*x4; 58 pc[3] = m4 = p2*x3 + p4*x4; 59 nz = bi[row+1] - diag_offset[row] - 1; 60 pv += 4; 61 for (j=0; j<nz; j++) { 62 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 63 x = rtmp + 4*pj[j]; 64 x[0] -= m1*x1 + m3*x2; 65 x[1] -= m2*x1 + m4*x2; 66 x[2] -= m1*x3 + m3*x4; 67 x[3] -= m2*x3 + m4*x4; 68 pv += 4; 69 } 70 PetscLogFlops(16*nz+12); 71 } 72 row = *ajtmp++; 73 } 74 /* finished row so stick it into b->a */ 75 pv = ba + 4*bi[i]; 76 pj = bj + bi[i]; 77 nz = bi[i+1] - bi[i]; 78 for (j=0; j<nz; j++) { 79 x = rtmp+4*pj[j]; 80 pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 81 pv += 4; 82 } 83 /* invert diagonal block */ 84 w = ba + 4*diag_offset[i]; 85 ierr = Kernel_A_gets_inverse_A_2(w);CHKERRQ(ierr); 86 } 87 88 ierr = PetscFree(rtmp);CHKERRQ(ierr); 89 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 90 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 91 C->factor = FACTOR_LU; 92 C->assembled = PETSC_TRUE; 93 PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */ 94 PetscFunctionReturn(0); 95 } 96 /* 97 Version for when blocks are 2 by 2 Using natural ordering 98 */ 99 #undef __FUNCT__ 100 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering" 101 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B) 102 { 103 Mat C = *B; 104 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 105 PetscErrorCode ierr; 106 PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 107 PetscInt *ajtmpold,*ajtmp,nz,row; 108 PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 109 MatScalar *pv,*v,*rtmp,*pc,*w,*x; 110 MatScalar p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4; 111 MatScalar *ba = b->a,*aa = a->a; 112 113 PetscFunctionBegin; 114 ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 115 116 for (i=0; i<n; i++) { 117 nz = bi[i+1] - bi[i]; 118 ajtmp = bj + bi[i]; 119 for (j=0; j<nz; j++) { 120 x = rtmp+4*ajtmp[j]; 121 x[0] = x[1] = x[2] = x[3] = 0.0; 122 } 123 /* load in initial (unfactored row) */ 124 nz = ai[i+1] - ai[i]; 125 ajtmpold = aj + ai[i]; 126 v = aa + 4*ai[i]; 127 for (j=0; j<nz; j++) { 128 x = rtmp+4*ajtmpold[j]; 129 x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 130 v += 4; 131 } 132 row = *ajtmp++; 133 while (row < i) { 134 pc = rtmp + 4*row; 135 p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 136 if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) { 137 pv = ba + 4*diag_offset[row]; 138 pj = bj + diag_offset[row] + 1; 139 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 140 pc[0] = m1 = p1*x1 + p3*x2; 141 pc[1] = m2 = p2*x1 + p4*x2; 142 pc[2] = m3 = p1*x3 + p3*x4; 143 pc[3] = m4 = p2*x3 + p4*x4; 144 nz = bi[row+1] - diag_offset[row] - 1; 145 pv += 4; 146 for (j=0; j<nz; j++) { 147 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 148 x = rtmp + 4*pj[j]; 149 x[0] -= m1*x1 + m3*x2; 150 x[1] -= m2*x1 + m4*x2; 151 x[2] -= m1*x3 + m3*x4; 152 x[3] -= m2*x3 + m4*x4; 153 pv += 4; 154 } 155 PetscLogFlops(16*nz+12); 156 } 157 row = *ajtmp++; 158 } 159 /* finished row so stick it into b->a */ 160 pv = ba + 4*bi[i]; 161 pj = bj + bi[i]; 162 nz = bi[i+1] - bi[i]; 163 for (j=0; j<nz; j++) { 164 x = rtmp+4*pj[j]; 165 pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 166 pv += 4; 167 } 168 /* invert diagonal block */ 169 w = ba + 4*diag_offset[i]; 170 ierr = Kernel_A_gets_inverse_A_2(w);CHKERRQ(ierr); 171 /*Kernel_A_gets_inverse_A(bs,w,v_pivots,v_work);*/ 172 } 173 174 ierr = PetscFree(rtmp);CHKERRQ(ierr); 175 C->factor = FACTOR_LU; 176 C->assembled = PETSC_TRUE; 177 PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */ 178 PetscFunctionReturn(0); 179 } 180 181 /* ----------------------------------------------------------- */ 182 /* 183 Version for when blocks are 1 by 1. 184 */ 185 #undef __FUNCT__ 186 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1" 187 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat A,MatFactorInfo *info,Mat *B) 188 { 189 Mat C = *B; 190 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 191 IS isrow = b->row,isicol = b->icol; 192 PetscErrorCode ierr; 193 PetscInt *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j; 194 PetscInt *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j; 195 PetscInt *diag_offset = b->diag,diag,*pj; 196 MatScalar *pv,*v,*rtmp,multiplier,*pc; 197 MatScalar *ba = b->a,*aa = a->a; 198 199 PetscFunctionBegin; 200 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 201 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 202 ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 203 204 for (i=0; i<n; i++) { 205 nz = bi[i+1] - bi[i]; 206 ajtmp = bj + bi[i]; 207 for (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0; 208 209 /* load in initial (unfactored row) */ 210 nz = ai[r[i]+1] - ai[r[i]]; 211 ajtmpold = aj + ai[r[i]]; 212 v = aa + ai[r[i]]; 213 for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] = v[j]; 214 215 row = *ajtmp++; 216 while (row < i) { 217 pc = rtmp + row; 218 if (*pc != 0.0) { 219 pv = ba + diag_offset[row]; 220 pj = bj + diag_offset[row] + 1; 221 multiplier = *pc * *pv++; 222 *pc = multiplier; 223 nz = bi[row+1] - diag_offset[row] - 1; 224 for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 225 PetscLogFlops(1+2*nz); 226 } 227 row = *ajtmp++; 228 } 229 /* finished row so stick it into b->a */ 230 pv = ba + bi[i]; 231 pj = bj + bi[i]; 232 nz = bi[i+1] - bi[i]; 233 for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];} 234 diag = diag_offset[i] - bi[i]; 235 /* check pivot entry for current row */ 236 if (pv[diag] == 0.0) { 237 SETERRQ(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot"); 238 } 239 pv[diag] = 1.0/pv[diag]; 240 } 241 242 ierr = PetscFree(rtmp);CHKERRQ(ierr); 243 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 244 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 245 C->factor = FACTOR_LU; 246 C->assembled = PETSC_TRUE; 247 PetscLogFlops(C->n); 248 PetscFunctionReturn(0); 249 } 250 251 252 /* ----------------------------------------------------------- */ 253 #undef __FUNCT__ 254 #define __FUNCT__ "MatLUFactor_SeqBAIJ" 255 PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,MatFactorInfo *info) 256 { 257 PetscErrorCode ierr; 258 Mat C; 259 260 PetscFunctionBegin; 261 ierr = MatLUFactorSymbolic(A,row,col,info,&C);CHKERRQ(ierr); 262 ierr = MatLUFactorNumeric(A,info,&C);CHKERRQ(ierr); 263 ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 264 PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol); 265 PetscFunctionReturn(0); 266 } 267 268 #include "src/mat/impls/sbaij/seq/sbaij.h" 269 #undef __FUNCT__ 270 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N" 271 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat A,MatFactorInfo *info,Mat *B) 272 { 273 PetscErrorCode ierr; 274 Mat C = *B; 275 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 276 Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 277 IS ip=b->row; 278 PetscInt *rip,i,j,mbs=a->mbs,bs=A->bs,*bi=b->i,*bj=b->j,*bcol; 279 PetscInt *ai=a->i,*aj=a->j; 280 PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 281 MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 282 PetscReal zeropivot,rs,shiftnz; 283 PetscTruth shiftpd; 284 ChShift_Ctx sctx; 285 PetscInt newshift; 286 287 PetscFunctionBegin; 288 if (bs > 1) { 289 if (!a->sbaijMat){ 290 ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr); 291 } 292 ierr = (a->sbaijMat)->ops->choleskyfactornumeric(a->sbaijMat,info,B);CHKERRQ(ierr); 293 ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr); 294 a->sbaijMat = PETSC_NULL; 295 PetscFunctionReturn(0); 296 } 297 298 /* initialization */ 299 shiftnz = info->shiftnz; 300 shiftpd = info->shiftpd; 301 zeropivot = info->zeropivot; 302 303 ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 304 nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 305 ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 306 jl = il + mbs; 307 rtmp = (MatScalar*)(jl + mbs); 308 309 sctx.shift_amount = 0; 310 sctx.nshift = 0; 311 do { 312 sctx.chshift = PETSC_FALSE; 313 for (i=0; i<mbs; i++) { 314 rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 315 } 316 317 for (k = 0; k<mbs; k++){ 318 bval = ba + bi[k]; 319 /* initialize k-th row by the perm[k]-th row of A */ 320 jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 321 for (j = jmin; j < jmax; j++){ 322 col = rip[aj[j]]; 323 if (col >= k){ /* only take upper triangular entry */ 324 rtmp[col] = aa[j]; 325 *bval++ = 0.0; /* for in-place factorization */ 326 } 327 } 328 329 /* shift the diagonal of the matrix */ 330 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 331 332 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 333 dk = rtmp[k]; 334 i = jl[k]; /* first row to be added to k_th row */ 335 336 while (i < k){ 337 nexti = jl[i]; /* next row to be added to k_th row */ 338 339 /* compute multiplier, update diag(k) and U(i,k) */ 340 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 341 uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 342 dk += uikdi*ba[ili]; 343 ba[ili] = uikdi; /* -U(i,k) */ 344 345 /* add multiple of row i to k-th row */ 346 jmin = ili + 1; jmax = bi[i+1]; 347 if (jmin < jmax){ 348 for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 349 /* update il and jl for row i */ 350 il[i] = jmin; 351 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 352 } 353 i = nexti; 354 } 355 356 /* shift the diagonals when zero pivot is detected */ 357 /* compute rs=sum of abs(off-diagonal) */ 358 rs = 0.0; 359 jmin = bi[k]+1; 360 nz = bi[k+1] - jmin; 361 if (nz){ 362 bcol = bj + jmin; 363 while (nz--){ 364 rs += PetscAbsScalar(rtmp[*bcol++]); 365 } 366 } 367 368 sctx.rs = rs; 369 sctx.pv = dk; 370 ierr = Mat_CholeskyCheckShift(info,&sctx,&newshift);CHKERRQ(ierr); 371 if (newshift == 1){ 372 break; /* sctx.shift_amount is updated */ 373 } else if (newshift == -1){ 374 SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs); 375 } 376 377 /* copy data into U(k,:) */ 378 ba[bi[k]] = 1.0/dk; /* U(k,k) */ 379 jmin = bi[k]+1; jmax = bi[k+1]; 380 if (jmin < jmax) { 381 for (j=jmin; j<jmax; j++){ 382 col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 383 } 384 /* add the k-th row into il and jl */ 385 il[k] = jmin; 386 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 387 } 388 } 389 } while (sctx.chshift); 390 ierr = PetscFree(il);CHKERRQ(ierr); 391 392 ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 393 C->factor = FACTOR_CHOLESKY; 394 C->assembled = PETSC_TRUE; 395 C->preallocated = PETSC_TRUE; 396 PetscLogFlops(C->m); 397 if (sctx.nshift){ 398 if (shiftnz) { 399 PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqBAIJ_1: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount); 400 } else if (shiftpd) { 401 PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqBAIJ_1: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount); 402 } 403 } 404 PetscFunctionReturn(0); 405 } 406 407 #undef __FUNCT__ 408 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering" 409 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *fact) 410 { 411 Mat C = *fact; 412 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 413 Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 414 PetscErrorCode ierr; 415 PetscInt i,j,am=a->mbs; 416 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 417 PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz; 418 MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 419 PetscReal zeropivot,rs,shiftnz; 420 PetscTruth shiftpd; 421 ChShift_Ctx sctx; 422 PetscInt newshift; 423 424 PetscFunctionBegin; 425 /* initialization */ 426 shiftnz = info->shiftnz; 427 shiftpd = info->shiftpd; 428 zeropivot = info->zeropivot; 429 430 nz = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar); 431 ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 432 jl = il + am; 433 rtmp = (MatScalar*)(jl + am); 434 435 sctx.shift_amount = 0; 436 sctx.nshift = 0; 437 do { 438 sctx.chshift = PETSC_FALSE; 439 for (i=0; i<am; i++) { 440 rtmp[i] = 0.0; jl[i] = am; il[0] = 0; 441 } 442 443 for (k = 0; k<am; k++){ 444 /* initialize k-th row with elements nonzero in row perm(k) of A */ 445 nz = ai[k+1] - ai[k]; 446 acol = aj + ai[k]; 447 aval = aa + ai[k]; 448 bval = ba + bi[k]; 449 while (nz -- ){ 450 if (*acol < k) { /* skip lower triangular entries */ 451 acol++; aval++; 452 } else { 453 rtmp[*acol++] = *aval++; 454 *bval++ = 0.0; /* for in-place factorization */ 455 } 456 } 457 458 /* shift the diagonal of the matrix */ 459 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 460 461 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 462 dk = rtmp[k]; 463 i = jl[k]; /* first row to be added to k_th row */ 464 465 while (i < k){ 466 nexti = jl[i]; /* next row to be added to k_th row */ 467 /* compute multiplier, update D(k) and U(i,k) */ 468 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 469 uikdi = - ba[ili]*ba[bi[i]]; 470 dk += uikdi*ba[ili]; 471 ba[ili] = uikdi; /* -U(i,k) */ 472 473 /* add multiple of row i to k-th row ... */ 474 jmin = ili + 1; 475 nz = bi[i+1] - jmin; 476 if (nz > 0){ 477 bcol = bj + jmin; 478 bval = ba + jmin; 479 while (nz --) rtmp[*bcol++] += uikdi*(*bval++); 480 /* update il and jl for i-th row */ 481 il[i] = jmin; 482 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 483 } 484 i = nexti; 485 } 486 487 /* shift the diagonals when zero pivot is detected */ 488 /* compute rs=sum of abs(off-diagonal) */ 489 rs = 0.0; 490 jmin = bi[k]+1; 491 nz = bi[k+1] - jmin; 492 if (nz){ 493 bcol = bj + jmin; 494 while (nz--){ 495 rs += PetscAbsScalar(rtmp[*bcol++]); 496 } 497 } 498 499 sctx.rs = rs; 500 sctx.pv = dk; 501 ierr = Mat_CholeskyCheckShift(info,&sctx,&newshift);CHKERRQ(ierr); 502 if (newshift == 1){ 503 break; /* sctx.shift_amount is updated */ 504 } else if (newshift == -1){ 505 SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs); 506 } 507 508 /* copy data into U(k,:) */ 509 ba[bi[k]] = 1.0/dk; 510 jmin = bi[k]+1; 511 nz = bi[k+1] - jmin; 512 if (nz){ 513 bcol = bj + jmin; 514 bval = ba + jmin; 515 while (nz--){ 516 *bval++ = rtmp[*bcol]; 517 rtmp[*bcol++] = 0.0; 518 } 519 /* add k-th row into il and jl */ 520 il[k] = jmin; 521 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 522 } 523 } 524 } while (sctx.chshift); 525 ierr = PetscFree(il);CHKERRQ(ierr); 526 527 C->factor = FACTOR_CHOLESKY; 528 C->assembled = PETSC_TRUE; 529 C->preallocated = PETSC_TRUE; 530 PetscLogFlops(C->m); 531 if (sctx.nshift){ 532 if (shiftnz) { 533 PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqBAIJ_1_NaturalOrdering: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount); 534 } else if (shiftpd) { 535 PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqBAIJ_1_NaturalOrdering: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount); 536 } 537 } 538 PetscFunctionReturn(0); 539 } 540 541 #include "petscbt.h" 542 #include "src/mat/utils/freespace.h" 543 #undef __FUNCT__ 544 #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ" 545 PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 546 { 547 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 548 Mat_SeqSBAIJ *b; 549 Mat B; 550 PetscErrorCode ierr; 551 PetscTruth perm_identity; 552 PetscInt reallocs=0,*rip,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->bs,*ui; 553 PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 554 PetscInt nlnk,*lnk,*lnk_lvl,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr; 555 PetscReal fill=info->fill,levels=info->levels; 556 FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 557 FreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 558 PetscBT lnkbt; 559 560 PetscFunctionBegin; 561 if (bs > 1){ 562 if (!a->sbaijMat){ 563 ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr); 564 } 565 ierr = MatICCFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr); 566 B = *fact; 567 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 568 PetscFunctionReturn(0); 569 } 570 571 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 572 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 573 574 /* special case that simply copies fill pattern */ 575 if (!levels && perm_identity) { 576 ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 577 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 578 for (i=0; i<am; i++) { 579 ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */ 580 } 581 ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr); 582 B = *fact; 583 ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 584 ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr); 585 586 b = (Mat_SeqSBAIJ*)B->data; 587 uj = b->j; 588 for (i=0; i<am; i++) { 589 aj = a->j + a->diag[i]; 590 for (j=0; j<ui[i]; j++){ 591 *uj++ = *aj++; 592 } 593 b->ilen[i] = ui[i]; 594 } 595 ierr = PetscFree(ui);CHKERRQ(ierr); 596 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 597 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 598 599 B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 600 B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 601 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 602 PetscFunctionReturn(0); 603 } 604 605 /* initialization */ 606 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 607 ui[0] = 0; 608 ierr = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr); 609 610 /* jl: linked list for storing indices of the pivot rows 611 il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 612 ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 613 il = jl + am; 614 uj_ptr = (PetscInt**)(il + am); 615 uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 616 for (i=0; i<am; i++){ 617 jl[i] = am; il[i] = 0; 618 } 619 620 /* create and initialize a linked list for storing column indices of the active row k */ 621 nlnk = am + 1; 622 ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 623 624 /* initial FreeSpace size is fill*(ai[am]+1) */ 625 ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 626 current_space = free_space; 627 ierr = GetMoreSpace((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 628 current_space_lvl = free_space_lvl; 629 630 for (k=0; k<am; k++){ /* for each active row k */ 631 /* initialize lnk by the column indices of row rip[k] of A */ 632 nzk = 0; 633 ncols = ai[rip[k]+1] - ai[rip[k]]; 634 ncols_upper = 0; 635 cols = cols_lvl + am; 636 for (j=0; j<ncols; j++){ 637 i = rip[*(aj + ai[rip[k]] + j)]; 638 if (i >= k){ /* only take upper triangular entry */ 639 cols[ncols_upper] = i; 640 cols_lvl[ncols_upper] = -1; /* initialize level for nonzero entries */ 641 ncols_upper++; 642 } 643 } 644 ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 645 nzk += nlnk; 646 647 /* update lnk by computing fill-in for each pivot row to be merged in */ 648 prow = jl[k]; /* 1st pivot row */ 649 650 while (prow < k){ 651 nextprow = jl[prow]; 652 653 /* merge prow into k-th row */ 654 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 655 jmax = ui[prow+1]; 656 ncols = jmax-jmin; 657 i = jmin - ui[prow]; 658 cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 659 for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j); 660 ierr = PetscIncompleteLLAdd(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 661 nzk += nlnk; 662 663 /* update il and jl for prow */ 664 if (jmin < jmax){ 665 il[prow] = jmin; 666 j = *cols; jl[prow] = jl[j]; jl[j] = prow; 667 } 668 prow = nextprow; 669 } 670 671 /* if free space is not available, make more free space */ 672 if (current_space->local_remaining<nzk) { 673 i = am - k + 1; /* num of unfactored rows */ 674 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 675 ierr = GetMoreSpace(i,¤t_space);CHKERRQ(ierr); 676 ierr = GetMoreSpace(i,¤t_space_lvl);CHKERRQ(ierr); 677 reallocs++; 678 } 679 680 /* copy data into free_space and free_space_lvl, then initialize lnk */ 681 ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 682 683 /* add the k-th row into il and jl */ 684 if (nzk-1 > 0){ 685 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 686 jl[k] = jl[i]; jl[i] = k; 687 il[k] = ui[k] + 1; 688 } 689 uj_ptr[k] = current_space->array; 690 uj_lvl_ptr[k] = current_space_lvl->array; 691 692 current_space->array += nzk; 693 current_space->local_used += nzk; 694 current_space->local_remaining -= nzk; 695 696 current_space_lvl->array += nzk; 697 current_space_lvl->local_used += nzk; 698 current_space_lvl->local_remaining -= nzk; 699 700 ui[k+1] = ui[k] + nzk; 701 } 702 703 if (ai[am] != 0) { 704 PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]); 705 PetscLogInfo(A,"MatICCFactorSymbolic_SeqBAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af); 706 PetscLogInfo(A,"MatICCFactorSymbolic_SeqBAIJ:Run with -pc_cholesky_fill %g or use \n",af); 707 PetscLogInfo(A,"MatICCFactorSymbolic_SeqBAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af); 708 } else { 709 PetscLogInfo(A,"MatICCFactorSymbolic_SeqBAIJ:Empty matrix.\n"); 710 } 711 712 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 713 ierr = PetscFree(jl);CHKERRQ(ierr); 714 ierr = PetscFree(cols_lvl);CHKERRQ(ierr); 715 716 /* destroy list of free space and other temporary array(s) */ 717 ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 718 ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 719 ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 720 ierr = DestroySpace(free_space_lvl);CHKERRQ(ierr); 721 722 /* put together the new matrix in MATSEQSBAIJ format */ 723 ierr = MatCreate(PETSC_COMM_SELF,am,am,am,am,fact);CHKERRQ(ierr); 724 B = *fact; 725 ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 726 ierr = MatSeqSBAIJSetPreallocation(B,1,0,PETSC_NULL);CHKERRQ(ierr); 727 728 b = (Mat_SeqSBAIJ*)B->data; 729 ierr = PetscFree(b->imax);CHKERRQ(ierr); 730 b->singlemalloc = PETSC_FALSE; 731 /* the next line frees the default space generated by the Create() */ 732 ierr = PetscFree(b->a);CHKERRQ(ierr); 733 ierr = PetscFree(b->ilen);CHKERRQ(ierr); 734 ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 735 b->j = uj; 736 b->i = ui; 737 b->diag = 0; 738 b->ilen = 0; 739 b->imax = 0; 740 b->row = perm; 741 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 742 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 743 b->icol = perm; 744 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 745 ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 746 PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar))); 747 b->maxnz = b->nz = ui[am]; 748 749 B->factor = FACTOR_CHOLESKY; 750 B->info.factor_mallocs = reallocs; 751 B->info.fill_ratio_given = fill; 752 if (ai[am] != 0) { 753 B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 754 } else { 755 B->info.fill_ratio_needed = 0.0; 756 } 757 if (perm_identity){ 758 B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 759 B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 760 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 761 } else { 762 (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 763 } 764 PetscFunctionReturn(0); 765 } 766 767 #undef __FUNCT__ 768 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ" 769 PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 770 { 771 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 772 Mat_SeqSBAIJ *b; 773 Mat B; 774 PetscErrorCode ierr; 775 PetscTruth perm_identity; 776 PetscReal fill = info->fill; 777 PetscInt *rip,*riip,i,mbs=a->mbs,bs=A->bs,*ai=a->i,*aj=a->j,reallocs=0,prow; 778 PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 779 PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 780 FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 781 PetscBT lnkbt; 782 IS iperm; 783 784 PetscFunctionBegin; 785 if (bs > 1) { /* convert to seqsbaij */ 786 if (!a->sbaijMat){ 787 ierr = MatConvert(A,MATSEQSBAIJ,&a->sbaijMat);CHKERRQ(ierr); 788 } 789 ierr = MatCholeskyFactorSymbolic(a->sbaijMat,perm,info,fact);CHKERRQ(ierr); 790 B = *fact; 791 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 792 PetscFunctionReturn(0); 793 } 794 795 /* check whether perm is the identity mapping */ 796 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 797 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 798 799 if (!perm_identity){ 800 /* check if perm is symmetric! */ 801 ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 802 ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 803 for (i=0; i<mbs; i++) { 804 if (rip[i] != riip[i]) SETERRQ(PETSC_ERR_ARG_INCOMP,"Non-symmetric permutation, must use symmetric permutation"); 805 } 806 ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr); 807 ierr = ISDestroy(iperm);CHKERRQ(ierr); 808 } 809 810 /* initialization */ 811 ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 812 ui[0] = 0; 813 814 /* jl: linked list for storing indices of the pivot rows 815 il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */ 816 ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 817 il = jl + mbs; 818 cols = il + mbs; 819 ui_ptr = (PetscInt**)(cols + mbs); 820 for (i=0; i<mbs; i++){ 821 jl[i] = mbs; il[i] = 0; 822 } 823 824 /* create and initialize a linked list for storing column indices of the active row k */ 825 nlnk = mbs + 1; 826 ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 827 828 /* initial FreeSpace size is fill*(ai[mbs]+1) */ 829 ierr = GetMoreSpace((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr); 830 current_space = free_space; 831 832 for (k=0; k<mbs; k++){ /* for each active row k */ 833 /* initialize lnk by the column indices of row rip[k] of A */ 834 nzk = 0; 835 ncols = ai[rip[k]+1] - ai[rip[k]]; 836 ncols_upper = 0; 837 for (j=0; j<ncols; j++){ 838 i = rip[*(aj + ai[rip[k]] + j)]; 839 if (i >= k){ /* only take upper triangular entry */ 840 cols[ncols_upper] = i; 841 ncols_upper++; 842 } 843 } 844 ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 845 nzk += nlnk; 846 847 /* update lnk by computing fill-in for each pivot row to be merged in */ 848 prow = jl[k]; /* 1st pivot row */ 849 850 while (prow < k){ 851 nextprow = jl[prow]; 852 /* merge prow into k-th row */ 853 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */ 854 jmax = ui[prow+1]; 855 ncols = jmax-jmin; 856 uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */ 857 ierr = PetscLLAdd(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 858 nzk += nlnk; 859 860 /* update il and jl for prow */ 861 if (jmin < jmax){ 862 il[prow] = jmin; 863 j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 864 } 865 prow = nextprow; 866 } 867 868 /* if free space is not available, make more free space */ 869 if (current_space->local_remaining<nzk) { 870 i = mbs - k + 1; /* num of unfactored rows */ 871 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 872 ierr = GetMoreSpace(i,¤t_space);CHKERRQ(ierr); 873 reallocs++; 874 } 875 876 /* copy data into free space, then initialize lnk */ 877 ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 878 879 /* add the k-th row into il and jl */ 880 if (nzk-1 > 0){ 881 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 882 jl[k] = jl[i]; jl[i] = k; 883 il[k] = ui[k] + 1; 884 } 885 ui_ptr[k] = current_space->array; 886 current_space->array += nzk; 887 current_space->local_used += nzk; 888 current_space->local_remaining -= nzk; 889 890 ui[k+1] = ui[k] + nzk; 891 } 892 893 if (ai[mbs] != 0) { 894 PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 895 PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Reallocs %D Fill ratio:given %g needed %g\n",reallocs,fill,af); 896 PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Run with -pc_cholesky_fill %g or use \n",af); 897 PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:PCCholeskySetFill(pc,%g) for best performance.\n",af); 898 } else { 899 PetscLogInfo(A,"MatCholeskyFactorSymbolic_SeqSBAIJ:Empty matrix.\n"); 900 } 901 902 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 903 ierr = PetscFree(jl);CHKERRQ(ierr); 904 905 /* destroy list of free space and other temporary array(s) */ 906 ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 907 ierr = MakeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 908 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 909 910 /* put together the new matrix in MATSEQSBAIJ format */ 911 ierr = MatCreate(PETSC_COMM_SELF,mbs,mbs,mbs,mbs,fact);CHKERRQ(ierr); 912 B = *fact; 913 ierr = MatSetType(B,MATSEQSBAIJ);CHKERRQ(ierr); 914 ierr = MatSeqSBAIJSetPreallocation(B,bs,0,PETSC_NULL);CHKERRQ(ierr); 915 916 b = (Mat_SeqSBAIJ*)B->data; 917 ierr = PetscFree(b->imax);CHKERRQ(ierr); 918 b->singlemalloc = PETSC_FALSE; 919 /* the next line frees the default space generated by the Create() */ 920 ierr = PetscFree(b->a);CHKERRQ(ierr); 921 ierr = PetscFree(b->ilen);CHKERRQ(ierr); 922 ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 923 b->j = uj; 924 b->i = ui; 925 b->diag = 0; 926 b->ilen = 0; 927 b->imax = 0; 928 b->row = perm; 929 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 930 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 931 b->icol = perm; 932 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 933 ierr = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 934 PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar))); 935 b->maxnz = b->nz = ui[mbs]; 936 937 B->factor = FACTOR_CHOLESKY; 938 B->info.factor_mallocs = reallocs; 939 B->info.fill_ratio_given = fill; 940 if (ai[mbs] != 0) { 941 B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 942 } else { 943 B->info.fill_ratio_needed = 0.0; 944 } 945 if (perm_identity){ 946 B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 947 B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 948 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 949 } else { 950 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 951 } 952 PetscFunctionReturn(0); 953 } 954