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