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