1 2 #include <../src/mat/impls/baij/seq/baij.h> 3 #include <../src/mat/impls/sbaij/seq/sbaij.h> 4 #include <petsc/private/kernels/blockinvert.h> 5 #include <petscis.h> 6 7 /* 8 input: 9 F -- numeric factor 10 output: 11 nneg, nzero, npos: matrix inertia 12 */ 13 14 #undef __FUNCT__ 15 #define __FUNCT__ "MatGetInertia_SeqSBAIJ" 16 PetscErrorCode MatGetInertia_SeqSBAIJ(Mat F,PetscInt *nneig,PetscInt *nzero,PetscInt *npos) 17 { 18 Mat_SeqSBAIJ *fact_ptr = (Mat_SeqSBAIJ*)F->data; 19 MatScalar *dd = fact_ptr->a; 20 PetscInt mbs =fact_ptr->mbs,bs=F->rmap->bs,i,nneig_tmp,npos_tmp,*fi = fact_ptr->diag; 21 22 PetscFunctionBegin; 23 if (bs != 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"No support for bs: %D >1 yet",bs); 24 nneig_tmp = 0; npos_tmp = 0; 25 for (i=0; i<mbs; i++) { 26 if (PetscRealPart(dd[*fi]) > 0.0) npos_tmp++; 27 else if (PetscRealPart(dd[*fi]) < 0.0) nneig_tmp++; 28 fi++; 29 } 30 if (nneig) *nneig = nneig_tmp; 31 if (npos) *npos = npos_tmp; 32 if (nzero) *nzero = mbs - nneig_tmp - npos_tmp; 33 PetscFunctionReturn(0); 34 } 35 36 /* 37 Symbolic U^T*D*U factorization for SBAIJ format. Modified from SSF of YSMP. 38 Use Modified Sparse Row (MSR) storage for u and ju. See page 85, "Iterative Methods ..." by Saad. 39 */ 40 #undef __FUNCT__ 41 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_MSR" 42 PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(Mat F,Mat A,IS perm,const MatFactorInfo *info) 43 { 44 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b; 45 PetscErrorCode ierr; 46 const PetscInt *rip,*ai,*aj; 47 PetscInt i,mbs = a->mbs,*jutmp,bs = A->rmap->bs,bs2=a->bs2; 48 PetscInt m,reallocs = 0,prow; 49 PetscInt *jl,*q,jmin,jmax,juidx,nzk,qm,*iu,*ju,k,j,vj,umax,maxadd; 50 PetscReal f = info->fill; 51 PetscBool perm_identity; 52 53 PetscFunctionBegin; 54 /* check whether perm is the identity mapping */ 55 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 56 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 57 58 if (perm_identity) { /* without permutation */ 59 a->permute = PETSC_FALSE; 60 61 ai = a->i; aj = a->j; 62 } else { /* non-trivial permutation */ 63 a->permute = PETSC_TRUE; 64 65 ierr = MatReorderingSeqSBAIJ(A,perm);CHKERRQ(ierr); 66 67 ai = a->inew; aj = a->jnew; 68 } 69 70 /* initialization */ 71 ierr = PetscMalloc1(mbs+1,&iu);CHKERRQ(ierr); 72 umax = (PetscInt)(f*ai[mbs] + 1); umax += mbs + 1; 73 ierr = PetscMalloc1(umax,&ju);CHKERRQ(ierr); 74 iu[0] = mbs+1; 75 juidx = mbs + 1; /* index for ju */ 76 /* jl linked list for pivot row -- linked list for col index */ 77 ierr = PetscMalloc2(mbs,&jl,mbs,&q);CHKERRQ(ierr); 78 for (i=0; i<mbs; i++) { 79 jl[i] = mbs; 80 q[i] = 0; 81 } 82 83 /* for each row k */ 84 for (k=0; k<mbs; k++) { 85 for (i=0; i<mbs; i++) q[i] = 0; /* to be removed! */ 86 nzk = 0; /* num. of nz blocks in k-th block row with diagonal block excluded */ 87 q[k] = mbs; 88 /* initialize nonzero structure of k-th row to row rip[k] of A */ 89 jmin = ai[rip[k]] +1; /* exclude diag[k] */ 90 jmax = ai[rip[k]+1]; 91 for (j=jmin; j<jmax; j++) { 92 vj = rip[aj[j]]; /* col. value */ 93 if (vj > k) { 94 qm = k; 95 do { 96 m = qm; qm = q[m]; 97 } while (qm < vj); 98 if (qm == vj) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Duplicate entry in A\n"); 99 nzk++; 100 q[m] = vj; 101 q[vj] = qm; 102 } /* if (vj > k) */ 103 } /* for (j=jmin; j<jmax; j++) */ 104 105 /* modify nonzero structure of k-th row by computing fill-in 106 for each row i to be merged in */ 107 prow = k; 108 prow = jl[prow]; /* next pivot row (== mbs for symbolic factorization) */ 109 110 while (prow < k) { 111 /* merge row prow into k-th row */ 112 jmin = iu[prow] + 1; jmax = iu[prow+1]; 113 qm = k; 114 for (j=jmin; j<jmax; j++) { 115 vj = ju[j]; 116 do { 117 m = qm; qm = q[m]; 118 } while (qm < vj); 119 if (qm != vj) { 120 nzk++; q[m] = vj; q[vj] = qm; qm = vj; 121 } 122 } 123 prow = jl[prow]; /* next pivot row */ 124 } 125 126 /* add k to row list for first nonzero element in k-th row */ 127 if (nzk > 0) { 128 i = q[k]; /* col value of first nonzero element in U(k, k+1:mbs-1) */ 129 jl[k] = jl[i]; jl[i] = k; 130 } 131 iu[k+1] = iu[k] + nzk; 132 133 /* allocate more space to ju if needed */ 134 if (iu[k+1] > umax) { 135 /* estimate how much additional space we will need */ 136 /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 137 /* just double the memory each time */ 138 maxadd = umax; 139 if (maxadd < nzk) maxadd = (mbs-k)*(nzk+1)/2; 140 umax += maxadd; 141 142 /* allocate a longer ju */ 143 ierr = PetscMalloc1(umax,&jutmp);CHKERRQ(ierr); 144 ierr = PetscMemcpy(jutmp,ju,iu[k]*sizeof(PetscInt));CHKERRQ(ierr); 145 ierr = PetscFree(ju);CHKERRQ(ierr); 146 ju = jutmp; 147 reallocs++; /* count how many times we realloc */ 148 } 149 150 /* save nonzero structure of k-th row in ju */ 151 i=k; 152 while (nzk--) { 153 i = q[i]; 154 ju[juidx++] = i; 155 } 156 } 157 158 #if defined(PETSC_USE_INFO) 159 if (ai[mbs] != 0) { 160 PetscReal af = ((PetscReal)iu[mbs])/((PetscReal)ai[mbs]); 161 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);CHKERRQ(ierr); 162 ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr); 163 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);CHKERRQ(ierr); 164 ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 165 } else { 166 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 167 } 168 #endif 169 170 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 171 ierr = PetscFree2(jl,q);CHKERRQ(ierr); 172 173 /* put together the new matrix */ 174 ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(F,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 175 176 /* ierr = PetscLogObjectParent((PetscObject)B,(PetscObject)iperm);CHKERRQ(ierr); */ 177 b = (Mat_SeqSBAIJ*)(F)->data; 178 b->singlemalloc = PETSC_FALSE; 179 b->free_a = PETSC_TRUE; 180 b->free_ij = PETSC_TRUE; 181 182 ierr = PetscMalloc1((iu[mbs]+1)*bs2,&b->a);CHKERRQ(ierr); 183 b->j = ju; 184 b->i = iu; 185 b->diag = 0; 186 b->ilen = 0; 187 b->imax = 0; 188 b->row = perm; 189 190 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 191 192 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 193 194 b->icol = perm; 195 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 196 ierr = PetscMalloc1(bs*mbs+bs,&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((PetscObject)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 PetscBool 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; 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=NULL,current_space=NULL; 232 PetscBT lnkbt; 233 234 PetscFunctionBegin; 235 if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n); 236 ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr); 237 if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i); 238 if (bs > 1) { 239 ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(fact,A,perm,info);CHKERRQ(ierr); 240 PetscFunctionReturn(0); 241 } 242 243 /* check whether perm is the identity mapping */ 244 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 245 if (!perm_identity) SETERRQ(PETSC_COMM_SELF,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 = PetscMalloc1(mbs+1,&ui);CHKERRQ(ierr); 251 ierr = PetscMalloc1(mbs+1,&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,&ui_ptr,mbs,&il,mbs,&jl,mbs,&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(PetscRealIntMultTruncate(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_COMM_SELF,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 = PetscIntMultTruncate(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 321 current_space->array += nzk; 322 current_space->local_used += nzk; 323 current_space->local_remaining -= nzk; 324 325 ui[k+1] = ui[k] + nzk; 326 } 327 328 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 329 ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr); 330 331 /* destroy list of free space and other temporary array(s) */ 332 ierr = PetscMalloc1(ui[mbs]+1,&uj);CHKERRQ(ierr); 333 ierr = PetscFreeSpaceContiguous_Cholesky(&free_space,uj,mbs,ui,udiag);CHKERRQ(ierr); /* store matrix factor */ 334 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 335 336 /* put together the new matrix in MATSEQSBAIJ format */ 337 ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 338 339 b = (Mat_SeqSBAIJ*)fact->data; 340 b->singlemalloc = PETSC_FALSE; 341 b->free_a = PETSC_TRUE; 342 b->free_ij = PETSC_TRUE; 343 344 ierr = PetscMalloc1(ui[mbs]+1,&b->a);CHKERRQ(ierr); 345 346 b->j = uj; 347 b->i = ui; 348 b->diag = udiag; 349 b->free_diag = PETSC_TRUE; 350 b->ilen = 0; 351 b->imax = 0; 352 b->row = perm; 353 b->icol = perm; 354 355 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 356 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 357 358 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 359 360 ierr = PetscMalloc1(mbs+1,&b->solve_work);CHKERRQ(ierr); 361 ierr = PetscLogObjectMemory((PetscObject)fact,ui[mbs]*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 362 363 b->maxnz = b->nz = ui[mbs]; 364 365 fact->info.factor_mallocs = reallocs; 366 fact->info.fill_ratio_given = fill; 367 if (ai[mbs] != 0) { 368 fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs]; 369 } else { 370 fact->info.fill_ratio_needed = 0.0; 371 } 372 #if defined(PETSC_USE_INFO) 373 if (ai[mbs] != 0) { 374 PetscReal af = fact->info.fill_ratio_needed; 375 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr); 376 ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr); 377 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr); 378 } else { 379 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 380 } 381 #endif 382 fact->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering; 383 PetscFunctionReturn(0); 384 } 385 386 #undef __FUNCT__ 387 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqSBAIJ_inplace" 388 PetscErrorCode MatCholeskyFactorSymbolic_SeqSBAIJ_inplace(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 389 { 390 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data; 391 Mat_SeqSBAIJ *b; 392 PetscErrorCode ierr; 393 PetscBool perm_identity,missing; 394 PetscReal fill = info->fill; 395 const PetscInt *rip,*ai,*aj; 396 PetscInt i,mbs=a->mbs,bs=A->rmap->bs,reallocs=0,prow,d; 397 PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 398 PetscInt nlnk,*lnk,ncols,*cols,*uj,**ui_ptr,*uj_ptr; 399 PetscFreeSpaceList free_space=NULL,current_space=NULL; 400 PetscBT lnkbt; 401 402 PetscFunctionBegin; 403 ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 404 if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 405 406 /* 407 This code originally uses Modified Sparse Row (MSR) storage 408 (see page 85, "Iterative Methods ..." by Saad) for the output matrix B - bad choise! 409 Then it is rewritten so the factor B takes seqsbaij format. However the associated 410 MatCholeskyFactorNumeric_() have not been modified for the cases of bs>1 or !perm_identity, 411 thus the original code in MSR format is still used for these cases. 412 The code below should replace MatCholeskyFactorSymbolic_SeqSBAIJ_MSR() whenever 413 MatCholeskyFactorNumeric_() is modified for using sbaij symbolic factor. 414 */ 415 if (bs > 1) { 416 ierr = MatCholeskyFactorSymbolic_SeqSBAIJ_MSR(fact,A,perm,info);CHKERRQ(ierr); 417 PetscFunctionReturn(0); 418 } 419 420 /* check whether perm is the identity mapping */ 421 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 422 423 if (perm_identity) { 424 a->permute = PETSC_FALSE; 425 426 ai = a->i; aj = a->j; 427 } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Matrix reordering is not supported for sbaij matrix. Use aij format"); 428 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 429 430 /* initialization */ 431 ierr = PetscMalloc1(mbs+1,&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,&ui_ptr,mbs,&il,mbs,&jl,mbs,&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(PetscRealIntMultTruncate(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 477 j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 478 } 479 prow = nextprow; 480 } 481 482 /* if free space is not available, make more free space */ 483 if (current_space->local_remaining<nzk) { 484 i = mbs - k + 1; /* num of unfactored rows */ 485 i = PetscMin(PetscIntMultTruncate(i,nzk), PetscIntMultTruncate(i,i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 486 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 487 reallocs++; 488 } 489 490 /* copy data into free space, then initialize lnk */ 491 ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 492 493 /* add the k-th row into il and jl */ 494 if (nzk-1 > 0) { 495 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 496 jl[k] = jl[i]; jl[i] = k; 497 il[k] = ui[k] + 1; 498 } 499 ui_ptr[k] = current_space->array; 500 501 current_space->array += nzk; 502 current_space->local_used += nzk; 503 current_space->local_remaining -= nzk; 504 505 ui[k+1] = ui[k] + nzk; 506 } 507 508 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 509 ierr = PetscFree4(ui_ptr,il,jl,cols);CHKERRQ(ierr); 510 511 /* destroy list of free space and other temporary array(s) */ 512 ierr = PetscMalloc1(ui[mbs]+1,&uj);CHKERRQ(ierr); 513 ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 514 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 515 516 /* put together the new matrix in MATSEQSBAIJ format */ 517 ierr = MatSeqSBAIJSetPreallocation_SeqSBAIJ(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 518 519 b = (Mat_SeqSBAIJ*)fact->data; 520 b->singlemalloc = PETSC_FALSE; 521 b->free_a = PETSC_TRUE; 522 b->free_ij = PETSC_TRUE; 523 524 ierr = PetscMalloc1(ui[mbs]+1,&b->a);CHKERRQ(ierr); 525 526 b->j = uj; 527 b->i = ui; 528 b->diag = 0; 529 b->ilen = 0; 530 b->imax = 0; 531 b->row = perm; 532 533 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 534 535 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 536 b->icol = perm; 537 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 538 ierr = PetscMalloc1(mbs+1,&b->solve_work);CHKERRQ(ierr); 539 ierr = PetscLogObjectMemory((PetscObject)fact,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 540 b->maxnz = b->nz = ui[mbs]; 541 542 fact->info.factor_mallocs = reallocs; 543 fact->info.fill_ratio_given = fill; 544 if (ai[mbs] != 0) { 545 fact->info.fill_ratio_needed = ((PetscReal)ui[mbs])/ai[mbs]; 546 } else { 547 fact->info.fill_ratio_needed = 0.0; 548 } 549 #if defined(PETSC_USE_INFO) 550 if (ai[mbs] != 0) { 551 PetscReal af = fact->info.fill_ratio_needed; 552 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)fill,(double)af);CHKERRQ(ierr); 553 ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr); 554 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g) for best performance.\n",(double)af);CHKERRQ(ierr); 555 } else { 556 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 557 } 558 #endif 559 ierr = MatSeqSBAIJSetNumericFactorization_inplace(fact,perm_identity);CHKERRQ(ierr); 560 PetscFunctionReturn(0); 561 } 562 563 #undef __FUNCT__ 564 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N" 565 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N(Mat C,Mat A,const MatFactorInfo *info) 566 { 567 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data; 568 IS perm = b->row; 569 PetscErrorCode ierr; 570 const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j; 571 PetscInt i,j; 572 PetscInt *a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 573 PetscInt bs =A->rmap->bs,bs2 = a->bs2; 574 MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 575 MatScalar *u,*diag,*rtmp,*rtmp_ptr; 576 MatScalar *work; 577 PetscInt *pivots; 578 579 PetscFunctionBegin; 580 /* initialization */ 581 ierr = PetscCalloc1(bs2*mbs,&rtmp);CHKERRQ(ierr); 582 ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr); 583 for (i=0; i<mbs; i++) { 584 jl[i] = mbs; il[0] = 0; 585 } 586 ierr = PetscMalloc3(bs2,&dk,bs2,&uik,bs,&work);CHKERRQ(ierr); 587 ierr = PetscMalloc1(bs,&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 = PetscMalloc1(bs2*ai[mbs],&aa);CHKERRQ(ierr); 597 ierr = PetscMemcpy(aa,a->a,bs2*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 598 ierr = PetscMalloc1(ai[mbs],&a2anew);CHKERRQ(ierr); 599 ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr); 600 601 for (i=0; i<mbs; i++) { 602 jmin = ai[i]; jmax = ai[i+1]; 603 for (j=jmin; j<jmax; j++) { 604 while (a2anew[j] != j) { 605 k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 606 for (k1=0; k1<bs2; k1++) { 607 dk[k1] = aa[k*bs2+k1]; 608 aa[k*bs2+k1] = aa[j*bs2+k1]; 609 aa[j*bs2+k1] = dk[k1]; 610 } 611 } 612 /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 613 if (i > aj[j]) { 614 /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 615 ap = aa + j*bs2; /* ptr to the beginning of j-th block of aa */ 616 for (k=0; k<bs2; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */ 617 for (k=0; k<bs; k++) { /* j-th block of aa <- dk^T */ 618 for (k1=0; k1<bs; k1++) *ap++ = dk[k + bs*k1]; 619 } 620 } 621 } 622 } 623 ierr = PetscFree(a2anew);CHKERRQ(ierr); 624 } 625 626 /* for each row k */ 627 for (k = 0; k<mbs; k++) { 628 629 /*initialize k-th row with elements nonzero in row perm(k) of A */ 630 jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 631 632 ap = aa + jmin*bs2; 633 for (j = jmin; j < jmax; j++) { 634 vj = perm_ptr[aj[j]]; /* block col. index */ 635 rtmp_ptr = rtmp + vj*bs2; 636 for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++; 637 } 638 639 /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 640 ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 641 i = jl[k]; /* first row to be added to k_th row */ 642 643 while (i < k) { 644 nexti = jl[i]; /* next row to be added to k_th row */ 645 646 /* compute multiplier */ 647 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 648 649 /* uik = -inv(Di)*U_bar(i,k) */ 650 diag = ba + i*bs2; 651 u = ba + ili*bs2; 652 ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 653 PetscKernel_A_gets_A_minus_B_times_C(bs,uik,diag,u); 654 655 /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 656 PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u); 657 ierr = PetscLogFlops(4.0*bs*bs2);CHKERRQ(ierr); 658 659 /* update -U(i,k) */ 660 ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 661 662 /* add multiple of row i to k-th row ... */ 663 jmin = ili + 1; jmax = bi[i+1]; 664 if (jmin < jmax) { 665 for (j=jmin; j<jmax; j++) { 666 /* rtmp += -U(i,k)^T * U_bar(i,j) */ 667 rtmp_ptr = rtmp + bj[j]*bs2; 668 u = ba + j*bs2; 669 PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u); 670 } 671 ierr = PetscLogFlops(2.0*bs*bs2*(jmax-jmin));CHKERRQ(ierr); 672 673 /* ... add i to row list for next nonzero entry */ 674 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 675 j = bj[jmin]; 676 jl[i] = jl[j]; jl[j] = i; /* update jl */ 677 } 678 i = nexti; 679 } 680 681 /* save nonzero entries in k-th row of U ... */ 682 683 /* invert diagonal block */ 684 diag = ba+k*bs2; 685 ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr); 686 ierr = PetscKernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr); 687 688 jmin = bi[k]; jmax = bi[k+1]; 689 if (jmin < jmax) { 690 for (j=jmin; j<jmax; j++) { 691 vj = bj[j]; /* block col. index of U */ 692 u = ba + j*bs2; 693 rtmp_ptr = rtmp + vj*bs2; 694 for (k1=0; k1<bs2; k1++) { 695 *u++ = *rtmp_ptr; 696 *rtmp_ptr++ = 0.0; 697 } 698 } 699 700 /* ... add k to row list for first nonzero entry in k-th row */ 701 il[k] = jmin; 702 i = bj[jmin]; 703 jl[k] = jl[i]; jl[i] = k; 704 } 705 } 706 707 ierr = PetscFree(rtmp);CHKERRQ(ierr); 708 ierr = PetscFree2(il,jl);CHKERRQ(ierr); 709 ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr); 710 ierr = PetscFree(pivots);CHKERRQ(ierr); 711 if (a->permute) { 712 ierr = PetscFree(aa);CHKERRQ(ierr); 713 } 714 715 ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 716 717 C->ops->solve = MatSolve_SeqSBAIJ_N_inplace; 718 C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_inplace; 719 C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_N_inplace; 720 C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_N_inplace; 721 722 C->assembled = PETSC_TRUE; 723 C->preallocated = PETSC_TRUE; 724 725 ierr = PetscLogFlops(1.3333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 726 PetscFunctionReturn(0); 727 } 728 729 #undef __FUNCT__ 730 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering" 731 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) 732 { 733 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data; 734 PetscErrorCode ierr; 735 PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 736 PetscInt *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 737 PetscInt bs =A->rmap->bs,bs2 = a->bs2; 738 MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 739 MatScalar *u,*diag,*rtmp,*rtmp_ptr; 740 MatScalar *work; 741 PetscInt *pivots; 742 743 PetscFunctionBegin; 744 ierr = PetscCalloc1(bs2*mbs,&rtmp);CHKERRQ(ierr); 745 ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr); 746 for (i=0; i<mbs; i++) { 747 jl[i] = mbs; il[0] = 0; 748 } 749 ierr = PetscMalloc3(bs2,&dk,bs2,&uik,bs,&work);CHKERRQ(ierr); 750 ierr = PetscMalloc1(bs,&pivots);CHKERRQ(ierr); 751 752 ai = a->i; aj = a->j; aa = a->a; 753 754 /* for each row k */ 755 for (k = 0; k<mbs; k++) { 756 757 /*initialize k-th row with elements nonzero in row k of A */ 758 jmin = ai[k]; jmax = ai[k+1]; 759 ap = aa + jmin*bs2; 760 for (j = jmin; j < jmax; j++) { 761 vj = aj[j]; /* block col. index */ 762 rtmp_ptr = rtmp + vj*bs2; 763 for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++; 764 } 765 766 /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 767 ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 768 i = jl[k]; /* first row to be added to k_th row */ 769 770 while (i < k) { 771 nexti = jl[i]; /* next row to be added to k_th row */ 772 773 /* compute multiplier */ 774 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 775 776 /* uik = -inv(Di)*U_bar(i,k) */ 777 diag = ba + i*bs2; 778 u = ba + ili*bs2; 779 ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 780 PetscKernel_A_gets_A_minus_B_times_C(bs,uik,diag,u); 781 782 /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 783 PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u); 784 ierr = PetscLogFlops(2.0*bs*bs2);CHKERRQ(ierr); 785 786 /* update -U(i,k) */ 787 ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 788 789 /* add multiple of row i to k-th row ... */ 790 jmin = ili + 1; jmax = bi[i+1]; 791 if (jmin < jmax) { 792 for (j=jmin; j<jmax; j++) { 793 /* rtmp += -U(i,k)^T * U_bar(i,j) */ 794 rtmp_ptr = rtmp + bj[j]*bs2; 795 u = ba + j*bs2; 796 PetscKernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u); 797 } 798 ierr = PetscLogFlops(2.0*bs*bs2*(jmax-jmin));CHKERRQ(ierr); 799 800 /* ... add i to row list for next nonzero entry */ 801 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 802 j = bj[jmin]; 803 jl[i] = jl[j]; jl[j] = i; /* update jl */ 804 } 805 i = nexti; 806 } 807 808 /* save nonzero entries in k-th row of U ... */ 809 810 /* invert diagonal block */ 811 diag = ba+k*bs2; 812 ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr); 813 ierr = PetscKernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr); 814 815 jmin = bi[k]; jmax = bi[k+1]; 816 if (jmin < jmax) { 817 for (j=jmin; j<jmax; j++) { 818 vj = bj[j]; /* block col. index of U */ 819 u = ba + j*bs2; 820 rtmp_ptr = rtmp + vj*bs2; 821 for (k1=0; k1<bs2; k1++) { 822 *u++ = *rtmp_ptr; 823 *rtmp_ptr++ = 0.0; 824 } 825 } 826 827 /* ... add k to row list for first nonzero entry in k-th row */ 828 il[k] = jmin; 829 i = bj[jmin]; 830 jl[k] = jl[i]; jl[i] = k; 831 } 832 } 833 834 ierr = PetscFree(rtmp);CHKERRQ(ierr); 835 ierr = PetscFree2(il,jl);CHKERRQ(ierr); 836 ierr = PetscFree3(dk,uik,work);CHKERRQ(ierr); 837 ierr = PetscFree(pivots);CHKERRQ(ierr); 838 839 C->ops->solve = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace; 840 C->ops->solvetranspose = MatSolve_SeqSBAIJ_N_NaturalOrdering_inplace; 841 C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace; 842 C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_N_NaturalOrdering_inplace; 843 C->assembled = PETSC_TRUE; 844 C->preallocated = PETSC_TRUE; 845 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 = PetscCalloc1(4*mbs,&rtmp);CHKERRQ(ierr); 879 ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr); 880 for (i=0; i<mbs; i++) { 881 jl[i] = mbs; il[0] = 0; 882 } 883 ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr); 884 885 /* check permutation */ 886 if (!a->permute) { 887 ai = a->i; aj = a->j; aa = a->a; 888 } else { 889 ai = a->inew; aj = a->jnew; 890 ierr = PetscMalloc1(4*ai[mbs],&aa);CHKERRQ(ierr); 891 ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 892 ierr = PetscMalloc1(ai[mbs],&a2anew);CHKERRQ(ierr); 893 ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr); 894 895 for (i=0; i<mbs; i++) { 896 jmin = ai[i]; jmax = ai[i+1]; 897 for (j=jmin; j<jmax; j++) { 898 while (a2anew[j] != j) { 899 k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 900 for (k1=0; k1<4; k1++) { 901 dk[k1] = aa[k*4+k1]; 902 aa[k*4+k1] = aa[j*4+k1]; 903 aa[j*4+k1] = dk[k1]; 904 } 905 } 906 /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 907 if (i > aj[j]) { 908 /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 909 ap = aa + j*4; /* ptr to the beginning of the block */ 910 dk[1] = ap[1]; /* swap ap[1] and ap[2] */ 911 ap[1] = ap[2]; 912 ap[2] = dk[1]; 913 } 914 } 915 } 916 ierr = PetscFree(a2anew);CHKERRQ(ierr); 917 } 918 919 /* for each row k */ 920 for (k = 0; k<mbs; k++) { 921 922 /*initialize k-th row with elements nonzero in row perm(k) of A */ 923 jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 924 ap = aa + jmin*4; 925 for (j = jmin; j < jmax; j++) { 926 vj = perm_ptr[aj[j]]; /* block col. index */ 927 rtmp_ptr = rtmp + vj*4; 928 for (i=0; i<4; i++) *rtmp_ptr++ = *ap++; 929 } 930 931 /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 932 ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr); 933 i = jl[k]; /* first row to be added to k_th row */ 934 935 while (i < k) { 936 nexti = jl[i]; /* next row to be added to k_th row */ 937 938 /* compute multiplier */ 939 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 940 941 /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */ 942 diag = ba + i*4; 943 u = ba + ili*4; 944 uik[0] = -(diag[0]*u[0] + diag[2]*u[1]); 945 uik[1] = -(diag[1]*u[0] + diag[3]*u[1]); 946 uik[2] = -(diag[0]*u[2] + diag[2]*u[3]); 947 uik[3] = -(diag[1]*u[2] + diag[3]*u[3]); 948 949 /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */ 950 dk[0] += uik[0]*u[0] + uik[1]*u[1]; 951 dk[1] += uik[2]*u[0] + uik[3]*u[1]; 952 dk[2] += uik[0]*u[2] + uik[1]*u[3]; 953 dk[3] += uik[2]*u[2] + uik[3]*u[3]; 954 955 ierr = PetscLogFlops(16.0*2.0);CHKERRQ(ierr); 956 957 /* update -U(i,k): ba[ili] = uik */ 958 ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr); 959 960 /* add multiple of row i to k-th row ... */ 961 jmin = ili + 1; jmax = bi[i+1]; 962 if (jmin < jmax) { 963 for (j=jmin; j<jmax; j++) { 964 /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */ 965 rtmp_ptr = rtmp + bj[j]*4; 966 u = ba + j*4; 967 rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1]; 968 rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1]; 969 rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3]; 970 rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3]; 971 } 972 ierr = PetscLogFlops(16.0*(jmax-jmin));CHKERRQ(ierr); 973 974 /* ... add i to row list for next nonzero entry */ 975 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 976 j = bj[jmin]; 977 jl[i] = jl[j]; jl[j] = i; /* update jl */ 978 } 979 i = nexti; 980 } 981 982 /* save nonzero entries in k-th row of U ... */ 983 984 /* invert diagonal block */ 985 diag = ba+k*4; 986 ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr); 987 ierr = PetscKernel_A_gets_inverse_A_2(diag,shift);CHKERRQ(ierr); 988 989 jmin = bi[k]; jmax = bi[k+1]; 990 if (jmin < jmax) { 991 for (j=jmin; j<jmax; j++) { 992 vj = bj[j]; /* block col. index of U */ 993 u = ba + j*4; 994 rtmp_ptr = rtmp + vj*4; 995 for (k1=0; k1<4; k1++) { 996 *u++ = *rtmp_ptr; 997 *rtmp_ptr++ = 0.0; 998 } 999 } 1000 1001 /* ... add k to row list for first nonzero entry in k-th row */ 1002 il[k] = jmin; 1003 i = bj[jmin]; 1004 jl[k] = jl[i]; jl[i] = k; 1005 } 1006 } 1007 1008 ierr = PetscFree(rtmp);CHKERRQ(ierr); 1009 ierr = PetscFree2(il,jl);CHKERRQ(ierr); 1010 if (a->permute) { 1011 ierr = PetscFree(aa);CHKERRQ(ierr); 1012 } 1013 ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 1014 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 1020 ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 1021 PetscFunctionReturn(0); 1022 } 1023 1024 /* 1025 Version for when blocks are 2 by 2 Using natural ordering 1026 */ 1027 #undef __FUNCT__ 1028 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering" 1029 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) 1030 { 1031 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data; 1032 PetscErrorCode ierr; 1033 PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 1034 PetscInt *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 1035 MatScalar *ba = b->a,*aa,*ap,dk[8],uik[8]; 1036 MatScalar *u,*diag,*rtmp,*rtmp_ptr; 1037 PetscReal shift = info->shiftamount; 1038 1039 PetscFunctionBegin; 1040 /* initialization */ 1041 /* il and jl record the first nonzero element in each row of the accessing 1042 window U(0:k, k:mbs-1). 1043 jl: list of rows to be added to uneliminated rows 1044 i>= k: jl(i) is the first row to be added to row i 1045 i< k: jl(i) is the row following row i in some list of rows 1046 jl(i) = mbs indicates the end of a list 1047 il(i): points to the first nonzero element in columns k,...,mbs-1 of 1048 row i of U */ 1049 ierr = PetscCalloc1(4*mbs,&rtmp);CHKERRQ(ierr); 1050 ierr = PetscMalloc2(mbs,&il,mbs,&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 = PetscKernel_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 1155 ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 1156 PetscFunctionReturn(0); 1157 } 1158 1159 /* 1160 Numeric U^T*D*U factorization for SBAIJ format. 1161 Version for blocks are 1 by 1. 1162 */ 1163 #undef __FUNCT__ 1164 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace" 1165 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_inplace(Mat C,Mat A,const MatFactorInfo *info) 1166 { 1167 Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ*)C->data; 1168 IS ip=b->row; 1169 PetscErrorCode ierr; 1170 const PetscInt *ai,*aj,*rip; 1171 PetscInt *a2anew,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol; 1172 PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 1173 MatScalar *rtmp,*ba=b->a,*bval,*aa,dk,uikdi; 1174 PetscReal rs; 1175 FactorShiftCtx sctx; 1176 1177 PetscFunctionBegin; 1178 /* MatPivotSetUp(): initialize shift context sctx */ 1179 ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr); 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 = PetscMalloc1(nz,&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,&rtmp,mbs,&il,mbs,&jl);CHKERRQ(ierr); 1205 1206 do { 1207 sctx.newshift = PETSC_FALSE; 1208 for (i=0; i<mbs; i++) { 1209 rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1210 } 1211 1212 for (k = 0; k<mbs; k++) { 1213 /*initialize k-th row by the perm[k]-th row of A */ 1214 jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 1215 bval = ba + bi[k]; 1216 for (j = jmin; j < jmax; j++) { 1217 col = rip[aj[j]]; 1218 rtmp[col] = aa[j]; 1219 *bval++ = 0.0; /* for in-place factorization */ 1220 } 1221 1222 /* shift the diagonal of the matrix */ 1223 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 1224 1225 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 1226 dk = rtmp[k]; 1227 i = jl[k]; /* first row to be added to k_th row */ 1228 1229 while (i < k) { 1230 nexti = jl[i]; /* next row to be added to k_th row */ 1231 1232 /* compute multiplier, update diag(k) and U(i,k) */ 1233 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1234 uikdi = -ba[ili]*ba[bi[i]]; /* diagonal(k) */ 1235 dk += uikdi*ba[ili]; 1236 ba[ili] = uikdi; /* -U(i,k) */ 1237 1238 /* add multiple of row i to k-th row */ 1239 jmin = ili + 1; jmax = bi[i+1]; 1240 if (jmin < jmax) { 1241 for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 1242 ierr = PetscLogFlops(2.0*(jmax-jmin));CHKERRQ(ierr); 1243 1244 /* update il and jl for row i */ 1245 il[i] = jmin; 1246 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 1247 } 1248 i = nexti; 1249 } 1250 1251 /* shift the diagonals when zero pivot is detected */ 1252 /* compute rs=sum of abs(off-diagonal) */ 1253 rs = 0.0; 1254 jmin = bi[k]+1; 1255 nz = bi[k+1] - jmin; 1256 if (nz) { 1257 bcol = bj + jmin; 1258 while (nz--) { 1259 rs += PetscAbsScalar(rtmp[*bcol]); 1260 bcol++; 1261 } 1262 } 1263 1264 sctx.rs = rs; 1265 sctx.pv = dk; 1266 ierr = MatPivotCheck(C,A,info,&sctx,k);CHKERRQ(ierr); 1267 if (sctx.newshift) break; /* sctx.shift_amount is updated */ 1268 dk = sctx.pv; 1269 1270 /* copy data into U(k,:) */ 1271 ba[bi[k]] = 1.0/dk; /* U(k,k) */ 1272 jmin = bi[k]+1; jmax = bi[k+1]; 1273 if (jmin < jmax) { 1274 for (j=jmin; j<jmax; j++) { 1275 col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 1276 } 1277 /* add the k-th row into il and jl */ 1278 il[k] = jmin; 1279 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 1280 } 1281 } 1282 } while (sctx.newshift); 1283 ierr = PetscFree3(rtmp,il,jl);CHKERRQ(ierr); 1284 if (a->permute) {ierr = PetscFree(aa);CHKERRQ(ierr);} 1285 1286 ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 1287 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 1296 ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 1297 if (sctx.nshift) { 1298 if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 1299 ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr); 1300 } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { 1301 ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr); 1302 } 1303 } 1304 PetscFunctionReturn(0); 1305 } 1306 1307 /* 1308 Version for when blocks are 1 by 1 Using natural ordering under new datastructure 1309 Modified from MatCholeskyFactorNumeric_SeqAIJ() 1310 */ 1311 #undef __FUNCT__ 1312 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering" 1313 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info) 1314 { 1315 Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data; 1316 Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)B->data; 1317 PetscErrorCode ierr; 1318 PetscInt i,j,mbs=A->rmap->n,*bi=b->i,*bj=b->j,*bdiag=b->diag,*bjtmp; 1319 PetscInt *ai=a->i,*aj=a->j,*ajtmp; 1320 PetscInt k,jmin,jmax,*c2r,*il,col,nexti,ili,nz; 1321 MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 1322 FactorShiftCtx sctx; 1323 PetscReal rs; 1324 MatScalar d,*v; 1325 1326 PetscFunctionBegin; 1327 ierr = PetscMalloc3(mbs,&rtmp,mbs,&il,mbs,&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 1335 ierr = PetscMemzero(rtmp,mbs*sizeof(MatScalar));CHKERRQ(ierr); 1336 1337 for (i=0; i<mbs; i++) { 1338 /* calculate sum(|aij|)-RealPart(aii), amt of shift needed for this row */ 1339 d = (aa)[a->diag[i]]; 1340 rtmp[i] += -PetscRealPart(d); /* diagonal entry */ 1341 ajtmp = aj + ai[i] + 1; /* exclude diagonal */ 1342 v = aa + ai[i] + 1; 1343 nz = ai[i+1] - ai[i] - 1; 1344 for (j=0; j<nz; j++) { 1345 rtmp[i] += PetscAbsScalar(v[j]); 1346 rtmp[ajtmp[j]] += PetscAbsScalar(v[j]); 1347 } 1348 if (PetscRealPart(rtmp[i]) > sctx.shift_top) sctx.shift_top = PetscRealPart(rtmp[i]); 1349 } 1350 sctx.shift_top *= 1.1; 1351 sctx.nshift_max = 5; 1352 sctx.shift_lo = 0.; 1353 sctx.shift_hi = 1.; 1354 } 1355 1356 /* allocate working arrays 1357 c2r: linked list, keep track of pivot rows for a given column. c2r[col]: head of the list for a given col 1358 il: for active k row, il[i] gives the index of the 1st nonzero entry in U[i,k:n-1] in bj and ba arrays 1359 */ 1360 do { 1361 sctx.newshift = PETSC_FALSE; 1362 1363 for (i=0; i<mbs; i++) c2r[i] = mbs; 1364 if (mbs) il[0] = 0; 1365 1366 for (k = 0; k<mbs; k++) { 1367 /* zero rtmp */ 1368 nz = bi[k+1] - bi[k]; 1369 bjtmp = bj + bi[k]; 1370 for (j=0; j<nz; j++) rtmp[bjtmp[j]] = 0.0; 1371 1372 /* load in initial unfactored row */ 1373 bval = ba + bi[k]; 1374 jmin = ai[k]; jmax = ai[k+1]; 1375 for (j = jmin; j < jmax; j++) { 1376 col = aj[j]; 1377 rtmp[col] = aa[j]; 1378 *bval++ = 0.0; /* for in-place factorization */ 1379 } 1380 /* shift the diagonal of the matrix: ZeropivotApply() */ 1381 rtmp[k] += sctx.shift_amount; /* shift the diagonal of the matrix */ 1382 1383 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 1384 dk = rtmp[k]; 1385 i = c2r[k]; /* first row to be added to k_th row */ 1386 1387 while (i < k) { 1388 nexti = c2r[i]; /* next row to be added to k_th row */ 1389 1390 /* compute multiplier, update diag(k) and U(i,k) */ 1391 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1392 uikdi = -ba[ili]*ba[bdiag[i]]; /* diagonal(k) */ 1393 dk += uikdi*ba[ili]; /* update diag[k] */ 1394 ba[ili] = uikdi; /* -U(i,k) */ 1395 1396 /* add multiple of row i to k-th row */ 1397 jmin = ili + 1; jmax = bi[i+1]; 1398 if (jmin < jmax) { 1399 for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 1400 /* update il and c2r for row i */ 1401 il[i] = jmin; 1402 j = bj[jmin]; c2r[i] = c2r[j]; c2r[j] = i; 1403 } 1404 i = nexti; 1405 } 1406 1407 /* copy data into U(k,:) */ 1408 rs = 0.0; 1409 jmin = bi[k]; jmax = bi[k+1]-1; 1410 if (jmin < jmax) { 1411 for (j=jmin; j<jmax; j++) { 1412 col = bj[j]; ba[j] = rtmp[col]; rs += PetscAbsScalar(ba[j]); 1413 } 1414 /* add the k-th row into il and c2r */ 1415 il[k] = jmin; 1416 i = bj[jmin]; c2r[k] = c2r[i]; c2r[i] = k; 1417 } 1418 1419 sctx.rs = rs; 1420 sctx.pv = dk; 1421 ierr = MatPivotCheck(B,A,info,&sctx,k);CHKERRQ(ierr); 1422 if (sctx.newshift) break; 1423 dk = sctx.pv; 1424 1425 ba[bdiag[k]] = 1.0/dk; /* U(k,k) */ 1426 } 1427 } while (sctx.newshift); 1428 1429 ierr = PetscFree3(rtmp,il,c2r);CHKERRQ(ierr); 1430 1431 B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1432 B->ops->solves = MatSolves_SeqSBAIJ_1; 1433 B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 1434 B->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering; 1435 B->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering; 1436 1437 B->assembled = PETSC_TRUE; 1438 B->preallocated = PETSC_TRUE; 1439 1440 ierr = PetscLogFlops(B->rmap->n);CHKERRQ(ierr); 1441 1442 /* MatPivotView() */ 1443 if (sctx.nshift) { 1444 if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { 1445 ierr = PetscInfo4(A,"number of shift_pd tries %D, shift_amount %g, diagonal shifted up by %e fraction top_value %e\n",sctx.nshift,(double)sctx.shift_amount,(double)sctx.shift_fraction,(double)sctx.shift_top);CHKERRQ(ierr); 1446 } else if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 1447 ierr = PetscInfo2(A,"number of shift_nz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr); 1448 } else if (info->shifttype == (PetscReal)MAT_SHIFT_INBLOCKS) { 1449 ierr = PetscInfo2(A,"number of shift_inblocks applied %D, each shift_amount %g\n",sctx.nshift,(double)info->shiftamount);CHKERRQ(ierr); 1450 } 1451 } 1452 PetscFunctionReturn(0); 1453 } 1454 1455 #undef __FUNCT__ 1456 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace" 1457 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info) 1458 { 1459 Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ*)C->data; 1460 PetscErrorCode ierr; 1461 PetscInt i,j,mbs = a->mbs; 1462 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 1463 PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz; 1464 MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 1465 PetscReal rs; 1466 FactorShiftCtx sctx; 1467 1468 PetscFunctionBegin; 1469 /* MatPivotSetUp(): initialize shift context sctx */ 1470 ierr = PetscMemzero(&sctx,sizeof(FactorShiftCtx));CHKERRQ(ierr); 1471 1472 /* initialization */ 1473 /* il and jl record the first nonzero element in each row of the accessing 1474 window U(0:k, k:mbs-1). 1475 jl: list of rows to be added to uneliminated rows 1476 i>= k: jl(i) is the first row to be added to row i 1477 i< k: jl(i) is the row following row i in some list of rows 1478 jl(i) = mbs indicates the end of a list 1479 il(i): points to the first nonzero element in U(i,k:mbs-1) 1480 */ 1481 ierr = PetscMalloc1(mbs,&rtmp);CHKERRQ(ierr); 1482 ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr); 1483 1484 do { 1485 sctx.newshift = PETSC_FALSE; 1486 for (i=0; i<mbs; i++) { 1487 rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1488 } 1489 1490 for (k = 0; k<mbs; k++) { 1491 /*initialize k-th row with elements nonzero in row perm(k) of A */ 1492 nz = ai[k+1] - ai[k]; 1493 acol = aj + ai[k]; 1494 aval = aa + ai[k]; 1495 bval = ba + bi[k]; 1496 while (nz--) { 1497 rtmp[*acol++] = *aval++; 1498 *bval++ = 0.0; /* for in-place factorization */ 1499 } 1500 1501 /* shift the diagonal of the matrix */ 1502 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 1503 1504 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 1505 dk = rtmp[k]; 1506 i = jl[k]; /* first row to be added to k_th row */ 1507 1508 while (i < k) { 1509 nexti = jl[i]; /* next row to be added to k_th row */ 1510 /* compute multiplier, update D(k) and U(i,k) */ 1511 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1512 uikdi = -ba[ili]*ba[bi[i]]; 1513 dk += uikdi*ba[ili]; 1514 ba[ili] = uikdi; /* -U(i,k) */ 1515 1516 /* add multiple of row i to k-th row ... */ 1517 jmin = ili + 1; 1518 nz = bi[i+1] - jmin; 1519 if (nz > 0) { 1520 bcol = bj + jmin; 1521 bval = ba + jmin; 1522 ierr = PetscLogFlops(2.0*nz);CHKERRQ(ierr); 1523 while (nz--) rtmp[*bcol++] += uikdi*(*bval++); 1524 1525 /* update il and jl for i-th row */ 1526 il[i] = jmin; 1527 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 1528 } 1529 i = nexti; 1530 } 1531 1532 /* shift the diagonals when zero pivot is detected */ 1533 /* compute rs=sum of abs(off-diagonal) */ 1534 rs = 0.0; 1535 jmin = bi[k]+1; 1536 nz = bi[k+1] - jmin; 1537 if (nz) { 1538 bcol = bj + jmin; 1539 while (nz--) { 1540 rs += PetscAbsScalar(rtmp[*bcol]); 1541 bcol++; 1542 } 1543 } 1544 1545 sctx.rs = rs; 1546 sctx.pv = dk; 1547 ierr = MatPivotCheck(C,A,info,&sctx,k);CHKERRQ(ierr); 1548 if (sctx.newshift) break; /* sctx.shift_amount is updated */ 1549 dk = sctx.pv; 1550 1551 /* copy data into U(k,:) */ 1552 ba[bi[k]] = 1.0/dk; 1553 jmin = bi[k]+1; 1554 nz = bi[k+1] - jmin; 1555 if (nz) { 1556 bcol = bj + jmin; 1557 bval = ba + jmin; 1558 while (nz--) { 1559 *bval++ = rtmp[*bcol]; 1560 rtmp[*bcol++] = 0.0; 1561 } 1562 /* add k-th row into il and jl */ 1563 il[k] = jmin; 1564 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 1565 } 1566 } /* end of for (k = 0; k<mbs; k++) */ 1567 } while (sctx.newshift); 1568 ierr = PetscFree(rtmp);CHKERRQ(ierr); 1569 ierr = PetscFree2(il,jl);CHKERRQ(ierr); 1570 1571 C->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 1572 C->ops->solves = MatSolves_SeqSBAIJ_1_inplace; 1573 C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 1574 C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 1575 C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_1_NaturalOrdering_inplace; 1576 1577 C->assembled = PETSC_TRUE; 1578 C->preallocated = PETSC_TRUE; 1579 1580 ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 1581 if (sctx.nshift) { 1582 if (info->shifttype == (PetscReal)MAT_SHIFT_NONZERO) { 1583 ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr); 1584 } else if (info->shifttype == (PetscReal)MAT_SHIFT_POSITIVE_DEFINITE) { 1585 ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,(double)sctx.shift_amount);CHKERRQ(ierr); 1586 } 1587 } 1588 PetscFunctionReturn(0); 1589 } 1590 1591 #undef __FUNCT__ 1592 #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ" 1593 PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,const MatFactorInfo *info) 1594 { 1595 PetscErrorCode ierr; 1596 Mat C; 1597 1598 PetscFunctionBegin; 1599 ierr = MatGetFactor(A,"petsc",MAT_FACTOR_CHOLESKY,&C);CHKERRQ(ierr); 1600 ierr = MatCholeskyFactorSymbolic(C,A,perm,info);CHKERRQ(ierr); 1601 ierr = MatCholeskyFactorNumeric(C,A,info);CHKERRQ(ierr); 1602 1603 A->ops->solve = C->ops->solve; 1604 A->ops->solvetranspose = C->ops->solvetranspose; 1605 1606 ierr = MatHeaderMerge(A,&C);CHKERRQ(ierr); 1607 PetscFunctionReturn(0); 1608 } 1609 1610 1611