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