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