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