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