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