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