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