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