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