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