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 rs; 562 FactorShiftCtx sctx; 563 564 PetscFunctionBegin; 565 if (bs > 1) {/* convert A to a SBAIJ matrix and apply Cholesky factorization from it */ 566 if (!a->sbaijMat){ 567 ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 568 } 569 ierr = (a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info);CHKERRQ(ierr); 570 ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr); 571 a->sbaijMat = PETSC_NULL; 572 PetscFunctionReturn(0); 573 } 574 575 /* MatPivotSetUp(): initialize shift context sctx */ 576 ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr); 577 578 ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 579 ierr = PetscMalloc3(mbs,MatScalar,&rtmp,mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr); 580 581 sctx.shift_amount = 0.; 582 sctx.nshift = 0; 583 do { 584 sctx.newshift = PETSC_FALSE; 585 for (i=0; i<mbs; i++) { 586 rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 587 } 588 589 for (k = 0; k<mbs; k++){ 590 bval = ba + bi[k]; 591 /* initialize k-th row by the perm[k]-th row of A */ 592 jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 593 for (j = jmin; j < jmax; j++){ 594 col = rip[aj[j]]; 595 if (col >= k){ /* only take upper triangular entry */ 596 rtmp[col] = aa[j]; 597 *bval++ = 0.0; /* for in-place factorization */ 598 } 599 } 600 601 /* shift the diagonal of the matrix */ 602 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 603 604 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 605 dk = rtmp[k]; 606 i = jl[k]; /* first row to be added to k_th row */ 607 608 while (i < k){ 609 nexti = jl[i]; /* next row to be added to k_th row */ 610 611 /* compute multiplier, update diag(k) and U(i,k) */ 612 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 613 uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 614 dk += uikdi*ba[ili]; 615 ba[ili] = uikdi; /* -U(i,k) */ 616 617 /* add multiple of row i to k-th row */ 618 jmin = ili + 1; jmax = bi[i+1]; 619 if (jmin < jmax){ 620 for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 621 /* update il and jl for row i */ 622 il[i] = jmin; 623 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 624 } 625 i = nexti; 626 } 627 628 /* shift the diagonals when zero pivot is detected */ 629 /* compute rs=sum of abs(off-diagonal) */ 630 rs = 0.0; 631 jmin = bi[k]+1; 632 nz = bi[k+1] - jmin; 633 if (nz){ 634 bcol = bj + jmin; 635 while (nz--){ 636 rs += PetscAbsScalar(rtmp[*bcol]); 637 bcol++; 638 } 639 } 640 641 sctx.rs = rs; 642 sctx.pv = dk; 643 ierr = MatPivotCheck(info,&sctx,k);CHKERRQ(ierr); 644 if (sctx.newshift) break; 645 dk = sctx.pv; 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.newshift); 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 rs; 688 FactorShiftCtx sctx; 689 690 PetscFunctionBegin; 691 /* MatPivotSetUp(): initialize shift context sctx */ 692 ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr); 693 694 ierr = PetscMalloc3(am,MatScalar,&rtmp,am,PetscInt,&il,am,PetscInt,&jl);CHKERRQ(ierr); 695 696 do { 697 sctx.newshift = PETSC_FALSE; 698 for (i=0; i<am; i++) { 699 rtmp[i] = 0.0; jl[i] = am; il[0] = 0; 700 } 701 702 for (k = 0; k<am; k++){ 703 /* initialize k-th row with elements nonzero in row perm(k) of A */ 704 nz = ai[k+1] - ai[k]; 705 acol = aj + ai[k]; 706 aval = aa + ai[k]; 707 bval = ba + bi[k]; 708 while (nz -- ){ 709 if (*acol < k) { /* skip lower triangular entries */ 710 acol++; aval++; 711 } else { 712 rtmp[*acol++] = *aval++; 713 *bval++ = 0.0; /* for in-place factorization */ 714 } 715 } 716 717 /* shift the diagonal of the matrix */ 718 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 719 720 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 721 dk = rtmp[k]; 722 i = jl[k]; /* first row to be added to k_th row */ 723 724 while (i < k){ 725 nexti = jl[i]; /* next row to be added to k_th row */ 726 /* compute multiplier, update D(k) and U(i,k) */ 727 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 728 uikdi = - ba[ili]*ba[bi[i]]; 729 dk += uikdi*ba[ili]; 730 ba[ili] = uikdi; /* -U(i,k) */ 731 732 /* add multiple of row i to k-th row ... */ 733 jmin = ili + 1; 734 nz = bi[i+1] - jmin; 735 if (nz > 0){ 736 bcol = bj + jmin; 737 bval = ba + jmin; 738 while (nz --) rtmp[*bcol++] += uikdi*(*bval++); 739 /* update il and jl for i-th row */ 740 il[i] = jmin; 741 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 742 } 743 i = nexti; 744 } 745 746 /* shift the diagonals when zero pivot is detected */ 747 /* compute rs=sum of abs(off-diagonal) */ 748 rs = 0.0; 749 jmin = bi[k]+1; 750 nz = bi[k+1] - jmin; 751 if (nz){ 752 bcol = bj + jmin; 753 while (nz--){ 754 rs += PetscAbsScalar(rtmp[*bcol]); 755 bcol++; 756 } 757 } 758 759 sctx.rs = rs; 760 sctx.pv = dk; 761 ierr = MatPivotCheck(info,&sctx,k);CHKERRQ(ierr); 762 if (sctx.newshift) break; /* sctx.shift_amount is updated */ 763 dk = sctx.pv; 764 765 /* copy data into U(k,:) */ 766 ba[bi[k]] = 1.0/dk; 767 jmin = bi[k]+1; 768 nz = bi[k+1] - jmin; 769 if (nz){ 770 bcol = bj + jmin; 771 bval = ba + jmin; 772 while (nz--){ 773 *bval++ = rtmp[*bcol]; 774 rtmp[*bcol++] = 0.0; 775 } 776 /* add k-th row into il and jl */ 777 il[k] = jmin; 778 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 779 } 780 } 781 } while (sctx.newshift); 782 ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr); 783 784 C->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 785 C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 786 C->assembled = PETSC_TRUE; 787 C->preallocated = PETSC_TRUE; 788 ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 789 if (sctx.nshift){ 790 if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 791 ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 792 } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { 793 ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 794 } 795 } 796 PetscFunctionReturn(0); 797 } 798 799 #include "petscbt.h" 800 #include "../src/mat/utils/freespace.h" 801 #undef __FUNCT__ 802 #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ" 803 PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 804 { 805 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 806 Mat_SeqSBAIJ *b; 807 Mat B; 808 PetscErrorCode ierr; 809 PetscTruth perm_identity; 810 PetscInt reallocs=0,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui; 811 const PetscInt *rip; 812 PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 813 PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr; 814 PetscReal fill=info->fill,levels=info->levels; 815 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 816 PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 817 PetscBT lnkbt; 818 819 PetscFunctionBegin; 820 if (bs > 1){ 821 if (!a->sbaijMat){ 822 ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 823 } 824 (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */ 825 ierr = MatICCFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr); 826 PetscFunctionReturn(0); 827 } 828 829 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 830 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 831 832 /* special case that simply copies fill pattern */ 833 if (!levels && perm_identity) { 834 ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 835 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 836 for (i=0; i<am; i++) { 837 ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */ 838 } 839 B = fact; 840 ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr); 841 842 843 b = (Mat_SeqSBAIJ*)B->data; 844 uj = b->j; 845 for (i=0; i<am; i++) { 846 aj = a->j + a->diag[i]; 847 for (j=0; j<ui[i]; j++){ 848 *uj++ = *aj++; 849 } 850 b->ilen[i] = ui[i]; 851 } 852 ierr = PetscFree(ui);CHKERRQ(ierr); 853 B->factortype = MAT_FACTOR_NONE; 854 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 855 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 856 B->factortype = MAT_FACTOR_ICC; 857 858 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 859 PetscFunctionReturn(0); 860 } 861 862 /* initialization */ 863 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 864 ui[0] = 0; 865 ierr = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr); 866 867 /* jl: linked list for storing indices of the pivot rows 868 il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 869 ierr = PetscMalloc4(am,PetscInt*,&uj_ptr,am,PetscInt*,&uj_lvl_ptr,am,PetscInt,&il,am,PetscInt,&jl);CHKERRQ(ierr); 870 for (i=0; i<am; i++){ 871 jl[i] = am; il[i] = 0; 872 } 873 874 /* create and initialize a linked list for storing column indices of the active row k */ 875 nlnk = am + 1; 876 ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 877 878 /* initial FreeSpace size is fill*(ai[am]+1) */ 879 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 880 current_space = free_space; 881 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 882 current_space_lvl = free_space_lvl; 883 884 for (k=0; k<am; k++){ /* for each active row k */ 885 /* initialize lnk by the column indices of row rip[k] of A */ 886 nzk = 0; 887 ncols = ai[rip[k]+1] - ai[rip[k]]; 888 ncols_upper = 0; 889 cols = cols_lvl + am; 890 for (j=0; j<ncols; j++){ 891 i = rip[*(aj + ai[rip[k]] + j)]; 892 if (i >= k){ /* only take upper triangular entry */ 893 cols[ncols_upper] = i; 894 cols_lvl[ncols_upper] = -1; /* initialize level for nonzero entries */ 895 ncols_upper++; 896 } 897 } 898 ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 899 nzk += nlnk; 900 901 /* update lnk by computing fill-in for each pivot row to be merged in */ 902 prow = jl[k]; /* 1st pivot row */ 903 904 while (prow < k){ 905 nextprow = jl[prow]; 906 907 /* merge prow into k-th row */ 908 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 909 jmax = ui[prow+1]; 910 ncols = jmax-jmin; 911 i = jmin - ui[prow]; 912 cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 913 for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j); 914 ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 915 nzk += nlnk; 916 917 /* update il and jl for prow */ 918 if (jmin < jmax){ 919 il[prow] = jmin; 920 j = *cols; jl[prow] = jl[j]; jl[j] = prow; 921 } 922 prow = nextprow; 923 } 924 925 /* if free space is not available, make more free space */ 926 if (current_space->local_remaining<nzk) { 927 i = am - k + 1; /* num of unfactored rows */ 928 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 929 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 930 ierr = PetscFreeSpaceGet(i,¤t_space_lvl);CHKERRQ(ierr); 931 reallocs++; 932 } 933 934 /* copy data into free_space and free_space_lvl, then initialize lnk */ 935 ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 936 937 /* add the k-th row into il and jl */ 938 if (nzk-1 > 0){ 939 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 940 jl[k] = jl[i]; jl[i] = k; 941 il[k] = ui[k] + 1; 942 } 943 uj_ptr[k] = current_space->array; 944 uj_lvl_ptr[k] = current_space_lvl->array; 945 946 current_space->array += nzk; 947 current_space->local_used += nzk; 948 current_space->local_remaining -= nzk; 949 950 current_space_lvl->array += nzk; 951 current_space_lvl->local_used += nzk; 952 current_space_lvl->local_remaining -= nzk; 953 954 ui[k+1] = ui[k] + nzk; 955 } 956 957 #if defined(PETSC_USE_INFO) 958 if (ai[am] != 0) { 959 PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]); 960 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 961 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 962 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 963 } else { 964 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 965 } 966 #endif 967 968 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 969 ierr = PetscFree4(uj_ptr,uj_lvl_ptr,il,jl);CHKERRQ(ierr); 970 ierr = PetscFree(cols_lvl);CHKERRQ(ierr); 971 972 /* destroy list of free space and other temporary array(s) */ 973 ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 974 ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 975 ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 976 ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 977 978 /* put together the new matrix in MATSEQSBAIJ format */ 979 B = fact; 980 ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 981 982 b = (Mat_SeqSBAIJ*)B->data; 983 b->singlemalloc = PETSC_FALSE; 984 b->free_a = PETSC_TRUE; 985 b->free_ij = PETSC_TRUE; 986 ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 987 b->j = uj; 988 b->i = ui; 989 b->diag = 0; 990 b->ilen = 0; 991 b->imax = 0; 992 b->row = perm; 993 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 994 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 995 b->icol = perm; 996 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 997 ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 998 ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 999 b->maxnz = b->nz = ui[am]; 1000 1001 B->info.factor_mallocs = reallocs; 1002 B->info.fill_ratio_given = fill; 1003 if (ai[am] != 0.) { 1004 B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 1005 } else { 1006 B->info.fill_ratio_needed = 0.0; 1007 } 1008 if (perm_identity){ 1009 B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 1010 B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 1011 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 1012 } else { 1013 (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 1014 } 1015 PetscFunctionReturn(0); 1016 } 1017 1018 #undef __FUNCT__ 1019 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ" 1020 PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 1021 { 1022 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 1023 Mat_SeqSBAIJ *b; 1024 Mat B; 1025 PetscErrorCode ierr; 1026 PetscTruth perm_identity; 1027 PetscReal fill = info->fill; 1028 const PetscInt *rip; 1029 PetscInt i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow; 1030 PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 1031 PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 1032 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1033 PetscBT lnkbt; 1034 1035 PetscFunctionBegin; 1036 if (bs > 1) { /* convert to seqsbaij */ 1037 if (!a->sbaijMat){ 1038 ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 1039 } 1040 (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */ 1041 ierr = MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr); 1042 PetscFunctionReturn(0); 1043 } 1044 1045 /* check whether perm is the identity mapping */ 1046 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 1047 if (!perm_identity) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported"); 1048 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 1049 1050 /* initialization */ 1051 ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 1052 ui[0] = 0; 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:mbs-1) */ 1056 ierr = PetscMalloc4(mbs,PetscInt*,&ui_ptr,mbs,PetscInt,&il,mbs,PetscInt,&jl,mbs,PetscInt,&cols);CHKERRQ(ierr); 1057 for (i=0; i<mbs; i++){ 1058 jl[i] = mbs; il[i] = 0; 1059 } 1060 1061 /* create and initialize a linked list for storing column indices of the active row k */ 1062 nlnk = mbs + 1; 1063 ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1064 1065 /* initial FreeSpace size is fill*(ai[mbs]+1) */ 1066 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr); 1067 current_space = free_space; 1068 1069 for (k=0; k<mbs; k++){ /* for each active row k */ 1070 /* initialize lnk by the column indices of row rip[k] of A */ 1071 nzk = 0; 1072 ncols = ai[rip[k]+1] - ai[rip[k]]; 1073 ncols_upper = 0; 1074 for (j=0; j<ncols; j++){ 1075 i = rip[*(aj + ai[rip[k]] + j)]; 1076 if (i >= k){ /* only take upper triangular entry */ 1077 cols[ncols_upper] = i; 1078 ncols_upper++; 1079 } 1080 } 1081 ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1082 nzk += nlnk; 1083 1084 /* update lnk by computing fill-in for each pivot row to be merged in */ 1085 prow = jl[k]; /* 1st pivot row */ 1086 1087 while (prow < k){ 1088 nextprow = jl[prow]; 1089 /* merge prow into k-th row */ 1090 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */ 1091 jmax = ui[prow+1]; 1092 ncols = jmax-jmin; 1093 uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */ 1094 ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1095 nzk += nlnk; 1096 1097 /* update il and jl for prow */ 1098 if (jmin < jmax){ 1099 il[prow] = jmin; 1100 j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 1101 } 1102 prow = nextprow; 1103 } 1104 1105 /* if free space is not available, make more free space */ 1106 if (current_space->local_remaining<nzk) { 1107 i = mbs - k + 1; /* num of unfactored rows */ 1108 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 1109 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 1110 reallocs++; 1111 } 1112 1113 /* copy data into free space, then initialize lnk */ 1114 ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 1115 1116 /* add the k-th row into il and jl */ 1117 if (nzk-1 > 0){ 1118 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 1119 jl[k] = jl[i]; jl[i] = k; 1120 il[k] = ui[k] + 1; 1121 } 1122 ui_ptr[k] = current_space->array; 1123 current_space->array += nzk; 1124 current_space->local_used += nzk; 1125 current_space->local_remaining -= nzk; 1126 1127 ui[k+1] = ui[k] + nzk; 1128 } 1129 1130 #if defined(PETSC_USE_INFO) 1131 if (ai[mbs] != 0.) { 1132 PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 1133 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 1134 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 1135 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 1136 } else { 1137 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 1138 } 1139 #endif 1140 1141 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 1142 ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr); 1143 1144 /* destroy list of free space and other temporary array(s) */ 1145 ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 1146 ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 1147 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1148 1149 /* put together the new matrix in MATSEQSBAIJ format */ 1150 B = fact; 1151 ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 1152 1153 b = (Mat_SeqSBAIJ*)B->data; 1154 b->singlemalloc = PETSC_FALSE; 1155 b->free_a = PETSC_TRUE; 1156 b->free_ij = PETSC_TRUE; 1157 ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 1158 b->j = uj; 1159 b->i = ui; 1160 b->diag = 0; 1161 b->ilen = 0; 1162 b->imax = 0; 1163 b->row = perm; 1164 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 1165 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1166 b->icol = perm; 1167 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1168 ierr = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1169 ierr = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 1170 b->maxnz = b->nz = ui[mbs]; 1171 1172 B->info.factor_mallocs = reallocs; 1173 B->info.fill_ratio_given = fill; 1174 if (ai[mbs] != 0.) { 1175 B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 1176 } else { 1177 B->info.fill_ratio_needed = 0.0; 1178 } 1179 if (perm_identity){ 1180 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 1181 } else { 1182 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 1183 } 1184 PetscFunctionReturn(0); 1185 } 1186 1187 #undef __FUNCT__ 1188 #define __FUNCT__ "MatSolve_SeqBAIJ_N_NaturalOrdering" 1189 PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering(Mat A,Vec bb,Vec xx) 1190 { 1191 Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 1192 PetscErrorCode ierr; 1193 const PetscInt *ai=a->i,*aj=a->j,*adiag=a->diag,*vi; 1194 PetscInt i,k,n=a->mbs; 1195 PetscInt nz,bs=A->rmap->bs,bs2=a->bs2; 1196 MatScalar *aa=a->a,*v; 1197 PetscScalar *x,*b,*s,*t,*ls; 1198 1199 PetscFunctionBegin; 1200 ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1201 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1202 t = a->solve_work; 1203 1204 /* forward solve the lower triangular */ 1205 ierr = PetscMemcpy(t,b,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy 1st block of b to t */ 1206 1207 for (i=1; i<n; i++) { 1208 v = aa + bs2*ai[i]; 1209 vi = aj + ai[i]; 1210 nz = ai[i+1] - ai[i]; 1211 s = t + bs*i; 1212 ierr = PetscMemcpy(s,b+bs*i,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy i_th block of b to t */ 1213 for(k=0;k<nz;k++){ 1214 Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[k]); 1215 v += bs2; 1216 } 1217 } 1218 1219 /* backward solve the upper triangular */ 1220 ls = a->solve_work + A->cmap->n; 1221 for (i=n-1; i>=0; i--){ 1222 v = aa + bs2*(adiag[i+1]+1); 1223 vi = aj + adiag[i+1]+1; 1224 nz = adiag[i] - adiag[i+1]-1; 1225 ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 1226 for(k=0;k<nz;k++){ 1227 Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[k]); 1228 v += bs2; 1229 } 1230 Kernel_w_gets_A_times_v(bs,ls,aa+bs2*adiag[i],t+i*bs); /* *inv(diagonal[i]) */ 1231 ierr = PetscMemcpy(x+i*bs,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 1232 } 1233 1234 ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 1235 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1236 ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr); 1237 PetscFunctionReturn(0); 1238 } 1239 1240 #undef __FUNCT__ 1241 #define __FUNCT__ "MatSolve_SeqBAIJ_N" 1242 PetscErrorCode MatSolve_SeqBAIJ_N(Mat A,Vec bb,Vec xx) 1243 { 1244 Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 1245 IS iscol=a->col,isrow=a->row; 1246 PetscErrorCode ierr; 1247 const PetscInt *r,*c,*rout,*cout,*ai=a->i,*aj=a->j,*adiag=a->diag,*vi; 1248 PetscInt i,m,n=a->mbs; 1249 PetscInt nz,bs=A->rmap->bs,bs2=a->bs2; 1250 MatScalar *aa=a->a,*v; 1251 PetscScalar *x,*b,*s,*t,*ls; 1252 1253 PetscFunctionBegin; 1254 ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1255 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1256 t = a->solve_work; 1257 1258 ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 1259 ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 1260 1261 /* forward solve the lower triangular */ 1262 ierr = PetscMemcpy(t,b+bs*r[0],bs*sizeof(PetscScalar));CHKERRQ(ierr); 1263 for (i=1; i<n; i++) { 1264 v = aa + bs2*ai[i]; 1265 vi = aj + ai[i]; 1266 nz = ai[i+1] - ai[i]; 1267 s = t + bs*i; 1268 ierr = PetscMemcpy(s,b+bs*r[i],bs*sizeof(PetscScalar));CHKERRQ(ierr); 1269 for(m=0;m<nz;m++){ 1270 Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*vi[m]); 1271 v += bs2; 1272 } 1273 } 1274 1275 /* backward solve the upper triangular */ 1276 ls = a->solve_work + A->cmap->n; 1277 for (i=n-1; i>=0; i--){ 1278 v = aa + bs2*(adiag[i+1]+1); 1279 vi = aj + adiag[i+1]+1; 1280 nz = adiag[i] - adiag[i+1] - 1; 1281 ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 1282 for(m=0;m<nz;m++){ 1283 Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*vi[m]); 1284 v += bs2; 1285 } 1286 Kernel_w_gets_A_times_v(bs,ls,v,t+i*bs); /* *inv(diagonal[i]) */ 1287 ierr = PetscMemcpy(x + bs*c[i],t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 1288 } 1289 ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 1290 ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 1291 ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 1292 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1293 ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr); 1294 PetscFunctionReturn(0); 1295 } 1296 1297 #undef __FUNCT__ 1298 #define __FUNCT__ "BlockAbs_privat" 1299 PetscErrorCode BlockAbs_private(PetscInt nbs,PetscInt bs2,PetscScalar *blockarray,PetscReal *absarray) 1300 { 1301 PetscErrorCode ierr; 1302 PetscInt i,j; 1303 PetscFunctionBegin; 1304 ierr = PetscMemzero(absarray,(nbs+1)*sizeof(PetscReal));CHKERRQ(ierr); 1305 for (i=0; i<nbs; i++){ 1306 for (j=0; j<bs2; j++){ 1307 if (absarray[i] < PetscAbsScalar(blockarray[i*nbs+j])) absarray[i] = PetscAbsScalar(blockarray[i*nbs+j]); 1308 } 1309 } 1310 PetscFunctionReturn(0); 1311 } 1312 1313 #undef __FUNCT__ 1314 #define __FUNCT__ "MatILUDTFactor_SeqBAIJ" 1315 /* 1316 This needs to be renamed and called by the regular MatILUFactor_SeqBAIJ when drop tolerance is used 1317 */ 1318 PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact) 1319 { 1320 Mat B = *fact; 1321 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b; 1322 IS isicol; 1323 PetscErrorCode ierr; 1324 const PetscInt *r,*ic; 1325 PetscInt i,mbs=a->mbs,bs=A->rmap->bs,bs2=a->bs2,*ai=a->i,*aj=a->j,*ajtmp,*adiag; 1326 PetscInt *bi,*bj,*bdiag; 1327 1328 PetscInt row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au; 1329 PetscInt nlnk,*lnk; 1330 PetscBT lnkbt; 1331 PetscTruth row_identity,icol_identity,both_identity; 1332 MatScalar *aatmp,*pv,*batmp,*ba,*rtmp,*pc,*multiplier,*vtmp; 1333 const PetscInt *ics; 1334 PetscInt j,nz,*pj,*bjtmp,k,ncut,*jtmp; 1335 1336 PetscReal dt=info->dt; /* shift=info->shiftamount; */ 1337 PetscInt nnz_max; 1338 PetscTruth missing; 1339 PetscReal *vtmp_abs; 1340 MatScalar *v_work; 1341 PetscInt *v_pivots; 1342 1343 PetscFunctionBegin; 1344 /* ------- symbolic factorization, can be reused ---------*/ 1345 ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr); 1346 if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i); 1347 adiag=a->diag; 1348 1349 ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 1350 1351 /* bdiag is location of diagonal in factor */ 1352 ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 1353 1354 /* allocate row pointers bi */ 1355 ierr = PetscMalloc((2*mbs+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 1356 1357 /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */ 1358 dtcount = (PetscInt)info->dtcount; 1359 if (dtcount > mbs-1) dtcount = mbs-1; 1360 nnz_max = ai[mbs]+2*mbs*dtcount +2; 1361 /* printf("MatILUDTFactor_SeqBAIJ, bs %d, ai[mbs] %d, nnz_max %d, dtcount %d\n",bs,ai[mbs],nnz_max,dtcount); */ 1362 ierr = PetscMalloc(nnz_max*sizeof(PetscInt),&bj);CHKERRQ(ierr); 1363 nnz_max = nnz_max*bs2; 1364 ierr = PetscMalloc(nnz_max*sizeof(MatScalar),&ba);CHKERRQ(ierr); 1365 1366 /* put together the new matrix */ 1367 ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 1368 ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr); 1369 b = (Mat_SeqBAIJ*)(B)->data; 1370 b->free_a = PETSC_TRUE; 1371 b->free_ij = PETSC_TRUE; 1372 b->singlemalloc = PETSC_FALSE; 1373 b->a = ba; 1374 b->j = bj; 1375 b->i = bi; 1376 b->diag = bdiag; 1377 b->ilen = 0; 1378 b->imax = 0; 1379 b->row = isrow; 1380 b->col = iscol; 1381 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1382 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 1383 b->icol = isicol; 1384 ierr = PetscMalloc((bs*(mbs+1))*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1385 1386 ierr = PetscLogObjectMemory(B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 1387 b->maxnz = nnz_max/bs2; 1388 1389 (B)->factortype = MAT_FACTOR_ILUDT; 1390 (B)->info.factor_mallocs = 0; 1391 (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)(ai[mbs]*bs2)); 1392 CHKMEMQ; 1393 /* ------- end of symbolic factorization ---------*/ 1394 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 1395 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 1396 ics = ic; 1397 1398 /* linked list for storing column indices of the active row */ 1399 nlnk = mbs + 1; 1400 ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1401 1402 /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */ 1403 ierr = PetscMalloc2(mbs,PetscInt,&im,mbs,PetscInt,&jtmp);CHKERRQ(ierr); 1404 /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */ 1405 ierr = PetscMalloc2(mbs*bs2,MatScalar,&rtmp,mbs*bs2,MatScalar,&vtmp);CHKERRQ(ierr); 1406 ierr = PetscMalloc((mbs+1)*sizeof(PetscReal),&vtmp_abs);CHKERRQ(ierr); 1407 ierr = PetscMalloc3(bs,MatScalar,&v_work,bs2,MatScalar,&multiplier,bs,PetscInt,&v_pivots);CHKERRQ(ierr); 1408 1409 bi[0] = 0; 1410 bdiag[0] = (nnz_max/bs2)-1; /* location of diagonal in factor B */ 1411 bi[2*mbs+1] = bdiag[0]+1; /* endof bj and ba array */ 1412 for (i=0; i<mbs; i++) { 1413 /* copy initial fill into linked list */ 1414 nzi = 0; /* nonzeros for active row i */ 1415 nzi = ai[r[i]+1] - ai[r[i]]; 1416 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); 1417 nzi_al = adiag[r[i]] - ai[r[i]]; 1418 nzi_au = ai[r[i]+1] - adiag[r[i]] -1; 1419 /* printf("row %d, nzi_al/au %d %d\n",i,nzi_al,nzi_au); */ 1420 1421 /* load in initial unfactored row */ 1422 ajtmp = aj + ai[r[i]]; 1423 ierr = PetscLLAddPerm(nzi,ajtmp,ic,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1424 ierr = PetscMemzero(rtmp,mbs*bs2*sizeof(PetscScalar));CHKERRQ(ierr); 1425 aatmp = a->a + bs2*ai[r[i]]; 1426 for (j=0; j<nzi; j++) { 1427 ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],aatmp+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1428 } 1429 1430 /* add pivot rows into linked list */ 1431 row = lnk[mbs]; 1432 while (row < i) { 1433 nzi_bl = bi[row+1] - bi[row] + 1; 1434 bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */ 1435 ierr = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr); 1436 nzi += nlnk; 1437 row = lnk[row]; 1438 } 1439 1440 /* copy data from lnk into jtmp, then initialize lnk */ 1441 ierr = PetscLLClean(mbs,mbs,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr); 1442 1443 /* numerical factorization */ 1444 bjtmp = jtmp; 1445 row = *bjtmp++; /* 1st pivot row */ 1446 1447 while (row < i) { 1448 pc = rtmp + bs2*row; 1449 pv = ba + bs2*bdiag[row]; /* inv(diag) of the pivot row */ 1450 Kernel_A_gets_A_times_B(bs,pc,pv,multiplier); /* pc= multiplier = pc*inv(diag[row]) */ 1451 ierr = BlockAbs_private(1,bs2,pc,vtmp_abs);CHKERRQ(ierr); 1452 if (vtmp_abs[0] > dt){ /* apply tolerance dropping rule */ 1453 pj = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */ 1454 pv = ba + bs2*(bdiag[row+1] + 1); 1455 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */ 1456 for (j=0; j<nz; j++){ 1457 Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); 1458 } 1459 /* ierr = PetscLogFlops(bslog*(nz+1.0)-bs);CHKERRQ(ierr); */ 1460 } 1461 row = *bjtmp++; 1462 } 1463 1464 /* copy sparse rtmp into contiguous vtmp; separate L and U part */ 1465 nzi_bl = 0; j = 0; 1466 while (jtmp[j] < i){ /* L-part. Note: jtmp is sorted */ 1467 ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1468 nzi_bl++; j++; 1469 } 1470 nzi_bu = nzi - nzi_bl -1; 1471 /* printf("nzi %d, nzi_bl %d, nzi_bu %d\n",nzi,nzi_bl,nzi_bu); */ 1472 1473 while (j < nzi){ /* U-part */ 1474 ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1475 /* 1476 printf(" col %d: ",jtmp[j]); 1477 for (j1=0; j1<bs2; j1++) printf(" %g",*(vtmp+bs2*j+j1)); 1478 printf(" \n"); 1479 */ 1480 j++; 1481 } 1482 1483 ierr = BlockAbs_private(nzi,bs2,vtmp,vtmp_abs);CHKERRQ(ierr); 1484 /* 1485 printf(" row %d, nzi %d, vtmp_abs\n",i,nzi); 1486 for (j1=0; j1<nzi; j1++) printf(" (%d %g),",jtmp[j1],vtmp_abs[j1]); 1487 printf(" \n"); 1488 */ 1489 bjtmp = bj + bi[i]; 1490 batmp = ba + bs2*bi[i]; 1491 /* apply level dropping rule to L part */ 1492 ncut = nzi_al + dtcount; 1493 if (ncut < nzi_bl){ 1494 ierr = PetscSortSplitReal(ncut,nzi_bl,vtmp_abs,jtmp);CHKERRQ(ierr); 1495 ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr); 1496 } else { 1497 ncut = nzi_bl; 1498 } 1499 for (j=0; j<ncut; j++){ 1500 bjtmp[j] = jtmp[j]; 1501 ierr = PetscMemcpy(batmp+bs2*j,rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1502 /* 1503 printf(" col %d: ",bjtmp[j]); 1504 for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*j+j1)); 1505 printf("\n"); 1506 */ 1507 } 1508 bi[i+1] = bi[i] + ncut; 1509 nzi = ncut + 1; 1510 1511 /* apply level dropping rule to U part */ 1512 ncut = nzi_au + dtcount; 1513 if (ncut < nzi_bu){ 1514 ierr = PetscSortSplitReal(ncut,nzi_bu,vtmp_abs+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr); 1515 ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr); 1516 } else { 1517 ncut = nzi_bu; 1518 } 1519 nzi += ncut; 1520 1521 /* mark bdiagonal */ 1522 bdiag[i+1] = bdiag[i] - (ncut + 1); 1523 bi[2*mbs - i] = bi[2*mbs - i +1] - (ncut + 1); 1524 1525 bjtmp = bj + bdiag[i]; 1526 batmp = ba + bs2*bdiag[i]; 1527 ierr = PetscMemcpy(batmp,rtmp+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1528 *bjtmp = i; 1529 /* 1530 printf(" diag %d: ",*bjtmp); 1531 for (j=0; j<bs2; j++){ 1532 printf(" %g,",batmp[j]); 1533 } 1534 printf("\n"); 1535 */ 1536 bjtmp = bj + bdiag[i+1]+1; 1537 batmp = ba + (bdiag[i+1]+1)*bs2; 1538 1539 for (k=0; k<ncut; k++){ 1540 bjtmp[k] = jtmp[nzi_bl+1+k]; 1541 ierr = PetscMemcpy(batmp+bs2*k,rtmp+bs2*bjtmp[k],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1542 /* 1543 printf(" col %d:",bjtmp[k]); 1544 for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*k+j1)); 1545 printf("\n"); 1546 */ 1547 } 1548 1549 im[i] = nzi; /* used by PetscLLAddSortedLU() */ 1550 1551 /* invert diagonal block for simplier triangular solves - add shift??? */ 1552 batmp = ba + bs2*bdiag[i]; 1553 ierr = Kernel_A_gets_inverse_A(bs,batmp,v_pivots,v_work);CHKERRQ(ierr); 1554 } /* for (i=0; i<mbs; i++) */ 1555 ierr = PetscFree3(v_work,multiplier,v_pivots);CHKERRQ(ierr); 1556 1557 /* printf("end of L %d, beginning of U %d\n",bi[mbs],bdiag[mbs]); */ 1558 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]); 1559 1560 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 1561 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 1562 1563 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1564 1565 ierr = PetscFree2(im,jtmp);CHKERRQ(ierr); 1566 ierr = PetscFree2(rtmp,vtmp);CHKERRQ(ierr); 1567 1568 ierr = PetscLogFlops(bs2*B->cmap->n);CHKERRQ(ierr); 1569 b->maxnz = b->nz = bi[mbs] + bdiag[0] - bdiag[mbs]; 1570 1571 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 1572 ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr); 1573 both_identity = (PetscTruth) (row_identity && icol_identity); 1574 if (row_identity && icol_identity) { 1575 B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering; 1576 } else { 1577 B->ops->solve = MatSolve_SeqBAIJ_N; 1578 } 1579 1580 B->ops->solveadd = 0; 1581 B->ops->solvetranspose = 0; 1582 B->ops->solvetransposeadd = 0; 1583 B->ops->matsolve = 0; 1584 B->assembled = PETSC_TRUE; 1585 B->preallocated = PETSC_TRUE; 1586 PetscFunctionReturn(0); 1587 } 1588