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 PetscTruth missing; 1392 IS isicol; 1393 1394 PetscFunctionBegin; 1395 /* printf("MatILUFactorSymbolic_SeqAIJ_ilu0_newdatastruct ...\n"); */ 1396 if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n); 1397 ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr); 1398 if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i); 1399 ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 1400 1401 ierr = MatDuplicateNoCreate_SeqAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_FALSE);CHKERRQ(ierr); 1402 b = (Mat_SeqAIJ*)(fact)->data; 1403 1404 /* allocate matrix arrays for new data structure */ 1405 ierr = PetscMalloc3(ai[n]+1,PetscScalar,&b->a,ai[n]+1,PetscInt,&b->j,n+1,PetscInt,&b->i);CHKERRQ(ierr); 1406 ierr = PetscLogObjectMemory(fact,ai[n]*(sizeof(PetscScalar)+sizeof(PetscInt))+(n+1)*sizeof(PetscInt));CHKERRQ(ierr); 1407 b->singlemalloc = PETSC_TRUE; 1408 if (!b->diag){ 1409 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&b->diag);CHKERRQ(ierr); 1410 ierr = PetscLogObjectMemory(fact,(n+1)*sizeof(PetscInt));CHKERRQ(ierr); 1411 } 1412 bdiag = b->diag; 1413 1414 if (n > 0) { 1415 ierr = PetscMemzero(b->a,(ai[n])*sizeof(MatScalar));CHKERRQ(ierr); 1416 } 1417 1418 /* set bi and bj with new data structure */ 1419 bi = b->i; 1420 bj = b->j; 1421 1422 /* L part */ 1423 bi[0] = 0; 1424 for (i=0; i<n; i++){ 1425 nz = adiag[i] - ai[i]; 1426 bi[i+1] = bi[i] + nz; 1427 aj = a->j + ai[i]; 1428 for (j=0; j<nz; j++){ 1429 *bj = aj[j]; bj++; 1430 } 1431 } 1432 1433 /* U part */ 1434 bi_temp = bi[n]; 1435 bdiag[n] = bi[n]-1; 1436 for (i=n-1; i>=0; i--){ 1437 nz = ai[i+1] - adiag[i] - 1; 1438 bi_temp = bi_temp + nz + 1; 1439 aj = a->j + adiag[i] + 1; 1440 for (j=0; j<nz; j++){ 1441 *bj = aj[j]; bj++; 1442 } 1443 /* diag[i] */ 1444 *bj = i; bj++; 1445 bdiag[i] = bi_temp - 1; 1446 } 1447 1448 fact->factor = MAT_FACTOR_ILU; 1449 fact->info.factor_mallocs = 0; 1450 fact->info.fill_ratio_given = info->fill; 1451 fact->info.fill_ratio_needed = 1.0; 1452 fact->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_newdatastruct; 1453 1454 b = (Mat_SeqAIJ*)(fact)->data; 1455 b->row = isrow; 1456 b->col = iscol; 1457 b->icol = isicol; 1458 ierr = PetscMalloc((fact->rmap->n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1459 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1460 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 1461 PetscFunctionReturn(0); 1462 } 1463 1464 #undef __FUNCT__ 1465 #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ_newdatastruct" 1466 PetscErrorCode MatILUFactorSymbolic_SeqAIJ_newdatastruct(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 1467 { 1468 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 1469 IS isicol; 1470 PetscErrorCode ierr; 1471 const PetscInt *r,*ic; 1472 PetscInt n=A->rmap->n,*ai=a->i,*aj=a->j; 1473 PetscInt *bi,*cols,nnz,*cols_lvl; 1474 PetscInt *bdiag,prow,fm,nzbd,reallocs=0,dcount=0; 1475 PetscInt i,levels,diagonal_fill; 1476 PetscTruth col_identity,row_identity; 1477 PetscReal f; 1478 PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL; 1479 PetscBT lnkbt; 1480 PetscInt nzi,*bj,**bj_ptr,**bjlvl_ptr; 1481 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1482 PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 1483 1484 PetscFunctionBegin; 1485 /* printf("MatILUFactorSymbolic_SeqAIJ_newdatastruct ...\n"); */ 1486 levels = (PetscInt)info->levels; 1487 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 1488 ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 1489 1490 if (!levels && row_identity && col_identity) { 1491 /* special case: ilu(0) with natural ordering */ 1492 ierr = MatILUFactorSymbolic_SeqAIJ_ilu0_newdatastruct(fact,A,isrow,iscol,info);CHKERRQ(ierr); 1493 PetscFunctionReturn(0); 1494 } 1495 1496 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); 1497 ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 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 f = info->fill; 1517 diagonal_fill = (PetscInt)info->diagonal_fill; 1518 ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 1519 current_space = free_space; 1520 ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr); 1521 current_space_lvl = free_space_lvl; 1522 1523 for (i=0; i<n; i++) { 1524 nzi = 0; 1525 /* copy current row into linked list */ 1526 nnz = ai[r[i]+1] - ai[r[i]]; 1527 if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 1528 cols = aj + ai[r[i]]; 1529 lnk[i] = -1; /* marker to indicate if diagonal exists */ 1530 ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 1531 nzi += nlnk; 1532 1533 /* make sure diagonal entry is included */ 1534 if (diagonal_fill && lnk[i] == -1) { 1535 fm = n; 1536 while (lnk[fm] < i) fm = lnk[fm]; 1537 lnk[i] = lnk[fm]; /* insert diagonal into linked list */ 1538 lnk[fm] = i; 1539 lnk_lvl[i] = 0; 1540 nzi++; dcount++; 1541 } 1542 1543 /* add pivot rows into the active row */ 1544 nzbd = 0; 1545 prow = lnk[n]; 1546 while (prow < i) { 1547 nnz = bdiag[prow]; 1548 cols = bj_ptr[prow] + nnz + 1; 1549 cols_lvl = bjlvl_ptr[prow] + nnz + 1; 1550 nnz = bi[prow+1] - bi[prow] - nnz - 1; 1551 ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr); 1552 nzi += nlnk; 1553 prow = lnk[prow]; 1554 nzbd++; 1555 } 1556 bdiag[i] = nzbd; 1557 bi[i+1] = bi[i] + nzi; 1558 1559 /* if free space is not available, make more free space */ 1560 if (current_space->local_remaining<nzi) { 1561 nnz = 2*nzi*(n - i); /* estimated and max additional space needed */ 1562 ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 1563 ierr = PetscFreeSpaceGet(nnz,¤t_space_lvl);CHKERRQ(ierr); 1564 reallocs++; 1565 } 1566 1567 /* copy data into free_space and free_space_lvl, then initialize lnk */ 1568 ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 1569 bj_ptr[i] = current_space->array; 1570 bjlvl_ptr[i] = current_space_lvl->array; 1571 1572 /* make sure the active row i has diagonal entry */ 1573 if (*(bj_ptr[i]+bdiag[i]) != i) { 1574 SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\ 1575 try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i); 1576 } 1577 1578 current_space->array += nzi; 1579 current_space->local_used += nzi; 1580 current_space->local_remaining -= nzi; 1581 current_space_lvl->array += nzi; 1582 current_space_lvl->local_used += nzi; 1583 current_space_lvl->local_remaining -= nzi; 1584 } 1585 1586 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 1587 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 1588 1589 /* destroy list of free space and other temporary arrays */ 1590 ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 1591 1592 /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */ 1593 ierr = PetscFreeSpaceContiguous_LU_v2(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr); 1594 1595 ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1596 ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 1597 ierr = PetscFree(bj_ptr);CHKERRQ(ierr); 1598 1599 #if defined(PETSC_USE_INFO) 1600 { 1601 PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 1602 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr); 1603 ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 1604 ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr); 1605 ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 1606 if (diagonal_fill) { 1607 ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr); 1608 } 1609 } 1610 #endif 1611 1612 /* put together the new matrix */ 1613 ierr = MatSeqAIJSetPreallocation_SeqAIJ(fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 1614 ierr = PetscLogObjectParent(fact,isicol);CHKERRQ(ierr); 1615 b = (Mat_SeqAIJ*)(fact)->data; 1616 b->free_a = PETSC_TRUE; 1617 b->free_ij = PETSC_TRUE; 1618 b->singlemalloc = PETSC_FALSE; 1619 ierr = PetscMalloc((bdiag[0]+1)*sizeof(PetscScalar),&b->a);CHKERRQ(ierr); 1620 b->j = bj; 1621 b->i = bi; 1622 b->diag = bdiag; 1623 b->ilen = 0; 1624 b->imax = 0; 1625 b->row = isrow; 1626 b->col = iscol; 1627 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1628 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 1629 b->icol = isicol; 1630 ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1631 /* In b structure: Free imax, ilen, old a, old j. 1632 Allocate bdiag, solve_work, new a, new j */ 1633 ierr = PetscLogObjectMemory(fact,(bdiag[0]+1)*(sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr); 1634 b->maxnz = b->nz = bdiag[0]+1; 1635 (fact)->info.factor_mallocs = reallocs; 1636 (fact)->info.fill_ratio_given = f; 1637 (fact)->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]); 1638 (fact)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ_newdatastruct; 1639 /* ierr = MatILUFactorSymbolic_SeqAIJ_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr); */ 1640 PetscFunctionReturn(0); 1641 } 1642 1643 #undef __FUNCT__ 1644 #define __FUNCT__ "MatILUFactorSymbolic_SeqAIJ" 1645 PetscErrorCode MatILUFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 1646 { 1647 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data,*b; 1648 IS isicol; 1649 PetscErrorCode ierr; 1650 const PetscInt *r,*ic; 1651 PetscInt n=A->rmap->n,*ai=a->i,*aj=a->j,d; 1652 PetscInt *bi,*cols,nnz,*cols_lvl; 1653 PetscInt *bdiag,prow,fm,nzbd,reallocs=0,dcount=0; 1654 PetscInt i,levels,diagonal_fill; 1655 PetscTruth col_identity,row_identity; 1656 PetscReal f; 1657 PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL; 1658 PetscBT lnkbt; 1659 PetscInt nzi,*bj,**bj_ptr,**bjlvl_ptr; 1660 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 1661 PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 1662 PetscTruth missing; 1663 PetscTruth newdatastruct=PETSC_FALSE; 1664 1665 PetscFunctionBegin; 1666 ierr = PetscOptionsGetTruth(PETSC_NULL,"-ilu_new",&newdatastruct,PETSC_NULL);CHKERRQ(ierr); 1667 if(newdatastruct){ 1668 ierr = MatILUFactorSymbolic_SeqAIJ_newdatastruct(fact,A,isrow,iscol,info);CHKERRQ(ierr); 1669 PetscFunctionReturn(0); 1670 } 1671 1672 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); 1673 f = info->fill; 1674 levels = (PetscInt)info->levels; 1675 diagonal_fill = (PetscInt)info->diagonal_fill; 1676 ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 1677 1678 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 1679 ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 1680 if (!levels && row_identity && col_identity) { /* special case: ilu(0) with natural ordering */ 1681 ierr = MatDuplicateNoCreate_SeqAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_TRUE);CHKERRQ(ierr); 1682 (fact)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ; 1683 1684 fact->factor = MAT_FACTOR_ILU; 1685 (fact)->info.factor_mallocs = 0; 1686 (fact)->info.fill_ratio_given = info->fill; 1687 (fact)->info.fill_ratio_needed = 1.0; 1688 b = (Mat_SeqAIJ*)(fact)->data; 1689 ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 1690 if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 1691 b->row = isrow; 1692 b->col = iscol; 1693 b->icol = isicol; 1694 ierr = PetscMalloc(((fact)->rmap->n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1695 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1696 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 1697 ierr = MatILUFactorSymbolic_SeqAIJ_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr); 1698 PetscFunctionReturn(0); 1699 } 1700 1701 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 1702 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 1703 1704 /* get new row pointers */ 1705 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 1706 bi[0] = 0; 1707 /* bdiag is location of diagonal in factor */ 1708 ierr = PetscMalloc((n+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 1709 bdiag[0] = 0; 1710 1711 ierr = PetscMalloc((2*n+1)*sizeof(PetscInt**),&bj_ptr);CHKERRQ(ierr); 1712 bjlvl_ptr = (PetscInt**)(bj_ptr + n); 1713 1714 /* create a linked list for storing column indices of the active row */ 1715 nlnk = n + 1; 1716 ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 1717 1718 /* initial FreeSpace size is f*(ai[n]+1) */ 1719 ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space);CHKERRQ(ierr); 1720 current_space = free_space; 1721 ierr = PetscFreeSpaceGet((PetscInt)(f*(ai[n]+1)),&free_space_lvl);CHKERRQ(ierr); 1722 current_space_lvl = free_space_lvl; 1723 1724 for (i=0; i<n; i++) { 1725 nzi = 0; 1726 /* copy current row into linked list */ 1727 nnz = ai[r[i]+1] - ai[r[i]]; 1728 if (!nnz) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 1729 cols = aj + ai[r[i]]; 1730 lnk[i] = -1; /* marker to indicate if diagonal exists */ 1731 ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 1732 nzi += nlnk; 1733 1734 /* make sure diagonal entry is included */ 1735 if (diagonal_fill && lnk[i] == -1) { 1736 fm = n; 1737 while (lnk[fm] < i) fm = lnk[fm]; 1738 lnk[i] = lnk[fm]; /* insert diagonal into linked list */ 1739 lnk[fm] = i; 1740 lnk_lvl[i] = 0; 1741 nzi++; dcount++; 1742 } 1743 1744 /* add pivot rows into the active row */ 1745 nzbd = 0; 1746 prow = lnk[n]; 1747 while (prow < i) { 1748 nnz = bdiag[prow]; 1749 cols = bj_ptr[prow] + nnz + 1; 1750 cols_lvl = bjlvl_ptr[prow] + nnz + 1; 1751 nnz = bi[prow+1] - bi[prow] - nnz - 1; 1752 ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr); 1753 nzi += nlnk; 1754 prow = lnk[prow]; 1755 nzbd++; 1756 } 1757 bdiag[i] = nzbd; 1758 bi[i+1] = bi[i] + nzi; 1759 1760 /* if free space is not available, make more free space */ 1761 if (current_space->local_remaining<nzi) { 1762 nnz = nzi*(n - i); /* estimated and max additional space needed */ 1763 ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 1764 ierr = PetscFreeSpaceGet(nnz,¤t_space_lvl);CHKERRQ(ierr); 1765 reallocs++; 1766 } 1767 1768 /* copy data into free_space and free_space_lvl, then initialize lnk */ 1769 ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 1770 bj_ptr[i] = current_space->array; 1771 bjlvl_ptr[i] = current_space_lvl->array; 1772 1773 /* make sure the active row i has diagonal entry */ 1774 if (*(bj_ptr[i]+bdiag[i]) != i) { 1775 SETERRQ1(PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\ 1776 try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i); 1777 } 1778 1779 current_space->array += nzi; 1780 current_space->local_used += nzi; 1781 current_space->local_remaining -= nzi; 1782 current_space_lvl->array += nzi; 1783 current_space_lvl->local_used += nzi; 1784 current_space_lvl->local_remaining -= nzi; 1785 } 1786 1787 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 1788 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 1789 1790 /* destroy list of free space and other temporary arrays */ 1791 ierr = PetscMalloc((bi[n]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 1792 ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); /* copy free_space -> bj */ 1793 ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1794 ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 1795 ierr = PetscFree(bj_ptr);CHKERRQ(ierr); 1796 1797 #if defined(PETSC_USE_INFO) 1798 { 1799 PetscReal af = ((PetscReal)bi[n])/((PetscReal)ai[n]); 1800 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr); 1801 ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 1802 ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%G);\n",af);CHKERRQ(ierr); 1803 ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 1804 if (diagonal_fill) { 1805 ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals",dcount);CHKERRQ(ierr); 1806 } 1807 } 1808 #endif 1809 1810 /* put together the new matrix */ 1811 ierr = MatSeqAIJSetPreallocation_SeqAIJ(fact,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 1812 ierr = PetscLogObjectParent(fact,isicol);CHKERRQ(ierr); 1813 b = (Mat_SeqAIJ*)(fact)->data; 1814 b->free_a = PETSC_TRUE; 1815 b->free_ij = PETSC_TRUE; 1816 b->singlemalloc = PETSC_FALSE; 1817 ierr = PetscMalloc( (bi[n] )*sizeof(PetscScalar),&b->a);CHKERRQ(ierr); 1818 b->j = bj; 1819 b->i = bi; 1820 for (i=0; i<n; i++) bdiag[i] += bi[i]; 1821 b->diag = bdiag; 1822 b->ilen = 0; 1823 b->imax = 0; 1824 b->row = isrow; 1825 b->col = iscol; 1826 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1827 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 1828 b->icol = isicol; 1829 ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1830 /* In b structure: Free imax, ilen, old a, old j. 1831 Allocate bdiag, solve_work, new a, new j */ 1832 ierr = PetscLogObjectMemory(fact,(bi[n]-n) * (sizeof(PetscInt)+sizeof(PetscScalar)));CHKERRQ(ierr); 1833 b->maxnz = b->nz = bi[n] ; 1834 (fact)->info.factor_mallocs = reallocs; 1835 (fact)->info.fill_ratio_given = f; 1836 (fact)->info.fill_ratio_needed = ((PetscReal)bi[n])/((PetscReal)ai[n]); 1837 (fact)->ops->lufactornumeric = MatLUFactorNumeric_SeqAIJ; 1838 ierr = MatILUFactorSymbolic_SeqAIJ_Inode(fact,A,isrow,iscol,info);CHKERRQ(ierr); 1839 PetscFunctionReturn(0); 1840 } 1841 1842 #undef __FUNCT__ 1843 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ_newdatastruct" 1844 PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ_newdatastruct(Mat B,Mat A,const MatFactorInfo *info) 1845 { 1846 Mat C = B; 1847 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 1848 Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 1849 IS ip=b->row,iip = b->icol; 1850 PetscErrorCode ierr; 1851 const PetscInt *rip,*riip; 1852 PetscInt i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bdiag=b->diag,*bjtmp; 1853 PetscInt *ai=a->i,*aj=a->j; 1854 PetscInt k,jmin,jmax,*c2r,*il,col,nexti,ili,nz; 1855 MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 1856 PetscReal zeropivot,shiftnz; 1857 PetscReal shiftpd; 1858 PetscTruth perm_identity; 1859 1860 PetscFunctionBegin; 1861 1862 shiftnz = info->shiftnz; 1863 shiftpd = info->shiftpd; 1864 zeropivot = info->zeropivot; 1865 1866 ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 1867 ierr = ISGetIndices(iip,&riip);CHKERRQ(ierr); 1868 1869 /* allocate working arrays 1870 c2r: linked list, keep track of pivot rows for a given column. c2r[col]: head of the list for a given col 1871 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 1872 */ 1873 nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 1874 ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 1875 c2r = il + mbs; 1876 rtmp = (MatScalar*)(c2r + mbs); 1877 1878 for (i=0; i<mbs; i++) { 1879 c2r[i] = mbs; 1880 } 1881 il[0] = 0; 1882 1883 for (k = 0; k<mbs; k++){ 1884 /* zero rtmp */ 1885 nz = bi[k+1] - bi[k]; 1886 bjtmp = bj + bi[k]; 1887 for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0; 1888 1889 /* load in initial unfactored row */ 1890 bval = ba + bi[k]; 1891 jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 1892 for (j = jmin; j < jmax; j++){ 1893 col = riip[aj[j]]; 1894 if (col >= k){ /* only take upper triangular entry */ 1895 rtmp[col] = aa[j]; 1896 *bval++ = 0.0; /* for in-place factorization */ 1897 } 1898 } 1899 /* shift the diagonal of the matrix */ 1900 1901 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 1902 dk = rtmp[k]; 1903 i = c2r[k]; /* first row to be added to k_th row */ 1904 1905 while (i < k){ 1906 nexti = c2r[i]; /* next row to be added to k_th row */ 1907 1908 /* compute multiplier, update diag(k) and U(i,k) */ 1909 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1910 uikdi = - ba[ili]*ba[bdiag[i]]; /* diagonal(k) */ 1911 dk += uikdi*ba[ili]; /* update diag[k] */ 1912 ba[ili] = uikdi; /* -U(i,k) */ 1913 1914 /* add multiple of row i to k-th row */ 1915 jmin = ili + 1; jmax = bi[i+1]; 1916 if (jmin < jmax){ 1917 for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 1918 /* update il and c2r for row i */ 1919 il[i] = jmin; 1920 j = bj[jmin]; c2r[i] = c2r[j]; c2r[j] = i; 1921 } 1922 i = nexti; 1923 } 1924 1925 /* copy data into U(k,:) */ 1926 ba[bdiag[k]] = 1.0/dk; /* U(k,k) */ 1927 jmin = bi[k]; jmax = bi[k+1]-1; 1928 if (jmin < jmax) { 1929 for (j=jmin; j<jmax; j++){ 1930 col = bj[j]; ba[j] = rtmp[col]; 1931 } 1932 /* add the k-th row into il and c2r */ 1933 il[k] = jmin; 1934 i = bj[jmin]; c2r[k] = c2r[i]; c2r[i] = k; 1935 } 1936 } 1937 1938 ierr = PetscFree(il);CHKERRQ(ierr); 1939 ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 1940 ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr); 1941 1942 ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr); 1943 if (perm_identity){ 1944 (B)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_newdatastruct; 1945 (B)->ops->solvetranspose = 0; 1946 (B)->ops->forwardsolve = 0; 1947 (B)->ops->backwardsolve = 0; 1948 } else { 1949 (B)->ops->solve = MatSolve_SeqSBAIJ_1_newdatastruct; 1950 (B)->ops->solvetranspose = 0; 1951 (B)->ops->forwardsolve = 0; 1952 (B)->ops->backwardsolve = 0; 1953 } 1954 1955 C->assembled = PETSC_TRUE; 1956 C->preallocated = PETSC_TRUE; 1957 ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr); 1958 PetscFunctionReturn(0); 1959 } 1960 1961 #undef __FUNCT__ 1962 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqAIJ" 1963 PetscErrorCode MatCholeskyFactorNumeric_SeqAIJ(Mat B,Mat A,const MatFactorInfo *info) 1964 { 1965 Mat C = B; 1966 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data; 1967 Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 1968 IS ip=b->row,iip = b->icol; 1969 PetscErrorCode ierr; 1970 const PetscInt *rip,*riip; 1971 PetscInt i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bcol; 1972 PetscInt *ai=a->i,*aj=a->j; 1973 PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 1974 MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 1975 PetscReal zeropivot,rs,shiftnz; 1976 PetscReal shiftpd; 1977 ChShift_Ctx sctx; 1978 PetscInt newshift; 1979 PetscTruth perm_identity; 1980 1981 PetscFunctionBegin; 1982 shiftnz = info->shiftnz; 1983 shiftpd = info->shiftpd; 1984 zeropivot = info->zeropivot; 1985 1986 ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 1987 ierr = ISGetIndices(iip,&riip);CHKERRQ(ierr); 1988 1989 /* initialization */ 1990 nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 1991 ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 1992 jl = il + mbs; 1993 rtmp = (MatScalar*)(jl + mbs); 1994 1995 sctx.shift_amount = 0; 1996 sctx.nshift = 0; 1997 do { 1998 sctx.chshift = PETSC_FALSE; 1999 for (i=0; i<mbs; i++) { 2000 rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 2001 } 2002 2003 for (k = 0; k<mbs; k++){ 2004 bval = ba + bi[k]; 2005 /* initialize k-th row by the perm[k]-th row of A */ 2006 jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 2007 for (j = jmin; j < jmax; j++){ 2008 col = riip[aj[j]]; 2009 if (col >= k){ /* only take upper triangular entry */ 2010 rtmp[col] = aa[j]; 2011 *bval++ = 0.0; /* for in-place factorization */ 2012 } 2013 } 2014 /* shift the diagonal of the matrix */ 2015 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 2016 2017 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 2018 dk = rtmp[k]; 2019 i = jl[k]; /* first row to be added to k_th row */ 2020 2021 while (i < k){ 2022 nexti = jl[i]; /* next row to be added to k_th row */ 2023 2024 /* compute multiplier, update diag(k) and U(i,k) */ 2025 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 2026 uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 2027 dk += uikdi*ba[ili]; 2028 ba[ili] = uikdi; /* -U(i,k) */ 2029 2030 /* add multiple of row i to k-th row */ 2031 jmin = ili + 1; jmax = bi[i+1]; 2032 if (jmin < jmax){ 2033 for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 2034 /* update il and jl for row i */ 2035 il[i] = jmin; 2036 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 2037 } 2038 i = nexti; 2039 } 2040 2041 /* shift the diagonals when zero pivot is detected */ 2042 /* compute rs=sum of abs(off-diagonal) */ 2043 rs = 0.0; 2044 jmin = bi[k]+1; 2045 nz = bi[k+1] - jmin; 2046 bcol = bj + jmin; 2047 for (j=0; j<nz; j++) { 2048 rs += PetscAbsScalar(rtmp[bcol[j]]); 2049 } 2050 2051 sctx.rs = rs; 2052 sctx.pv = dk; 2053 ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 2054 2055 if (newshift == 1) { 2056 if (!sctx.shift_amount) { 2057 sctx.shift_amount = 1e-5; 2058 } 2059 break; 2060 } 2061 2062 /* copy data into U(k,:) */ 2063 ba[bi[k]] = 1.0/dk; /* U(k,k) */ 2064 jmin = bi[k]+1; jmax = bi[k+1]; 2065 if (jmin < jmax) { 2066 for (j=jmin; j<jmax; j++){ 2067 col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 2068 } 2069 /* add the k-th row into il and jl */ 2070 il[k] = jmin; 2071 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 2072 } 2073 } 2074 } while (sctx.chshift); 2075 ierr = PetscFree(il);CHKERRQ(ierr); 2076 2077 ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 2078 ierr = ISRestoreIndices(iip,&riip);CHKERRQ(ierr); 2079 2080 ierr = ISIdentity(ip,&perm_identity);CHKERRQ(ierr); 2081 if (perm_identity){ 2082 (B)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 2083 (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 2084 (B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering; 2085 (B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering; 2086 } else { 2087 (B)->ops->solve = MatSolve_SeqSBAIJ_1; 2088 (B)->ops->solvetranspose = MatSolve_SeqSBAIJ_1; 2089 (B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1; 2090 (B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1; 2091 } 2092 2093 C->assembled = PETSC_TRUE; 2094 C->preallocated = PETSC_TRUE; 2095 ierr = PetscLogFlops(C->rmap->n);CHKERRQ(ierr); 2096 if (sctx.nshift){ 2097 if (shiftnz) { 2098 ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 2099 } else if (shiftpd) { 2100 ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 2101 } 2102 } 2103 PetscFunctionReturn(0); 2104 } 2105 2106 /* 2107 icc() under revised new data structure. 2108 Factored arrays bj and ba are stored as 2109 U(0,:),...,U(i,:),U(n-1,:) 2110 2111 ui=fact->i is an array of size n+1, in which 2112 ui+ 2113 ui[i]: points to 1st entry of U(i,:),i=0,...,n-1 2114 ui[n]: points to U(n-1,n-1)+1 2115 2116 udiag=fact->diag is an array of size n,in which 2117 udiag[i]: points to diagonal of U(i,:), i=0,...,n-1 2118 2119 U(i,:) contains udiag[i] as its last entry, i.e., 2120 U(i,:) = (u[i,i+1],...,u[i,n-1],diag[i]) 2121 */ 2122 2123 #undef __FUNCT__ 2124 #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ_newdatastruct" 2125 PetscErrorCode MatICCFactorSymbolic_SeqAIJ_newdatastruct(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 2126 { 2127 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 2128 Mat_SeqSBAIJ *b; 2129 PetscErrorCode ierr; 2130 PetscTruth perm_identity,missing; 2131 PetscInt reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui,*udiag; 2132 const PetscInt *rip,*riip; 2133 PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 2134 PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL,d; 2135 PetscInt ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr; 2136 PetscReal fill=info->fill,levels=info->levels; 2137 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 2138 PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 2139 PetscBT lnkbt; 2140 IS iperm; 2141 2142 PetscFunctionBegin; 2143 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); 2144 ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 2145 if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 2146 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 2147 ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 2148 2149 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 2150 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&udiag);CHKERRQ(ierr); 2151 ui[0] = 0; 2152 2153 /* ICC(0) without matrix ordering: simply rearrange column indices */ 2154 if (!levels && perm_identity) { 2155 2156 for (i=0; i<am; i++) { 2157 ui[i+1] = ui[i] + ai[i+1] - a->diag[i]; 2158 udiag[i] = ui[i+1] - 1; /* points to the last entry of U(i,:) */ 2159 } 2160 ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 2161 cols = uj; 2162 for (i=0; i<am; i++) { 2163 aj = a->j + a->diag[i] + 1; /* 1st entry of U(i,:) without diagonal */ 2164 ncols = ui[i+1] - ui[i] - 1; 2165 for (j=0; j<ncols; j++) *cols++ = aj[j]; 2166 *cols++ = i; /* diagoanl is located as the last entry of U(i,:) */ 2167 } 2168 } else { /* case: levels>0 || (levels=0 && !perm_identity) */ 2169 if (levels) SETERRQ(PETSC_ERR_ARG_WRONG,"not done yet"); 2170 ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 2171 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 2172 2173 /* initialization */ 2174 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr); 2175 2176 /* jl: linked list for storing indices of the pivot rows 2177 il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 2178 ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 2179 il = jl + am; 2180 uj_ptr = (PetscInt**)(il + am); 2181 uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 2182 for (i=0; i<am; i++){ 2183 jl[i] = am; il[i] = 0; 2184 } 2185 2186 /* create and initialize a linked list for storing column indices of the active row k */ 2187 nlnk = am + 1; 2188 ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 2189 2190 /* initial FreeSpace size is fill*(ai[am]+1) */ 2191 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 2192 current_space = free_space; 2193 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 2194 current_space_lvl = free_space_lvl; 2195 2196 for (k=0; k<am; k++){ /* for each active row k */ 2197 /* initialize lnk by the column indices of row rip[k] of A */ 2198 nzk = 0; 2199 ncols = ai[rip[k]+1] - ai[rip[k]]; 2200 if (!ncols) SETERRQ2(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k); 2201 ncols_upper = 0; 2202 for (j=0; j<ncols; j++){ 2203 i = *(aj + ai[rip[k]] + j); /* unpermuted column index */ 2204 if (riip[i] >= k){ /* only take upper triangular entry */ 2205 ajtmp[ncols_upper] = i; 2206 ncols_upper++; 2207 } 2208 } 2209 ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 2210 nzk += nlnk; 2211 2212 /* update lnk by computing fill-in for each pivot row to be merged in */ 2213 prow = jl[k]; /* 1st pivot row */ 2214 2215 while (prow < k){ 2216 nextprow = jl[prow]; 2217 2218 /* merge prow into k-th row */ 2219 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 2220 jmax = ui[prow+1]; 2221 ncols = jmax-jmin; 2222 i = jmin - ui[prow]; 2223 cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 2224 uj = uj_lvl_ptr[prow] + i; /* levels of cols */ 2225 j = *(uj - 1); 2226 ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr); 2227 nzk += nlnk; 2228 2229 /* update il and jl for prow */ 2230 if (jmin < jmax){ 2231 il[prow] = jmin; 2232 j = *cols; jl[prow] = jl[j]; jl[j] = prow; 2233 } 2234 prow = nextprow; 2235 } 2236 2237 /* if free space is not available, make more free space */ 2238 if (current_space->local_remaining<nzk) { 2239 i = am - k + 1; /* num of unfactored rows */ 2240 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 2241 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 2242 ierr = PetscFreeSpaceGet(i,¤t_space_lvl);CHKERRQ(ierr); 2243 reallocs++; 2244 } 2245 2246 /* copy data into free_space and free_space_lvl, then initialize lnk */ 2247 if (nzk == 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k); 2248 ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 2249 2250 /* add the k-th row into il and jl */ 2251 if (nzk > 1){ 2252 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 2253 jl[k] = jl[i]; jl[i] = k; 2254 il[k] = ui[k] + 1; 2255 } 2256 uj_ptr[k] = current_space->array; 2257 uj_lvl_ptr[k] = current_space_lvl->array; 2258 2259 current_space->array += nzk; 2260 current_space->local_used += nzk; 2261 current_space->local_remaining -= nzk; 2262 2263 current_space_lvl->array += nzk; 2264 current_space_lvl->local_used += nzk; 2265 current_space_lvl->local_remaining -= nzk; 2266 2267 ui[k+1] = ui[k] + nzk; 2268 } 2269 2270 #if defined(PETSC_USE_INFO) 2271 if (ai[am] != 0) { 2272 PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]); 2273 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 2274 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 2275 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 2276 } else { 2277 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 2278 } 2279 #endif 2280 2281 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 2282 ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr); 2283 ierr = PetscFree(jl);CHKERRQ(ierr); 2284 ierr = PetscFree(ajtmp);CHKERRQ(ierr); 2285 2286 /* destroy list of free space and other temporary array(s) */ 2287 ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 2288 ierr = PetscFreeSpaceContiguous_Cholesky(&free_space,uj,am,ui,udiag);CHKERRQ(ierr); 2289 ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 2290 ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 2291 2292 } /* end of case: levels>0 || (levels=0 && !perm_identity) */ 2293 2294 /* put together the new matrix in MATSEQSBAIJ format */ 2295 2296 b = (Mat_SeqSBAIJ*)(fact)->data; 2297 b->singlemalloc = PETSC_FALSE; 2298 ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 2299 b->j = uj; 2300 b->i = ui; 2301 b->diag = udiag; 2302 b->free_diag = PETSC_TRUE; 2303 b->ilen = 0; 2304 b->imax = 0; 2305 b->row = perm; 2306 b->col = perm; 2307 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 2308 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 2309 b->icol = iperm; 2310 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 2311 ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 2312 ierr = PetscLogObjectMemory((fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 2313 b->maxnz = b->nz = ui[am]; 2314 b->free_a = PETSC_TRUE; 2315 b->free_ij = PETSC_TRUE; 2316 2317 (fact)->info.factor_mallocs = reallocs; 2318 (fact)->info.fill_ratio_given = fill; 2319 if (ai[am] != 0) { 2320 (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 2321 } else { 2322 (fact)->info.fill_ratio_needed = 0.0; 2323 } 2324 (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ_newdatastruct; 2325 PetscFunctionReturn(0); 2326 } 2327 2328 #undef __FUNCT__ 2329 #define __FUNCT__ "MatICCFactorSymbolic_SeqAIJ" 2330 PetscErrorCode MatICCFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 2331 { 2332 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 2333 Mat_SeqSBAIJ *b; 2334 PetscErrorCode ierr; 2335 PetscTruth perm_identity,missing; 2336 PetscInt reallocs=0,i,*ai=a->i,*aj=a->j,am=A->rmap->n,*ui,*udiag; 2337 const PetscInt *rip,*riip; 2338 PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 2339 PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL,d; 2340 PetscInt ncols,ncols_upper,*cols,*ajtmp,*uj,**uj_ptr,**uj_lvl_ptr; 2341 PetscReal fill=info->fill,levels=info->levels; 2342 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 2343 PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 2344 PetscBT lnkbt; 2345 IS iperm; 2346 PetscTruth newdatastruct=PETSC_FALSE; 2347 2348 PetscFunctionBegin; 2349 ierr = PetscOptionsGetTruth(PETSC_NULL,"-icc_new",&newdatastruct,PETSC_NULL);CHKERRQ(ierr); 2350 if(newdatastruct){ 2351 ierr = MatICCFactorSymbolic_SeqAIJ_newdatastruct(fact,A,perm,info);CHKERRQ(ierr); 2352 PetscFunctionReturn(0); 2353 } 2354 2355 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); 2356 ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 2357 if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 2358 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 2359 ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 2360 2361 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 2362 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&udiag);CHKERRQ(ierr); 2363 ui[0] = 0; 2364 2365 /* ICC(0) without matrix ordering: simply copies fill pattern */ 2366 if (!levels && perm_identity) { 2367 2368 for (i=0; i<am; i++) { 2369 ui[i+1] = ui[i] + ai[i+1] - a->diag[i]; 2370 udiag[i] = ui[i]; 2371 } 2372 ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 2373 cols = uj; 2374 for (i=0; i<am; i++) { 2375 aj = a->j + a->diag[i]; 2376 ncols = ui[i+1] - ui[i]; 2377 for (j=0; j<ncols; j++) *cols++ = *aj++; 2378 } 2379 } else { /* case: levels>0 || (levels=0 && !perm_identity) */ 2380 ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 2381 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 2382 2383 /* initialization */ 2384 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ajtmp);CHKERRQ(ierr); 2385 2386 /* jl: linked list for storing indices of the pivot rows 2387 il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 2388 ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 2389 il = jl + am; 2390 uj_ptr = (PetscInt**)(il + am); 2391 uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 2392 for (i=0; i<am; i++){ 2393 jl[i] = am; il[i] = 0; 2394 } 2395 2396 /* create and initialize a linked list for storing column indices of the active row k */ 2397 nlnk = am + 1; 2398 ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 2399 2400 /* initial FreeSpace size is fill*(ai[am]+1) */ 2401 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 2402 current_space = free_space; 2403 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 2404 current_space_lvl = free_space_lvl; 2405 2406 for (k=0; k<am; k++){ /* for each active row k */ 2407 /* initialize lnk by the column indices of row rip[k] of A */ 2408 nzk = 0; 2409 ncols = ai[rip[k]+1] - ai[rip[k]]; 2410 if (!ncols) SETERRQ2(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k); 2411 ncols_upper = 0; 2412 for (j=0; j<ncols; j++){ 2413 i = *(aj + ai[rip[k]] + j); /* unpermuted column index */ 2414 if (riip[i] >= k){ /* only take upper triangular entry */ 2415 ajtmp[ncols_upper] = i; 2416 ncols_upper++; 2417 } 2418 } 2419 ierr = PetscIncompleteLLInit(ncols_upper,ajtmp,am,riip,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 2420 nzk += nlnk; 2421 2422 /* update lnk by computing fill-in for each pivot row to be merged in */ 2423 prow = jl[k]; /* 1st pivot row */ 2424 2425 while (prow < k){ 2426 nextprow = jl[prow]; 2427 2428 /* merge prow into k-th row */ 2429 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 2430 jmax = ui[prow+1]; 2431 ncols = jmax-jmin; 2432 i = jmin - ui[prow]; 2433 cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 2434 uj = uj_lvl_ptr[prow] + i; /* levels of cols */ 2435 j = *(uj - 1); 2436 ierr = PetscICCLLAddSorted(ncols,cols,levels,uj,am,nlnk,lnk,lnk_lvl,lnkbt,j);CHKERRQ(ierr); 2437 nzk += nlnk; 2438 2439 /* update il and jl for prow */ 2440 if (jmin < jmax){ 2441 il[prow] = jmin; 2442 j = *cols; jl[prow] = jl[j]; jl[j] = prow; 2443 } 2444 prow = nextprow; 2445 } 2446 2447 /* if free space is not available, make more free space */ 2448 if (current_space->local_remaining<nzk) { 2449 i = am - k + 1; /* num of unfactored rows */ 2450 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 2451 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 2452 ierr = PetscFreeSpaceGet(i,¤t_space_lvl);CHKERRQ(ierr); 2453 reallocs++; 2454 } 2455 2456 /* copy data into free_space and free_space_lvl, then initialize lnk */ 2457 if (nzk == 0) SETERRQ1(PETSC_ERR_ARG_WRONG,"Empty row %D in ICC matrix factor",k); 2458 ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 2459 2460 /* add the k-th row into il and jl */ 2461 if (nzk > 1){ 2462 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 2463 jl[k] = jl[i]; jl[i] = k; 2464 il[k] = ui[k] + 1; 2465 } 2466 uj_ptr[k] = current_space->array; 2467 uj_lvl_ptr[k] = current_space_lvl->array; 2468 2469 current_space->array += nzk; 2470 current_space->local_used += nzk; 2471 current_space->local_remaining -= nzk; 2472 2473 current_space_lvl->array += nzk; 2474 current_space_lvl->local_used += nzk; 2475 current_space_lvl->local_remaining -= nzk; 2476 2477 ui[k+1] = ui[k] + nzk; 2478 } 2479 2480 #if defined(PETSC_USE_INFO) 2481 if (ai[am] != 0) { 2482 PetscReal af = (PetscReal)ui[am]/((PetscReal)ai[am]); 2483 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 2484 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 2485 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 2486 } else { 2487 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 2488 } 2489 #endif 2490 2491 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 2492 ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr); 2493 ierr = PetscFree(jl);CHKERRQ(ierr); 2494 ierr = PetscFree(ajtmp);CHKERRQ(ierr); 2495 2496 /* destroy list of free space and other temporary array(s) */ 2497 ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 2498 ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 2499 ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 2500 ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 2501 2502 } /* end of case: levels>0 || (levels=0 && !perm_identity) */ 2503 2504 /* put together the new matrix in MATSEQSBAIJ format */ 2505 2506 b = (Mat_SeqSBAIJ*)(fact)->data; 2507 b->singlemalloc = PETSC_FALSE; 2508 ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 2509 b->j = uj; 2510 b->i = ui; 2511 b->diag = udiag; 2512 b->free_diag = PETSC_TRUE; 2513 b->ilen = 0; 2514 b->imax = 0; 2515 b->row = perm; 2516 b->col = perm; 2517 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 2518 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 2519 b->icol = iperm; 2520 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 2521 ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 2522 ierr = PetscLogObjectMemory((fact),(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 2523 b->maxnz = b->nz = ui[am]; 2524 b->free_a = PETSC_TRUE; 2525 b->free_ij = PETSC_TRUE; 2526 2527 (fact)->info.factor_mallocs = reallocs; 2528 (fact)->info.fill_ratio_given = fill; 2529 if (ai[am] != 0) { 2530 (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 2531 } else { 2532 (fact)->info.fill_ratio_needed = 0.0; 2533 } 2534 (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 2535 PetscFunctionReturn(0); 2536 } 2537 2538 #undef __FUNCT__ 2539 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqAIJ" 2540 PetscErrorCode MatCholeskyFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 2541 { 2542 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 2543 Mat_SeqSBAIJ *b; 2544 PetscErrorCode ierr; 2545 PetscTruth perm_identity; 2546 PetscReal fill = info->fill; 2547 const PetscInt *rip,*riip; 2548 PetscInt i,am=A->rmap->n,*ai=a->i,*aj=a->j,reallocs=0,prow; 2549 PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 2550 PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 2551 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 2552 PetscBT lnkbt; 2553 IS iperm; 2554 2555 PetscFunctionBegin; 2556 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); 2557 /* check whether perm is the identity mapping */ 2558 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 2559 ierr = ISInvertPermutation(perm,PETSC_DECIDE,&iperm);CHKERRQ(ierr); 2560 ierr = ISGetIndices(iperm,&riip);CHKERRQ(ierr); 2561 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 2562 2563 /* initialization */ 2564 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 2565 ui[0] = 0; 2566 2567 /* jl: linked list for storing indices of the pivot rows 2568 il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 2569 ierr = PetscMalloc((3*am+1)*sizeof(PetscInt)+am*sizeof(PetscInt**),&jl);CHKERRQ(ierr); 2570 il = jl + am; 2571 cols = il + am; 2572 ui_ptr = (PetscInt**)(cols + am); 2573 for (i=0; i<am; i++){ 2574 jl[i] = am; il[i] = 0; 2575 } 2576 2577 /* create and initialize a linked list for storing column indices of the active row k */ 2578 nlnk = am + 1; 2579 ierr = PetscLLCreate(am,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 2580 2581 /* initial FreeSpace size is fill*(ai[am]+1) */ 2582 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 2583 current_space = free_space; 2584 2585 for (k=0; k<am; k++){ /* for each active row k */ 2586 /* initialize lnk by the column indices of row rip[k] of A */ 2587 nzk = 0; 2588 ncols = ai[rip[k]+1] - ai[rip[k]]; 2589 if (!ncols) SETERRQ2(PETSC_ERR_MAT_CH_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",rip[k],k); 2590 ncols_upper = 0; 2591 for (j=0; j<ncols; j++){ 2592 i = riip[*(aj + ai[rip[k]] + j)]; 2593 if (i >= k){ /* only take upper triangular entry */ 2594 cols[ncols_upper] = i; 2595 ncols_upper++; 2596 } 2597 } 2598 ierr = PetscLLAdd(ncols_upper,cols,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 2599 nzk += nlnk; 2600 2601 /* update lnk by computing fill-in for each pivot row to be merged in */ 2602 prow = jl[k]; /* 1st pivot row */ 2603 2604 while (prow < k){ 2605 nextprow = jl[prow]; 2606 /* merge prow into k-th row */ 2607 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 2608 jmax = ui[prow+1]; 2609 ncols = jmax-jmin; 2610 uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 2611 ierr = PetscLLAddSorted(ncols,uj_ptr,am,nlnk,lnk,lnkbt);CHKERRQ(ierr); 2612 nzk += nlnk; 2613 2614 /* update il and jl for prow */ 2615 if (jmin < jmax){ 2616 il[prow] = jmin; 2617 j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 2618 } 2619 prow = nextprow; 2620 } 2621 2622 /* if free space is not available, make more free space */ 2623 if (current_space->local_remaining<nzk) { 2624 i = am - k + 1; /* num of unfactored rows */ 2625 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 2626 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 2627 reallocs++; 2628 } 2629 2630 /* copy data into free space, then initialize lnk */ 2631 ierr = PetscLLClean(am,am,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 2632 2633 /* add the k-th row into il and jl */ 2634 if (nzk-1 > 0){ 2635 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 2636 jl[k] = jl[i]; jl[i] = k; 2637 il[k] = ui[k] + 1; 2638 } 2639 ui_ptr[k] = current_space->array; 2640 current_space->array += nzk; 2641 current_space->local_used += nzk; 2642 current_space->local_remaining -= nzk; 2643 2644 ui[k+1] = ui[k] + nzk; 2645 } 2646 2647 #if defined(PETSC_USE_INFO) 2648 if (ai[am] != 0) { 2649 PetscReal af = (PetscReal)(ui[am])/((PetscReal)ai[am]); 2650 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 2651 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 2652 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 2653 } else { 2654 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 2655 } 2656 #endif 2657 2658 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 2659 ierr = ISRestoreIndices(iperm,&riip);CHKERRQ(ierr); 2660 ierr = PetscFree(jl);CHKERRQ(ierr); 2661 2662 /* destroy list of free space and other temporary array(s) */ 2663 ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 2664 ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 2665 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 2666 2667 /* put together the new matrix in MATSEQSBAIJ format */ 2668 2669 b = (Mat_SeqSBAIJ*)(fact)->data; 2670 b->singlemalloc = PETSC_FALSE; 2671 b->free_a = PETSC_TRUE; 2672 b->free_ij = PETSC_TRUE; 2673 ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 2674 b->j = uj; 2675 b->i = ui; 2676 b->diag = 0; 2677 b->ilen = 0; 2678 b->imax = 0; 2679 b->row = perm; 2680 b->col = perm; 2681 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 2682 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 2683 b->icol = iperm; 2684 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 2685 ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 2686 ierr = PetscLogObjectMemory(fact,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 2687 b->maxnz = b->nz = ui[am]; 2688 2689 (fact)->info.factor_mallocs = reallocs; 2690 (fact)->info.fill_ratio_given = fill; 2691 if (ai[am] != 0) { 2692 (fact)->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 2693 } else { 2694 (fact)->info.fill_ratio_needed = 0.0; 2695 } 2696 (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqAIJ; 2697 PetscFunctionReturn(0); 2698 } 2699 2700 #undef __FUNCT__ 2701 #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering_newdatastruct" 2702 PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_newdatastruct(Mat A,Vec bb,Vec xx) 2703 { 2704 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 2705 PetscErrorCode ierr; 2706 PetscInt n = A->rmap->n; 2707 const PetscInt *ai = a->i,*aj = a->j,*vi; 2708 PetscScalar *x,sum; 2709 const PetscScalar *b; 2710 const MatScalar *aa = a->a,*v; 2711 PetscInt i,nz; 2712 2713 PetscFunctionBegin; 2714 if (!n) PetscFunctionReturn(0); 2715 2716 ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 2717 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2718 2719 /* forward solve the lower triangular */ 2720 x[0] = b[0]; 2721 v = aa; 2722 vi = aj; 2723 for (i=1; i<n; i++) { 2724 nz = ai[i+1] - ai[i]; 2725 sum = b[i]; 2726 PetscSparseDenseMinusDot(sum,x,v,vi,nz); 2727 v += nz; 2728 vi += nz; 2729 x[i] = sum; 2730 } 2731 2732 /* backward solve the upper triangular */ 2733 v = aa + ai[n+1]; 2734 vi = aj + ai[n+1]; 2735 for (i=n-1; i>=0; i--){ 2736 nz = ai[2*n-i +1] - ai[2*n-i]-1; 2737 sum = x[i]; 2738 PetscSparseDenseMinusDot(sum,x,v,vi,nz); 2739 v += nz; 2740 vi += nz; vi++; 2741 x[i] = *v++ *sum; /* x[i]=aa[adiag[i]]*sum; v++; */ 2742 } 2743 2744 ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr); 2745 ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 2746 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2747 PetscFunctionReturn(0); 2748 } 2749 2750 #undef __FUNCT__ 2751 #define __FUNCT__ "MatSolve_SeqAIJ_NaturalOrdering_newdatastruct_v2" 2752 PetscErrorCode MatSolve_SeqAIJ_NaturalOrdering_newdatastruct_v2(Mat A,Vec bb,Vec xx) 2753 { 2754 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 2755 PetscErrorCode ierr; 2756 PetscInt n = A->rmap->n; 2757 const PetscInt *ai = a->i,*aj = a->j,*adiag = a->diag,*vi; 2758 PetscScalar *x,sum; 2759 const PetscScalar *b; 2760 const MatScalar *aa = a->a,*v; 2761 PetscInt i,nz; 2762 2763 PetscFunctionBegin; 2764 if (!n) PetscFunctionReturn(0); 2765 2766 ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 2767 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2768 2769 /* forward solve the lower triangular */ 2770 x[0] = b[0]; 2771 v = aa; 2772 vi = aj; 2773 for (i=1; i<n; i++) { 2774 nz = ai[i+1] - ai[i]; 2775 sum = b[i]; 2776 PetscSparseDenseMinusDot(sum,x,v,vi,nz); 2777 v += nz; 2778 vi += nz; 2779 x[i] = sum; 2780 } 2781 2782 /* backward solve the upper triangular */ 2783 /* v = aa + ai[n+1]; 2784 vi = aj + ai[n+1]; */ 2785 v = aa + adiag[n-1]; 2786 vi = aj + adiag[n-1]; 2787 for (i=n-1; i>=0; i--){ 2788 nz = adiag[i] - adiag[i+1]-1; 2789 sum = x[i]; 2790 PetscSparseDenseMinusDot(sum,x,v,vi,nz); 2791 v += nz; 2792 vi += nz; vi++; 2793 x[i] = *v++ *sum; /* x[i]=aa[adiag[i]]*sum; v++; */ 2794 } 2795 2796 ierr = PetscLogFlops(2.0*a->nz - A->cmap->n);CHKERRQ(ierr); 2797 ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 2798 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2799 PetscFunctionReturn(0); 2800 } 2801 2802 #undef __FUNCT__ 2803 #define __FUNCT__ "MatSolve_SeqAIJ_newdatastruct" 2804 PetscErrorCode MatSolve_SeqAIJ_newdatastruct(Mat A,Vec bb,Vec xx) 2805 { 2806 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 2807 IS iscol = a->col,isrow = a->row; 2808 PetscErrorCode ierr; 2809 PetscInt i,n=A->rmap->n,*vi,*ai=a->i,*aj=a->j,nz,k; 2810 const PetscInt *rout,*cout,*r,*c; 2811 PetscScalar *x,*tmp,*tmps,sum; 2812 const PetscScalar *b; 2813 const MatScalar *aa = a->a,*v; 2814 2815 PetscFunctionBegin; 2816 if (!n) PetscFunctionReturn(0); 2817 2818 ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 2819 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2820 tmp = a->solve_work; 2821 2822 ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 2823 ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 2824 2825 /* forward solve the lower triangular */ 2826 tmp[0] = b[*r++]; 2827 tmps = tmp; 2828 v = aa; 2829 vi = aj; 2830 for (i=1; i<n; i++) { 2831 nz = ai[i+1] - ai[i]; 2832 sum = b[*r++]; 2833 PetscSparseDenseMinusDot(sum,tmps,v,vi,nz); 2834 tmp[i] = sum; 2835 v += nz; vi += nz; 2836 } 2837 2838 /* backward solve the upper triangular */ 2839 k = n+1; 2840 v = aa + ai[k]; /* 1st entry of U(n-1,:) */ 2841 vi = aj + ai[k]; 2842 for (i=n-1; i>=0; i--){ 2843 k = 2*n-i; 2844 nz = ai[k +1] - ai[k] - 1; 2845 sum = tmp[i]; 2846 PetscSparseDenseMinusDot(sum,tmps,v,vi,nz); 2847 x[c[i]] = tmp[i] = sum*v[nz]; /* v[nz] = aa[adiag[i]] */ 2848 v += nz+1; vi += nz+1; 2849 } 2850 2851 ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 2852 ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 2853 ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 2854 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2855 ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr); 2856 PetscFunctionReturn(0); 2857 } 2858 2859 #undef __FUNCT__ 2860 #define __FUNCT__ "MatSolve_SeqAIJ_newdatastruct_v2" 2861 PetscErrorCode MatSolve_SeqAIJ_newdatastruct_v2(Mat A,Vec bb,Vec xx) 2862 { 2863 Mat_SeqAIJ *a = (Mat_SeqAIJ*)A->data; 2864 IS iscol = a->col,isrow = a->row; 2865 PetscErrorCode ierr; 2866 PetscInt i,n=A->rmap->n,*vi,*ai=a->i,*aj=a->j,*adiag = a->diag,nz; 2867 const PetscInt *rout,*cout,*r,*c; 2868 PetscScalar *x,*tmp,sum; 2869 const PetscScalar *b; 2870 const MatScalar *aa = a->a,*v; 2871 2872 PetscFunctionBegin; 2873 if (!n) PetscFunctionReturn(0); 2874 2875 ierr = VecGetArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 2876 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 2877 tmp = a->solve_work; 2878 2879 ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 2880 ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout; 2881 2882 /* forward solve the lower triangular */ 2883 /* tmp[0] = b[*r++]; */ 2884 tmp[0] = b[r[0]]; 2885 v = aa; 2886 vi = aj; 2887 for (i=1; i<n; i++) { 2888 nz = ai[i+1] - ai[i]; 2889 sum = b[r[i]]; 2890 PetscSparseDenseMinusDot(sum,tmp,v,vi,nz); 2891 tmp[i] = sum; 2892 v += nz; vi += nz; 2893 } 2894 2895 /* backward solve the upper triangular */ 2896 /* v = aa + ai[k]; *//* 1st entry of U(n-1,:) */ 2897 /* vi = aj + ai[k]; */ 2898 v = aa + adiag[n-1]; 2899 vi = aj + adiag[n-1]; 2900 for (i=n-1; i>=0; i--){ 2901 /* nz = ai[k +1] - ai[k] - 1;*/ 2902 nz = adiag[i]-adiag[i+1]-1; 2903 sum = tmp[i]; 2904 PetscSparseDenseMinusDot(sum,tmp,v,vi,nz); 2905 x[c[i]] = tmp[i] = sum*v[nz]; /* v[nz] = aa[adiag[i]] */ 2906 v += nz+1; vi += nz+1; 2907 } 2908 2909 ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 2910 ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 2911 ierr = VecRestoreArray(bb,(PetscScalar**)&b);CHKERRQ(ierr); 2912 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 2913 ierr = PetscLogFlops(2*a->nz - A->cmap->n);CHKERRQ(ierr); 2914 PetscFunctionReturn(0); 2915 } 2916 2917 #undef __FUNCT__ 2918 #define __FUNCT__ "MatILUDTFactor_SeqAIJ" 2919 PetscErrorCode MatILUDTFactor_SeqAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact) 2920 { 2921 Mat B = *fact; 2922 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b; 2923 IS isicol; 2924 PetscErrorCode ierr; 2925 const PetscInt *r,*ic; 2926 PetscInt i,n=A->rmap->n,*ai=a->i,*aj=a->j,*ajtmp,*adiag; 2927 PetscInt *bi,*bj,*bdiag,*bdiag_rev; 2928 PetscInt row,nzi,nzi_bl,nzi_bu,*im,nzi_al,nzi_au; 2929 PetscInt nlnk,*lnk; 2930 PetscBT lnkbt; 2931 PetscTruth row_identity,icol_identity,both_identity; 2932 MatScalar *aatmp,*pv,*batmp,*ba,*rtmp,*pc,multiplier,*vtmp,diag_tmp; 2933 const PetscInt *ics; 2934 PetscInt j,nz,*pj,*bjtmp,k,ncut,*jtmp; 2935 PetscReal dt=info->dt,dtcol=info->dtcol,shift=info->shiftinblocks; 2936 PetscInt dtcount=(PetscInt)info->dtcount,nnz_max; 2937 PetscTruth missing; 2938 2939 PetscFunctionBegin; 2940 2941 if (dt == PETSC_DEFAULT) dt = 0.005; 2942 if (dtcol == PETSC_DEFAULT) dtcol = 0.01; /* XXX unused! */ 2943 if (dtcount == PETSC_DEFAULT) dtcount = (PetscInt)(1.5*a->rmax); 2944 2945 /* ------- symbolic factorization, can be reused ---------*/ 2946 ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr); 2947 if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i); 2948 adiag=a->diag; 2949 2950 ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 2951 2952 /* bdiag is location of diagonal in factor */ 2953 ierr = PetscMalloc((2*n+2)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 2954 bdiag_rev = bdiag + n+1; 2955 2956 /* allocate row pointers bi */ 2957 ierr = PetscMalloc((2*n+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 2958 2959 /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */ 2960 if (dtcount > n-1) dtcount = n-1; /* diagonal is excluded */ 2961 nnz_max = ai[n]+2*n*dtcount+2; 2962 2963 ierr = PetscMalloc((nnz_max+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 2964 ierr = PetscMalloc((nnz_max+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr); 2965 2966 /* put together the new matrix */ 2967 ierr = MatSeqAIJSetPreallocation_SeqAIJ(B,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 2968 ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr); 2969 b = (Mat_SeqAIJ*)(B)->data; 2970 b->free_a = PETSC_TRUE; 2971 b->free_ij = PETSC_TRUE; 2972 b->singlemalloc = PETSC_FALSE; 2973 b->a = ba; 2974 b->j = bj; 2975 b->i = bi; 2976 b->diag = bdiag; 2977 b->ilen = 0; 2978 b->imax = 0; 2979 b->row = isrow; 2980 b->col = iscol; 2981 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 2982 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 2983 b->icol = isicol; 2984 ierr = PetscMalloc((n+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 2985 2986 ierr = PetscLogObjectMemory(B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 2987 b->maxnz = nnz_max; 2988 2989 (B)->factor = MAT_FACTOR_ILUDT; 2990 (B)->info.factor_mallocs = 0; 2991 (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)ai[n]); 2992 CHKMEMQ; 2993 /* ------- end of symbolic factorization ---------*/ 2994 2995 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 2996 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 2997 ics = ic; 2998 2999 /* linked list for storing column indices of the active row */ 3000 nlnk = n + 1; 3001 ierr = PetscLLCreate(n,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 3002 3003 /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */ 3004 ierr = PetscMalloc((2*n+1)*sizeof(PetscInt),&im);CHKERRQ(ierr); 3005 jtmp = im + n; 3006 /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */ 3007 ierr = PetscMalloc((2*n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 3008 ierr = PetscMemzero(rtmp,(n+1)*sizeof(PetscScalar));CHKERRQ(ierr); 3009 vtmp = rtmp + n; 3010 3011 bi[0] = 0; 3012 bdiag[0] = nnz_max-1; /* location of diag[0] in factor B */ 3013 bdiag_rev[n] = bdiag[0]; 3014 bi[2*n+1] = bdiag[0]+1; /* endof bj and ba array */ 3015 for (i=0; i<n; i++) { 3016 /* copy initial fill into linked list */ 3017 nzi = 0; /* nonzeros for active row i */ 3018 nzi = ai[r[i]+1] - ai[r[i]]; 3019 if (!nzi) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 3020 nzi_al = adiag[r[i]] - ai[r[i]]; 3021 nzi_au = ai[r[i]+1] - adiag[r[i]] -1; 3022 ajtmp = aj + ai[r[i]]; 3023 ierr = PetscLLAddPerm(nzi,ajtmp,ic,n,nlnk,lnk,lnkbt);CHKERRQ(ierr); 3024 3025 /* load in initial (unfactored row) */ 3026 aatmp = a->a + ai[r[i]]; 3027 for (j=0; j<nzi; j++) { 3028 rtmp[ics[*ajtmp++]] = *aatmp++; 3029 } 3030 3031 /* add pivot rows into linked list */ 3032 row = lnk[n]; 3033 while (row < i ) { 3034 nzi_bl = bi[row+1] - bi[row] + 1; 3035 bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */ 3036 ierr = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr); 3037 nzi += nlnk; 3038 row = lnk[row]; 3039 } 3040 3041 /* copy data from lnk into jtmp, then initialize lnk */ 3042 ierr = PetscLLClean(n,n,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr); 3043 3044 /* numerical factorization */ 3045 bjtmp = jtmp; 3046 row = *bjtmp++; /* 1st pivot row */ 3047 while ( row < i ) { 3048 pc = rtmp + row; 3049 pv = ba + bdiag[row]; /* 1./(diag of the pivot row) */ 3050 multiplier = (*pc) * (*pv); 3051 *pc = multiplier; 3052 if (PetscAbsScalar(*pc) > dt){ /* apply tolerance dropping rule */ 3053 pj = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */ 3054 pv = ba + bdiag[row+1] + 1; 3055 /* if (multiplier < -1.0 or multiplier >1.0) printf("row/prow %d, %d, multiplier %g\n",i,row,multiplier); */ 3056 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */ 3057 for (j=0; j<nz; j++) rtmp[*pj++] -= multiplier * (*pv++); 3058 ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); 3059 } 3060 row = *bjtmp++; 3061 } 3062 3063 /* copy sparse rtmp into contiguous vtmp; separate L and U part */ 3064 diag_tmp = rtmp[i]; /* save diagonal value - may not needed?? */ 3065 nzi_bl = 0; j = 0; 3066 while (jtmp[j] < i){ /* Note: jtmp is sorted */ 3067 vtmp[j] = rtmp[jtmp[j]]; rtmp[jtmp[j]]=0.0; 3068 nzi_bl++; j++; 3069 } 3070 nzi_bu = nzi - nzi_bl -1; 3071 while (j < nzi){ 3072 vtmp[j] = rtmp[jtmp[j]]; rtmp[jtmp[j]]=0.0; 3073 j++; 3074 } 3075 3076 bjtmp = bj + bi[i]; 3077 batmp = ba + bi[i]; 3078 /* apply level dropping rule to L part */ 3079 ncut = nzi_al + dtcount; 3080 if (ncut < nzi_bl){ 3081 ierr = PetscSortSplit(ncut,nzi_bl,vtmp,jtmp);CHKERRQ(ierr); 3082 ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr); 3083 } else { 3084 ncut = nzi_bl; 3085 } 3086 for (j=0; j<ncut; j++){ 3087 bjtmp[j] = jtmp[j]; 3088 batmp[j] = vtmp[j]; 3089 /* printf(" (%d,%g),",bjtmp[j],batmp[j]); */ 3090 } 3091 bi[i+1] = bi[i] + ncut; 3092 nzi = ncut + 1; 3093 3094 /* apply level dropping rule to U part */ 3095 ncut = nzi_au + dtcount; 3096 if (ncut < nzi_bu){ 3097 ierr = PetscSortSplit(ncut,nzi_bu,vtmp+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr); 3098 ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr); 3099 } else { 3100 ncut = nzi_bu; 3101 } 3102 nzi += ncut; 3103 3104 /* mark bdiagonal */ 3105 bdiag[i+1] = bdiag[i] - (ncut + 1); 3106 bdiag_rev[n-i-1] = bdiag[i+1]; 3107 bi[2*n - i] = bi[2*n - i +1] - (ncut + 1); 3108 bjtmp = bj + bdiag[i]; 3109 batmp = ba + bdiag[i]; 3110 *bjtmp = i; 3111 *batmp = diag_tmp; /* rtmp[i]; */ 3112 if (*batmp == 0.0) { 3113 *batmp = dt+shift; 3114 /* printf(" row %d add shift %g\n",i,shift); */ 3115 } 3116 *batmp = 1.0/(*batmp); /* invert diagonal entries for simplier triangular solves */ 3117 /* printf(" (%d,%g),",*bjtmp,*batmp); */ 3118 3119 bjtmp = bj + bdiag[i+1]+1; 3120 batmp = ba + bdiag[i+1]+1; 3121 for (k=0; k<ncut; k++){ 3122 bjtmp[k] = jtmp[nzi_bl+1+k]; 3123 batmp[k] = vtmp[nzi_bl+1+k]; 3124 /* printf(" (%d,%g),",bjtmp[k],batmp[k]); */ 3125 } 3126 /* printf("\n"); */ 3127 3128 im[i] = nzi; /* used by PetscLLAddSortedLU() */ 3129 /* 3130 printf("row %d: bi %d, bdiag %d\n",i,bi[i],bdiag[i]); 3131 printf(" ----------------------------\n"); 3132 */ 3133 } /* for (i=0; i<n; i++) */ 3134 /* printf("end of L %d, beginning of U %d\n",bi[n],bdiag[n]); */ 3135 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]); 3136 3137 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 3138 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 3139 3140 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 3141 ierr = PetscFree(im);CHKERRQ(ierr); 3142 ierr = PetscFree(rtmp);CHKERRQ(ierr); 3143 3144 ierr = PetscLogFlops(B->cmap->n);CHKERRQ(ierr); 3145 b->maxnz = b->nz = bi[n] + bdiag[0] - bdiag[n]; 3146 3147 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 3148 ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr); 3149 both_identity = (PetscTruth) (row_identity && icol_identity); 3150 if (row_identity && icol_identity) { 3151 B->ops->solve = MatSolve_SeqAIJ_NaturalOrdering_newdatastruct; 3152 } else { 3153 B->ops->solve = MatSolve_SeqAIJ_newdatastruct; 3154 } 3155 3156 B->ops->lufactorsymbolic = MatILUDTFactorSymbolic_SeqAIJ; 3157 B->ops->lufactornumeric = MatILUDTFactorNumeric_SeqAIJ; 3158 B->ops->solveadd = 0; 3159 B->ops->solvetranspose = 0; 3160 B->ops->solvetransposeadd = 0; 3161 B->ops->matsolve = 0; 3162 B->assembled = PETSC_TRUE; 3163 B->preallocated = PETSC_TRUE; 3164 PetscFunctionReturn(0); 3165 } 3166 3167 /* a wraper of MatILUDTFactor_SeqAIJ() */ 3168 #undef __FUNCT__ 3169 #define __FUNCT__ "MatILUDTFactorSymbolic_SeqAIJ" 3170 PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorSymbolic_SeqAIJ(Mat fact,Mat A,IS row,IS col,const MatFactorInfo *info) 3171 { 3172 PetscErrorCode ierr; 3173 3174 PetscFunctionBegin; 3175 ierr = MatILUDTFactor_SeqAIJ(A,row,col,info,&fact);CHKERRQ(ierr); 3176 3177 fact->ops->lufactornumeric = MatILUDTFactorNumeric_SeqAIJ; 3178 PetscFunctionReturn(0); 3179 } 3180 3181 /* 3182 same as MatLUFactorNumeric_SeqAIJ(), except using contiguous array matrix factors 3183 - intend to replace existing MatLUFactorNumeric_SeqAIJ() 3184 */ 3185 #undef __FUNCT__ 3186 #define __FUNCT__ "MatILUDTFactorNumeric_SeqAIJ" 3187 PetscErrorCode PETSCMAT_DLLEXPORT MatILUDTFactorNumeric_SeqAIJ(Mat fact,Mat A,const MatFactorInfo *info) 3188 { 3189 Mat C=fact; 3190 Mat_SeqAIJ *a=(Mat_SeqAIJ*)A->data,*b=(Mat_SeqAIJ *)C->data; 3191 IS isrow = b->row,isicol = b->icol; 3192 PetscErrorCode ierr; 3193 const PetscInt *r,*ic,*ics; 3194 PetscInt i,j,k,n=A->rmap->n,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 3195 PetscInt *ajtmp,*bjtmp,nz,nzl,nzu,row,*bdiag = b->diag,*pj; 3196 MatScalar *rtmp,*pc,multiplier,*v,*pv,*aa=a->a; 3197 PetscReal dt=info->dt,shift=info->shiftinblocks; 3198 PetscTruth row_identity, col_identity; 3199 3200 PetscFunctionBegin; 3201 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 3202 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 3203 ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 3204 ics = ic; 3205 3206 for (i=0; i<n; i++){ 3207 /* initialize rtmp array */ 3208 nzl = bi[i+1] - bi[i]; /* num of nozeros in L(i,:) */ 3209 bjtmp = bj + bi[i]; 3210 for (j=0; j<nzl; j++) rtmp[*bjtmp++] = 0.0; 3211 rtmp[i] = 0.0; 3212 nzu = bdiag[i] - bdiag[i+1]; /* num of nozeros in U(i,:) */ 3213 bjtmp = bj + bdiag[i+1] + 1; 3214 for (j=0; j<nzu; j++) rtmp[*bjtmp++] = 0.0; 3215 3216 /* load in initial unfactored row of A */ 3217 /* printf("row %d\n",i); */ 3218 nz = ai[r[i]+1] - ai[r[i]]; 3219 ajtmp = aj + ai[r[i]]; 3220 v = aa + ai[r[i]]; 3221 for (j=0; j<nz; j++) { 3222 rtmp[ics[*ajtmp++]] = v[j]; 3223 /* printf(" (%d,%g),",ics[ajtmp[j]],rtmp[ics[ajtmp[j]]]); */ 3224 } 3225 /* printf("\n"); */ 3226 3227 /* numerical factorization */ 3228 bjtmp = bj + bi[i]; /* point to 1st entry of L(i,:) */ 3229 nzl = bi[i+1] - bi[i]; /* num of entries in L(i,:) */ 3230 k = 0; 3231 while (k < nzl){ 3232 row = *bjtmp++; 3233 /* printf(" prow %d\n",row); */ 3234 pc = rtmp + row; 3235 pv = b->a + bdiag[row]; /* 1./(diag of the pivot row) */ 3236 multiplier = (*pc) * (*pv); 3237 *pc = multiplier; 3238 if (PetscAbsScalar(multiplier) > dt){ 3239 pj = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */ 3240 pv = b->a + bdiag[row+1] + 1; 3241 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */ 3242 for (j=0; j<nz; j++) rtmp[*pj++] -= multiplier * (*pv++); 3243 /* ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); */ 3244 } 3245 k++; 3246 } 3247 3248 /* finished row so stick it into b->a */ 3249 /* L-part */ 3250 pv = b->a + bi[i] ; 3251 pj = bj + bi[i] ; 3252 nzl = bi[i+1] - bi[i]; 3253 for (j=0; j<nzl; j++) { 3254 pv[j] = rtmp[pj[j]]; 3255 /* printf(" (%d,%g),",pj[j],pv[j]); */ 3256 } 3257 3258 /* diagonal: invert diagonal entries for simplier triangular solves */ 3259 if (rtmp[i] == 0.0) rtmp[i] = dt+shift; 3260 b->a[bdiag[i]] = 1.0/rtmp[i]; 3261 /* printf(" (%d,%g),",i,b->a[bdiag[i]]); */ 3262 3263 /* U-part */ 3264 pv = b->a + bdiag[i+1] + 1; 3265 pj = bj + bdiag[i+1] + 1; 3266 nzu = bdiag[i] - bdiag[i+1] - 1; 3267 for (j=0; j<nzu; j++) { 3268 pv[j] = rtmp[pj[j]]; 3269 /* printf(" (%d,%g),",pj[j],pv[j]); */ 3270 } 3271 /* printf("\n"); */ 3272 } 3273 3274 ierr = PetscFree(rtmp);CHKERRQ(ierr); 3275 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 3276 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 3277 3278 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 3279 ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr); 3280 if (row_identity && col_identity) { 3281 C->ops->solve = MatSolve_SeqAIJ_NaturalOrdering_newdatastruct; 3282 } else { 3283 C->ops->solve = MatSolve_SeqAIJ_newdatastruct; 3284 } 3285 C->ops->solveadd = 0; 3286 C->ops->solvetranspose = 0; 3287 C->ops->solvetransposeadd = 0; 3288 C->ops->matsolve = 0; 3289 C->assembled = PETSC_TRUE; 3290 C->preallocated = PETSC_TRUE; 3291 ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr); 3292 PetscFunctionReturn(0); 3293 } 3294