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