1 #define PETSCMAT_DLL 2 3 #include "src/mat/impls/baij/seq/baij.h" 4 #include "src/mat/impls/sbaij/seq/sbaij.h" 5 #include "src/inline/ilu.h" 6 #include "include/petscis.h" 7 8 #if !defined(PETSC_USE_COMPLEX) 9 /* 10 input: 11 F -- numeric factor 12 output: 13 nneg, nzero, npos: matrix inertia 14 */ 15 16 #undef __FUNCT__ 17 #define __FUNCT__ "MatGetInertia_SeqSBAIJ" 18 PetscErrorCode MatGetInertia_SeqSBAIJ(Mat F,PetscInt *nneig,PetscInt *nzero,PetscInt *npos) 19 { 20 Mat_SeqSBAIJ *fact_ptr = (Mat_SeqSBAIJ*)F->data; 21 MatScalar *dd = fact_ptr->a; 22 PetscInt mbs=fact_ptr->mbs,bs=F->rmap.bs,i,nneig_tmp,npos_tmp,*fi = fact_ptr->i; 23 24 PetscFunctionBegin; 25 if (bs != 1) SETERRQ1(PETSC_ERR_SUP,"No support for bs: %D >1 yet",bs); 26 nneig_tmp = 0; npos_tmp = 0; 27 for (i=0; i<mbs; i++){ 28 if (PetscRealPart(dd[*fi]) > 0.0){ 29 npos_tmp++; 30 } else if (PetscRealPart(dd[*fi]) < 0.0){ 31 nneig_tmp++; 32 } 33 fi++; 34 } 35 if (nneig) *nneig = nneig_tmp; 36 if (npos) *npos = npos_tmp; 37 if (nzero) *nzero = mbs - nneig_tmp - npos_tmp; 38 39 PetscFunctionReturn(0); 40 } 41 #endif /* !defined(PETSC_USE_COMPLEX) */ 42 43 /* 44 Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP. 45 Use Modified Sparse Row (MSR) storage for u and ju. See page 85, "Iterative Methods ..." by Saad. 46 */ 47 #undef __FUNCT__ 48 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_MSR" 49 PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat A,IS perm,MatFactorInfo *info,Mat *B) 50 { 51 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b; 52 PetscErrorCode ierr; 53 PetscInt *rip,i,mbs = a->mbs,*ai,*aj; 54 PetscInt *jutmp,bs = A->rmap.bs,bs2=a->bs2; 55 PetscInt m,reallocs = 0,prow; 56 PetscInt *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd; 57 PetscReal f = info->fill; 58 PetscTruth perm_identity; 59 60 PetscFunctionBegin; 61 /* check whether perm is the identity mapping */ 62 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 63 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 64 65 if (perm_identity){ /* without permutation */ 66 a->permute = PETSC_FALSE; 67 ai = a->i; aj = a->j; 68 } else { /* non-trivial permutation */ 69 a->permute = PETSC_TRUE; 70 ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr); 71 ai = a->inew; aj = a->jnew; 72 } 73 74 /* initialization */ 75 ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&iu);CHKERRQ(ierr); 76 umax = (PetscInt)(f*ai[mbs] + 1); umax += mbs + 1; 77 ierr = PetscMalloc(umax*sizeof(PetscInt),&ju);CHKERRQ(ierr); 78 iu[0] = mbs+1; 79 juidx = mbs + 1; /* index for ju */ 80 ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&jl);CHKERRQ(ierr); /* linked list for pivot row */ 81 q = jl + mbs; /* linked list for col index */ 82 for (i=0; i<mbs; i++){ 83 jl[i] = mbs; 84 q[i] = 0; 85 } 86 87 /* for each row k */ 88 for (k=0; k<mbs; k++){ 89 for (i=0; i<mbs; i++) q[i] = 0; /* to be removed! */ 90 nzk = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */ 91 q[k] = mbs; 92 /* initialize nonzero structure of k-th row to row rip[k] of A */ 93 jmin = ai[rip[k]] +1; /* exclude diag[k] */ 94 jmax = ai[rip[k]+1]; 95 for (j=jmin; j<jmax; j++){ 96 vj = rip[aj[j]]; /* col. value */ 97 if(vj > k){ 98 qm = k; 99 do { 100 m = qm; qm = q[m]; 101 } while(qm < vj); 102 if (qm == vj) { 103 SETERRQ(PETSC_ERR_PLIB,"Duplicate entry in A\n"); 104 } 105 nzk++; 106 q[m] = vj; 107 q[vj] = qm; 108 } /* if(vj > k) */ 109 } /* for (j=jmin; j<jmax; j++) */ 110 111 /* modify nonzero structure of k-th row by computing fill-in 112 for each row i to be merged in */ 113 prow = k; 114 prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */ 115 116 while (prow < k){ 117 /* merge row prow into k-th row */ 118 jmin = iu[prow] + 1; jmax = iu[prow+1]; 119 qm = k; 120 for (j=jmin; j<jmax; j++){ 121 vj = ju[j]; 122 do { 123 m = qm; qm = q[m]; 124 } while (qm < vj); 125 if (qm != vj){ 126 nzk++; q[m] = vj; q[vj] = qm; qm = vj; 127 } 128 } 129 prow = jl[prow]; /* next pivot row */ 130 } 131 132 /* add k to row list for first nonzero element in k-th row */ 133 if (nzk > 0){ 134 i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */ 135 jl[k] = jl[i]; jl[i] = k; 136 } 137 iu[k+1] = iu[k] + nzk; 138 139 /* allocate more space to ju if needed */ 140 if (iu[k+1] > umax) { 141 /* estimate how much additional space we will need */ 142 /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 143 /* just double the memory each time */ 144 maxadd = umax; 145 if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2; 146 umax += maxadd; 147 148 /* allocate a longer ju */ 149 ierr = PetscMalloc(umax*sizeof(PetscInt),&jutmp);CHKERRQ(ierr); 150 ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(PetscInt));CHKERRQ(ierr); 151 ierr = PetscFree(ju);CHKERRQ(ierr); 152 ju = jutmp; 153 reallocs++; /* count how many times we realloc */ 154 } 155 156 /* save nonzero structure of k-th row in ju */ 157 i=k; 158 while (nzk --) { 159 i = q[i]; 160 ju[juidx++] = i; 161 } 162 } 163 164 #if defined(PETSC_USE_INFO) 165 if (ai[mbs] != 0) { 166 PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); 167 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,f,af);CHKERRQ(ierr); 168 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 169 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G);\n",af);CHKERRQ(ierr); 170 ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 171 } else { 172 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 173 } 174 #endif 175 176 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 177 ierr = PetscFree(jl);CHKERRQ(ierr); 178 179 /* put together the new matrix */ 180 ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 181 182 /* ierr = PetscLogObjectParent(*B,iperm);CHKERRQ(ierr); */ 183 b = (Mat_SeqSBAIJ*)(*B)->data; 184 b->singlemalloc = PETSC_FALSE; 185 b->free_a = PETSC_TRUE; 186 b->free_ij = PETSC_TRUE; 187 ierr = PetscMalloc((iu[mbs]+1)*sizeof(MatScalar)*bs2,&b->a);CHKERRQ(ierr); 188 b->j = ju; 189 b->i = iu; 190 b->diag = 0; 191 b->ilen = 0; 192 b->imax = 0; 193 b->row = perm; 194 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 195 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 196 b->icol = perm; 197 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 198 ierr = PetscMalloc((bs*mbs+bs)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 199 /* In b structure: Free imax, ilen, old a, old j. 200 Allocate idnew, solve_work, new a, new j */ 201 ierr = PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 202 b->maxnz = b->nz = iu[mbs]; 203 204 (*B)->factor = MAT_FACTOR_CHOLESKY; 205 (*B)->info.factor_mallocs = reallocs; 206 (*B)->info.fill_ratio_given = f; 207 if (ai[mbs] != 0) { 208 (*B)->info.fill_ratio_needed = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); 209 } else { 210 (*B)->info.fill_ratio_needed = 0.0; 211 } 212 213 if (perm_identity){ 214 switch (bs) { 215 case 1: 216 (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 217 (*B)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 218 (*B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering; 219 (*B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering; 220 ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=1\n");CHKERRQ(ierr); 221 break; 222 case 2: 223 (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 224 (*B)->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 225 (*B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering; 226 (*B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering; 227 ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=2\n");CHKERRQ(ierr); 228 break; 229 case 3: 230 (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 231 (*B)->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 232 (*B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_3_NaturalOrdering; 233 (*B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering; 234 ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=3\n");CHKERRQ(ierr); 235 break; 236 case 4: 237 (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 238 (*B)->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 239 (*B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_4_NaturalOrdering; 240 (*B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_4_NaturalOrdering; 241 ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=4\n");CHKERRQ(ierr); 242 break; 243 case 5: 244 (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 245 (*B)->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 246 (*B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_5_NaturalOrdering; 247 (*B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_5_NaturalOrdering; 248 ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=5\n");CHKERRQ(ierr); 249 break; 250 case 6: 251 (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 252 (*B)->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 253 (*B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_6_NaturalOrdering; 254 (*B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_6_NaturalOrdering; 255 ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=6\n");CHKERRQ(ierr); 256 break; 257 case 7: 258 (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 259 (*B)->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 260 (*B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_7_NaturalOrdering; 261 (*B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_7_NaturalOrdering; 262 ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=7\n");CHKERRQ(ierr); 263 break; 264 default: 265 (*B)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering; 266 (*B)->ops->solve = MatSolve_SeqSBAIJ_N_NaturalOrdering; 267 (*B)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering; 268 (*B)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering; 269 ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS>7\n");CHKERRQ(ierr); 270 break; 271 } 272 } 273 PetscFunctionReturn(0); 274 } 275 /* 276 Symbolic U^T*D*U factorization for SBAIJ format. 277 */ 278 #include "petscbt.h" 279 #include "src/mat/utils/freespace.h" 280 #undef __FUNCT__ 281 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ" 282 PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info,Mat *fact) 283 { 284 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 285 Mat_SeqSBAIJ *b; 286 PetscErrorCode ierr; 287 PetscTruth perm_identity,missing; 288 PetscReal fill = info->fill; 289 PetscInt *rip,i,mbs=a->mbs,bs=A->rmap.bs,*ai,*aj,reallocs=0,prow,d; 290 PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 291 PetscInt nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr; 292 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 293 PetscBT lnkbt; 294 295 PetscFunctionBegin; 296 ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 297 if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 298 299 /* 300 This code originally uses Modified Sparse Row (MSR) storage 301 (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choise! 302 Then it is rewritten so the factor B takes seqsbaij format. However the associated 303 MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity, 304 thus the original code in MSR format is still used for these cases. 305 The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever 306 MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor. 307 */ 308 if (bs > 1){ 309 ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(A,perm,info,fact);CHKERRQ(ierr); 310 PetscFunctionReturn(0); 311 } 312 313 /* check whether perm is the identity mapping */ 314 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 315 316 if (perm_identity){ 317 a->permute = PETSC_FALSE; 318 ai = a->i; aj = a->j; 319 } else { 320 SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format"); 321 /* There are bugs for reordeing. Needs further work. 322 MatReordering for sbaij cannot be efficient. User should use aij formt! */ 323 a->permute = PETSC_TRUE; 324 ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr); 325 ai = a->inew; aj = a->jnew; 326 } 327 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 328 329 /* initialization */ 330 ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 331 ui[0] = 0; 332 333 /* jl: linked list for storing indices of the pivot rows 334 il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */ 335 ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 336 il = jl + mbs; 337 cols = il + mbs; 338 ui_ptr = (PetscInt**)(cols + mbs); 339 340 for (i=0; i<mbs; i++){ 341 jl[i] = mbs; il[i] = 0; 342 } 343 344 /* create and initialize a linked list for storing column indices of the active row k */ 345 nlnk = mbs + 1; 346 ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 347 348 /* initial FreeSpace size is fill*(ai[mbs]+1) */ 349 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr); 350 current_space = free_space; 351 352 for (k=0; k<mbs; k++){ /* for each active row k */ 353 /* initialize lnk by the column indices of row rip[k] of A */ 354 nzk = 0; 355 ncols = ai[rip[k]+1] - ai[rip[k]]; 356 for (j=0; j<ncols; j++){ 357 i = *(aj + ai[rip[k]] + j); 358 cols[j] = rip[i]; 359 } 360 ierr = PetscLLAdd(ncols,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 361 nzk += nlnk; 362 363 /* update lnk by computing fill-in for each pivot row to be merged in */ 364 prow = jl[k]; /* 1st pivot row */ 365 366 while (prow < k){ 367 nextprow = jl[prow]; 368 /* merge prow into k-th row */ 369 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */ 370 jmax = ui[prow+1]; 371 ncols = jmax-jmin; 372 uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */ 373 ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 374 nzk += nlnk; 375 376 /* update il and jl for prow */ 377 if (jmin < jmax){ 378 il[prow] = jmin; 379 j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 380 } 381 prow = nextprow; 382 } 383 384 /* if free space is not available, make more free space */ 385 if (current_space->local_remaining<nzk) { 386 i = mbs - k + 1; /* num of unfactored rows */ 387 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 388 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 389 reallocs++; 390 } 391 392 /* copy data into free space, then initialize lnk */ 393 ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 394 395 /* add the k-th row into il and jl */ 396 if (nzk-1 > 0){ 397 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 398 jl[k] = jl[i]; jl[i] = k; 399 il[k] = ui[k] + 1; 400 } 401 ui_ptr[k] = current_space->array; 402 current_space->array += nzk; 403 current_space->local_used += nzk; 404 current_space->local_remaining -= nzk; 405 406 ui[k+1] = ui[k] + nzk; 407 } 408 409 #if defined(PETSC_USE_INFO) 410 if (ai[mbs] != 0) { 411 PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 412 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 413 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 414 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 415 } else { 416 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 417 } 418 #endif 419 420 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 421 ierr = PetscFree(jl);CHKERRQ(ierr); 422 423 /* destroy list of free space and other temporary array(s) */ 424 ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 425 ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 426 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 427 428 /* put together the new matrix in MATSEQSBAIJ format */ 429 ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(*fact,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 430 431 b = (Mat_SeqSBAIJ*)(*fact)->data; 432 b->singlemalloc = PETSC_FALSE; 433 b->free_a = PETSC_TRUE; 434 b->free_ij = PETSC_TRUE; 435 ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 436 b->j = uj; 437 b->i = ui; 438 b->diag = 0; 439 b->ilen = 0; 440 b->imax = 0; 441 b->row = perm; 442 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 443 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 444 b->icol = perm; 445 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 446 ierr = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 447 ierr = PetscLogObjectMemory(*fact,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 448 b->maxnz = b->nz = ui[mbs]; 449 450 (*fact)->factor = MAT_FACTOR_CHOLESKY; 451 (*fact)->info.factor_mallocs = reallocs; 452 (*fact)->info.fill_ratio_given = fill; 453 if (ai[mbs] != 0) { 454 (*fact)->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 455 } else { 456 (*fact)->info.fill_ratio_needed = 0.0; 457 } 458 459 if (perm_identity){ 460 switch (bs) { 461 case 1: 462 (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 463 (*fact)->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 464 (*fact)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering; 465 (*fact)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering; 466 ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=1\n");CHKERRQ(ierr); 467 break; 468 case 2: 469 (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering; 470 (*fact)->ops->solve = MatSolve_SeqSBAIJ_2_NaturalOrdering; 471 (*fact)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_2_NaturalOrdering; 472 (*fact)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_2_NaturalOrdering; 473 ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=2\n");CHKERRQ(ierr); 474 break; 475 case 3: 476 (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering; 477 (*fact)->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering; 478 (*fact)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_3_NaturalOrdering; 479 (*fact)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering; 480 ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=3\n");CHKERRQ(ierr); 481 break; 482 case 4: 483 (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering; 484 (*fact)->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering; 485 (*fact)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_4_NaturalOrdering; 486 (*fact)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_4_NaturalOrdering; 487 ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=4\n");CHKERRQ(ierr); 488 break; 489 case 5: 490 (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_5_NaturalOrdering; 491 (*fact)->ops->solve = MatSolve_SeqSBAIJ_5_NaturalOrdering; 492 (*fact)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_5_NaturalOrdering; 493 (*fact)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_5_NaturalOrdering; 494 ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=5\n");CHKERRQ(ierr); 495 break; 496 case 6: 497 (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_6_NaturalOrdering; 498 (*fact)->ops->solve = MatSolve_SeqSBAIJ_6_NaturalOrdering; 499 (*fact)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_6_NaturalOrdering; 500 (*fact)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_6_NaturalOrdering; 501 ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=6\n");CHKERRQ(ierr); 502 break; 503 case 7: 504 (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_7_NaturalOrdering; 505 (*fact)->ops->solve = MatSolve_SeqSBAIJ_7_NaturalOrdering; 506 (*fact)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_7_NaturalOrdering; 507 (*fact)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_7_NaturalOrdering; 508 ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS=7\n");CHKERRQ(ierr); 509 break; 510 default: 511 (*fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering; 512 (*fact)->ops->solve = MatSolve_SeqSBAIJ_N_NaturalOrdering; 513 (*fact)->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering; 514 (*fact)->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering; 515 ierr = PetscInfo(A,"Using special in-place natural ordering factor and solve BS>7\n");CHKERRQ(ierr); 516 break; 517 } 518 } 519 PetscFunctionReturn(0); 520 } 521 #undef __FUNCT__ 522 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N" 523 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat A,MatFactorInfo *info,Mat *B) 524 { 525 Mat C = *B; 526 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 527 IS perm = b->row; 528 PetscErrorCode ierr; 529 PetscInt *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 530 PetscInt *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 531 PetscInt bs=A->rmap.bs,bs2 = a->bs2,bslog = 0; 532 MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 533 MatScalar *u,*diag,*rtmp,*rtmp_ptr; 534 MatScalar *work; 535 PetscInt *pivots; 536 537 PetscFunctionBegin; 538 /* initialization */ 539 ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 540 ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr); 541 ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 542 jl = il + mbs; 543 for (i=0; i<mbs; i++) { 544 jl[i] = mbs; il[0] = 0; 545 } 546 ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr); 547 uik = dk + bs2; 548 work = uik + bs2; 549 ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr); 550 551 ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr); 552 553 /* check permutation */ 554 if (!a->permute){ 555 ai = a->i; aj = a->j; aa = a->a; 556 } else { 557 ai = a->inew; aj = a->jnew; 558 ierr = PetscMalloc(bs2*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 559 ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 560 ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr); 561 ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr); 562 563 /* flops in while loop */ 564 bslog = 2*bs*bs2; 565 566 for (i=0; i<mbs; i++){ 567 jmin = ai[i]; jmax = ai[i+1]; 568 for (j=jmin; j<jmax; j++){ 569 while (a2anew[j] != j){ 570 k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 571 for (k1=0; k1<bs2; k1++){ 572 dk[k1] = aa[k*bs2+k1]; 573 aa[k*bs2+k1] = aa[j*bs2+k1]; 574 aa[j*bs2+k1] = dk[k1]; 575 } 576 } 577 /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 578 if (i > aj[j]){ 579 /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 580 ap = aa + j*bs2; /* ptr to the beginning of j-th block of aa */ 581 for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */ 582 for (k=0; k<bs; k++){ /* j-th block of aa <- dk^T */ 583 for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1]; 584 } 585 } 586 } 587 } 588 ierr = PetscFree(a2anew);CHKERRQ(ierr); 589 } 590 591 /* for each row k */ 592 for (k = 0; k<mbs; k++){ 593 594 /*initialize k-th row with elements nonzero in row perm(k) of A */ 595 jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 596 597 ap = aa + jmin*bs2; 598 for (j = jmin; j < jmax; j++){ 599 vj = perm_ptr[aj[j]]; /* block col. index */ 600 rtmp_ptr = rtmp + vj*bs2; 601 for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++; 602 } 603 604 /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 605 ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 606 i = jl[k]; /* first row to be added to k_th row */ 607 608 while (i < k){ 609 nexti = jl[i]; /* next row to be added to k_th row */ 610 611 /* compute multiplier */ 612 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 613 614 /* uik = -inv(Di)*U_bar(i,k) */ 615 diag = ba + i*bs2; 616 u = ba + ili*bs2; 617 ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 618 Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u); 619 620 /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 621 Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u); 622 ierr = PetscLogFlops(bslog*2);CHKERRQ(ierr); 623 624 /* update -U(i,k) */ 625 ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 626 627 /* add multiple of row i to k-th row ... */ 628 jmin = ili + 1; jmax = bi[i+1]; 629 if (jmin < jmax){ 630 for (j=jmin; j<jmax; j++) { 631 /* rtmp += -U(i,k)^T * U_bar(i,j) */ 632 rtmp_ptr = rtmp + bj[j]*bs2; 633 u = ba + j*bs2; 634 Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u); 635 } 636 ierr = PetscLogFlops(bslog*(jmax-jmin));CHKERRQ(ierr); 637 638 /* ... add i to row list for next nonzero entry */ 639 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 640 j = bj[jmin]; 641 jl[i] = jl[j]; jl[j] = i; /* update jl */ 642 } 643 i = nexti; 644 } 645 646 /* save nonzero entries in k-th row of U ... */ 647 648 /* invert diagonal block */ 649 diag = ba+k*bs2; 650 ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr); 651 ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr); 652 653 jmin = bi[k]; jmax = bi[k+1]; 654 if (jmin < jmax) { 655 for (j=jmin; j<jmax; j++){ 656 vj = bj[j]; /* block col. index of U */ 657 u = ba + j*bs2; 658 rtmp_ptr = rtmp + vj*bs2; 659 for (k1=0; k1<bs2; k1++){ 660 *u++ = *rtmp_ptr; 661 *rtmp_ptr++ = 0.0; 662 } 663 } 664 665 /* ... add k to row list for first nonzero entry in k-th row */ 666 il[k] = jmin; 667 i = bj[jmin]; 668 jl[k] = jl[i]; jl[i] = k; 669 } 670 } 671 672 ierr = PetscFree(rtmp);CHKERRQ(ierr); 673 ierr = PetscFree(il);CHKERRQ(ierr); 674 ierr = PetscFree(dk);CHKERRQ(ierr); 675 ierr = PetscFree(pivots);CHKERRQ(ierr); 676 if (a->permute){ 677 ierr = PetscFree(aa);CHKERRQ(ierr); 678 } 679 680 ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 681 C->factor = MAT_FACTOR_CHOLESKY; 682 C->assembled = PETSC_TRUE; 683 C->preallocated = PETSC_TRUE; 684 ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 685 PetscFunctionReturn(0); 686 } 687 688 #undef __FUNCT__ 689 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering" 690 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B) 691 { 692 Mat C = *B; 693 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 694 PetscErrorCode ierr; 695 PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 696 PetscInt *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 697 PetscInt bs=A->rmap.bs,bs2 = a->bs2,bslog; 698 MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 699 MatScalar *u,*diag,*rtmp,*rtmp_ptr; 700 MatScalar *work; 701 PetscInt *pivots; 702 703 PetscFunctionBegin; 704 ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 705 ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr); 706 ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 707 jl = il + mbs; 708 for (i=0; i<mbs; i++) { 709 jl[i] = mbs; il[0] = 0; 710 } 711 ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr); 712 uik = dk + bs2; 713 work = uik + bs2; 714 ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr); 715 716 ai = a->i; aj = a->j; aa = a->a; 717 718 /* flops in while loop */ 719 bslog = 2*bs*bs2; 720 721 /* for each row k */ 722 for (k = 0; k<mbs; k++){ 723 724 /*initialize k-th row with elements nonzero in row k of A */ 725 jmin = ai[k]; jmax = ai[k+1]; 726 ap = aa + jmin*bs2; 727 for (j = jmin; j < jmax; j++){ 728 vj = aj[j]; /* block col. index */ 729 rtmp_ptr = rtmp + vj*bs2; 730 for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++; 731 } 732 733 /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 734 ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 735 i = jl[k]; /* first row to be added to k_th row */ 736 737 while (i < k){ 738 nexti = jl[i]; /* next row to be added to k_th row */ 739 740 /* compute multiplier */ 741 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 742 743 /* uik = -inv(Di)*U_bar(i,k) */ 744 diag = ba + i*bs2; 745 u = ba + ili*bs2; 746 ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 747 Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u); 748 749 /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 750 Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u); 751 ierr = PetscLogFlops(bslog*2);CHKERRQ(ierr); 752 753 /* update -U(i,k) */ 754 ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 755 756 /* add multiple of row i to k-th row ... */ 757 jmin = ili + 1; jmax = bi[i+1]; 758 if (jmin < jmax){ 759 for (j=jmin; j<jmax; j++) { 760 /* rtmp += -U(i,k)^T * U_bar(i,j) */ 761 rtmp_ptr = rtmp + bj[j]*bs2; 762 u = ba + j*bs2; 763 Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u); 764 } 765 ierr = PetscLogFlops(bslog*(jmax-jmin));CHKERRQ(ierr); 766 767 /* ... add i to row list for next nonzero entry */ 768 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 769 j = bj[jmin]; 770 jl[i] = jl[j]; jl[j] = i; /* update jl */ 771 } 772 i = nexti; 773 } 774 775 /* save nonzero entries in k-th row of U ... */ 776 777 /* invert diagonal block */ 778 diag = ba+k*bs2; 779 ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr); 780 ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr); 781 782 jmin = bi[k]; jmax = bi[k+1]; 783 if (jmin < jmax) { 784 for (j=jmin; j<jmax; j++){ 785 vj = bj[j]; /* block col. index of U */ 786 u = ba + j*bs2; 787 rtmp_ptr = rtmp + vj*bs2; 788 for (k1=0; k1<bs2; k1++){ 789 *u++ = *rtmp_ptr; 790 *rtmp_ptr++ = 0.0; 791 } 792 } 793 794 /* ... add k to row list for first nonzero entry in k-th row */ 795 il[k] = jmin; 796 i = bj[jmin]; 797 jl[k] = jl[i]; jl[i] = k; 798 } 799 } 800 801 ierr = PetscFree(rtmp);CHKERRQ(ierr); 802 ierr = PetscFree(il);CHKERRQ(ierr); 803 ierr = PetscFree(dk);CHKERRQ(ierr); 804 ierr = PetscFree(pivots);CHKERRQ(ierr); 805 806 C->factor = MAT_FACTOR_CHOLESKY; 807 C->assembled = PETSC_TRUE; 808 C->preallocated = PETSC_TRUE; 809 ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 810 PetscFunctionReturn(0); 811 } 812 813 /* 814 Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP. 815 Version for blocks 2 by 2. 816 */ 817 #undef __FUNCT__ 818 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2" 819 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat A,MatFactorInfo *info,Mat *B) 820 { 821 Mat C = *B; 822 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 823 IS perm = b->row; 824 PetscErrorCode ierr; 825 PetscInt *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 826 PetscInt *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 827 MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 828 MatScalar *u,*diag,*rtmp,*rtmp_ptr; 829 PetscReal shift = info->shiftinblocks; 830 831 PetscFunctionBegin; 832 /* initialization */ 833 /* il and jl record the first nonzero element in each row of the accessing 834 window U(0:k, k:mbs-1). 835 jl: list of rows to be added to uneliminated rows 836 i>= k: jl(i) is the first row to be added to row i 837 i< k: jl(i) is the row following row i in some list of rows 838 jl(i) = mbs indicates the end of a list 839 il(i): points to the first nonzero element in columns k,...,mbs-1 of 840 row i of U */ 841 ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 842 ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); 843 ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 844 jl = il + mbs; 845 for (i=0; i<mbs; i++) { 846 jl[i] = mbs; il[0] = 0; 847 } 848 ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr); 849 uik = dk + 4; 850 ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr); 851 852 /* check permutation */ 853 if (!a->permute){ 854 ai = a->i; aj = a->j; aa = a->a; 855 } else { 856 ai = a->inew; aj = a->jnew; 857 ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 858 ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 859 ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr); 860 ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr); 861 862 for (i=0; i<mbs; i++){ 863 jmin = ai[i]; jmax = ai[i+1]; 864 for (j=jmin; j<jmax; j++){ 865 while (a2anew[j] != j){ 866 k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 867 for (k1=0; k1<4; k1++){ 868 dk[k1] = aa[k*4+k1]; 869 aa[k*4+k1] = aa[j*4+k1]; 870 aa[j*4+k1] = dk[k1]; 871 } 872 } 873 /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 874 if (i > aj[j]){ 875 /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 876 ap = aa + j*4; /* ptr to the beginning of the block */ 877 dk[1] = ap[1]; /* swap ap[1] and ap[2] */ 878 ap[1] = ap[2]; 879 ap[2] = dk[1]; 880 } 881 } 882 } 883 ierr = PetscFree(a2anew);CHKERRQ(ierr); 884 } 885 886 /* for each row k */ 887 for (k = 0; k<mbs; k++){ 888 889 /*initialize k-th row with elements nonzero in row perm(k) of A */ 890 jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 891 ap = aa + jmin*4; 892 for (j = jmin; j < jmax; j++){ 893 vj = perm_ptr[aj[j]]; /* block col. index */ 894 rtmp_ptr = rtmp + vj*4; 895 for (i=0; i<4; i++) *rtmp_ptr++ = *ap++; 896 } 897 898 /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 899 ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr); 900 i = jl[k]; /* first row to be added to k_th row */ 901 902 while (i < k){ 903 nexti = jl[i]; /* next row to be added to k_th row */ 904 905 /* compute multiplier */ 906 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 907 908 /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */ 909 diag = ba + i*4; 910 u = ba + ili*4; 911 uik[0] = -(diag[0]*u[0] + diag[2]*u[1]); 912 uik[1] = -(diag[1]*u[0] + diag[3]*u[1]); 913 uik[2] = -(diag[0]*u[2] + diag[2]*u[3]); 914 uik[3] = -(diag[1]*u[2] + diag[3]*u[3]); 915 916 /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */ 917 dk[0] += uik[0]*u[0] + uik[1]*u[1]; 918 dk[1] += uik[2]*u[0] + uik[3]*u[1]; 919 dk[2] += uik[0]*u[2] + uik[1]*u[3]; 920 dk[3] += uik[2]*u[2] + uik[3]*u[3]; 921 922 ierr = PetscLogFlops(16*2);CHKERRQ(ierr); 923 924 /* update -U(i,k): ba[ili] = uik */ 925 ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr); 926 927 /* add multiple of row i to k-th row ... */ 928 jmin = ili + 1; jmax = bi[i+1]; 929 if (jmin < jmax){ 930 for (j=jmin; j<jmax; j++) { 931 /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */ 932 rtmp_ptr = rtmp + bj[j]*4; 933 u = ba + j*4; 934 rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1]; 935 rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1]; 936 rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3]; 937 rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3]; 938 } 939 ierr = PetscLogFlops(16*(jmax-jmin));CHKERRQ(ierr); 940 941 /* ... add i to row list for next nonzero entry */ 942 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 943 j = bj[jmin]; 944 jl[i] = jl[j]; jl[j] = i; /* update jl */ 945 } 946 i = nexti; 947 } 948 949 /* save nonzero entries in k-th row of U ... */ 950 951 /* invert diagonal block */ 952 diag = ba+k*4; 953 ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr); 954 ierr = Kernel_A_gets_inverse_A_2(diag,shift);CHKERRQ(ierr); 955 956 jmin = bi[k]; jmax = bi[k+1]; 957 if (jmin < jmax) { 958 for (j=jmin; j<jmax; j++){ 959 vj = bj[j]; /* block col. index of U */ 960 u = ba + j*4; 961 rtmp_ptr = rtmp + vj*4; 962 for (k1=0; k1<4; k1++){ 963 *u++ = *rtmp_ptr; 964 *rtmp_ptr++ = 0.0; 965 } 966 } 967 968 /* ... add k to row list for first nonzero entry in k-th row */ 969 il[k] = jmin; 970 i = bj[jmin]; 971 jl[k] = jl[i]; jl[i] = k; 972 } 973 } 974 975 ierr = PetscFree(rtmp);CHKERRQ(ierr); 976 ierr = PetscFree(il);CHKERRQ(ierr); 977 ierr = PetscFree(dk);CHKERRQ(ierr); 978 if (a->permute) { 979 ierr = PetscFree(aa);CHKERRQ(ierr); 980 } 981 ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 982 C->factor = MAT_FACTOR_CHOLESKY; 983 C->assembled = PETSC_TRUE; 984 C->preallocated = PETSC_TRUE; 985 ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 986 PetscFunctionReturn(0); 987 } 988 989 /* 990 Version for when blocks are 2 by 2 Using natural ordering 991 */ 992 #undef __FUNCT__ 993 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering" 994 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B) 995 { 996 Mat C = *B; 997 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 998 PetscErrorCode ierr; 999 PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 1000 PetscInt *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 1001 MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 1002 MatScalar *u,*diag,*rtmp,*rtmp_ptr; 1003 PetscReal shift = info->shiftinblocks; 1004 1005 PetscFunctionBegin; 1006 /* initialization */ 1007 /* il and jl record the first nonzero element in each row of the accessing 1008 window U(0:k, k:mbs-1). 1009 jl: list of rows to be added to uneliminated rows 1010 i>= k: jl(i) is the first row to be added to row i 1011 i< k: jl(i) is the row following row i in some list of rows 1012 jl(i) = mbs indicates the end of a list 1013 il(i): points to the first nonzero element in columns k,...,mbs-1 of 1014 row i of U */ 1015 ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 1016 ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); 1017 ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 1018 jl = il + mbs; 1019 for (i=0; i<mbs; i++) { 1020 jl[i] = mbs; il[0] = 0; 1021 } 1022 ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr); 1023 uik = dk + 4; 1024 1025 ai = a->i; aj = a->j; aa = a->a; 1026 1027 /* for each row k */ 1028 for (k = 0; k<mbs; k++){ 1029 1030 /*initialize k-th row with elements nonzero in row k of A */ 1031 jmin = ai[k]; jmax = ai[k+1]; 1032 ap = aa + jmin*4; 1033 for (j = jmin; j < jmax; j++){ 1034 vj = aj[j]; /* block col. index */ 1035 rtmp_ptr = rtmp + vj*4; 1036 for (i=0; i<4; i++) *rtmp_ptr++ = *ap++; 1037 } 1038 1039 /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 1040 ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr); 1041 i = jl[k]; /* first row to be added to k_th row */ 1042 1043 while (i < k){ 1044 nexti = jl[i]; /* next row to be added to k_th row */ 1045 1046 /* compute multiplier */ 1047 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1048 1049 /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */ 1050 diag = ba + i*4; 1051 u = ba + ili*4; 1052 uik[0] = -(diag[0]*u[0] + diag[2]*u[1]); 1053 uik[1] = -(diag[1]*u[0] + diag[3]*u[1]); 1054 uik[2] = -(diag[0]*u[2] + diag[2]*u[3]); 1055 uik[3] = -(diag[1]*u[2] + diag[3]*u[3]); 1056 1057 /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */ 1058 dk[0] += uik[0]*u[0] + uik[1]*u[1]; 1059 dk[1] += uik[2]*u[0] + uik[3]*u[1]; 1060 dk[2] += uik[0]*u[2] + uik[1]*u[3]; 1061 dk[3] += uik[2]*u[2] + uik[3]*u[3]; 1062 1063 ierr = PetscLogFlops(16*2);CHKERRQ(ierr); 1064 1065 /* update -U(i,k): ba[ili] = uik */ 1066 ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr); 1067 1068 /* add multiple of row i to k-th row ... */ 1069 jmin = ili + 1; jmax = bi[i+1]; 1070 if (jmin < jmax){ 1071 for (j=jmin; j<jmax; j++) { 1072 /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */ 1073 rtmp_ptr = rtmp + bj[j]*4; 1074 u = ba + j*4; 1075 rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1]; 1076 rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1]; 1077 rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3]; 1078 rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3]; 1079 } 1080 ierr = PetscLogFlops(16*(jmax-jmin));CHKERRQ(ierr); 1081 1082 /* ... add i to row list for next nonzero entry */ 1083 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 1084 j = bj[jmin]; 1085 jl[i] = jl[j]; jl[j] = i; /* update jl */ 1086 } 1087 i = nexti; 1088 } 1089 1090 /* save nonzero entries in k-th row of U ... */ 1091 1092 /* invert diagonal block */ 1093 diag = ba+k*4; 1094 ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr); 1095 ierr = Kernel_A_gets_inverse_A_2(diag,shift);CHKERRQ(ierr); 1096 1097 jmin = bi[k]; jmax = bi[k+1]; 1098 if (jmin < jmax) { 1099 for (j=jmin; j<jmax; j++){ 1100 vj = bj[j]; /* block col. index of U */ 1101 u = ba + j*4; 1102 rtmp_ptr = rtmp + vj*4; 1103 for (k1=0; k1<4; k1++){ 1104 *u++ = *rtmp_ptr; 1105 *rtmp_ptr++ = 0.0; 1106 } 1107 } 1108 1109 /* ... add k to row list for first nonzero entry in k-th row */ 1110 il[k] = jmin; 1111 i = bj[jmin]; 1112 jl[k] = jl[i]; jl[i] = k; 1113 } 1114 } 1115 1116 ierr = PetscFree(rtmp);CHKERRQ(ierr); 1117 ierr = PetscFree(il);CHKERRQ(ierr); 1118 ierr = PetscFree(dk);CHKERRQ(ierr); 1119 1120 C->factor = MAT_FACTOR_CHOLESKY; 1121 C->assembled = PETSC_TRUE; 1122 C->preallocated = PETSC_TRUE; 1123 ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 1124 PetscFunctionReturn(0); 1125 } 1126 1127 /* 1128 Numeric U^T*D*U factorization for SBAIJ format. 1129 Version for blocks are 1 by 1. 1130 */ 1131 #undef __FUNCT__ 1132 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1" 1133 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat A,MatFactorInfo *info,Mat *B) 1134 { 1135 Mat C = *B; 1136 Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data; 1137 IS ip=b->row; 1138 PetscErrorCode ierr; 1139 PetscInt *rip,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol; 1140 PetscInt *ai,*aj,*a2anew; 1141 PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 1142 MatScalar *rtmp,*ba=b->a,*bval,*aa,dk,uikdi; 1143 PetscReal zeropivot,rs,shiftnz; 1144 PetscReal shiftpd; 1145 ChShift_Ctx sctx; 1146 PetscInt newshift; 1147 1148 PetscFunctionBegin; 1149 /* initialization */ 1150 shiftnz = info->shiftnz; 1151 shiftpd = info->shiftpd; 1152 zeropivot = info->zeropivot; 1153 1154 ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 1155 if (!a->permute){ 1156 ai = a->i; aj = a->j; aa = a->a; 1157 } else { 1158 ai = a->inew; aj = a->jnew; 1159 nz = ai[mbs]; 1160 ierr = PetscMalloc(nz*sizeof(MatScalar),&aa);CHKERRQ(ierr); 1161 a2anew = a->a2anew; 1162 bval = a->a; 1163 for (j=0; j<nz; j++){ 1164 aa[a2anew[j]] = *(bval++); 1165 } 1166 } 1167 1168 /* initialization */ 1169 /* il and jl record the first nonzero element in each row of the accessing 1170 window U(0:k, k:mbs-1). 1171 jl: list of rows to be added to uneliminated rows 1172 i>= k: jl(i) is the first row to be added to row i 1173 i< k: jl(i) is the row following row i in some list of rows 1174 jl(i) = mbs indicates the end of a list 1175 il(i): points to the first nonzero element in columns k,...,mbs-1 of 1176 row i of U */ 1177 nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 1178 ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 1179 jl = il + mbs; 1180 rtmp = (MatScalar*)(jl + mbs); 1181 1182 sctx.shift_amount = 0; 1183 sctx.nshift = 0; 1184 do { 1185 sctx.chshift = PETSC_FALSE; 1186 for (i=0; i<mbs; i++) { 1187 rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1188 } 1189 1190 for (k = 0; k<mbs; k++){ 1191 /*initialize k-th row by the perm[k]-th row of A */ 1192 jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 1193 bval = ba + bi[k]; 1194 for (j = jmin; j < jmax; j++){ 1195 col = rip[aj[j]]; 1196 rtmp[col] = aa[j]; 1197 *bval++ = 0.0; /* for in-place factorization */ 1198 } 1199 1200 /* shift the diagonal of the matrix */ 1201 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 1202 1203 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 1204 dk = rtmp[k]; 1205 i = jl[k]; /* first row to be added to k_th row */ 1206 1207 while (i < k){ 1208 nexti = jl[i]; /* next row to be added to k_th row */ 1209 1210 /* compute multiplier, update diag(k) and U(i,k) */ 1211 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1212 uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 1213 dk += uikdi*ba[ili]; 1214 ba[ili] = uikdi; /* -U(i,k) */ 1215 1216 /* add multiple of row i to k-th row */ 1217 jmin = ili + 1; jmax = bi[i+1]; 1218 if (jmin < jmax){ 1219 for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 1220 ierr = PetscLogFlops(2*(jmax-jmin));CHKERRQ(ierr); 1221 1222 /* update il and jl for row i */ 1223 il[i] = jmin; 1224 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 1225 } 1226 i = nexti; 1227 } 1228 1229 /* shift the diagonals when zero pivot is detected */ 1230 /* compute rs=sum of abs(off-diagonal) */ 1231 rs = 0.0; 1232 jmin = bi[k]+1; 1233 nz = bi[k+1] - jmin; 1234 if (nz){ 1235 bcol = bj + jmin; 1236 while (nz--){ 1237 rs += PetscAbsScalar(rtmp[*bcol]); 1238 bcol++; 1239 } 1240 } 1241 1242 sctx.rs = rs; 1243 sctx.pv = dk; 1244 ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 1245 if (newshift == 1) break; /* sctx.shift_amount is updated */ 1246 1247 /* copy data into U(k,:) */ 1248 ba[bi[k]] = 1.0/dk; /* U(k,k) */ 1249 jmin = bi[k]+1; jmax = bi[k+1]; 1250 if (jmin < jmax) { 1251 for (j=jmin; j<jmax; j++){ 1252 col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 1253 } 1254 /* add the k-th row into il and jl */ 1255 il[k] = jmin; 1256 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 1257 } 1258 } 1259 } while (sctx.chshift); 1260 ierr = PetscFree(il);CHKERRQ(ierr); 1261 if (a->permute){ierr = PetscFree(aa);CHKERRQ(ierr);} 1262 1263 ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 1264 C->factor = MAT_FACTOR_CHOLESKY; 1265 C->assembled = PETSC_TRUE; 1266 C->preallocated = PETSC_TRUE; 1267 ierr = PetscLogFlops(C->rmap.N);CHKERRQ(ierr); 1268 if (sctx.nshift){ 1269 if (shiftnz) { 1270 ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 1271 } else if (shiftpd) { 1272 ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 1273 } 1274 } 1275 PetscFunctionReturn(0); 1276 } 1277 1278 /* 1279 Version for when blocks are 1 by 1 Using natural ordering 1280 */ 1281 #undef __FUNCT__ 1282 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering" 1283 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B) 1284 { 1285 Mat C = *B; 1286 Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data; 1287 PetscErrorCode ierr; 1288 PetscInt i,j,mbs = a->mbs; 1289 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 1290 PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz; 1291 MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 1292 PetscReal zeropivot,rs,shiftnz; 1293 PetscReal shiftpd; 1294 ChShift_Ctx sctx; 1295 PetscInt newshift; 1296 1297 PetscFunctionBegin; 1298 /* initialization */ 1299 shiftnz = info->shiftnz; 1300 shiftpd = info->shiftpd; 1301 zeropivot = info->zeropivot; 1302 1303 /* il and jl record the first nonzero element in each row of the accessing 1304 window U(0:k, k:mbs-1). 1305 jl: list of rows to be added to uneliminated rows 1306 i>= k: jl(i) is the first row to be added to row i 1307 i< k: jl(i) is the row following row i in some list of rows 1308 jl(i) = mbs indicates the end of a list 1309 il(i): points to the first nonzero element in U(i,k:mbs-1) 1310 */ 1311 ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 1312 ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 1313 jl = il + mbs; 1314 1315 sctx.shift_amount = 0; 1316 sctx.nshift = 0; 1317 do { 1318 sctx.chshift = PETSC_FALSE; 1319 for (i=0; i<mbs; i++) { 1320 rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1321 } 1322 1323 for (k = 0; k<mbs; k++){ 1324 /*initialize k-th row with elements nonzero in row perm(k) of A */ 1325 nz = ai[k+1] - ai[k]; 1326 acol = aj + ai[k]; 1327 aval = aa + ai[k]; 1328 bval = ba + bi[k]; 1329 while (nz -- ){ 1330 rtmp[*acol++] = *aval++; 1331 *bval++ = 0.0; /* for in-place factorization */ 1332 } 1333 1334 /* shift the diagonal of the matrix */ 1335 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 1336 1337 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 1338 dk = rtmp[k]; 1339 i = jl[k]; /* first row to be added to k_th row */ 1340 1341 while (i < k){ 1342 nexti = jl[i]; /* next row to be added to k_th row */ 1343 /* compute multiplier, update D(k) and U(i,k) */ 1344 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1345 uikdi = - ba[ili]*ba[bi[i]]; 1346 dk += uikdi*ba[ili]; 1347 ba[ili] = uikdi; /* -U(i,k) */ 1348 1349 /* add multiple of row i to k-th row ... */ 1350 jmin = ili + 1; 1351 nz = bi[i+1] - jmin; 1352 if (nz > 0){ 1353 bcol = bj + jmin; 1354 bval = ba + jmin; 1355 ierr = PetscLogFlops(2*nz);CHKERRQ(ierr); 1356 while (nz --) rtmp[*bcol++] += uikdi*(*bval++); 1357 1358 /* update il and jl for i-th row */ 1359 il[i] = jmin; 1360 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 1361 } 1362 i = nexti; 1363 } 1364 1365 /* shift the diagonals when zero pivot is detected */ 1366 /* compute rs=sum of abs(off-diagonal) */ 1367 rs = 0.0; 1368 jmin = bi[k]+1; 1369 nz = bi[k+1] - jmin; 1370 if (nz){ 1371 bcol = bj + jmin; 1372 while (nz--){ 1373 rs += PetscAbsScalar(rtmp[*bcol]); 1374 bcol++; 1375 } 1376 } 1377 1378 sctx.rs = rs; 1379 sctx.pv = dk; 1380 ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 1381 if (newshift == 1) break; /* sctx.shift_amount is updated */ 1382 1383 /* copy data into U(k,:) */ 1384 ba[bi[k]] = 1.0/dk; 1385 jmin = bi[k]+1; 1386 nz = bi[k+1] - jmin; 1387 if (nz){ 1388 bcol = bj + jmin; 1389 bval = ba + jmin; 1390 while (nz--){ 1391 *bval++ = rtmp[*bcol]; 1392 rtmp[*bcol++] = 0.0; 1393 } 1394 /* add k-th row into il and jl */ 1395 il[k] = jmin; 1396 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 1397 } 1398 } /* end of for (k = 0; k<mbs; k++) */ 1399 } while (sctx.chshift); 1400 ierr = PetscFree(rtmp);CHKERRQ(ierr); 1401 ierr = PetscFree(il);CHKERRQ(ierr); 1402 1403 C->factor = MAT_FACTOR_CHOLESKY; 1404 C->assembled = PETSC_TRUE; 1405 C->preallocated = PETSC_TRUE; 1406 ierr = PetscLogFlops(C->rmap.N);CHKERRQ(ierr); 1407 if (sctx.nshift){ 1408 if (shiftnz) { 1409 ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 1410 } else if (shiftpd) { 1411 ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 1412 } 1413 } 1414 PetscFunctionReturn(0); 1415 } 1416 1417 #undef __FUNCT__ 1418 #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ" 1419 PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info) 1420 { 1421 PetscErrorCode ierr; 1422 Mat C; 1423 1424 PetscFunctionBegin; 1425 ierr = MatCholeskyFactorSymbolic(A,perm,info,&C);CHKERRQ(ierr); 1426 ierr = MatCholeskyFactorNumeric(A,info,&C);CHKERRQ(ierr); 1427 ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 1428 PetscFunctionReturn(0); 1429 } 1430 1431 1432