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