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