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