1 2 /* 3 Factorization code for BAIJ format. 4 */ 5 #include <../src/mat/impls/baij/seq/baij.h> 6 #include <../src/mat/blockinvert.h> 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2" 10 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat B,Mat A,const MatFactorInfo *info) 11 { 12 Mat C=B; 13 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 14 IS isrow = b->row,isicol = b->icol; 15 PetscErrorCode ierr; 16 const PetscInt *r,*ic; 17 PetscInt i,j,k,nz,nzL,row,*pj; 18 const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,bs2=a->bs2; 19 const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag; 20 MatScalar *rtmp,*pc,*mwork,*pv; 21 MatScalar *aa=a->a,*v; 22 PetscInt flg; 23 PetscReal shift = info->shiftamount; 24 25 PetscFunctionBegin; 26 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 27 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 28 29 /* generate work space needed by the factorization */ 30 ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr); 31 ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 32 33 for (i=0; i<n; i++){ 34 /* zero rtmp */ 35 /* L part */ 36 nz = bi[i+1] - bi[i]; 37 bjtmp = bj + bi[i]; 38 for (j=0; j<nz; j++){ 39 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 40 } 41 42 /* U part */ 43 nz = bdiag[i] - bdiag[i+1]; 44 bjtmp = bj + bdiag[i+1]+1; 45 for (j=0; j<nz; j++){ 46 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 47 } 48 49 /* load in initial (unfactored row) */ 50 nz = ai[r[i]+1] - ai[r[i]]; 51 ajtmp = aj + ai[r[i]]; 52 v = aa + bs2*ai[r[i]]; 53 for (j=0; j<nz; j++) { 54 ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 55 } 56 57 /* elimination */ 58 bjtmp = bj + bi[i]; 59 nzL = bi[i+1] - bi[i]; 60 for(k=0;k < nzL;k++) { 61 row = bjtmp[k]; 62 pc = rtmp + bs2*row; 63 for (flg=0,j=0; j<bs2; j++) { if (pc[j] != (PetscScalar)0.0) { flg = 1; break; }} 64 if (flg) { 65 pv = b->a + bs2*bdiag[row]; 66 /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 67 ierr = Kernel_A_gets_A_times_B_2(pc,pv,mwork);CHKERRQ(ierr); 68 69 pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */ 70 pv = b->a + bs2*(bdiag[row+1]+1); 71 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */ 72 for (j=0; j<nz; j++) { 73 /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 74 /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 75 v = rtmp + 4*pj[j]; 76 ierr = Kernel_A_gets_A_minus_B_times_C_2(v,pc,pv);CHKERRQ(ierr); 77 pv += 4; 78 } 79 ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 80 } 81 } 82 83 /* finished row so stick it into b->a */ 84 /* L part */ 85 pv = b->a + bs2*bi[i] ; 86 pj = b->j + bi[i] ; 87 nz = bi[i+1] - bi[i]; 88 for (j=0; j<nz; j++) { 89 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 90 } 91 92 /* Mark diagonal and invert diagonal for simplier triangular solves */ 93 pv = b->a + bs2*bdiag[i]; 94 pj = b->j + bdiag[i]; 95 ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 96 /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 97 ierr = Kernel_A_gets_inverse_A_2(pv,shift);CHKERRQ(ierr); 98 99 /* U part */ 100 pv = b->a + bs2*(bdiag[i+1]+1); 101 pj = b->j + bdiag[i+1]+1; 102 nz = bdiag[i] - bdiag[i+1] - 1; 103 for (j=0; j<nz; j++){ 104 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 105 } 106 } 107 108 ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 109 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 110 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 111 C->ops->solve = MatSolve_SeqBAIJ_2; 112 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2; 113 114 C->assembled = PETSC_TRUE; 115 ierr = PetscLogFlops(1.333333333333*2*2*2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 116 PetscFunctionReturn(0); 117 } 118 119 #undef __FUNCT__ 120 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering" 121 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info) 122 { 123 Mat C=B; 124 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 125 PetscErrorCode ierr; 126 PetscInt i,j,k,nz,nzL,row,*pj; 127 const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,bs2=a->bs2; 128 const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag; 129 MatScalar *rtmp,*pc,*mwork,*pv; 130 MatScalar *aa=a->a,*v; 131 PetscInt flg; 132 PetscReal shift = info->shiftamount; 133 134 PetscFunctionBegin; 135 /* generate work space needed by the factorization */ 136 ierr = PetscMalloc2(bs2*n,MatScalar,&rtmp,bs2,MatScalar,&mwork);CHKERRQ(ierr); 137 ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 138 139 for (i=0; i<n; i++){ 140 /* zero rtmp */ 141 /* L part */ 142 nz = bi[i+1] - bi[i]; 143 bjtmp = bj + bi[i]; 144 for (j=0; j<nz; j++){ 145 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 146 } 147 148 /* U part */ 149 nz = bdiag[i] - bdiag[i+1]; 150 bjtmp = bj + bdiag[i+1]+1; 151 for (j=0; j<nz; j++){ 152 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 153 } 154 155 /* load in initial (unfactored row) */ 156 nz = ai[i+1] - ai[i]; 157 ajtmp = aj + ai[i]; 158 v = aa + bs2*ai[i]; 159 for (j=0; j<nz; j++) { 160 ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 161 } 162 163 /* elimination */ 164 bjtmp = bj + bi[i]; 165 nzL = bi[i+1] - bi[i]; 166 for(k=0;k < nzL;k++) { 167 row = bjtmp[k]; 168 pc = rtmp + bs2*row; 169 for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=(PetscScalar)0.0) { flg = 1; break; }} 170 if (flg) { 171 pv = b->a + bs2*bdiag[row]; 172 /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 173 ierr = Kernel_A_gets_A_times_B_2(pc,pv,mwork);CHKERRQ(ierr); 174 175 pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */ 176 pv = b->a + bs2*(bdiag[row+1]+1); 177 nz = bdiag[row]-bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */ 178 for (j=0; j<nz; j++) { 179 /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 180 /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 181 v = rtmp + 4*pj[j]; 182 ierr = Kernel_A_gets_A_minus_B_times_C_2(v,pc,pv);CHKERRQ(ierr); 183 pv += 4; 184 } 185 ierr = PetscLogFlops(16*nz+12);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 186 } 187 } 188 189 /* finished row so stick it into b->a */ 190 /* L part */ 191 pv = b->a + bs2*bi[i] ; 192 pj = b->j + bi[i] ; 193 nz = bi[i+1] - bi[i]; 194 for (j=0; j<nz; j++) { 195 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 196 } 197 198 /* Mark diagonal and invert diagonal for simplier triangular solves */ 199 pv = b->a + bs2*bdiag[i]; 200 pj = b->j + bdiag[i]; 201 ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 202 /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 203 ierr = Kernel_A_gets_inverse_A_2(pv,shift);CHKERRQ(ierr); 204 205 /* U part */ 206 /* 207 pv = b->a + bs2*bi[2*n-i]; 208 pj = b->j + bi[2*n-i]; 209 nz = bi[2*n-i+1] - bi[2*n-i] - 1; 210 */ 211 pv = b->a + bs2*(bdiag[i+1]+1); 212 pj = b->j + bdiag[i+1]+1; 213 nz = bdiag[i] - bdiag[i+1] - 1; 214 for (j=0; j<nz; j++){ 215 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 216 } 217 } 218 ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 219 220 C->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering; 221 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering; 222 C->assembled = PETSC_TRUE; 223 ierr = PetscLogFlops(1.333333333333*2*2*2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 224 PetscFunctionReturn(0); 225 } 226 227 #undef __FUNCT__ 228 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_inplace" 229 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_inplace(Mat B,Mat A,const MatFactorInfo *info) 230 { 231 Mat C = B; 232 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 233 IS isrow = b->row,isicol = b->icol; 234 PetscErrorCode ierr; 235 const PetscInt *r,*ic; 236 PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 237 PetscInt *ajtmpold,*ajtmp,nz,row; 238 PetscInt *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj; 239 MatScalar *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4; 240 MatScalar p1,p2,p3,p4; 241 MatScalar *ba = b->a,*aa = a->a; 242 PetscReal shift = info->shiftamount; 243 244 PetscFunctionBegin; 245 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 246 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 247 ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 248 249 for (i=0; i<n; i++) { 250 nz = bi[i+1] - bi[i]; 251 ajtmp = bj + bi[i]; 252 for (j=0; j<nz; j++) { 253 x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0; 254 } 255 /* load in initial (unfactored row) */ 256 idx = r[i]; 257 nz = ai[idx+1] - ai[idx]; 258 ajtmpold = aj + ai[idx]; 259 v = aa + 4*ai[idx]; 260 for (j=0; j<nz; j++) { 261 x = rtmp+4*ic[ajtmpold[j]]; 262 x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 263 v += 4; 264 } 265 row = *ajtmp++; 266 while (row < i) { 267 pc = rtmp + 4*row; 268 p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 269 if (p1 != (PetscScalar)0.0 || p2 != (PetscScalar)0.0 || p3 != (PetscScalar)0.0 || p4 != (PetscScalar)0.0) { 270 pv = ba + 4*diag_offset[row]; 271 pj = bj + diag_offset[row] + 1; 272 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 273 pc[0] = m1 = p1*x1 + p3*x2; 274 pc[1] = m2 = p2*x1 + p4*x2; 275 pc[2] = m3 = p1*x3 + p3*x4; 276 pc[3] = m4 = p2*x3 + p4*x4; 277 nz = bi[row+1] - diag_offset[row] - 1; 278 pv += 4; 279 for (j=0; j<nz; j++) { 280 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 281 x = rtmp + 4*pj[j]; 282 x[0] -= m1*x1 + m3*x2; 283 x[1] -= m2*x1 + m4*x2; 284 x[2] -= m1*x3 + m3*x4; 285 x[3] -= m2*x3 + m4*x4; 286 pv += 4; 287 } 288 ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr); 289 } 290 row = *ajtmp++; 291 } 292 /* finished row so stick it into b->a */ 293 pv = ba + 4*bi[i]; 294 pj = bj + bi[i]; 295 nz = bi[i+1] - bi[i]; 296 for (j=0; j<nz; j++) { 297 x = rtmp+4*pj[j]; 298 pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 299 pv += 4; 300 } 301 /* invert diagonal block */ 302 w = ba + 4*diag_offset[i]; 303 ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr); 304 } 305 306 ierr = PetscFree(rtmp);CHKERRQ(ierr); 307 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 308 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 309 C->ops->solve = MatSolve_SeqBAIJ_2_inplace; 310 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_inplace; 311 C->assembled = PETSC_TRUE; 312 ierr = PetscLogFlops(1.333333333333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 313 PetscFunctionReturn(0); 314 } 315 /* 316 Version for when blocks are 2 by 2 Using natural ordering 317 */ 318 #undef __FUNCT__ 319 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace" 320 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info) 321 { 322 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 323 PetscErrorCode ierr; 324 PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 325 PetscInt *ajtmpold,*ajtmp,nz,row; 326 PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 327 MatScalar *pv,*v,*rtmp,*pc,*w,*x; 328 MatScalar p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4; 329 MatScalar *ba = b->a,*aa = a->a; 330 PetscReal shift = info->shiftamount; 331 332 PetscFunctionBegin; 333 ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 334 for (i=0; i<n; i++) { 335 nz = bi[i+1] - bi[i]; 336 ajtmp = bj + bi[i]; 337 for (j=0; j<nz; j++) { 338 x = rtmp+4*ajtmp[j]; 339 x[0] = x[1] = x[2] = x[3] = 0.0; 340 } 341 /* load in initial (unfactored row) */ 342 nz = ai[i+1] - ai[i]; 343 ajtmpold = aj + ai[i]; 344 v = aa + 4*ai[i]; 345 for (j=0; j<nz; j++) { 346 x = rtmp+4*ajtmpold[j]; 347 x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 348 v += 4; 349 } 350 row = *ajtmp++; 351 while (row < i) { 352 pc = rtmp + 4*row; 353 p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 354 if (p1 != (PetscScalar)0.0 || p2 != (PetscScalar)0.0 || p3 != (PetscScalar)0.0 || p4 != (PetscScalar)0.0) { 355 pv = ba + 4*diag_offset[row]; 356 pj = bj + diag_offset[row] + 1; 357 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 358 pc[0] = m1 = p1*x1 + p3*x2; 359 pc[1] = m2 = p2*x1 + p4*x2; 360 pc[2] = m3 = p1*x3 + p3*x4; 361 pc[3] = m4 = p2*x3 + p4*x4; 362 nz = bi[row+1] - diag_offset[row] - 1; 363 pv += 4; 364 for (j=0; j<nz; j++) { 365 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 366 x = rtmp + 4*pj[j]; 367 x[0] -= m1*x1 + m3*x2; 368 x[1] -= m2*x1 + m4*x2; 369 x[2] -= m1*x3 + m3*x4; 370 x[3] -= m2*x3 + m4*x4; 371 pv += 4; 372 } 373 ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr); 374 } 375 row = *ajtmp++; 376 } 377 /* finished row so stick it into b->a */ 378 pv = ba + 4*bi[i]; 379 pj = bj + bi[i]; 380 nz = bi[i+1] - bi[i]; 381 for (j=0; j<nz; j++) { 382 x = rtmp+4*pj[j]; 383 pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 384 /* 385 printf(" col %d:",pj[j]); 386 PetscInt j1; 387 for (j1=0; j1<4; j1++) printf(" %g,",*(pv+j1)); 388 printf("\n"); 389 */ 390 pv += 4; 391 } 392 /* invert diagonal block */ 393 w = ba + 4*diag_offset[i]; 394 /* 395 printf(" \n%d -th: diag: ",i); 396 for (j=0; j<4; j++){ 397 printf(" %g,",w[j]); 398 } 399 printf("\n----------------------------\n"); 400 */ 401 ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr); 402 } 403 404 ierr = PetscFree(rtmp);CHKERRQ(ierr); 405 C->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering_inplace; 406 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering_inplace; 407 C->assembled = PETSC_TRUE; 408 ierr = PetscLogFlops(1.333333333333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 409 PetscFunctionReturn(0); 410 } 411 412 /* ----------------------------------------------------------- */ 413 /* 414 Version for when blocks are 1 by 1. 415 */ 416 #undef __FUNCT__ 417 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1" 418 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat B,Mat A,const MatFactorInfo *info) 419 { 420 Mat C=B; 421 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 422 IS isrow = b->row,isicol = b->icol; 423 PetscErrorCode ierr; 424 const PetscInt *r,*ic,*ics; 425 const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*bdiag=b->diag; 426 PetscInt i,j,k,nz,nzL,row,*pj; 427 const PetscInt *ajtmp,*bjtmp; 428 MatScalar *rtmp,*pc,multiplier,*pv; 429 const MatScalar *aa=a->a,*v; 430 PetscBool row_identity,col_identity; 431 FactorShiftCtx sctx; 432 const PetscInt *ddiag; 433 PetscReal rs; 434 MatScalar d; 435 436 PetscFunctionBegin; 437 /* MatPivotSetUp(): initialize shift context sctx */ 438 ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr); 439 440 if (info->shifttype == (PetscReal) MAT_SHIFT_POSITIVE_DEFINITE) { /* set sctx.shift_top=max{rs} */ 441 ddiag = a->diag; 442 sctx.shift_top = info->zeropivot; 443 for (i=0; i<n; i++) { 444 /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */ 445 d = (aa)[ddiag[i]]; 446 rs = -PetscAbsScalar(d) - PetscRealPart(d); 447 v = aa+ai[i]; 448 nz = ai[i+1] - ai[i]; 449 for (j=0; j<nz; j++) 450 rs += PetscAbsScalar(v[j]); 451 if (rs>sctx.shift_top) sctx.shift_top = rs; 452 } 453 sctx.shift_top *= 1.1; 454 sctx.nshift_max = 5; 455 sctx.shift_lo = 0.; 456 sctx.shift_hi = 1.; 457 } 458 459 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 460 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 461 ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 462 ics = ic; 463 464 do { 465 sctx.newshift = PETSC_FALSE; 466 for (i=0; i<n; i++){ 467 /* zero rtmp */ 468 /* L part */ 469 nz = bi[i+1] - bi[i]; 470 bjtmp = bj + bi[i]; 471 for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0; 472 473 /* U part */ 474 nz = bdiag[i]-bdiag[i+1]; 475 bjtmp = bj + bdiag[i+1]+1; 476 for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0; 477 478 /* load in initial (unfactored row) */ 479 nz = ai[r[i]+1] - ai[r[i]]; 480 ajtmp = aj + ai[r[i]]; 481 v = aa + ai[r[i]]; 482 for (j=0; j<nz; j++) { 483 rtmp[ics[ajtmp[j]]] = v[j]; 484 } 485 /* ZeropivotApply() */ 486 rtmp[i] += sctx.shift_amount; /* shift the diagonal of the matrix */ 487 488 /* elimination */ 489 bjtmp = bj + bi[i]; 490 row = *bjtmp++; 491 nzL = bi[i+1] - bi[i]; 492 for(k=0; k < nzL;k++) { 493 pc = rtmp + row; 494 if (*pc != (PetscScalar)0.0) { 495 pv = b->a + bdiag[row]; 496 multiplier = *pc * (*pv); 497 *pc = multiplier; 498 pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */ 499 pv = b->a + bdiag[row+1]+1; 500 nz = bdiag[row]-bdiag[row+1]-1; /* num of entries in U(row,:) excluding diag */ 501 for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 502 ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); 503 } 504 row = *bjtmp++; 505 } 506 507 /* finished row so stick it into b->a */ 508 rs = 0.0; 509 /* L part */ 510 pv = b->a + bi[i] ; 511 pj = b->j + bi[i] ; 512 nz = bi[i+1] - bi[i]; 513 for (j=0; j<nz; j++) { 514 pv[j] = rtmp[pj[j]]; rs += PetscAbsScalar(pv[j]); 515 } 516 517 /* U part */ 518 pv = b->a + bdiag[i+1]+1; 519 pj = b->j + bdiag[i+1]+1; 520 nz = bdiag[i] - bdiag[i+1]-1; 521 for (j=0; j<nz; j++) { 522 pv[j] = rtmp[pj[j]]; rs += PetscAbsScalar(pv[j]); 523 } 524 525 sctx.rs = rs; 526 sctx.pv = rtmp[i]; 527 ierr = MatPivotCheck(A,info,&sctx,i);CHKERRQ(ierr); 528 if(sctx.newshift) break; /* break for-loop */ 529 rtmp[i] = sctx.pv; /* sctx.pv might be updated in the case of MAT_SHIFT_INBLOCKS */ 530 531 /* Mark diagonal and invert diagonal for simplier triangular solves */ 532 pv = b->a + bdiag[i]; 533 *pv = (PetscScalar)1.0/rtmp[i]; 534 535 } /* endof for (i=0; i<n; i++){ */ 536 537 /* MatPivotRefine() */ 538 if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE && !sctx.newshift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max){ 539 /* 540 * if no shift in this attempt & shifting & started shifting & can refine, 541 * then try lower shift 542 */ 543 sctx.shift_hi = sctx.shift_fraction; 544 sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.; 545 sctx.shift_amount = sctx.shift_fraction * sctx.shift_top; 546 sctx.newshift = PETSC_TRUE; 547 sctx.nshift++; 548 } 549 } while (sctx.newshift); 550 551 ierr = PetscFree(rtmp);CHKERRQ(ierr); 552 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 553 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 554 555 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 556 ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr); 557 if (row_identity && col_identity) { 558 C->ops->solve = MatSolve_SeqBAIJ_1_NaturalOrdering; 559 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering; 560 } else { 561 C->ops->solve = MatSolve_SeqBAIJ_1; 562 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1; 563 } 564 C->assembled = PETSC_TRUE; 565 ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr); 566 567 /* MatShiftView(A,info,&sctx) */ 568 if (sctx.nshift){ 569 if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { 570 ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,sctx.shift_fraction,sctx.shift_top);CHKERRQ(ierr); 571 } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 572 ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 573 } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS){ 574 ierr = PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %G\n",sctx.nshift,info->shiftamount);CHKERRQ(ierr); 575 } 576 } 577 PetscFunctionReturn(0); 578 } 579 580 #undef __FUNCT__ 581 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1_inplace" 582 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info) 583 { 584 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 585 IS isrow = b->row,isicol = b->icol; 586 PetscErrorCode ierr; 587 const PetscInt *r,*ic; 588 PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 589 PetscInt *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j; 590 PetscInt *diag_offset = b->diag,diag,*pj; 591 MatScalar *pv,*v,*rtmp,multiplier,*pc; 592 MatScalar *ba = b->a,*aa = a->a; 593 PetscBool row_identity, col_identity; 594 595 PetscFunctionBegin; 596 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 597 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 598 ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 599 600 for (i=0; i<n; i++) { 601 nz = bi[i+1] - bi[i]; 602 ajtmp = bj + bi[i]; 603 for (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0; 604 605 /* load in initial (unfactored row) */ 606 nz = ai[r[i]+1] - ai[r[i]]; 607 ajtmpold = aj + ai[r[i]]; 608 v = aa + ai[r[i]]; 609 for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] = v[j]; 610 611 row = *ajtmp++; 612 while (row < i) { 613 pc = rtmp + row; 614 if (*pc != 0.0) { 615 pv = ba + diag_offset[row]; 616 pj = bj + diag_offset[row] + 1; 617 multiplier = *pc * *pv++; 618 *pc = multiplier; 619 nz = bi[row+1] - diag_offset[row] - 1; 620 for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 621 ierr = PetscLogFlops(1.0+2.0*nz);CHKERRQ(ierr); 622 } 623 row = *ajtmp++; 624 } 625 /* finished row so stick it into b->a */ 626 pv = ba + bi[i]; 627 pj = bj + bi[i]; 628 nz = bi[i+1] - bi[i]; 629 for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];} 630 diag = diag_offset[i] - bi[i]; 631 /* check pivot entry for current row */ 632 if (pv[diag] == 0.0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot: row in original ordering %D in permuted ordering %D",r[i],i); 633 pv[diag] = 1.0/pv[diag]; 634 } 635 636 ierr = PetscFree(rtmp);CHKERRQ(ierr); 637 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 638 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 639 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 640 ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr); 641 if (row_identity && col_identity) { 642 C->ops->solve = MatSolve_SeqBAIJ_1_NaturalOrdering_inplace; 643 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering_inplace; 644 } else { 645 C->ops->solve = MatSolve_SeqBAIJ_1_inplace; 646 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_inplace; 647 } 648 C->assembled = PETSC_TRUE; 649 ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr); 650 PetscFunctionReturn(0); 651 } 652 653 EXTERN_C_BEGIN 654 #undef __FUNCT__ 655 #define __FUNCT__ "MatGetFactor_seqbaij_petsc" 656 PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,MatFactorType ftype,Mat *B) 657 { 658 PetscInt n = A->rmap->n; 659 PetscErrorCode ierr; 660 661 PetscFunctionBegin; 662 ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr); 663 ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 664 if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) { 665 ierr = MatSetType(*B,MATSEQBAIJ);CHKERRQ(ierr); 666 (*B)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqBAIJ; 667 (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ; 668 } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 669 ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr); 670 ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 671 (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqBAIJ; 672 (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ; 673 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Factor type not supported"); 674 (*B)->factortype = ftype; 675 PetscFunctionReturn(0); 676 } 677 EXTERN_C_END 678 679 EXTERN_C_BEGIN 680 #undef __FUNCT__ 681 #define __FUNCT__ "MatGetFactorAvailable_seqbaij_petsc" 682 PetscErrorCode MatGetFactorAvailable_seqbaij_petsc(Mat A,MatFactorType ftype,PetscBool *flg) 683 { 684 PetscFunctionBegin; 685 *flg = PETSC_TRUE; 686 PetscFunctionReturn(0); 687 } 688 EXTERN_C_END 689 690 /* ----------------------------------------------------------- */ 691 #undef __FUNCT__ 692 #define __FUNCT__ "MatLUFactor_SeqBAIJ" 693 PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,const MatFactorInfo *info) 694 { 695 PetscErrorCode ierr; 696 Mat C; 697 698 PetscFunctionBegin; 699 ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr); 700 ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr); 701 ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr); 702 A->ops->solve = C->ops->solve; 703 A->ops->solvetranspose = C->ops->solvetranspose; 704 ierr = MatHeaderMerge(A,C);CHKERRQ(ierr); 705 ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr); 706 PetscFunctionReturn(0); 707 } 708 709 #include <../src/mat/impls/sbaij/seq/sbaij.h> 710 #undef __FUNCT__ 711 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N" 712 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C,Mat A,const MatFactorInfo *info) 713 { 714 PetscErrorCode ierr; 715 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 716 Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 717 IS ip=b->row; 718 const PetscInt *rip; 719 PetscInt i,j,mbs=a->mbs,bs=A->rmap->bs,*bi=b->i,*bj=b->j,*bcol; 720 PetscInt *ai=a->i,*aj=a->j; 721 PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 722 MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 723 PetscReal rs; 724 FactorShiftCtx sctx; 725 726 PetscFunctionBegin; 727 if (bs > 1) {/* convert A to a SBAIJ matrix and apply Cholesky factorization from it */ 728 if (!a->sbaijMat){ 729 ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 730 } 731 ierr = (a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info);CHKERRQ(ierr); 732 ierr = MatDestroy(&a->sbaijMat);CHKERRQ(ierr); 733 PetscFunctionReturn(0); 734 } 735 736 /* MatPivotSetUp(): initialize shift context sctx */ 737 ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr); 738 739 ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 740 ierr = PetscMalloc3(mbs,MatScalar,&rtmp,mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr); 741 742 sctx.shift_amount = 0.; 743 sctx.nshift = 0; 744 do { 745 sctx.newshift = PETSC_FALSE; 746 for (i=0; i<mbs; i++) { 747 rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 748 } 749 750 for (k = 0; k<mbs; k++){ 751 bval = ba + bi[k]; 752 /* initialize k-th row by the perm[k]-th row of A */ 753 jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 754 for (j = jmin; j < jmax; j++){ 755 col = rip[aj[j]]; 756 if (col >= k){ /* only take upper triangular entry */ 757 rtmp[col] = aa[j]; 758 *bval++ = 0.0; /* for in-place factorization */ 759 } 760 } 761 762 /* shift the diagonal of the matrix */ 763 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 764 765 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 766 dk = rtmp[k]; 767 i = jl[k]; /* first row to be added to k_th row */ 768 769 while (i < k){ 770 nexti = jl[i]; /* next row to be added to k_th row */ 771 772 /* compute multiplier, update diag(k) and U(i,k) */ 773 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 774 uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 775 dk += uikdi*ba[ili]; 776 ba[ili] = uikdi; /* -U(i,k) */ 777 778 /* add multiple of row i to k-th row */ 779 jmin = ili + 1; jmax = bi[i+1]; 780 if (jmin < jmax){ 781 for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 782 /* update il and jl for row i */ 783 il[i] = jmin; 784 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 785 } 786 i = nexti; 787 } 788 789 /* shift the diagonals when zero pivot is detected */ 790 /* compute rs=sum of abs(off-diagonal) */ 791 rs = 0.0; 792 jmin = bi[k]+1; 793 nz = bi[k+1] - jmin; 794 if (nz){ 795 bcol = bj + jmin; 796 while (nz--){ 797 rs += PetscAbsScalar(rtmp[*bcol]); 798 bcol++; 799 } 800 } 801 802 sctx.rs = rs; 803 sctx.pv = dk; 804 ierr = MatPivotCheck(A,info,&sctx,k);CHKERRQ(ierr); 805 if (sctx.newshift) break; 806 dk = sctx.pv; 807 808 /* copy data into U(k,:) */ 809 ba[bi[k]] = 1.0/dk; /* U(k,k) */ 810 jmin = bi[k]+1; jmax = bi[k+1]; 811 if (jmin < jmax) { 812 for (j=jmin; j<jmax; j++){ 813 col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 814 } 815 /* add the k-th row into il and jl */ 816 il[k] = jmin; 817 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 818 } 819 } 820 } while (sctx.newshift); 821 ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr); 822 823 ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 824 C->assembled = PETSC_TRUE; 825 C->preallocated = PETSC_TRUE; 826 ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 827 if (sctx.nshift){ 828 if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { 829 ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 830 } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 831 ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 832 } 833 } 834 PetscFunctionReturn(0); 835 } 836 837 #undef __FUNCT__ 838 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering" 839 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) 840 { 841 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 842 Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 843 PetscErrorCode ierr; 844 PetscInt i,j,am=a->mbs; 845 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 846 PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz; 847 MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 848 PetscReal rs; 849 FactorShiftCtx sctx; 850 851 PetscFunctionBegin; 852 /* MatPivotSetUp(): initialize shift context sctx */ 853 ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr); 854 855 ierr = PetscMalloc3(am,MatScalar,&rtmp,am,PetscInt,&il,am,PetscInt,&jl);CHKERRQ(ierr); 856 857 do { 858 sctx.newshift = PETSC_FALSE; 859 for (i=0; i<am; i++) { 860 rtmp[i] = 0.0; jl[i] = am; il[0] = 0; 861 } 862 863 for (k = 0; k<am; k++){ 864 /* initialize k-th row with elements nonzero in row perm(k) of A */ 865 nz = ai[k+1] - ai[k]; 866 acol = aj + ai[k]; 867 aval = aa + ai[k]; 868 bval = ba + bi[k]; 869 while (nz -- ){ 870 if (*acol < k) { /* skip lower triangular entries */ 871 acol++; aval++; 872 } else { 873 rtmp[*acol++] = *aval++; 874 *bval++ = 0.0; /* for in-place factorization */ 875 } 876 } 877 878 /* shift the diagonal of the matrix */ 879 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 880 881 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 882 dk = rtmp[k]; 883 i = jl[k]; /* first row to be added to k_th row */ 884 885 while (i < k){ 886 nexti = jl[i]; /* next row to be added to k_th row */ 887 /* compute multiplier, update D(k) and U(i,k) */ 888 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 889 uikdi = - ba[ili]*ba[bi[i]]; 890 dk += uikdi*ba[ili]; 891 ba[ili] = uikdi; /* -U(i,k) */ 892 893 /* add multiple of row i to k-th row ... */ 894 jmin = ili + 1; 895 nz = bi[i+1] - jmin; 896 if (nz > 0){ 897 bcol = bj + jmin; 898 bval = ba + jmin; 899 while (nz --) rtmp[*bcol++] += uikdi*(*bval++); 900 /* update il and jl for i-th row */ 901 il[i] = jmin; 902 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 903 } 904 i = nexti; 905 } 906 907 /* shift the diagonals when zero pivot is detected */ 908 /* compute rs=sum of abs(off-diagonal) */ 909 rs = 0.0; 910 jmin = bi[k]+1; 911 nz = bi[k+1] - jmin; 912 if (nz){ 913 bcol = bj + jmin; 914 while (nz--){ 915 rs += PetscAbsScalar(rtmp[*bcol]); 916 bcol++; 917 } 918 } 919 920 sctx.rs = rs; 921 sctx.pv = dk; 922 ierr = MatPivotCheck(A,info,&sctx,k);CHKERRQ(ierr); 923 if (sctx.newshift) break; /* sctx.shift_amount is updated */ 924 dk = sctx.pv; 925 926 /* copy data into U(k,:) */ 927 ba[bi[k]] = 1.0/dk; 928 jmin = bi[k]+1; 929 nz = bi[k+1] - jmin; 930 if (nz){ 931 bcol = bj + jmin; 932 bval = ba + jmin; 933 while (nz--){ 934 *bval++ = rtmp[*bcol]; 935 rtmp[*bcol++] = 0.0; 936 } 937 /* add k-th row into il and jl */ 938 il[k] = jmin; 939 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 940 } 941 } 942 } while (sctx.newshift); 943 ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr); 944 945 C->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 946 C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 947 C->assembled = PETSC_TRUE; 948 C->preallocated = PETSC_TRUE; 949 ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 950 if (sctx.nshift){ 951 if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 952 ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 953 } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { 954 ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 955 } 956 } 957 PetscFunctionReturn(0); 958 } 959 960 #include <petscbt.h> 961 #include <../src/mat/utils/freespace.h> 962 #undef __FUNCT__ 963 #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ" 964 PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 965 { 966 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 967 Mat_SeqSBAIJ *b; 968 Mat B; 969 PetscErrorCode ierr; 970 PetscBool perm_identity; 971 PetscInt reallocs=0,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui; 972 const PetscInt *rip; 973 PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 974 PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr; 975 PetscReal fill=info->fill,levels=info->levels; 976 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 977 PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 978 PetscBT lnkbt; 979 980 PetscFunctionBegin; 981 if (bs > 1){ 982 if (!a->sbaijMat){ 983 ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 984 } 985 (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */ 986 ierr = MatICCFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr); 987 PetscFunctionReturn(0); 988 } 989 990 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 991 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 992 993 /* special case that simply copies fill pattern */ 994 if (!levels && perm_identity) { 995 ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 996 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 997 for (i=0; i<am; i++) { 998 ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */ 999 } 1000 B = fact; 1001 ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr); 1002 1003 1004 b = (Mat_SeqSBAIJ*)B->data; 1005 uj = b->j; 1006 for (i=0; i<am; i++) { 1007 aj = a->j + a->diag[i]; 1008 for (j=0; j<ui[i]; j++){ 1009 *uj++ = *aj++; 1010 } 1011 b->ilen[i] = ui[i]; 1012 } 1013 ierr = PetscFree(ui);CHKERRQ(ierr); 1014 B->factortype = MAT_FACTOR_NONE; 1015 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1016 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1017 B->factortype = MAT_FACTOR_ICC; 1018 1019 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 1020 PetscFunctionReturn(0); 1021 } 1022 1023 /* initialization */ 1024 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 1025 ui[0] = 0; 1026 ierr = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr); 1027 1028 /* jl: linked list for storing indices of the pivot rows 1029 il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 1030 ierr = PetscMalloc4(am,PetscInt*,&uj_ptr,am,PetscInt*,&uj_lvl_ptr,am,PetscInt,&il,am,PetscInt,&jl);CHKERRQ(ierr); 1031 for (i=0; i<am; i++){ 1032 jl[i] = am; il[i] = 0; 1033 } 1034 1035 /* create and initialize a linked list for storing column indices of the active row k */ 1036 nlnk = am + 1; 1037 ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 1038 1039 /* initial FreeSpace size is fill*(ai[am]+am)/2 */ 1040 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+am)/2),&free_space);CHKERRQ(ierr); 1041 current_space = free_space; 1042 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+am)/2),&free_space_lvl);CHKERRQ(ierr); 1043 current_space_lvl = free_space_lvl; 1044 1045 for (k=0; k<am; k++){ /* for each active row k */ 1046 /* initialize lnk by the column indices of row rip[k] of A */ 1047 nzk = 0; 1048 ncols = ai[rip[k]+1] - ai[rip[k]]; 1049 ncols_upper = 0; 1050 cols = cols_lvl + am; 1051 for (j=0; j<ncols; j++){ 1052 i = rip[*(aj + ai[rip[k]] + j)]; 1053 if (i >= k){ /* only take upper triangular entry */ 1054 cols[ncols_upper] = i; 1055 cols_lvl[ncols_upper] = -1; /* initialize level for nonzero entries */ 1056 ncols_upper++; 1057 } 1058 } 1059 ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 1060 nzk += nlnk; 1061 1062 /* update lnk by computing fill-in for each pivot row to be merged in */ 1063 prow = jl[k]; /* 1st pivot row */ 1064 1065 while (prow < k){ 1066 nextprow = jl[prow]; 1067 1068 /* merge prow into k-th row */ 1069 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 1070 jmax = ui[prow+1]; 1071 ncols = jmax-jmin; 1072 i = jmin - ui[prow]; 1073 cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 1074 for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j); 1075 ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 1076 nzk += nlnk; 1077 1078 /* update il and jl for prow */ 1079 if (jmin < jmax){ 1080 il[prow] = jmin; 1081 j = *cols; jl[prow] = jl[j]; jl[j] = prow; 1082 } 1083 prow = nextprow; 1084 } 1085 1086 /* if free space is not available, make more free space */ 1087 if (current_space->local_remaining<nzk) { 1088 i = am - k + 1; /* num of unfactored rows */ 1089 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 1090 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 1091 ierr = PetscFreeSpaceGet(i,¤t_space_lvl);CHKERRQ(ierr); 1092 reallocs++; 1093 } 1094 1095 /* copy data into free_space and free_space_lvl, then initialize lnk */ 1096 ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 1097 1098 /* add the k-th row into il and jl */ 1099 if (nzk-1 > 0){ 1100 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 1101 jl[k] = jl[i]; jl[i] = k; 1102 il[k] = ui[k] + 1; 1103 } 1104 uj_ptr[k] = current_space->array; 1105 uj_lvl_ptr[k] = current_space_lvl->array; 1106 1107 current_space->array += nzk; 1108 current_space->local_used += nzk; 1109 current_space->local_remaining -= nzk; 1110 1111 current_space_lvl->array += nzk; 1112 current_space_lvl->local_used += nzk; 1113 current_space_lvl->local_remaining -= nzk; 1114 1115 ui[k+1] = ui[k] + nzk; 1116 } 1117 1118 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 1119 ierr = PetscFree4(uj_ptr,uj_lvl_ptr,il,jl);CHKERRQ(ierr); 1120 ierr = PetscFree(cols_lvl);CHKERRQ(ierr); 1121 1122 /* copy free_space into uj and free free_space; set uj in new datastructure; */ 1123 ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 1124 ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 1125 ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1126 ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 1127 1128 /* put together the new matrix in MATSEQSBAIJ format */ 1129 B = fact; 1130 ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 1131 1132 b = (Mat_SeqSBAIJ*)B->data; 1133 b->singlemalloc = PETSC_FALSE; 1134 b->free_a = PETSC_TRUE; 1135 b->free_ij = PETSC_TRUE; 1136 ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 1137 b->j = uj; 1138 b->i = ui; 1139 b->diag = 0; 1140 b->ilen = 0; 1141 b->imax = 0; 1142 b->row = perm; 1143 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 1144 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1145 b->icol = perm; 1146 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1147 ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1148 ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 1149 b->maxnz = b->nz = ui[am]; 1150 1151 B->info.factor_mallocs = reallocs; 1152 B->info.fill_ratio_given = fill; 1153 if (ai[am] != 0.) { 1154 /* nonzeros in lower triangular part of A (includign diagonals)= (ai[am]+am)/2 */ 1155 B->info.fill_ratio_needed = ((PetscReal)2*ui[am])/(ai[am]+am); 1156 } else { 1157 B->info.fill_ratio_needed = 0.0; 1158 } 1159 #if defined(PETSC_USE_INFO) 1160 if (ai[am] != 0) { 1161 PetscReal af = B->info.fill_ratio_needed; 1162 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 1163 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 1164 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 1165 } else { 1166 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 1167 } 1168 #endif 1169 if (perm_identity){ 1170 B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 1171 B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 1172 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 1173 } else { 1174 (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 1175 } 1176 PetscFunctionReturn(0); 1177 } 1178 1179 #undef __FUNCT__ 1180 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ" 1181 PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 1182 { 1183 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1184 Mat_SeqSBAIJ *b; 1185 Mat B; 1186 PetscErrorCode ierr; 1187 PetscBool perm_identity; 1188 PetscReal fill = info->fill; 1189 const PetscInt *rip; 1190 PetscInt i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow; 1191 PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 1192 PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 1193 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1194 PetscBT lnkbt; 1195 1196 PetscFunctionBegin; 1197 if (bs > 1) { /* convert to seqsbaij */ 1198 if (!a->sbaijMat){ 1199 ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 1200 } 1201 (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */ 1202 ierr = MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr); 1203 PetscFunctionReturn(0); 1204 } 1205 1206 /* check whether perm is the identity mapping */ 1207 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 1208 if (!perm_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported"); 1209 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 1210 1211 /* initialization */ 1212 ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 1213 ui[0] = 0; 1214 1215 /* jl: linked list for storing indices of the pivot rows 1216 il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */ 1217 ierr = PetscMalloc4(mbs,PetscInt*,&ui_ptr,mbs,PetscInt,&il,mbs,PetscInt,&jl,mbs,PetscInt,&cols);CHKERRQ(ierr); 1218 for (i=0; i<mbs; i++){ 1219 jl[i] = mbs; il[i] = 0; 1220 } 1221 1222 /* create and initialize a linked list for storing column indices of the active row k */ 1223 nlnk = mbs + 1; 1224 ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1225 1226 /* initial FreeSpace size is fill* (ai[mbs]+mbs)/2 */ 1227 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+mbs)/2),&free_space);CHKERRQ(ierr); 1228 current_space = free_space; 1229 1230 for (k=0; k<mbs; k++){ /* for each active row k */ 1231 /* initialize lnk by the column indices of row rip[k] of A */ 1232 nzk = 0; 1233 ncols = ai[rip[k]+1] - ai[rip[k]]; 1234 ncols_upper = 0; 1235 for (j=0; j<ncols; j++){ 1236 i = rip[*(aj + ai[rip[k]] + j)]; 1237 if (i >= k){ /* only take upper triangular entry */ 1238 cols[ncols_upper] = i; 1239 ncols_upper++; 1240 } 1241 } 1242 ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1243 nzk += nlnk; 1244 1245 /* update lnk by computing fill-in for each pivot row to be merged in */ 1246 prow = jl[k]; /* 1st pivot row */ 1247 1248 while (prow < k){ 1249 nextprow = jl[prow]; 1250 /* merge prow into k-th row */ 1251 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */ 1252 jmax = ui[prow+1]; 1253 ncols = jmax-jmin; 1254 uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */ 1255 ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1256 nzk += nlnk; 1257 1258 /* update il and jl for prow */ 1259 if (jmin < jmax){ 1260 il[prow] = jmin; 1261 j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 1262 } 1263 prow = nextprow; 1264 } 1265 1266 /* if free space is not available, make more free space */ 1267 if (current_space->local_remaining<nzk) { 1268 i = mbs - k + 1; /* num of unfactored rows */ 1269 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 1270 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 1271 reallocs++; 1272 } 1273 1274 /* copy data into free space, then initialize lnk */ 1275 ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 1276 1277 /* add the k-th row into il and jl */ 1278 if (nzk-1 > 0){ 1279 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 1280 jl[k] = jl[i]; jl[i] = k; 1281 il[k] = ui[k] + 1; 1282 } 1283 ui_ptr[k] = current_space->array; 1284 current_space->array += nzk; 1285 current_space->local_used += nzk; 1286 current_space->local_remaining -= nzk; 1287 1288 ui[k+1] = ui[k] + nzk; 1289 } 1290 1291 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 1292 ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr); 1293 1294 /* copy free_space into uj and free free_space; set uj in new datastructure; */ 1295 ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 1296 ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 1297 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1298 1299 /* put together the new matrix in MATSEQSBAIJ format */ 1300 B = fact; 1301 ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 1302 1303 b = (Mat_SeqSBAIJ*)B->data; 1304 b->singlemalloc = PETSC_FALSE; 1305 b->free_a = PETSC_TRUE; 1306 b->free_ij = PETSC_TRUE; 1307 ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 1308 b->j = uj; 1309 b->i = ui; 1310 b->diag = 0; 1311 b->ilen = 0; 1312 b->imax = 0; 1313 b->row = perm; 1314 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 1315 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1316 b->icol = perm; 1317 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1318 ierr = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1319 ierr = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 1320 b->maxnz = b->nz = ui[mbs]; 1321 1322 B->info.factor_mallocs = reallocs; 1323 B->info.fill_ratio_given = fill; 1324 if (ai[mbs] != 0.) { 1325 /* nonzeros in lower triangular part of A = (ai[mbs]+mbs)/2 */ 1326 B->info.fill_ratio_needed = ((PetscReal)2*ui[mbs])/(ai[mbs]+mbs); 1327 } else { 1328 B->info.fill_ratio_needed = 0.0; 1329 } 1330 #if defined(PETSC_USE_INFO) 1331 if (ai[mbs] != 0.) { 1332 PetscReal af = B->info.fill_ratio_needed; 1333 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 1334 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 1335 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 1336 } else { 1337 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 1338 } 1339 #endif 1340 if (perm_identity){ 1341 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 1342 } else { 1343 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 1344 } 1345 PetscFunctionReturn(0); 1346 } 1347 1348 #undef __FUNCT__ 1349 #define __FUNCT__ "MatSolve_SeqBAIJ_N_NaturalOrdering" 1350 PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering(Mat A,Vec bb,Vec xx) 1351 { 1352 Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 1353 PetscErrorCode ierr; 1354 const PetscInt *ai=a->i,*aj=a->j,*adiag=a->diag,*vi; 1355 PetscInt i,k,n=a->mbs; 1356 PetscInt nz,bs=A->rmap->bs,bs2=a->bs2; 1357 MatScalar *aa=a->a,*v; 1358 PetscScalar *x,*b,*s,*t,*ls; 1359 1360 PetscFunctionBegin; 1361 ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1362 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1363 t = a->solve_work; 1364 1365 /* forward solve the lower triangular */ 1366 ierr = PetscMemcpy(t,b,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy 1st block of b to t */ 1367 1368 for (i=1; i<n; i++) { 1369 v = aa + bs2*ai[i]; 1370 vi = aj + ai[i]; 1371 nz = ai[i+1] - ai[i]; 1372 s = t + bs*i; 1373 ierr = PetscMemcpy(s,b+bs*i,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy i_th block of b to t */ 1374 for(k=0;k<nz;k++){ 1375 Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[k]); 1376 v += bs2; 1377 } 1378 } 1379 1380 /* backward solve the upper triangular */ 1381 ls = a->solve_work + A->cmap->n; 1382 for (i=n-1; i>=0; i--){ 1383 v = aa + bs2*(adiag[i+1]+1); 1384 vi = aj + adiag[i+1]+1; 1385 nz = adiag[i] - adiag[i+1]-1; 1386 ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 1387 for(k=0;k<nz;k++){ 1388 Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[k]); 1389 v += bs2; 1390 } 1391 Kernel_w_gets_A_times_v(bs,ls,aa+bs2*adiag[i],t+i*bs); /* *inv(diagonal[i]) */ 1392 ierr = PetscMemcpy(x+i*bs,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 1393 } 1394 1395 ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 1396 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1397 ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr); 1398 PetscFunctionReturn(0); 1399 } 1400 1401 #undef __FUNCT__ 1402 #define __FUNCT__ "MatSolve_SeqBAIJ_N" 1403 PetscErrorCode MatSolve_SeqBAIJ_N(Mat A,Vec bb,Vec xx) 1404 { 1405 Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 1406 IS iscol=a->col,isrow=a->row; 1407 PetscErrorCode ierr; 1408 const PetscInt *r,*c,*rout,*cout,*ai=a->i,*aj=a->j,*adiag=a->diag,*vi; 1409 PetscInt i,m,n=a->mbs; 1410 PetscInt nz,bs=A->rmap->bs,bs2=a->bs2; 1411 MatScalar *aa=a->a,*v; 1412 PetscScalar *x,*b,*s,*t,*ls; 1413 1414 PetscFunctionBegin; 1415 ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1416 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1417 t = a->solve_work; 1418 1419 ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 1420 ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 1421 1422 /* forward solve the lower triangular */ 1423 ierr = PetscMemcpy(t,b+bs*r[0],bs*sizeof(PetscScalar));CHKERRQ(ierr); 1424 for (i=1; i<n; i++) { 1425 v = aa + bs2*ai[i]; 1426 vi = aj + ai[i]; 1427 nz = ai[i+1] - ai[i]; 1428 s = t + bs*i; 1429 ierr = PetscMemcpy(s,b+bs*r[i],bs*sizeof(PetscScalar));CHKERRQ(ierr); 1430 for(m=0;m<nz;m++){ 1431 Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[m]); 1432 v += bs2; 1433 } 1434 } 1435 1436 /* backward solve the upper triangular */ 1437 ls = a->solve_work + A->cmap->n; 1438 for (i=n-1; i>=0; i--){ 1439 v = aa + bs2*(adiag[i+1]+1); 1440 vi = aj + adiag[i+1]+1; 1441 nz = adiag[i] - adiag[i+1] - 1; 1442 ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 1443 for(m=0;m<nz;m++){ 1444 Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[m]); 1445 v += bs2; 1446 } 1447 Kernel_w_gets_A_times_v(bs,ls,v,t+i*bs); /* *inv(diagonal[i]) */ 1448 ierr = PetscMemcpy(x + bs*c[i],t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 1449 } 1450 ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 1451 ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 1452 ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 1453 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1454 ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr); 1455 PetscFunctionReturn(0); 1456 } 1457 1458 #undef __FUNCT__ 1459 #define __FUNCT__ "BlockAbs_privat" 1460 PetscErrorCode BlockAbs_private(PetscInt nbs,PetscInt bs2,PetscScalar *blockarray,PetscReal *absarray) 1461 { 1462 PetscErrorCode ierr; 1463 PetscInt i,j; 1464 PetscFunctionBegin; 1465 ierr = PetscMemzero(absarray,(nbs+1)*sizeof(PetscReal));CHKERRQ(ierr); 1466 for (i=0; i<nbs; i++){ 1467 for (j=0; j<bs2; j++){ 1468 if (absarray[i] < PetscAbsScalar(blockarray[i*nbs+j])) absarray[i] = PetscAbsScalar(blockarray[i*nbs+j]); 1469 } 1470 } 1471 PetscFunctionReturn(0); 1472 } 1473 1474 #undef __FUNCT__ 1475 #define __FUNCT__ "MatILUDTFactor_SeqBAIJ" 1476 /* 1477 This needs to be renamed and called by the regular MatILUFactor_SeqBAIJ when drop tolerance is used 1478 */ 1479 PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact) 1480 { 1481 Mat B = *fact; 1482 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b; 1483 IS isicol; 1484 PetscErrorCode ierr; 1485 const PetscInt *r,*ic; 1486 PetscInt i,mbs=a->mbs,bs=A->rmap->bs,bs2=a->bs2,*ai=a->i,*aj=a->j,*ajtmp,*adiag; 1487 PetscInt *bi,*bj,*bdiag; 1488 1489 PetscInt row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au; 1490 PetscInt nlnk,*lnk; 1491 PetscBT lnkbt; 1492 PetscBool row_identity,icol_identity; 1493 MatScalar *aatmp,*pv,*batmp,*ba,*rtmp,*pc,*multiplier,*vtmp; 1494 PetscInt j,nz,*pj,*bjtmp,k,ncut,*jtmp; 1495 1496 PetscReal dt=info->dt; /* shift=info->shiftamount; */ 1497 PetscInt nnz_max; 1498 PetscBool missing; 1499 PetscReal *vtmp_abs; 1500 MatScalar *v_work; 1501 PetscInt *v_pivots; 1502 1503 PetscFunctionBegin; 1504 /* ------- symbolic factorization, can be reused ---------*/ 1505 ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr); 1506 if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i); 1507 adiag=a->diag; 1508 1509 ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 1510 1511 /* bdiag is location of diagonal in factor */ 1512 ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 1513 1514 /* allocate row pointers bi */ 1515 ierr = PetscMalloc((2*mbs+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 1516 1517 /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */ 1518 dtcount = (PetscInt)info->dtcount; 1519 if (dtcount > mbs-1) dtcount = mbs-1; 1520 nnz_max = ai[mbs]+2*mbs*dtcount +2; 1521 /* printf("MatILUDTFactor_SeqBAIJ, bs %d, ai[mbs] %d, nnz_max %d, dtcount %d\n",bs,ai[mbs],nnz_max,dtcount); */ 1522 ierr = PetscMalloc(nnz_max*sizeof(PetscInt),&bj);CHKERRQ(ierr); 1523 nnz_max = nnz_max*bs2; 1524 ierr = PetscMalloc(nnz_max*sizeof(MatScalar),&ba);CHKERRQ(ierr); 1525 1526 /* put together the new matrix */ 1527 ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 1528 ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr); 1529 b = (Mat_SeqBAIJ*)(B)->data; 1530 b->free_a = PETSC_TRUE; 1531 b->free_ij = PETSC_TRUE; 1532 b->singlemalloc = PETSC_FALSE; 1533 b->a = ba; 1534 b->j = bj; 1535 b->i = bi; 1536 b->diag = bdiag; 1537 b->ilen = 0; 1538 b->imax = 0; 1539 b->row = isrow; 1540 b->col = iscol; 1541 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1542 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 1543 b->icol = isicol; 1544 ierr = PetscMalloc((bs*(mbs+1))*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1545 1546 ierr = PetscLogObjectMemory(B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 1547 b->maxnz = nnz_max/bs2; 1548 1549 (B)->factortype = MAT_FACTOR_ILUDT; 1550 (B)->info.factor_mallocs = 0; 1551 (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)(ai[mbs]*bs2)); 1552 CHKMEMQ; 1553 /* ------- end of symbolic factorization ---------*/ 1554 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 1555 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 1556 1557 /* linked list for storing column indices of the active row */ 1558 nlnk = mbs + 1; 1559 ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1560 1561 /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */ 1562 ierr = PetscMalloc2(mbs,PetscInt,&im,mbs,PetscInt,&jtmp);CHKERRQ(ierr); 1563 /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */ 1564 ierr = PetscMalloc2(mbs*bs2,MatScalar,&rtmp,mbs*bs2,MatScalar,&vtmp);CHKERRQ(ierr); 1565 ierr = PetscMalloc((mbs+1)*sizeof(PetscReal),&vtmp_abs);CHKERRQ(ierr); 1566 ierr = PetscMalloc3(bs,MatScalar,&v_work,bs2,MatScalar,&multiplier,bs,PetscInt,&v_pivots);CHKERRQ(ierr); 1567 1568 bi[0] = 0; 1569 bdiag[0] = (nnz_max/bs2)-1; /* location of diagonal in factor B */ 1570 bi[2*mbs+1] = bdiag[0]+1; /* endof bj and ba array */ 1571 for (i=0; i<mbs; i++) { 1572 /* copy initial fill into linked list */ 1573 nzi = 0; /* nonzeros for active row i */ 1574 nzi = ai[r[i]+1] - ai[r[i]]; 1575 if (!nzi) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 1576 nzi_al = adiag[r[i]] - ai[r[i]]; 1577 nzi_au = ai[r[i]+1] - adiag[r[i]] -1; 1578 /* printf("row %d, nzi_al/au %d %d\n",i,nzi_al,nzi_au); */ 1579 1580 /* load in initial unfactored row */ 1581 ajtmp = aj + ai[r[i]]; 1582 ierr = PetscLLAddPerm(nzi,ajtmp,ic,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1583 ierr = PetscMemzero(rtmp,mbs*bs2*sizeof(PetscScalar));CHKERRQ(ierr); 1584 aatmp = a->a + bs2*ai[r[i]]; 1585 for (j=0; j<nzi; j++) { 1586 ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],aatmp+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1587 } 1588 1589 /* add pivot rows into linked list */ 1590 row = lnk[mbs]; 1591 while (row < i) { 1592 nzi_bl = bi[row+1] - bi[row] + 1; 1593 bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */ 1594 ierr = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr); 1595 nzi += nlnk; 1596 row = lnk[row]; 1597 } 1598 1599 /* copy data from lnk into jtmp, then initialize lnk */ 1600 ierr = PetscLLClean(mbs,mbs,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr); 1601 1602 /* numerical factorization */ 1603 bjtmp = jtmp; 1604 row = *bjtmp++; /* 1st pivot row */ 1605 1606 while (row < i) { 1607 pc = rtmp + bs2*row; 1608 pv = ba + bs2*bdiag[row]; /* inv(diag) of the pivot row */ 1609 Kernel_A_gets_A_times_B(bs,pc,pv,multiplier); /* pc= multiplier = pc*inv(diag[row]) */ 1610 ierr = BlockAbs_private(1,bs2,pc,vtmp_abs);CHKERRQ(ierr); 1611 if (vtmp_abs[0] > dt){ /* apply tolerance dropping rule */ 1612 pj = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */ 1613 pv = ba + bs2*(bdiag[row+1] + 1); 1614 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */ 1615 for (j=0; j<nz; j++){ 1616 Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); 1617 } 1618 /* ierr = PetscLogFlops(bslog*(nz+1.0)-bs);CHKERRQ(ierr); */ 1619 } 1620 row = *bjtmp++; 1621 } 1622 1623 /* copy sparse rtmp into contiguous vtmp; separate L and U part */ 1624 nzi_bl = 0; j = 0; 1625 while (jtmp[j] < i){ /* L-part. Note: jtmp is sorted */ 1626 ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1627 nzi_bl++; j++; 1628 } 1629 nzi_bu = nzi - nzi_bl -1; 1630 /* printf("nzi %d, nzi_bl %d, nzi_bu %d\n",nzi,nzi_bl,nzi_bu); */ 1631 1632 while (j < nzi){ /* U-part */ 1633 ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1634 /* 1635 printf(" col %d: ",jtmp[j]); 1636 for (j1=0; j1<bs2; j1++) printf(" %g",*(vtmp+bs2*j+j1)); 1637 printf(" \n"); 1638 */ 1639 j++; 1640 } 1641 1642 ierr = BlockAbs_private(nzi,bs2,vtmp,vtmp_abs);CHKERRQ(ierr); 1643 /* 1644 printf(" row %d, nzi %d, vtmp_abs\n",i,nzi); 1645 for (j1=0; j1<nzi; j1++) printf(" (%d %g),",jtmp[j1],vtmp_abs[j1]); 1646 printf(" \n"); 1647 */ 1648 bjtmp = bj + bi[i]; 1649 batmp = ba + bs2*bi[i]; 1650 /* apply level dropping rule to L part */ 1651 ncut = nzi_al + dtcount; 1652 if (ncut < nzi_bl){ 1653 ierr = PetscSortSplitReal(ncut,nzi_bl,vtmp_abs,jtmp);CHKERRQ(ierr); 1654 ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr); 1655 } else { 1656 ncut = nzi_bl; 1657 } 1658 for (j=0; j<ncut; j++){ 1659 bjtmp[j] = jtmp[j]; 1660 ierr = PetscMemcpy(batmp+bs2*j,rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1661 /* 1662 printf(" col %d: ",bjtmp[j]); 1663 for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*j+j1)); 1664 printf("\n"); 1665 */ 1666 } 1667 bi[i+1] = bi[i] + ncut; 1668 nzi = ncut + 1; 1669 1670 /* apply level dropping rule to U part */ 1671 ncut = nzi_au + dtcount; 1672 if (ncut < nzi_bu){ 1673 ierr = PetscSortSplitReal(ncut,nzi_bu,vtmp_abs+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr); 1674 ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr); 1675 } else { 1676 ncut = nzi_bu; 1677 } 1678 nzi += ncut; 1679 1680 /* mark bdiagonal */ 1681 bdiag[i+1] = bdiag[i] - (ncut + 1); 1682 bi[2*mbs - i] = bi[2*mbs - i +1] - (ncut + 1); 1683 1684 bjtmp = bj + bdiag[i]; 1685 batmp = ba + bs2*bdiag[i]; 1686 ierr = PetscMemcpy(batmp,rtmp+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1687 *bjtmp = i; 1688 /* 1689 printf(" diag %d: ",*bjtmp); 1690 for (j=0; j<bs2; j++){ 1691 printf(" %g,",batmp[j]); 1692 } 1693 printf("\n"); 1694 */ 1695 bjtmp = bj + bdiag[i+1]+1; 1696 batmp = ba + (bdiag[i+1]+1)*bs2; 1697 1698 for (k=0; k<ncut; k++){ 1699 bjtmp[k] = jtmp[nzi_bl+1+k]; 1700 ierr = PetscMemcpy(batmp+bs2*k,rtmp+bs2*bjtmp[k],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1701 /* 1702 printf(" col %d:",bjtmp[k]); 1703 for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*k+j1)); 1704 printf("\n"); 1705 */ 1706 } 1707 1708 im[i] = nzi; /* used by PetscLLAddSortedLU() */ 1709 1710 /* invert diagonal block for simplier triangular solves - add shift??? */ 1711 batmp = ba + bs2*bdiag[i]; 1712 ierr = Kernel_A_gets_inverse_A(bs,batmp,v_pivots,v_work);CHKERRQ(ierr); 1713 } /* for (i=0; i<mbs; i++) */ 1714 ierr = PetscFree3(v_work,multiplier,v_pivots);CHKERRQ(ierr); 1715 1716 /* printf("end of L %d, beginning of U %d\n",bi[mbs],bdiag[mbs]); */ 1717 if (bi[mbs] >= bdiag[mbs]) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"end of L array %d cannot >= the beginning of U array %d",bi[mbs],bdiag[mbs]); 1718 1719 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 1720 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 1721 1722 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1723 1724 ierr = PetscFree2(im,jtmp);CHKERRQ(ierr); 1725 ierr = PetscFree2(rtmp,vtmp);CHKERRQ(ierr); 1726 1727 ierr = PetscLogFlops(bs2*B->cmap->n);CHKERRQ(ierr); 1728 b->maxnz = b->nz = bi[mbs] + bdiag[0] - bdiag[mbs]; 1729 1730 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 1731 ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr); 1732 if (row_identity && icol_identity) { 1733 B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering; 1734 } else { 1735 B->ops->solve = MatSolve_SeqBAIJ_N; 1736 } 1737 1738 B->ops->solveadd = 0; 1739 B->ops->solvetranspose = 0; 1740 B->ops->solvetransposeadd = 0; 1741 B->ops->matsolve = 0; 1742 B->assembled = PETSC_TRUE; 1743 B->preallocated = PETSC_TRUE; 1744 PetscFunctionReturn(0); 1745 } 1746