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