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