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