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