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