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