1 #define PETSCMAT_DLL 2 3 #include "../src/mat/impls/aij/seq/aij.h" 4 #include "../src/inline/dot.h" 5 #define PETSC_USE_WHILE_KERNELS 6 #include "../src/inline/spops.h" 7 #include "petscbt.h" 8 #include "../src/mat/utils/freespace.h" 9 10 EXTERN_C_BEGIN 11 #undef __FUNCT__ 12 #define __FUNCT__ "MatOrdering_Flow_SeqAIJ" 13 /* 14 Computes an ordering to get most of the large numerical values in the lower triangular part of the matrix 15 */ 16 PetscErrorCode MatOrdering_Flow_SeqAIJ(Mat mat,const MatOrderingType type,IS *irow,IS *icol) 17 { 18 Mat_SeqAIJ *a = (Mat_SeqAIJ*)mat->data; 19 PetscErrorCode ierr; 20 PetscInt i,j,jj,k, kk,n = mat->rmap->n, current = 0, newcurrent = 0,*order; 21 const PetscInt *ai = a->i, *aj = a->j; 22 const PetscScalar *aa = a->a; 23 PetscTruth *done; 24 PetscReal best,past = 0,future; 25 26 PetscFunctionBegin; 27 /* pick initial row */ 28 best = -1; 29 for (i=0; i<n; i++) { 30 future = 0; 31 for (j=ai[i]; j<ai[i+1]; j++) { 32 if (aj[j] != i) future += PetscAbsScalar(aa[j]); else past = PetscAbsScalar(aa[j]); 33 } 34 if (!future) future = 1.e-10; /* if there is zero in the upper diagonal part want to rank this row high */ 35 if (past/future > best) { 36 best = past/future; 37 current = i; 38 } 39 } 40 41 ierr = PetscMalloc(n*sizeof(PetscTruth),&done);CHKERRQ(ierr); 42 ierr = PetscMalloc(n*sizeof(PetscInt),&order);CHKERRQ(ierr); 43 ierr = PetscMemzero(done,n*sizeof(PetscTruth));CHKERRQ(ierr); 44 order[0] = current; 45 for (i=0; i<n-1; i++) { 46 done[current] = PETSC_TRUE; 47 best = -1; 48 /* loop over all neighbors of current pivot */ 49 for (j=ai[current]; j<ai[current+1]; j++) { 50 jj = aj[j]; 51 if (done[jj]) continue; 52 /* loop over columns of potential next row computing weights for below and above diagonal */ 53 past = future = 0.0; 54 for (k=ai[jj]; k<ai[jj+1]; k++) { 55 kk = aj[k]; 56 if (done[kk]) past += PetscAbsScalar(aa[k]); 57 else if (kk != jj) future += PetscAbsScalar(aa[k]); 58 } 59 if (!future) future = 1.e-10; /* if there is zero in the upper diagonal part want to rank this row high */ 60 if (past/future > best) { 61 best = past/future; 62 newcurrent = jj; 63 } 64 } 65 if (best == -1) { /* no neighbors to select from so select best of all that remain */ 66 best = -1; 67 for (k=0; k<n; k++) { 68 if (done[k]) continue; 69 future = 0; 70 past = 0; 71 for (j=ai[k]; j<ai[k+1]; j++) { 72 kk = aj[j]; 73 if (done[kk]) past += PetscAbsScalar(aa[j]); 74 else if (kk != k) future += PetscAbsScalar(aa[j]); 75 } 76 if (!future) future = 1.e-10; /* if there is zero in the upper diagonal part want to rank this row high */ 77 if (past/future > best) { 78 best = past/future; 79 newcurrent = k; 80 } 81 } 82 } 83 if (current == newcurrent) SETERRQ(PETSC_ERR_PLIB,"newcurrent cannot be current"); 84 current = newcurrent; 85 order[i+1] = current; 86 } 87 ierr = ISCreateGeneral(PETSC_COMM_SELF,n,order,irow);CHKERRQ(ierr); 88 *icol = *irow; 89 ierr = PetscObjectReference((PetscObject)*irow);CHKERRQ(ierr); 90 ierr = PetscFree(done);CHKERRQ(ierr); 91 ierr = PetscFree(order);CHKERRQ(ierr); 92 PetscFunctionReturn(0); 93 } 94 EXTERN_C_END 95 96 EXTERN_C_BEGIN 97 #undef __FUNCT__ 98 #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc" 99 PetscErrorCode MatGetFactorAvailable_seqaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg) 100 { 101 PetscFunctionBegin; 102 *flg = PETSC_TRUE; 103 PetscFunctionReturn(0); 104 } 105 EXTERN_C_END 106 107 EXTERN_C_BEGIN 108 #undef __FUNCT__ 109 #define __FUNCT__ "MatGetFactor_seqaij_petsc" 110 PetscErrorCode MatGetFactor_seqaij_petsc(Mat A,MatFactorType ftype,Mat *B) 111 { 112 PetscInt n = A->rmap->n; 113 PetscErrorCode ierr; 114 115 PetscFunctionBegin; 116 ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr); 117 ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 118 if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT){ 119 ierr = MatSetType(*B,MATSEQAIJ);CHKERRQ(ierr); 120 (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqAIJ; 121 (*B)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqAIJ; 122 (*B)->ops->iludtfactor = MatILUDTFactor_SeqAIJ; 123 } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 124 ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr); 125 ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 126 (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqAIJ; 127 (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqAIJ; 128 } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported"); 129 (*B)->factor = ftype; 130 PetscFunctionReturn(0); 131 } 132 EXTERN_C_END 133 134 #undef __FUNCT__ 135 #define __FUNCT__ "MatLUFactorSymbolic_SeqAIJ" 136 PetscErrorCode MatLUFactorSymbolic_SeqAIJ(Mat B,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 137 { 138 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 139 IS isicol; 140 PetscErrorCode ierr; 141 const PetscInt *r,*ic; 142 PetscInt i,n=A->rmap->n,*ai=a->i,*aj=a->j; 143 PetscInt *bi,*bj,*ajtmp; 144 PetscInt *bdiag,row,nnz,nzi,reallocs=0,nzbd,*im; 145 PetscReal f; 146 PetscInt nlnk,*lnk,k,**bi_ptr; 147 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 148 PetscBT lnkbt; 149 150 PetscFunctionBegin; 151 if (A->rmap->N != A->cmap->N) SETERRQ(PETSC_ERR_ARG_WRONG,"matrix must be square"); 152 ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 153 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 154 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 155 156 /* get new row pointers */ 157 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 158 bi[0] = 0; 159 160 /* bdiag is location of diagonal in factor */ 161 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 162 bdiag[0] = 0; 163 164 /* linked list for storing column indices of the active row */ 165 nlnk = n + 1; 166 ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 167 168 ierr = PetscMalloc2(n+1,PetscInt**,&bi_ptr,n+1,PetscInt,&im);CHKERRQ(ierr); 169 170 /* initial FreeSpace size is f*(ai[n]+1) */ 171 f = info->fill; 172 ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 173 current_space = free_space; 174 175 for (i=0; i<n; i++) { 176 /* copy previous fill into linked list */ 177 nzi = 0; 178 nnz = ai[r[i]+1] - ai[r[i]]; 179 if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 180 ajtmp = aj + ai[r[i]]; 181 ierr = PetscLLAddPerm(nnz,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 182 nzi += nlnk; 183 184 /* add pivot rows into linked list */ 185 row = lnk[n]; 186 while (row < i) { 187 nzbd = bdiag[row] - bi[row] + 1; /* num of entries in the row with column index <= row */ 188 ajtmp = bi_ptr[row] + nzbd; /* points to the entry next to the diagonal */ 189 ierr = PetscLLAddSortedLU(ajtmp,row,nlnk,lnk,lnkbt,i,nzbd,im);CHKERRQ(ierr); 190 nzi += nlnk; 191 row = lnk[row]; 192 } 193 bi[i+1] = bi[i] + nzi; 194 im[i] = nzi; 195 196 /* mark bdiag */ 197 nzbd = 0; 198 nnz = nzi; 199 k = lnk[n]; 200 while (nnz-- && k < i){ 201 nzbd++; 202 k = lnk[k]; 203 } 204 bdiag[i] = bi[i] + nzbd; 205 206 /* if free space is not available, make more free space */ 207 if (current_space->local_remaining<nzi) { 208 nnz = (n - i)*nzi; /* estimated and max additional space needed */ 209 ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 210 reallocs++; 211 } 212 213 /* copy data into free space, then initialize lnk */ 214 ierr = PetscLLClean(n,n,nzi,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 215 bi_ptr[i] = current_space->array; 216 current_space->array += nzi; 217 current_space->local_used += nzi; 218 current_space->local_remaining -= nzi; 219 } 220 #if defined(PETSC_USE_INFO) 221 if (ai[n] != 0) { 222 PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 223 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr); 224 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 225 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr); 226 ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 227 } else { 228 ierr = PetscInfo(A,"Empty matrix\n");CHKERRQ(ierr); 229 } 230 #endif 231 232 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 233 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 234 235 /* destroy list of free space and other temporary array(s) */ 236 ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 237 ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 238 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 239 ierr = PetscFree2(bi_ptr,im);CHKERRQ(ierr); 240 241 /* put together the new matrix */ 242 ierr = MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 243 ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr); 244 b = (Mat_SeqAIJ*)(B)->data; 245 b->free_a = PETSC_TRUE; 246 b->free_ij = PETSC_TRUE; 247 b->singlemalloc = PETSC_FALSE; 248 ierr = PetscMalloc((bi[n]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr); 249 b->j = bj; 250 b->i = bi; 251 b->diag = bdiag; 252 b->ilen = 0; 253 b->imax = 0; 254 b->row = isrow; 255 b->col = iscol; 256 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 257 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 258 b->icol = isicol; 259 ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 260 261 /* In b structure: Free imax, ilen, old a, old j. Allocate solve_work, new a, new j */ 262 ierr = PetscLogObjectMemory(B,(bi[n]-n)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr); 263 b->maxnz = b->nz = bi[n] ; 264 265 (B)->factor = MAT_FACTOR_LU; 266 (B)->info.factor_mallocs = reallocs; 267 (B)->info.fill_ratio_given = f; 268 269 if (ai[n]) { 270 (B)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 271 } else { 272 (B)->info.fill_ratio_needed = 0.0; 273 } 274 (B)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ; 275 (B)->ops->solve = MatSolve_SeqAIJ; 276 (B)->ops->solvetranspose = MatSolveTranspose_SeqAIJ; 277 /* switch to inodes if appropriate */ 278 ierr = MatLUFactorSymbolic_Inode(B,A,isrow,iscol,info);CHKERRQ(ierr); 279 PetscFunctionReturn(0); 280 } 281 282 /* 283 Trouble in factorization, should we dump the original matrix? 284 */ 285 #undef __FUNCT__ 286 #define __FUNCT__ "MatFactorDumpMatrix" 287 PetscErrorCode MatFactorDumpMatrix(Mat A) 288 { 289 PetscErrorCode ierr; 290 PetscTruth flg = PETSC_FALSE; 291 292 PetscFunctionBegin; 293 ierr = PetscOptionsGetTruth(PETSC_NULL,"-mat_factor_dump_on_error",&flg,PETSC_NULL);CHKERRQ(ierr); 294 if (flg) { 295 PetscViewer viewer; 296 char filename[PETSC_MAX_PATH_LEN]; 297 298 ierr = PetscSNPrintf(filename,PETSC_MAX_PATH_LEN,"matrix_factor_error.%d",PetscGlobalRank);CHKERRQ(ierr); 299 ierr = PetscViewerBinaryOpen(((PetscObject)A)->comm,filename,FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); 300 ierr = MatView(A,viewer);CHKERRQ(ierr); 301 ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 302 } 303 PetscFunctionReturn(0); 304 } 305 306 extern PetscErrorCode MatSolve_Inode(Mat,Vec,Vec); 307 308 /* ----------------------------------------------------------- */ 309 extern PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_iludt(Mat,Vec,Vec); 310 extern PetscErrorCode MatSolve_SeqAIJ_iludt(Mat,Vec,Vec); 311 312 #undef __FUNCT__ 313 #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_newdatastruct" 314 PetscErrorCode MatLUFactorNumeric_SeqAIJ_newdatastruct(Mat B,Mat A,const MatFactorInfo *info) 315 { 316 Mat C=B; 317 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data; 318 IS isrow = b->row,isicol = b->icol; 319 PetscErrorCode ierr; 320 const PetscInt *r,*ic,*ics; 321 PetscInt i,j,k,n=A->rmap->n,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 322 PetscInt *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj; 323 MatScalar *rtmp,*pc,multiplier,*v,*pv,*aa=a->a; 324 PetscReal shift=info->shiftinblocks; 325 PetscTruth row_identity, col_identity; 326 327 PetscFunctionBegin; 328 /* printf("MatLUFactorNumeric_SeqAIJ_newdatastruct is called ...\n"); */ 329 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 330 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 331 ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 332 ics = ic; 333 334 for (i=0; i<n; i++){ 335 /* zero rtmp */ 336 /* L part */ 337 nz = bi[i+1] - bi[i]; 338 bjtmp = bj + bi[i]; 339 for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0; 340 341 /* U part */ 342 /* nz = bi[2*n-i] - bi[2*n-i - 1] - 1; */ 343 nz = bdiag[i] - bdiag[i+1]; 344 bjtmp = bj + bdiag[i+1]+1; /* bi[2*n-i - 1]+1; */ 345 for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0; 346 347 /* load in initial (unfactored row) */ 348 nz = ai[r[i]+1] - ai[r[i]]; 349 ajtmp = aj + ai[r[i]]; 350 v = aa + ai[r[i]]; 351 for (j=0; j<nz; j++) { 352 rtmp[ics[ajtmp[j]]] = v[j]; 353 } 354 if (rtmp[ics[r[i]]] == 0.0){ 355 rtmp[ics[r[i]]] += shift; /* shift the diagonal of the matrix */ 356 /* printf("row %d, shift %g\n",i,shift); */ 357 } 358 359 /* elimination */ 360 bjtmp = bj + bi[i]; 361 row = *bjtmp++; 362 nzL = bi[i+1] - bi[i]; 363 k = 0; 364 while (k < nzL) { 365 pc = rtmp + row; 366 if (*pc != 0.0) { 367 pv = b->a + bdiag[row]; 368 multiplier = *pc * (*pv); 369 *pc = multiplier; 370 pj = b->j + bdiag[row+1] + 1; 371 pv = b->a + bdiag[row+1] + 1; 372 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diag */ 373 for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 374 ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); 375 } 376 row = *bjtmp++; k++; 377 } 378 379 /* finished row so stick it into b->a */ 380 /* L part */ 381 pv = b->a + bi[i] ; 382 pj = b->j + bi[i] ; 383 nz = bi[i+1] - bi[i]; 384 for (j=0; j<nz; j++) { 385 pv[j] = rtmp[pj[j]]; 386 } 387 388 /* Mark diagonal and invert diagonal for simplier triangular solves */ 389 pv = b->a + bdiag[i]; 390 pj = b->j + bdiag[i]; 391 /* if (*pj != i)SETERRQ2(PETSC_ERR_SUP,"row %d != *pj %d",i,*pj) */ 392 *pv = 1.0/rtmp[*pj]; 393 394 /* U part */ 395 pv = b->a + bdiag[i+1]+1; 396 pj = b->j + bdiag[i+1]+1; 397 nz = bdiag[i] - bdiag[i+1]-1; 398 for (j=0; j<nz; j++) { 399 pv[j] = rtmp[pj[j]]; 400 } 401 402 } 403 ierr = PetscFree(rtmp);CHKERRQ(ierr); 404 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 405 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 406 407 408 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 409 ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr); 410 if (row_identity && col_identity) { 411 C->ops->solve = MatSolve_SeqAIJ_NaturalOrdering_iludt; 412 } else { 413 C->ops->solve = MatSolve_SeqAIJ_iludt; 414 } 415 416 C->ops->solveadd = 0; 417 C->ops->solvetranspose = 0; 418 C->ops->solvetransposeadd = 0; 419 C->ops->matsolve = 0; 420 C->assembled = PETSC_TRUE; 421 C->preallocated = PETSC_TRUE; 422 ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr); 423 PetscFunctionReturn(0); 424 } 425 426 #undef __FUNCT__ 427 #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ" 428 PetscErrorCode MatLUFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info) 429 { 430 Mat C=B; 431 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data; 432 IS isrow = b->row,isicol = b->icol; 433 PetscErrorCode ierr; 434 const PetscInt *r,*ic,*ics; 435 PetscInt i,j,n=A->rmap->n,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 436 PetscInt *ajtmp,*bjtmp,nz,row,*diag_offset = b->diag,diag,*pj; 437 MatScalar *rtmp,*pc,multiplier,*v,*pv,d,*aa=a->a; 438 PetscReal rs=0.0; 439 LUShift_Ctx sctx; 440 PetscInt newshift,*ddiag; 441 442 PetscFunctionBegin; 443 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 444 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 445 ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 446 ics = ic; 447 448 sctx.shift_top = 0; 449 sctx.nshift_max = 0; 450 sctx.shift_lo = 0; 451 sctx.shift_hi = 0; 452 sctx.shift_fraction = 0; 453 454 /* if both shift schemes are chosen by user, only use info->shiftpd */ 455 if (info->shiftpd) { /* set sctx.shift_top=max{rs} */ 456 ddiag = a->diag; 457 sctx.shift_top = info->zeropivot; 458 for (i=0; i<n; i++) { 459 /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */ 460 d = (aa)[ddiag[i]]; 461 rs = -PetscAbsScalar(d) - PetscRealPart(d); 462 v = aa+ai[i]; 463 nz = ai[i+1] - ai[i]; 464 for (j=0; j<nz; j++) 465 rs += PetscAbsScalar(v[j]); 466 if (rs>sctx.shift_top) sctx.shift_top = rs; 467 } 468 sctx.shift_top *= 1.1; 469 sctx.nshift_max = 5; 470 sctx.shift_lo = 0.; 471 sctx.shift_hi = 1.; 472 } 473 474 sctx.shift_amount = 0.0; 475 sctx.nshift = 0; 476 do { 477 sctx.lushift = PETSC_FALSE; 478 for (i=0; i<n; i++){ 479 nz = bi[i+1] - bi[i]; 480 bjtmp = bj + bi[i]; 481 for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0; 482 483 /* load in initial (unfactored row) */ 484 nz = ai[r[i]+1] - ai[r[i]]; 485 ajtmp = aj + ai[r[i]]; 486 v = aa + ai[r[i]]; 487 for (j=0; j<nz; j++) { 488 rtmp[ics[ajtmp[j]]] = v[j]; 489 } 490 rtmp[ics[r[i]]] += sctx.shift_amount; /* shift the diagonal of the matrix */ 491 if (sctx.shift_amount > 0.0) printf("row %d, shift %g\n",i,sctx.shift_amount); 492 493 row = *bjtmp++; 494 while (row < i) { 495 pc = rtmp + row; 496 if (*pc != 0.0) { 497 pv = b->a + diag_offset[row]; 498 pj = b->j + diag_offset[row] + 1; 499 multiplier = *pc / *pv++; 500 *pc = multiplier; 501 nz = bi[row+1] - diag_offset[row] - 1; 502 for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 503 ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); 504 } 505 row = *bjtmp++; 506 } 507 /* finished row so stick it into b->a */ 508 pv = b->a + bi[i] ; 509 pj = b->j + bi[i] ; 510 nz = bi[i+1] - bi[i]; 511 diag = diag_offset[i] - bi[i]; 512 rs = -PetscAbsScalar(pv[diag]); 513 for (j=0; j<nz; j++) { 514 pv[j] = rtmp[pj[j]]; 515 rs += PetscAbsScalar(pv[j]); 516 } 517 518 /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */ 519 sctx.rs = rs; 520 sctx.pv = pv[diag]; 521 ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr); 522 if (newshift == 1) break; 523 } 524 525 if (info->shiftpd && !sctx.lushift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) { 526 /* 527 * if no shift in this attempt & shifting & started shifting & can refine, 528 * then try lower shift 529 */ 530 sctx.shift_hi = sctx.shift_fraction; 531 sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.; 532 sctx.shift_amount = sctx.shift_fraction * sctx.shift_top; 533 sctx.lushift = PETSC_TRUE; 534 sctx.nshift++; 535 } 536 } while (sctx.lushift); 537 538 /* invert diagonal entries for simplier triangular solves */ 539 for (i=0; i<n; i++) { 540 b->a[diag_offset[i]] = 1.0/b->a[diag_offset[i]]; 541 } 542 ierr = PetscFree(rtmp);CHKERRQ(ierr); 543 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 544 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 545 if (b->inode.use) { 546 C->ops->solve = MatSolve_Inode; 547 } else { 548 PetscTruth row_identity, col_identity; 549 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 550 ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr); 551 if (row_identity && col_identity) { 552 C->ops->solve = MatSolve_SeqAIJ_NaturalOrdering; 553 } else { 554 C->ops->solve = MatSolve_SeqAIJ; 555 } 556 } 557 C->ops->solveadd = MatSolveAdd_SeqAIJ; 558 C->ops->solvetranspose = MatSolveTranspose_SeqAIJ; 559 C->ops->solvetransposeadd = MatSolveTransposeAdd_SeqAIJ; 560 C->ops->matsolve = MatMatSolve_SeqAIJ; 561 C->assembled = PETSC_TRUE; 562 C->preallocated = PETSC_TRUE; 563 ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr); 564 if (sctx.nshift){ 565 if (info->shiftpd) { 566 ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,sctx.shift_fraction,sctx.shift_top);CHKERRQ(ierr); 567 } else if (info->shiftnz) { 568 ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 569 } 570 } 571 PetscFunctionReturn(0); 572 } 573 574 /* 575 This routine implements inplace ILU(0) with row or/and column permutations. 576 Input: 577 A - original matrix 578 Output; 579 A - a->i (rowptr) is same as original rowptr, but factored i-the row is stored in rowperm[i] 580 a->j (col index) is permuted by the inverse of colperm, then sorted 581 a->a reordered accordingly with a->j 582 a->diag (ptr to diagonal elements) is updated. 583 */ 584 #undef __FUNCT__ 585 #define __FUNCT__ "MatLUFactorNumeric_SeqAIJ_InplaceWithPerm" 586 PetscErrorCode MatLUFactorNumeric_SeqAIJ_InplaceWithPerm(Mat B,Mat A,const MatFactorInfo *info) 587 { 588 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 589 IS isrow = a->row,isicol = a->icol; 590 PetscErrorCode ierr; 591 const PetscInt *r,*ic,*ics; 592 PetscInt i,j,n=A->rmap->n,*ai=a->i,*aj=a->j; 593 PetscInt *ajtmp,nz,row; 594 PetscInt *diag = a->diag,nbdiag,*pj; 595 PetscScalar *rtmp,*pc,multiplier,d; 596 MatScalar *v,*pv; 597 PetscReal rs; 598 LUShift_Ctx sctx; 599 PetscInt newshift; 600 601 PetscFunctionBegin; 602 if (A != B) SETERRQ(PETSC_ERR_ARG_INCOMP,"input and output matrix must have same address"); 603 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 604 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 605 ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&rtmp);CHKERRQ(ierr); 606 ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr); 607 ics = ic; 608 609 sctx.shift_top = 0; 610 sctx.nshift_max = 0; 611 sctx.shift_lo = 0; 612 sctx.shift_hi = 0; 613 sctx.shift_fraction = 0; 614 615 /* if both shift schemes are chosen by user, only use info->shiftpd */ 616 if (info->shiftpd) { /* set sctx.shift_top=max{rs} */ 617 sctx.shift_top = 0; 618 for (i=0; i<n; i++) { 619 /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */ 620 d = (a->a)[diag[i]]; 621 rs = -PetscAbsScalar(d) - PetscRealPart(d); 622 v = a->a+ai[i]; 623 nz = ai[i+1] - ai[i]; 624 for (j=0; j<nz; j++) 625 rs += PetscAbsScalar(v[j]); 626 if (rs>sctx.shift_top) sctx.shift_top = rs; 627 } 628 if (sctx.shift_top < info->zeropivot) sctx.shift_top = info->zeropivot; 629 sctx.shift_top *= 1.1; 630 sctx.nshift_max = 5; 631 sctx.shift_lo = 0.; 632 sctx.shift_hi = 1.; 633 } 634 635 sctx.shift_amount = 0; 636 sctx.nshift = 0; 637 do { 638 sctx.lushift = PETSC_FALSE; 639 for (i=0; i<n; i++){ 640 /* load in initial unfactored row */ 641 nz = ai[r[i]+1] - ai[r[i]]; 642 ajtmp = aj + ai[r[i]]; 643 v = a->a + ai[r[i]]; 644 /* sort permuted ajtmp and values v accordingly */ 645 for (j=0; j<nz; j++) ajtmp[j] = ics[ajtmp[j]]; 646 ierr = PetscSortIntWithScalarArray(nz,ajtmp,v);CHKERRQ(ierr); 647 648 diag[r[i]] = ai[r[i]]; 649 for (j=0; j<nz; j++) { 650 rtmp[ajtmp[j]] = v[j]; 651 if (ajtmp[j] < i) diag[r[i]]++; /* update a->diag */ 652 } 653 rtmp[r[i]] += sctx.shift_amount; /* shift the diagonal of the matrix */ 654 655 row = *ajtmp++; 656 while (row < i) { 657 pc = rtmp + row; 658 if (*pc != 0.0) { 659 pv = a->a + diag[r[row]]; 660 pj = aj + diag[r[row]] + 1; 661 662 multiplier = *pc / *pv++; 663 *pc = multiplier; 664 nz = ai[r[row]+1] - diag[r[row]] - 1; 665 for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 666 ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); 667 } 668 row = *ajtmp++; 669 } 670 /* finished row so overwrite it onto a->a */ 671 pv = a->a + ai[r[i]] ; 672 pj = aj + ai[r[i]] ; 673 nz = ai[r[i]+1] - ai[r[i]]; 674 nbdiag = diag[r[i]] - ai[r[i]]; /* num of entries before the diagonal */ 675 676 rs = 0.0; 677 for (j=0; j<nz; j++) { 678 pv[j] = rtmp[pj[j]]; 679 if (j != nbdiag) rs += PetscAbsScalar(pv[j]); 680 } 681 682 /* 9/13/02 Victor Eijkhout suggested scaling zeropivot by rs for matrices with funny scalings */ 683 sctx.rs = rs; 684 sctx.pv = pv[nbdiag]; 685 ierr = MatLUCheckShift_inline(info,sctx,i,newshift);CHKERRQ(ierr); 686 if (newshift == 1) break; 687 } 688 689 if (info->shiftpd && !sctx.lushift && sctx.shift_fraction>0 && sctx.nshift<sctx.nshift_max) { 690 /* 691 * if no shift in this attempt & shifting & started shifting & can refine, 692 * then try lower shift 693 */ 694 sctx.shift_hi = sctx.shift_fraction; 695 sctx.shift_fraction = (sctx.shift_hi+sctx.shift_lo)/2.; 696 sctx.shift_amount = sctx.shift_fraction * sctx.shift_top; 697 sctx.lushift = PETSC_TRUE; 698 sctx.nshift++; 699 } 700 } while (sctx.lushift); 701 702 /* invert diagonal entries for simplier triangular solves */ 703 for (i=0; i<n; i++) { 704 a->a[diag[r[i]]] = 1.0/a->a[diag[r[i]]]; 705 } 706 707 ierr = PetscFree(rtmp);CHKERRQ(ierr); 708 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 709 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 710 A->ops->solve = MatSolve_SeqAIJ_InplaceWithPerm; 711 A->ops->solveadd = MatSolveAdd_SeqAIJ; 712 A->ops->solvetranspose = MatSolveTranspose_SeqAIJ; 713 A->ops->solvetransposeadd = MatSolveTransposeAdd_SeqAIJ; 714 A->assembled = PETSC_TRUE; 715 A->preallocated = PETSC_TRUE; 716 ierr = PetscLogFlops(A->cmap->n);CHKERRQ(ierr); 717 if (sctx.nshift){ 718 if (info->shiftpd) { 719 ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %G, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,sctx.shift_amount,sctx.shift_fraction,sctx.shift_top);CHKERRQ(ierr); 720 } else if (info->shiftnz) { 721 ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 722 } 723 } 724 PetscFunctionReturn(0); 725 } 726 727 /* ----------------------------------------------------------- */ 728 #undef __FUNCT__ 729 #define __FUNCT__ "MatLUFactor_SeqAIJ" 730 PetscErrorCode MatLUFactor_SeqAIJ(Mat A,IS row,IS col,const MatFactorInfo *info) 731 { 732 PetscErrorCode ierr; 733 Mat C; 734 735 PetscFunctionBegin; 736 ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr); 737 ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr); 738 ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr); 739 A->ops->solve = C->ops->solve; 740 A->ops->solvetranspose = C->ops->solvetranspose; 741 ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 742 ierr = PetscLogObjectParent(A,((Mat_SeqAIJ*)(A->data))->icol);CHKERRQ(ierr); 743 PetscFunctionReturn(0); 744 } 745 /* ----------------------------------------------------------- */ 746 747 748 #undef __FUNCT__ 749 #define __FUNCT__ "MatSolve_SeqAIJ" 750 PetscErrorCode MatSolve_SeqAIJ(Mat A,Vec bb,Vec xx) 751 { 752 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 753 IS iscol = a->col,isrow = a->row; 754 PetscErrorCode ierr; 755 PetscInt i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j; 756 PetscInt nz; 757 const PetscInt *rout,*cout,*r,*c; 758 PetscScalar *x,*tmp,*tmps,sum; 759 const PetscScalar *b; 760 const MatScalar *aa = a->a,*v; 761 762 PetscFunctionBegin; 763 if (!n) PetscFunctionReturn(0); 764 765 ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 766 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 767 tmp = a->solve_work; 768 769 ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 770 ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 771 772 /* forward solve the lower triangular */ 773 tmp[0] = b[*r++]; 774 tmps = tmp; 775 for (i=1; i<n; i++) { 776 v = aa + ai[i] ; 777 vi = aj + ai[i] ; 778 nz = a->diag[i] - ai[i]; 779 sum = b[*r++]; 780 PetscSparseDenseMinusDot(sum,tmps,v,vi,nz); 781 tmp[i] = sum; 782 } 783 784 /* backward solve the upper triangular */ 785 for (i=n-1; i>=0; i--){ 786 v = aa + a->diag[i] + 1; 787 vi = aj + a->diag[i] + 1; 788 nz = ai[i+1] - a->diag[i] - 1; 789 sum = tmp[i]; 790 PetscSparseDenseMinusDot(sum,tmps,v,vi,nz); 791 x[*c--] = tmp[i] = sum*aa[a->diag[i]]; 792 } 793 794 ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 795 ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 796 ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 797 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 798 ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr); 799 PetscFunctionReturn(0); 800 } 801 802 #undef __FUNCT__ 803 #define __FUNCT__ "MatMatSolve_SeqAIJ" 804 PetscErrorCode MatMatSolve_SeqAIJ(Mat A,Mat B,Mat X) 805 { 806 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 807 IS iscol = a->col,isrow = a->row; 808 PetscErrorCode ierr; 809 PetscInt i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j; 810 PetscInt nz,neq; 811 const PetscInt *rout,*cout,*r,*c; 812 PetscScalar *x,*b,*tmp,*tmps,sum; 813 const MatScalar *aa = a->a,*v; 814 PetscTruth bisdense,xisdense; 815 816 PetscFunctionBegin; 817 if (!n) PetscFunctionReturn(0); 818 819 ierr = PetscTypeCompare((PetscObject)B,MATSEQDENSE,&bisdense);CHKERRQ(ierr); 820 if (!bisdense) SETERRQ(PETSC_ERR_ARG_INCOMP,"B matrix must be a SeqDense matrix"); 821 ierr = PetscTypeCompare((PetscObject)X,MATSEQDENSE,&xisdense);CHKERRQ(ierr); 822 if (!xisdense) SETERRQ(PETSC_ERR_ARG_INCOMP,"X matrix must be a SeqDense matrix"); 823 824 ierr = MatGetArray(B,&b);CHKERRQ(ierr); 825 ierr = MatGetArray(X,&x);CHKERRQ(ierr); 826 827 tmp = a->solve_work; 828 ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 829 ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 830 831 for (neq=0; neq<B->cmap->n; neq++){ 832 /* forward solve the lower triangular */ 833 tmp[0] = b[r[0]]; 834 tmps = tmp; 835 for (i=1; i<n; i++) { 836 v = aa + ai[i] ; 837 vi = aj + ai[i] ; 838 nz = a->diag[i] - ai[i]; 839 sum = b[r[i]]; 840 PetscSparseDenseMinusDot(sum,tmps,v,vi,nz); 841 tmp[i] = sum; 842 } 843 /* backward solve the upper triangular */ 844 for (i=n-1; i>=0; i--){ 845 v = aa + a->diag[i] + 1; 846 vi = aj + a->diag[i] + 1; 847 nz = ai[i+1] - a->diag[i] - 1; 848 sum = tmp[i]; 849 PetscSparseDenseMinusDot(sum,tmps,v,vi,nz); 850 x[c[i]] = tmp[i] = sum*aa[a->diag[i]]; 851 } 852 853 b += n; 854 x += n; 855 } 856 ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 857 ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 858 ierr = MatRestoreArray(B,&b);CHKERRQ(ierr); 859 ierr = MatRestoreArray(X,&x);CHKERRQ(ierr); 860 ierr = PetscLogFlops(B->cmap->n*(2.0*a->nz - n));CHKERRQ(ierr); 861 PetscFunctionReturn(0); 862 } 863 864 #undef __FUNCT__ 865 #define __FUNCT__ "MatSolve_SeqAIJ_InplaceWithPerm" 866 PetscErrorCode MatSolve_SeqAIJ_InplaceWithPerm(Mat A,Vec bb,Vec xx) 867 { 868 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 869 IS iscol = a->col,isrow = a->row; 870 PetscErrorCode ierr; 871 const PetscInt *r,*c,*rout,*cout; 872 PetscInt i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j; 873 PetscInt nz,row; 874 PetscScalar *x,*b,*tmp,*tmps,sum; 875 const MatScalar *aa = a->a,*v; 876 877 PetscFunctionBegin; 878 if (!n) PetscFunctionReturn(0); 879 880 ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 881 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 882 tmp = a->solve_work; 883 884 ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 885 ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 886 887 /* forward solve the lower triangular */ 888 tmp[0] = b[*r++]; 889 tmps = tmp; 890 for (row=1; row<n; row++) { 891 i = rout[row]; /* permuted row */ 892 v = aa + ai[i] ; 893 vi = aj + ai[i] ; 894 nz = a->diag[i] - ai[i]; 895 sum = b[*r++]; 896 PetscSparseDenseMinusDot(sum,tmps,v,vi,nz); 897 tmp[row] = sum; 898 } 899 900 /* backward solve the upper triangular */ 901 for (row=n-1; row>=0; row--){ 902 i = rout[row]; /* permuted row */ 903 v = aa + a->diag[i] + 1; 904 vi = aj + a->diag[i] + 1; 905 nz = ai[i+1] - a->diag[i] - 1; 906 sum = tmp[row]; 907 PetscSparseDenseMinusDot(sum,tmps,v,vi,nz); 908 x[*c--] = tmp[row] = sum*aa[a->diag[i]]; 909 } 910 911 ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 912 ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 913 ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 914 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 915 ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr); 916 PetscFunctionReturn(0); 917 } 918 919 /* ----------------------------------------------------------- */ 920 #undef __FUNCT__ 921 #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering" 922 PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering(Mat A,Vec bb,Vec xx) 923 { 924 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 925 PetscErrorCode ierr; 926 PetscInt n = A->rmap->n; 927 const PetscInt *ai = a->i,*aj = a->j,*adiag = a->diag,*vi; 928 PetscScalar *x; 929 const PetscScalar *b; 930 const MatScalar *aa = a->a; 931 #if !defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 932 PetscInt adiag_i,i,nz,ai_i; 933 const MatScalar *v; 934 PetscScalar sum; 935 #endif 936 937 PetscFunctionBegin; 938 if (!n) PetscFunctionReturn(0); 939 940 ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 941 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 942 943 #if defined(PETSC_USE_FORTRAN_KERNEL_SOLVEAIJ) 944 fortransolveaij_(&n,x,ai,aj,adiag,aa,b); 945 #else 946 /* forward solve the lower triangular */ 947 x[0] = b[0]; 948 for (i=1; i<n; i++) { 949 ai_i = ai[i]; 950 v = aa + ai_i; 951 vi = aj + ai_i; 952 nz = adiag[i] - ai_i; 953 sum = b[i]; 954 PetscSparseDenseMinusDot(sum,x,v,vi,nz); 955 x[i] = sum; 956 } 957 958 /* backward solve the upper triangular */ 959 for (i=n-1; i>=0; i--){ 960 adiag_i = adiag[i]; 961 v = aa + adiag_i + 1; 962 vi = aj + adiag_i + 1; 963 nz = ai[i+1] - adiag_i - 1; 964 sum = x[i]; 965 PetscSparseDenseMinusDot(sum,x,v,vi,nz); 966 x[i] = sum*aa[adiag_i]; 967 } 968 #endif 969 ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr); 970 ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 971 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 972 PetscFunctionReturn(0); 973 } 974 975 #undef __FUNCT__ 976 #define __FUNCT__ "MatSolveAdd_SeqAIJ" 977 PetscErrorCode MatSolveAdd_SeqAIJ(Mat A,Vec bb,Vec yy,Vec xx) 978 { 979 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 980 IS iscol = a->col,isrow = a->row; 981 PetscErrorCode ierr; 982 PetscInt i, n = A->rmap->n,*vi,*ai = a->i,*aj = a->j; 983 PetscInt nz; 984 const PetscInt *rout,*cout,*r,*c; 985 PetscScalar *x,*b,*tmp,sum; 986 const MatScalar *aa = a->a,*v; 987 988 PetscFunctionBegin; 989 if (yy != xx) {ierr = VecCopy(yy,xx);CHKERRQ(ierr);} 990 991 ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 992 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 993 tmp = a->solve_work; 994 995 ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 996 ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 997 998 /* forward solve the lower triangular */ 999 tmp[0] = b[*r++]; 1000 for (i=1; i<n; i++) { 1001 v = aa + ai[i] ; 1002 vi = aj + ai[i] ; 1003 nz = a->diag[i] - ai[i]; 1004 sum = b[*r++]; 1005 while (nz--) sum -= *v++ * tmp[*vi++ ]; 1006 tmp[i] = sum; 1007 } 1008 1009 /* backward solve the upper triangular */ 1010 for (i=n-1; i>=0; i--){ 1011 v = aa + a->diag[i] + 1; 1012 vi = aj + a->diag[i] + 1; 1013 nz = ai[i+1] - a->diag[i] - 1; 1014 sum = tmp[i]; 1015 while (nz--) sum -= *v++ * tmp[*vi++ ]; 1016 tmp[i] = sum*aa[a->diag[i]]; 1017 x[*c--] += tmp[i]; 1018 } 1019 1020 ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 1021 ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 1022 ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 1023 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1024 ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 1025 1026 PetscFunctionReturn(0); 1027 } 1028 /* -------------------------------------------------------------------*/ 1029 #undef __FUNCT__ 1030 #define __FUNCT__ "MatSolveTranspose_SeqAIJ" 1031 PetscErrorCode MatSolveTranspose_SeqAIJ(Mat A,Vec bb,Vec xx) 1032 { 1033 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1034 IS iscol = a->col,isrow = a->row; 1035 PetscErrorCode ierr; 1036 const PetscInt *rout,*cout,*r,*c; 1037 PetscInt i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j; 1038 PetscInt nz,*diag = a->diag; 1039 PetscScalar *x,*b,*tmp,s1; 1040 const MatScalar *aa = a->a,*v; 1041 1042 PetscFunctionBegin; 1043 ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1044 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1045 tmp = a->solve_work; 1046 1047 ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 1048 ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 1049 1050 /* copy the b into temp work space according to permutation */ 1051 for (i=0; i<n; i++) tmp[i] = b[c[i]]; 1052 1053 /* forward solve the U^T */ 1054 for (i=0; i<n; i++) { 1055 v = aa + diag[i] ; 1056 vi = aj + diag[i] + 1; 1057 nz = ai[i+1] - diag[i] - 1; 1058 s1 = tmp[i]; 1059 s1 *= (*v++); /* multiply by inverse of diagonal entry */ 1060 while (nz--) { 1061 tmp[*vi++ ] -= (*v++)*s1; 1062 } 1063 tmp[i] = s1; 1064 } 1065 1066 /* backward solve the L^T */ 1067 for (i=n-1; i>=0; i--){ 1068 v = aa + diag[i] - 1 ; 1069 vi = aj + diag[i] - 1 ; 1070 nz = diag[i] - ai[i]; 1071 s1 = tmp[i]; 1072 while (nz--) { 1073 tmp[*vi-- ] -= (*v--)*s1; 1074 } 1075 } 1076 1077 /* copy tmp into x according to permutation */ 1078 for (i=0; i<n; i++) x[r[i]] = tmp[i]; 1079 1080 ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 1081 ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 1082 ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 1083 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1084 1085 ierr = PetscLogFlops(2.0*a->nz-A->cmap->n);CHKERRQ(ierr); 1086 PetscFunctionReturn(0); 1087 } 1088 1089 #undef __FUNCT__ 1090 #define __FUNCT__ "MatSolveTransposeAdd_SeqAIJ" 1091 PetscErrorCode MatSolveTransposeAdd_SeqAIJ(Mat A,Vec bb,Vec zz,Vec xx) 1092 { 1093 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1094 IS iscol = a->col,isrow = a->row; 1095 PetscErrorCode ierr; 1096 const PetscInt *r,*c,*rout,*cout; 1097 PetscInt i,n = A->rmap->n,*vi,*ai = a->i,*aj = a->j; 1098 PetscInt nz,*diag = a->diag; 1099 PetscScalar *x,*b,*tmp; 1100 const MatScalar *aa = a->a,*v; 1101 1102 PetscFunctionBegin; 1103 if (zz != xx) {ierr = VecCopy(zz,xx);CHKERRQ(ierr);} 1104 1105 ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1106 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1107 tmp = a->solve_work; 1108 1109 ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 1110 ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 1111 1112 /* copy the b into temp work space according to permutation */ 1113 for (i=0; i<n; i++) tmp[i] = b[c[i]]; 1114 1115 /* forward solve the U^T */ 1116 for (i=0; i<n; i++) { 1117 v = aa + diag[i] ; 1118 vi = aj + diag[i] + 1; 1119 nz = ai[i+1] - diag[i] - 1; 1120 tmp[i] *= *v++; 1121 while (nz--) { 1122 tmp[*vi++ ] -= (*v++)*tmp[i]; 1123 } 1124 } 1125 1126 /* backward solve the L^T */ 1127 for (i=n-1; i>=0; i--){ 1128 v = aa + diag[i] - 1 ; 1129 vi = aj + diag[i] - 1 ; 1130 nz = diag[i] - ai[i]; 1131 while (nz--) { 1132 tmp[*vi-- ] -= (*v--)*tmp[i]; 1133 } 1134 } 1135 1136 /* copy tmp into x according to permutation */ 1137 for (i=0; i<n; i++) x[r[i]] += tmp[i]; 1138 1139 ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 1140 ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 1141 ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 1142 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1143 1144 ierr = PetscLogFlops(2.0*a->nz);CHKERRQ(ierr); 1145 PetscFunctionReturn(0); 1146 } 1147 /* ----------------------------------------------------------------*/ 1148 EXTERN PetscErrorCode Mat_CheckInode(Mat,PetscTruth); 1149 EXTERN PetscErrorCode MatDuplicateNoCreate_SeqAIJ(Mat,Mat,MatDuplicateOption); 1150 1151 /* ilu(0) with natural ordering under new data structure */ 1152 #undef __FUNCT__ 1153 #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ_newdatastruct" 1154 PetscErrorCode MatILUFactorSymbolic_SeqAIJ_newdatastruct(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 1155 { 1156 1157 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 1158 PetscErrorCode ierr; 1159 PetscInt n=A->rmap->n,*ai=a->i,*aj,*adiag=a->diag; 1160 PetscInt i,j,nz,*bi,*bj,*bdiag,*bdiag_rev; 1161 1162 PetscFunctionBegin; 1163 ierr = MatDuplicateNoCreate_SeqAIJ(fact,A,MAT_DO_NOT_COPY_VALUES);CHKERRQ(ierr); 1164 b = (Mat_SeqAIJ*)(fact)->data; 1165 1166 /* replace matrix arrays with single allocations, then reset values */ 1167 ierr = PetscFree3(b->a,b->j,b->i);CHKERRQ(ierr); 1168 ierr = PetscFree(b->diag);CHKERRQ(ierr); 1169 1170 ierr = PetscMalloc((2*n+2)*sizeof(PetscInt),&b->i);CHKERRQ(ierr); 1171 ierr = PetscMalloc((ai[n]+2)*sizeof(PetscInt),&b->j);CHKERRQ(ierr); 1172 ierr = PetscMalloc((ai[n]+2)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr); 1173 b->singlemalloc = PETSC_FALSE; 1174 ierr = PetscMalloc((2*n+2)*sizeof(PetscInt),&b->diag);CHKERRQ(ierr); 1175 bdiag = b->diag; 1176 bdiag_rev = bdiag + n+1; 1177 1178 if (n > 0) { 1179 ierr = PetscMemzero(b->a,(ai[n])*sizeof(PetscScalar));CHKERRQ(ierr); 1180 } 1181 1182 /* set bi and bj with new data structure */ 1183 bi = b->i; 1184 bj = b->j; 1185 1186 /* L part */ 1187 bi[0] = 0; 1188 for (i=0; i<n; i++){ 1189 nz = adiag[i] - ai[i]; 1190 bi[i+1] = bi[i] + nz; 1191 aj = a->j + ai[i]; 1192 for (j=0; j<nz; j++){ 1193 *bj = aj[j]; bj++; 1194 } 1195 } 1196 1197 bj++; /* diag[n] */ 1198 bdiag[n] = bi[n]; 1199 1200 /* U part */ 1201 for (i=n-1; i>=0; i--){ 1202 nz = ai[i+1] - adiag[i] - 1; 1203 bi[n+ (n-i)] = bi[n+(n-i)-1] + nz + 1; 1204 aj = a->j + adiag[i] + 1; 1205 for (j=0; j<nz; j++){ 1206 *bj = aj[j]; bj++; 1207 } 1208 /* diag[i] */ 1209 *bj = i; bj++; 1210 bdiag[i] = bi[n+ (n-i)] ; 1211 bdiag_rev[n-i] = bdiag[i]; 1212 } 1213 bdiag_rev[0] = bdiag[n]; 1214 PetscFunctionReturn(0); 1215 } 1216 1217 #undef __FUNCT__ 1218 #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ" 1219 PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 1220 { 1221 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 1222 IS isicol; 1223 PetscErrorCode ierr; 1224 const PetscInt *r,*ic; 1225 PetscInt n=A->rmap->n,*ai=a->i,*aj=a->j,d; 1226 PetscInt *bi,*cols,nnz,*cols_lvl; 1227 PetscInt *bdiag,prow,fm,nzbd,reallocs=0,dcount=0; 1228 PetscInt i,levels,diagonal_fill; 1229 PetscTruth col_identity,row_identity; 1230 PetscReal f; 1231 PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL; 1232 PetscBT lnkbt; 1233 PetscInt nzi,*bj,**bj_ptr,**bjlvl_ptr; 1234 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1235 PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 1236 PetscTruth missing; 1237 1238 PetscFunctionBegin; 1239 if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n); 1240 f = info->fill; 1241 levels = (PetscInt)info->levels; 1242 diagonal_fill = (PetscInt)info->diagonal_fill; 1243 ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 1244 1245 /* special case that simply copies fill pattern */ 1246 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 1247 ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 1248 if (!levels && row_identity && col_identity) { 1249 1250 PetscTruth newdatastruct=PETSC_FALSE; 1251 ierr = PetscOptionsGetTruth(PETSC_NULL,"-ilu_new",&newdatastruct,PETSC_NULL);CHKERRQ(ierr); 1252 if (newdatastruct){ 1253 ierr = MatILUFactorSymbolic_SeqAIJ_newdatastruct(fact,A,isrow,iscol,info);CHKERRQ(ierr); 1254 (fact)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_newdatastruct; 1255 } else { 1256 ierr = MatDuplicateNoCreate_SeqAIJ(fact,A,MAT_DO_NOT_COPY_VALUES);CHKERRQ(ierr); 1257 (fact)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ; 1258 } 1259 1260 fact->factor = MAT_FACTOR_ILU; 1261 (fact)->info.factor_mallocs = 0; 1262 (fact)->info.fill_ratio_given = info->fill; 1263 (fact)->info.fill_ratio_needed = 1.0; 1264 b = (Mat_SeqAIJ*)(fact)->data; 1265 ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 1266 if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 1267 b->row = isrow; 1268 b->col = iscol; 1269 b->icol = isicol; 1270 ierr = PetscMalloc(((fact)->rmap->n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1271 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1272 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 1273 ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr); 1274 PetscFunctionReturn(0); 1275 } 1276 1277 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 1278 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 1279 1280 /* get new row pointers */ 1281 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 1282 bi[0] = 0; 1283 /* bdiag is location of diagonal in factor */ 1284 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 1285 bdiag[0] = 0; 1286 1287 ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr); 1288 bjlvl_ptr = (PetscInt**)(bj_ptr + n); 1289 1290 /* create a linked list for storing column indices of the active row */ 1291 nlnk = n + 1; 1292 ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 1293 1294 /* initial FreeSpace size is f*(ai[n]+1) */ 1295 ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 1296 current_space = free_space; 1297 ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr); 1298 current_space_lvl = free_space_lvl; 1299 1300 for (i=0; i<n; i++) { 1301 nzi = 0; 1302 /* copy current row into linked list */ 1303 nnz = ai[r[i]+1] - ai[r[i]]; 1304 if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 1305 cols = aj + ai[r[i]]; 1306 lnk[i] = -1; /* marker to indicate if diagonal exists */ 1307 ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 1308 nzi += nlnk; 1309 1310 /* make sure diagonal entry is included */ 1311 if (diagonal_fill && lnk[i] == -1) { 1312 fm = n; 1313 while (lnk[fm] < i) fm = lnk[fm]; 1314 lnk[i] = lnk[fm]; /* insert diagonal into linked list */ 1315 lnk[fm] = i; 1316 lnk_lvl[i] = 0; 1317 nzi++; dcount++; 1318 } 1319 1320 /* add pivot rows into the active row */ 1321 nzbd = 0; 1322 prow = lnk[n]; 1323 while (prow < i) { 1324 nnz = bdiag[prow]; 1325 cols = bj_ptr[prow] + nnz + 1; 1326 cols_lvl = bjlvl_ptr[prow] + nnz + 1; 1327 nnz = bi[prow+1] - bi[prow] - nnz - 1; 1328 ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr); 1329 nzi += nlnk; 1330 prow = lnk[prow]; 1331 nzbd++; 1332 } 1333 bdiag[i] = nzbd; 1334 bi[i+1] = bi[i] + nzi; 1335 1336 /* if free space is not available, make more free space */ 1337 if (current_space->local_remaining<nzi) { 1338 nnz = nzi*(n - i); /* estimated and max additional space needed */ 1339 ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 1340 ierr = PetscFreeSpaceGet(nnz,¤t_space_lvl);CHKERRQ(ierr); 1341 reallocs++; 1342 } 1343 1344 /* copy data into free_space and free_space_lvl, then initialize lnk */ 1345 ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 1346 bj_ptr[i] = current_space->array; 1347 bjlvl_ptr[i] = current_space_lvl->array; 1348 1349 /* make sure the active row i has diagonal entry */ 1350 if (*(bj_ptr[i]+bdiag[i]) != i) { 1351 SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\ 1352 try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i); 1353 } 1354 1355 current_space->array += nzi; 1356 current_space->local_used += nzi; 1357 current_space->local_remaining -= nzi; 1358 current_space_lvl->array += nzi; 1359 current_space_lvl->local_used += nzi; 1360 current_space_lvl->local_remaining -= nzi; 1361 } 1362 1363 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 1364 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 1365 1366 /* destroy list of free space and other temporary arrays */ 1367 ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 1368 ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 1369 ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1370 ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 1371 ierr = PetscFree(bj_ptr);CHKERRQ(ierr); 1372 1373 #if defined(PETSC_USE_INFO) 1374 { 1375 PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 1376 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr); 1377 ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 1378 ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr); 1379 ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 1380 if (diagonal_fill) { 1381 ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr); 1382 } 1383 } 1384 #endif 1385 1386 /* put together the new matrix */ 1387 ierr = MatSeqAIJSetPreallocation_SeqAIJ(fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 1388 ierr = PetscLogObjectParent(fact,isicol);CHKERRQ(ierr); 1389 b = (Mat_SeqAIJ*)(fact)->data; 1390 b->free_a = PETSC_TRUE; 1391 b->free_ij = PETSC_TRUE; 1392 b->singlemalloc = PETSC_FALSE; 1393 ierr = PetscMalloc( (bi[n] )*sizeof(PetscScalar),&b->a);CHKERRQ(ierr); 1394 b->j = bj; 1395 b->i = bi; 1396 for (i=0; i<n; i++) bdiag[i] += bi[i]; 1397 b->diag = bdiag; 1398 b->ilen = 0; 1399 b->imax = 0; 1400 b->row = isrow; 1401 b->col = iscol; 1402 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1403 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 1404 b->icol = isicol; 1405 ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1406 /* In b structure: Free imax, ilen, old a, old j. 1407 Allocate bdiag, solve_work, new a, new j */ 1408 ierr = PetscLogObjectMemory(fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr); 1409 b->maxnz = b->nz = bi[n] ; 1410 (fact)->info.factor_mallocs = reallocs; 1411 (fact)->info.fill_ratio_given = f; 1412 (fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 1413 (fact)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ; 1414 ierr = MatILUFactorSymbolic_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr); 1415 PetscFunctionReturn(0); 1416 } 1417 1418 #include "../src/mat/impls/sbaij/seq/sbaij.h" 1419 #undef __FUNCT__ 1420 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ" 1421 PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info) 1422 { 1423 Mat C = B; 1424 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 1425 Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 1426 IS ip=b->row,iip = b->icol; 1427 PetscErrorCode ierr; 1428 const PetscInt *rip,*riip; 1429 PetscInt i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bcol; 1430 PetscInt *ai=a->i,*aj=a->j; 1431 PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 1432 MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 1433 PetscReal zeropivot,rs,shiftnz; 1434 PetscReal shiftpd; 1435 ChShift_Ctx sctx; 1436 PetscInt newshift; 1437 PetscTruth perm_identity; 1438 1439 PetscFunctionBegin; 1440 1441 shiftnz = info->shiftnz; 1442 shiftpd = info->shiftpd; 1443 zeropivot = info->zeropivot; 1444 1445 ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 1446 ierr = ISGetIndices(iip,&riip);CHKERRQ(ierr); 1447 1448 /* initialization */ 1449 nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 1450 ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 1451 jl = il + mbs; 1452 rtmp = (MatScalar*)(jl + mbs); 1453 1454 sctx.shift_amount = 0; 1455 sctx.nshift = 0; 1456 do { 1457 sctx.chshift = PETSC_FALSE; 1458 for (i=0; i<mbs; i++) { 1459 rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1460 } 1461 1462 for (k = 0; k<mbs; k++){ 1463 bval = ba + bi[k]; 1464 /* initialize k-th row by the perm[k]-th row of A */ 1465 jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 1466 for (j = jmin; j < jmax; j++){ 1467 col = riip[aj[j]]; 1468 if (col >= k){ /* only take upper triangular entry */ 1469 rtmp[col] = aa[j]; 1470 *bval++ = 0.0; /* for in-place factorization */ 1471 } 1472 } 1473 /* shift the diagonal of the matrix */ 1474 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 1475 1476 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 1477 dk = rtmp[k]; 1478 i = jl[k]; /* first row to be added to k_th row */ 1479 1480 while (i < k){ 1481 nexti = jl[i]; /* next row to be added to k_th row */ 1482 1483 /* compute multiplier, update diag(k) and U(i,k) */ 1484 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1485 uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 1486 dk += uikdi*ba[ili]; 1487 ba[ili] = uikdi; /* -U(i,k) */ 1488 1489 /* add multiple of row i to k-th row */ 1490 jmin = ili + 1; jmax = bi[i+1]; 1491 if (jmin < jmax){ 1492 for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 1493 /* update il and jl for row i */ 1494 il[i] = jmin; 1495 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 1496 } 1497 i = nexti; 1498 } 1499 1500 /* shift the diagonals when zero pivot is detected */ 1501 /* compute rs=sum of abs(off-diagonal) */ 1502 rs = 0.0; 1503 jmin = bi[k]+1; 1504 nz = bi[k+1] - jmin; 1505 bcol = bj + jmin; 1506 while (nz--){ 1507 rs += PetscAbsScalar(rtmp[*bcol]); 1508 bcol++; 1509 } 1510 1511 sctx.rs = rs; 1512 sctx.pv = dk; 1513 ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 1514 1515 if (newshift == 1) { 1516 if (!sctx.shift_amount) { 1517 sctx.shift_amount = 1e-5; 1518 } 1519 break; 1520 } 1521 1522 /* copy data into U(k,:) */ 1523 ba[bi[k]] = 1.0/dk; /* U(k,k) */ 1524 jmin = bi[k]+1; jmax = bi[k+1]; 1525 if (jmin < jmax) { 1526 for (j=jmin; j<jmax; j++){ 1527 col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 1528 } 1529 /* add the k-th row into il and jl */ 1530 il[k] = jmin; 1531 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 1532 } 1533 } 1534 } while (sctx.chshift); 1535 ierr = PetscFree(il);CHKERRQ(ierr); 1536 1537 ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 1538 ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr); 1539 1540 ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr); 1541 if (perm_identity){ 1542 (B)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1543 (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1544 (B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering; 1545 (B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering; 1546 } else { 1547 (B)->ops->solve = MatSolve_SeqSBAIJ_1; 1548 (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1; 1549 (B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1; 1550 (B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1; 1551 } 1552 1553 C->assembled = PETSC_TRUE; 1554 C->preallocated = PETSC_TRUE; 1555 ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr); 1556 if (sctx.nshift){ 1557 if (shiftnz) { 1558 ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 1559 } else if (shiftpd) { 1560 ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 1561 } 1562 } 1563 PetscFunctionReturn(0); 1564 } 1565 1566 #undef __FUNCT__ 1567 #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ" 1568 PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 1569 { 1570 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1571 Mat_SeqSBAIJ *b; 1572 PetscErrorCode ierr; 1573 PetscTruth perm_identity,missing; 1574 PetscInt reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui; 1575 const PetscInt *rip,*riip; 1576 PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 1577 PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL,d; 1578 PetscInt ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr; 1579 PetscReal fill=info->fill,levels=info->levels; 1580 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1581 PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 1582 PetscBT lnkbt; 1583 IS iperm; 1584 1585 PetscFunctionBegin; 1586 if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n); 1587 ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 1588 if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 1589 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 1590 ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 1591 1592 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 1593 ui[0] = 0; 1594 1595 /* ICC(0) without matrix ordering: simply copies fill pattern */ 1596 if (!levels && perm_identity) { 1597 1598 for (i=0; i<am; i++) { 1599 ui[i+1] = ui[i] + ai[i+1] - a->diag[i]; 1600 } 1601 ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 1602 cols = uj; 1603 for (i=0; i<am; i++) { 1604 aj = a->j + a->diag[i]; 1605 ncols = ui[i+1] - ui[i]; 1606 for (j=0; j<ncols; j++) *cols++ = *aj++; 1607 } 1608 } else { /* case: levels>0 || (levels=0 && !perm_identity) */ 1609 ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 1610 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 1611 1612 /* initialization */ 1613 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr); 1614 1615 /* jl: linked list for storing indices of the pivot rows 1616 il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 1617 ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 1618 il = jl + am; 1619 uj_ptr = (PetscInt**)(il + am); 1620 uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 1621 for (i=0; i<am; i++){ 1622 jl[i] = am; il[i] = 0; 1623 } 1624 1625 /* create and initialize a linked list for storing column indices of the active row k */ 1626 nlnk = am + 1; 1627 ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 1628 1629 /* initial FreeSpace size is fill*(ai[am]+1) */ 1630 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 1631 current_space = free_space; 1632 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 1633 current_space_lvl = free_space_lvl; 1634 1635 for (k=0; k<am; k++){ /* for each active row k */ 1636 /* initialize lnk by the column indices of row rip[k] of A */ 1637 nzk = 0; 1638 ncols = ai[rip[k]+1] - ai[rip[k]]; 1639 if (!ncols) SETERRQ2(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k); 1640 ncols_upper = 0; 1641 for (j=0; j<ncols; j++){ 1642 i = *(aj + ai[rip[k]] + j); /* unpermuted column index */ 1643 if (riip[i] >= k){ /* only take upper triangular entry */ 1644 ajtmp[ncols_upper] = i; 1645 ncols_upper++; 1646 } 1647 } 1648 ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 1649 nzk += nlnk; 1650 1651 /* update lnk by computing fill-in for each pivot row to be merged in */ 1652 prow = jl[k]; /* 1st pivot row */ 1653 1654 while (prow < k){ 1655 nextprow = jl[prow]; 1656 1657 /* merge prow into k-th row */ 1658 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 1659 jmax = ui[prow+1]; 1660 ncols = jmax-jmin; 1661 i = jmin - ui[prow]; 1662 cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 1663 uj = uj_lvl_ptr[prow] + i; /* levels of cols */ 1664 j = *(uj - 1); 1665 ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr); 1666 nzk += nlnk; 1667 1668 /* update il and jl for prow */ 1669 if (jmin < jmax){ 1670 il[prow] = jmin; 1671 j = *cols; jl[prow] = jl[j]; jl[j] = prow; 1672 } 1673 prow = nextprow; 1674 } 1675 1676 /* if free space is not available, make more free space */ 1677 if (current_space->local_remaining<nzk) { 1678 i = am - k + 1; /* num of unfactored rows */ 1679 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 1680 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 1681 ierr = PetscFreeSpaceGet(i,¤t_space_lvl);CHKERRQ(ierr); 1682 reallocs++; 1683 } 1684 1685 /* copy data into free_space and free_space_lvl, then initialize lnk */ 1686 if (nzk == 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k); 1687 ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 1688 1689 /* add the k-th row into il and jl */ 1690 if (nzk > 1){ 1691 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 1692 jl[k] = jl[i]; jl[i] = k; 1693 il[k] = ui[k] + 1; 1694 } 1695 uj_ptr[k] = current_space->array; 1696 uj_lvl_ptr[k] = current_space_lvl->array; 1697 1698 current_space->array += nzk; 1699 current_space->local_used += nzk; 1700 current_space->local_remaining -= nzk; 1701 1702 current_space_lvl->array += nzk; 1703 current_space_lvl->local_used += nzk; 1704 current_space_lvl->local_remaining -= nzk; 1705 1706 ui[k+1] = ui[k] + nzk; 1707 } 1708 1709 #if defined(PETSC_USE_INFO) 1710 if (ai[am] != 0) { 1711 PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]); 1712 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 1713 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 1714 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 1715 } else { 1716 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 1717 } 1718 #endif 1719 1720 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 1721 ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr); 1722 ierr = PetscFree(jl);CHKERRQ(ierr); 1723 ierr = PetscFree(ajtmp);CHKERRQ(ierr); 1724 1725 /* destroy list of free space and other temporary array(s) */ 1726 ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 1727 ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 1728 ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1729 ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 1730 1731 } /* end of case: levels>0 || (levels=0 && !perm_identity) */ 1732 1733 /* put together the new matrix in MATSEQSBAIJ format */ 1734 1735 b = (Mat_SeqSBAIJ*)(fact)->data; 1736 b->singlemalloc = PETSC_FALSE; 1737 ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 1738 b->j = uj; 1739 b->i = ui; 1740 b->diag = 0; 1741 b->ilen = 0; 1742 b->imax = 0; 1743 b->row = perm; 1744 b->col = perm; 1745 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1746 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1747 b->icol = iperm; 1748 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 1749 ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1750 ierr = PetscLogObjectMemory((fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 1751 b->maxnz = b->nz = ui[am]; 1752 b->free_a = PETSC_TRUE; 1753 b->free_ij = PETSC_TRUE; 1754 1755 (fact)->info.factor_mallocs = reallocs; 1756 (fact)->info.fill_ratio_given = fill; 1757 if (ai[am] != 0) { 1758 (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 1759 } else { 1760 (fact)->info.fill_ratio_needed = 0.0; 1761 } 1762 (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 1763 PetscFunctionReturn(0); 1764 } 1765 1766 #undef __FUNCT__ 1767 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ" 1768 PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 1769 { 1770 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1771 Mat_SeqSBAIJ *b; 1772 PetscErrorCode ierr; 1773 PetscTruth perm_identity; 1774 PetscReal fill = info->fill; 1775 const PetscInt *rip,*riip; 1776 PetscInt i,am=A->rmap->n,*ai=a->i,*aj=a->j,reallocs=0,prow; 1777 PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 1778 PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 1779 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1780 PetscBT lnkbt; 1781 IS iperm; 1782 1783 PetscFunctionBegin; 1784 if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n); 1785 /* check whether perm is the identity mapping */ 1786 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 1787 ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 1788 ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 1789 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 1790 1791 /* initialization */ 1792 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 1793 ui[0] = 0; 1794 1795 /* jl: linked list for storing indices of the pivot rows 1796 il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 1797 ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 1798 il = jl + am; 1799 cols = il + am; 1800 ui_ptr = (PetscInt**)(cols + am); 1801 for (i=0; i<am; i++){ 1802 jl[i] = am; il[i] = 0; 1803 } 1804 1805 /* create and initialize a linked list for storing column indices of the active row k */ 1806 nlnk = am + 1; 1807 ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1808 1809 /* initial FreeSpace size is fill*(ai[am]+1) */ 1810 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 1811 current_space = free_space; 1812 1813 for (k=0; k<am; k++){ /* for each active row k */ 1814 /* initialize lnk by the column indices of row rip[k] of A */ 1815 nzk = 0; 1816 ncols = ai[rip[k]+1] - ai[rip[k]]; 1817 if (!ncols) SETERRQ2(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k); 1818 ncols_upper = 0; 1819 for (j=0; j<ncols; j++){ 1820 i = riip[*(aj + ai[rip[k]] + j)]; 1821 if (i >= k){ /* only take upper triangular entry */ 1822 cols[ncols_upper] = i; 1823 ncols_upper++; 1824 } 1825 } 1826 ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1827 nzk += nlnk; 1828 1829 /* update lnk by computing fill-in for each pivot row to be merged in */ 1830 prow = jl[k]; /* 1st pivot row */ 1831 1832 while (prow < k){ 1833 nextprow = jl[prow]; 1834 /* merge prow into k-th row */ 1835 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 1836 jmax = ui[prow+1]; 1837 ncols = jmax-jmin; 1838 uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 1839 ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1840 nzk += nlnk; 1841 1842 /* update il and jl for prow */ 1843 if (jmin < jmax){ 1844 il[prow] = jmin; 1845 j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 1846 } 1847 prow = nextprow; 1848 } 1849 1850 /* if free space is not available, make more free space */ 1851 if (current_space->local_remaining<nzk) { 1852 i = am - k + 1; /* num of unfactored rows */ 1853 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 1854 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 1855 reallocs++; 1856 } 1857 1858 /* copy data into free space, then initialize lnk */ 1859 ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 1860 1861 /* add the k-th row into il and jl */ 1862 if (nzk-1 > 0){ 1863 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 1864 jl[k] = jl[i]; jl[i] = k; 1865 il[k] = ui[k] + 1; 1866 } 1867 ui_ptr[k] = current_space->array; 1868 current_space->array += nzk; 1869 current_space->local_used += nzk; 1870 current_space->local_remaining -= nzk; 1871 1872 ui[k+1] = ui[k] + nzk; 1873 } 1874 1875 #if defined(PETSC_USE_INFO) 1876 if (ai[am] != 0) { 1877 PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]); 1878 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 1879 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 1880 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 1881 } else { 1882 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 1883 } 1884 #endif 1885 1886 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 1887 ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr); 1888 ierr = PetscFree(jl);CHKERRQ(ierr); 1889 1890 /* destroy list of free space and other temporary array(s) */ 1891 ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 1892 ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 1893 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1894 1895 /* put together the new matrix in MATSEQSBAIJ format */ 1896 1897 b = (Mat_SeqSBAIJ*)(fact)->data; 1898 b->singlemalloc = PETSC_FALSE; 1899 b->free_a = PETSC_TRUE; 1900 b->free_ij = PETSC_TRUE; 1901 ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 1902 b->j = uj; 1903 b->i = ui; 1904 b->diag = 0; 1905 b->ilen = 0; 1906 b->imax = 0; 1907 b->row = perm; 1908 b->col = perm; 1909 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1910 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 1911 b->icol = iperm; 1912 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 1913 ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1914 ierr = PetscLogObjectMemory(fact,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 1915 b->maxnz = b->nz = ui[am]; 1916 1917 (fact)->info.factor_mallocs = reallocs; 1918 (fact)->info.fill_ratio_given = fill; 1919 if (ai[am] != 0) { 1920 (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 1921 } else { 1922 (fact)->info.fill_ratio_needed = 0.0; 1923 } 1924 (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 1925 PetscFunctionReturn(0); 1926 } 1927 1928 #undef __FUNCT__ 1929 #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering_iludt" 1930 PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_iludt(Mat A,Vec bb,Vec xx) 1931 { 1932 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1933 PetscErrorCode ierr; 1934 PetscInt n = A->rmap->n; 1935 const PetscInt *ai = a->i,*aj = a->j,*vi,*adiag_rev=a->diag+n+1; 1936 PetscScalar *x,sum; 1937 const PetscScalar *b; 1938 const MatScalar *aa = a->a,*v; 1939 PetscInt i,nz,k; 1940 1941 PetscFunctionBegin; 1942 if (!n) PetscFunctionReturn(0); 1943 1944 ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 1945 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1946 1947 /* forward solve the lower triangular */ 1948 x[0] = b[0]; 1949 v = aa; 1950 vi = aj; 1951 for (i=1; i<n; i++) { 1952 nz = ai[i+1] - ai[i]; 1953 sum = b[i]; 1954 PetscSparseDenseMinusDot(sum,x,v,vi,nz); 1955 /* while (nz--) sum -= *v++ * x[*vi++];*/ 1956 v += nz; 1957 vi += nz; 1958 x[i] = sum; 1959 } 1960 1961 /* backward solve the upper triangular */ 1962 v = aa + adiag_rev[0] + 1; 1963 vi = aj + adiag_rev[0] + 1; 1964 k = 0; /* i + k = n */ 1965 for (i=n-1; i>=0; i--){ 1966 /* nz = adiag[i] - adiag[i+1] - 1; */ 1967 nz = adiag_rev[k+1] - adiag_rev[k] - 1; 1968 sum = x[i]; 1969 PetscSparseDenseMinusDot(sum,x,v,vi,nz); 1970 /* while (nz--) sum -= *v++ * x[*vi++]; */ 1971 v += nz; 1972 vi += nz; 1973 x[i] = *v++ *sum; /* x[i]=aa[adiag[i]]*sum; v++; */ 1974 vi++; k++; 1975 } 1976 1977 ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr); 1978 ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 1979 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1980 PetscFunctionReturn(0); 1981 } 1982 1983 #undef __FUNCT__ 1984 #define __FUNCT__ "MatSolve_SeqAIJ_iludt" 1985 PetscErrorCode MatSolve_SeqAIJ_iludt(Mat A,Vec bb,Vec xx) 1986 { 1987 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 1988 IS iscol = a->col,isrow = a->row; 1989 PetscErrorCode ierr; 1990 PetscInt i,n=A->rmap->n,*vi,*ai = a->i,*aj = a->j,*adiag=a->diag; 1991 PetscInt nz; 1992 const PetscInt *rout,*cout,*r,*c; 1993 PetscScalar *x,*tmp,*tmps; 1994 const PetscScalar *b; 1995 const MatScalar *aa = a->a,*v; 1996 1997 PetscFunctionBegin; 1998 if (!n) PetscFunctionReturn(0); 1999 2000 ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 2001 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2002 tmp = a->solve_work; 2003 2004 ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 2005 ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 2006 2007 /* forward solve the lower triangular */ 2008 tmp[0] = b[*r++]; 2009 tmps = tmp; 2010 v = aa; 2011 vi = aj; 2012 for (i=1; i<n; i++) { 2013 nz = ai[i+1] - ai[i]; 2014 tmp[i] = b[*r++]; 2015 PetscSparseDenseMinusDot(tmp[i],tmps,v,vi,nz); 2016 v += nz; vi += nz; 2017 } 2018 2019 /* backward solve the upper triangular */ 2020 v = aa + adiag[n] + 1; 2021 vi = aj + adiag[n] + 1; 2022 for (i=n-1; i>=0; i--){ 2023 nz = adiag[i] - adiag[i+1] - 1; 2024 PetscSparseDenseMinusDot(tmp[i],tmps,v,vi,nz); 2025 x[*c--] = tmp[i] = tmp[i]*aa[adiag[i]]; 2026 v += nz+1; vi += nz+1; 2027 } 2028 2029 ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 2030 ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 2031 ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 2032 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2033 ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr); 2034 PetscFunctionReturn(0); 2035 } 2036 2037 #undef __FUNCT__ 2038 #define __FUNCT__ "MatILUDTFactor_SeqAIJ" 2039 PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact) 2040 { 2041 Mat B = *fact; 2042 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b; 2043 IS isicol; 2044 PetscErrorCode ierr; 2045 const PetscInt *r,*ic; 2046 PetscInt i,n=A->rmap->n,*ai=a->i,*aj=a->j,*ajtmp,*adiag; 2047 PetscInt *bi,*bj,*bdiag,*bdiag_rev; 2048 PetscInt row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au; 2049 PetscInt nlnk,*lnk; 2050 PetscBT lnkbt; 2051 PetscTruth row_identity,icol_identity,both_identity; 2052 MatScalar *aatmp,*pv,*batmp,*ba,*rtmp,*pc,multiplier,*vtmp,diag_tmp; 2053 const PetscInt *ics; 2054 PetscInt j,nz,*pj,*bjtmp,k,ncut,*jtmp; 2055 PetscReal dt=info->dt,shift=info->shiftinblocks; 2056 PetscInt nnz_max; 2057 PetscTruth missing; 2058 2059 PetscFunctionBegin; 2060 printf("MatILUDTFactor_SeqAIJ is callled ...\n"); 2061 /* ------- symbolic factorization, can be reused ---------*/ 2062 ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr); 2063 if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i); 2064 adiag=a->diag; 2065 2066 ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 2067 2068 /* bdiag is location of diagonal in factor */ 2069 ierr = PetscMalloc((2*n+2)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 2070 bdiag_rev = bdiag + n+1; 2071 2072 /* allocate row pointers bi */ 2073 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 2074 2075 /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */ 2076 dtcount = (PetscInt)info->dtcount; 2077 if (dtcount > n-1) dtcount = n-1; /* diagonal is excluded */ 2078 nnz_max = ai[n]+2*n*dtcount+2; 2079 /* if (nnz_max > n*n) nnz_max = n*n+1; n*n may cause overflow ?? */ 2080 ierr = PetscMalloc(nnz_max*sizeof(PetscInt),&bj);CHKERRQ(ierr); 2081 ierr = PetscMalloc(nnz_max*sizeof(MatScalar),&ba);CHKERRQ(ierr); 2082 2083 /* put together the new matrix */ 2084 ierr = MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 2085 ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr); 2086 b = (Mat_SeqAIJ*)(B)->data; 2087 b->free_a = PETSC_TRUE; 2088 b->free_ij = PETSC_TRUE; 2089 b->singlemalloc = PETSC_FALSE; 2090 b->a = ba; 2091 b->j = bj; 2092 b->i = bi; 2093 b->diag = bdiag; 2094 b->ilen = 0; 2095 b->imax = 0; 2096 b->row = isrow; 2097 b->col = iscol; 2098 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 2099 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 2100 b->icol = isicol; 2101 ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 2102 2103 ierr = PetscLogObjectMemory(B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 2104 b->maxnz = nnz_max; 2105 2106 (B)->factor = MAT_FACTOR_ILUDT; 2107 (B)->info.factor_mallocs = 0; 2108 (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)ai[n]); 2109 CHKMEMQ; 2110 /* ------- end of symbolic factorization ---------*/ 2111 2112 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 2113 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 2114 ics = ic; 2115 2116 /* linked list for storing column indices of the active row */ 2117 nlnk = n + 1; 2118 ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 2119 2120 /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */ 2121 ierr = PetscMalloc((2*n+1)*sizeof(PetscInt),&im);CHKERRQ(ierr); 2122 jtmp = im + n; 2123 /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */ 2124 ierr = PetscMalloc((2*n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 2125 ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr); 2126 vtmp = rtmp + n; 2127 2128 bi[0] = 0; 2129 bdiag[0] = nnz_max-1; /* location of diagonal in factor B */ 2130 bdiag_rev[n] = bdiag[0]; 2131 for (i=0; i<n; i++) { 2132 /* copy initial fill into linked list */ 2133 nzi = 0; /* nonzeros for active row i */ 2134 nzi = ai[r[i]+1] - ai[r[i]]; 2135 if (!nzi) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 2136 nzi_al = adiag[r[i]] - ai[r[i]]; 2137 nzi_au = ai[r[i]+1] - adiag[r[i]] -1; 2138 ajtmp = aj + ai[r[i]]; 2139 ierr = PetscLLAddPerm(nzi,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 2140 2141 /* load in initial (unfactored row) */ 2142 aatmp = a->a + ai[r[i]]; 2143 for (j=0; j<nzi; j++) { 2144 rtmp[ics[*ajtmp++]] = *aatmp++; 2145 } 2146 2147 /* add pivot rows into linked list */ 2148 row = lnk[n]; 2149 while (row < i ) { 2150 nzi_bl = bi[row+1] - bi[row] + 1; 2151 bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */ 2152 ierr = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr); 2153 nzi += nlnk; 2154 row = lnk[row]; 2155 } 2156 2157 /* copy data from lnk into jtmp, then initialize lnk */ 2158 ierr = PetscLLClean(n,n,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr); 2159 2160 /* numerical factorization */ 2161 bjtmp = jtmp; 2162 row = *bjtmp++; /* 1st pivot row */ 2163 while ( row < i ) { 2164 pc = rtmp + row; 2165 pv = ba + bdiag[row]; /* 1./(diag of the pivot row) */ 2166 multiplier = (*pc) * (*pv); 2167 *pc = multiplier; 2168 if (PetscAbsScalar(*pc) > dt){ /* apply tolerance dropping rule */ 2169 pj = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */ 2170 pv = ba + bdiag[row+1] + 1; 2171 /* if (multiplier < -1.0 or multiplier >1.0) printf("row/prow %d, %d, multiplier %g\n",i,row,multiplier); */ 2172 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */ 2173 for (j=0; j<nz; j++) rtmp[*pj++] -= multiplier * (*pv++); 2174 ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); 2175 } 2176 row = *bjtmp++; 2177 } 2178 2179 /* copy sparse rtmp into contiguous vtmp; separate L and U part */ 2180 diag_tmp = rtmp[i]; /* save diagonal value - may not needed?? */ 2181 nzi_bl = 0; j = 0; 2182 while (jtmp[j] < i){ /* Note: jtmp is sorted */ 2183 vtmp[j] = rtmp[jtmp[j]]; rtmp[jtmp[j]]=0.0; 2184 nzi_bl++; j++; 2185 } 2186 nzi_bu = nzi - nzi_bl -1; 2187 while (j < nzi){ 2188 vtmp[j] = rtmp[jtmp[j]]; rtmp[jtmp[j]]=0.0; 2189 j++; 2190 } 2191 2192 bjtmp = bj + bi[i]; 2193 batmp = ba + bi[i]; 2194 /* apply level dropping rule to L part */ 2195 ncut = nzi_al + dtcount; 2196 if (ncut < nzi_bl){ 2197 ierr = PetscSortSplit(ncut,nzi_bl,vtmp,jtmp);CHKERRQ(ierr); 2198 ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr); 2199 } else { 2200 ncut = nzi_bl; 2201 } 2202 for (j=0; j<ncut; j++){ 2203 bjtmp[j] = jtmp[j]; 2204 batmp[j] = vtmp[j]; 2205 /* printf(" (%d,%g),",bjtmp[j],batmp[j]); */ 2206 } 2207 bi[i+1] = bi[i] + ncut; 2208 nzi = ncut + 1; 2209 2210 /* apply level dropping rule to U part */ 2211 ncut = nzi_au + dtcount; 2212 if (ncut < nzi_bu){ 2213 ierr = PetscSortSplit(ncut,nzi_bu,vtmp+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr); 2214 ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr); 2215 } else { 2216 ncut = nzi_bu; 2217 } 2218 nzi += ncut; 2219 2220 /* mark bdiagonal */ 2221 bdiag[i+1] = bdiag[i] - (ncut + 1); 2222 bdiag_rev[n-i-1] = bdiag[i+1]; 2223 bjtmp = bj + bdiag[i]; 2224 batmp = ba + bdiag[i]; 2225 *bjtmp = i; 2226 *batmp = diag_tmp; /* rtmp[i]; */ 2227 if (*batmp == 0.0) { 2228 *batmp = dt+shift; 2229 printf(" row %d add shift %g\n",i,shift); 2230 } 2231 *batmp = 1.0/(*batmp); /* invert diagonal entries for simplier triangular solves */ 2232 /* printf(" (%d,%g),",*bjtmp,*batmp); */ 2233 2234 bjtmp = bj + bdiag[i+1]+1; 2235 batmp = ba + bdiag[i+1]+1; 2236 for (k=0; k<ncut; k++){ 2237 bjtmp[k] = jtmp[nzi_bl+1+k]; 2238 batmp[k] = vtmp[nzi_bl+1+k]; 2239 /* printf(" (%d,%g),",bjtmp[k],batmp[k]); */ 2240 } 2241 /* printf("\n"); */ 2242 2243 im[i] = nzi; /* used by PetscLLAddSortedLU() */ 2244 /* 2245 printf("row %d: bi %d, bdiag %d\n",i,bi[i],bdiag[i]); 2246 printf(" ----------------------------\n"); 2247 */ 2248 } /* for (i=0; i<n; i++) */ 2249 /* printf("end of L %d, beginning of U %d\n",bi[n],bdiag[n]); */ 2250 if (bi[n] >= bdiag[n]) SETERRQ2(PETSC_ERR_ARG_SIZ,"end of L array %d cannot >= the beginning of U array %d",bi[n],bdiag[n]); 2251 2252 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 2253 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 2254 2255 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 2256 ierr = PetscFree(im);CHKERRQ(ierr); 2257 ierr = PetscFree(rtmp);CHKERRQ(ierr); 2258 2259 ierr = PetscLogFlops(B->cmap->n);CHKERRQ(ierr); 2260 b->maxnz = b->nz = bi[n] + bdiag[0] - bdiag[n]; 2261 2262 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 2263 ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr); 2264 both_identity = (PetscTruth) (row_identity && icol_identity); 2265 if (row_identity && icol_identity) { 2266 B->ops->solve = MatSolve_SeqAIJ_NaturalOrdering_iludt; 2267 } else { 2268 B->ops->solve = MatSolve_SeqAIJ_iludt; 2269 } 2270 2271 B->ops->lufactorsymbolic = MatILUDTFactorSymbolic_SeqAIJ; 2272 B->ops->lufactornumeric = MatILUDTFactorNumeric_SeqAIJ; 2273 B->ops->solveadd = 0; 2274 B->ops->solvetranspose = 0; 2275 B->ops->solvetransposeadd = 0; 2276 B->ops->matsolve = 0; 2277 B->assembled = PETSC_TRUE; 2278 B->preallocated = PETSC_TRUE; 2279 PetscFunctionReturn(0); 2280 } 2281 2282 /* a wraper of MatILUDTFactor_SeqAIJ() */ 2283 #undef __FUNCT__ 2284 #define __FUNCT__ "MatILUDTFactorSymbolic_SeqAIJ" 2285 PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info) 2286 { 2287 PetscErrorCode ierr; 2288 2289 PetscFunctionBegin; 2290 printf("MatILUDTFactorSymbolic_SeqAIJ is called...\n"); 2291 ierr = MatILUDTFactor_SeqAIJ(A,row,col,info,&fact);CHKERRQ(ierr); 2292 2293 fact->ops->lufactornumeric = MatILUDTFactorNumeric_SeqAIJ; 2294 PetscFunctionReturn(0); 2295 } 2296 2297 /* 2298 same as MatLUFactorNumeric_SeqAIJ(), except using contiguous array matrix factors 2299 - intend to replace existing MatLUFactorNumeric_SeqAIJ() 2300 */ 2301 #undef __FUNCT__ 2302 #define __FUNCT__ "MatILUDTFactorNumeric_SeqAIJ" 2303 PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorNumeric_SeqAIJ(Mat fact,Mat A,const MatFactorInfo *info) 2304 { 2305 Mat C=fact; 2306 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data; 2307 IS isrow = b->row,isicol = b->icol; 2308 PetscErrorCode ierr; 2309 const PetscInt *r,*ic,*ics; 2310 PetscInt i,j,k,n=A->rmap->n,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 2311 PetscInt *ajtmp,*bjtmp,nz,nzl,nzu,row,*bdiag = b->diag,*pj; 2312 MatScalar *rtmp,*pc,multiplier,*v,*pv,*aa=a->a; 2313 PetscReal dt=info->dt,shift=info->shiftinblocks; 2314 PetscTruth row_identity, col_identity; 2315 2316 PetscFunctionBegin; 2317 printf("MatILUDTFactorNumeric_SeqAIJ is called...\n"); 2318 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 2319 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 2320 ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 2321 ics = ic; 2322 2323 for (i=0; i<n; i++){ 2324 /* initialize rtmp array */ 2325 nzl = bi[i+1] - bi[i]; /* num of nozeros in L(i,:) */ 2326 bjtmp = bj + bi[i]; 2327 for (j=0; j<nzl; j++) rtmp[*bjtmp++] = 0.0; 2328 rtmp[i] = 0.0; 2329 nzu = bdiag[i] - bdiag[i+1]; /* num of nozeros in U(i,:) */ 2330 bjtmp = bj + bdiag[i+1] + 1; 2331 for (j=0; j<nzu; j++) rtmp[*bjtmp++] = 0.0; 2332 2333 /* load in initial unfactored row of A */ 2334 /* printf("row %d\n",i); */ 2335 nz = ai[r[i]+1] - ai[r[i]]; 2336 ajtmp = aj + ai[r[i]]; 2337 v = aa + ai[r[i]]; 2338 for (j=0; j<nz; j++) { 2339 rtmp[ics[*ajtmp++]] = v[j]; 2340 /* printf(" (%d,%g),",ics[ajtmp[j]],rtmp[ics[ajtmp[j]]]); */ 2341 } 2342 /* printf("\n"); */ 2343 2344 /* numerical factorization */ 2345 bjtmp = bj + bi[i]; /* point to 1st entry of L(i,:) */ 2346 nzl = bi[i+1] - bi[i]; /* num of entries in L(i,:) */ 2347 k = 0; 2348 while (k < nzl){ 2349 row = *bjtmp++; 2350 /* printf(" prow %d\n",row); */ 2351 pc = rtmp + row; 2352 pv = b->a + bdiag[row]; /* 1./(diag of the pivot row) */ 2353 multiplier = (*pc) * (*pv); 2354 *pc = multiplier; 2355 if (PetscAbsScalar(multiplier) > dt){ 2356 pj = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */ 2357 pv = b->a + bdiag[row+1] + 1; 2358 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */ 2359 for (j=0; j<nz; j++) rtmp[*pj++] -= multiplier * (*pv++); 2360 /* ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); */ 2361 } 2362 k++; 2363 } 2364 2365 /* finished row so stick it into b->a */ 2366 /* L-part */ 2367 pv = b->a + bi[i] ; 2368 pj = bj + bi[i] ; 2369 nzl = bi[i+1] - bi[i]; 2370 for (j=0; j<nzl; j++) { 2371 pv[j] = rtmp[pj[j]]; 2372 /* printf(" (%d,%g),",pj[j],pv[j]); */ 2373 } 2374 2375 /* diagonal: invert diagonal entries for simplier triangular solves */ 2376 if (rtmp[i] == 0.0) rtmp[i] = dt+shift; 2377 b->a[bdiag[i]] = 1.0/rtmp[i]; 2378 /* printf(" (%d,%g),",i,b->a[bdiag[i]]); */ 2379 2380 /* U-part */ 2381 pv = b->a + bdiag[i+1] + 1; 2382 pj = bj + bdiag[i+1] + 1; 2383 nzu = bdiag[i] - bdiag[i+1] - 1; 2384 for (j=0; j<nzu; j++) { 2385 pv[j] = rtmp[pj[j]]; 2386 /* printf(" (%d,%g),",pj[j],pv[j]); */ 2387 } 2388 /* printf("\n"); */ 2389 } 2390 2391 ierr = PetscFree(rtmp);CHKERRQ(ierr); 2392 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 2393 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 2394 2395 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 2396 ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr); 2397 if (row_identity && col_identity) { 2398 C->ops->solve = MatSolve_SeqAIJ_NaturalOrdering_iludt; 2399 } else { 2400 C->ops->solve = MatSolve_SeqAIJ_iludt; 2401 } 2402 C->ops->solveadd = 0; 2403 C->ops->solvetranspose = 0; 2404 C->ops->solvetransposeadd = 0; 2405 C->ops->matsolve = 0; 2406 C->assembled = PETSC_TRUE; 2407 C->preallocated = PETSC_TRUE; 2408 ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr); 2409 PetscFunctionReturn(0); 2410 } 2411