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 ierr = PetscLogObjectMemory(*B,(iu[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 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 = PetscLLAddSorted(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 ierr = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 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 ierr = PetscMalloc(bs2*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 667 ierr = PetscMemzero(rtmp,bs2*mbs*sizeof(MatScalar));CHKERRQ(ierr); 668 ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 669 jl = il + mbs; 670 for (i=0; i<mbs; i++) { 671 jl[i] = mbs; il[0] = 0; 672 } 673 ierr = PetscMalloc((2*bs2+bs)*sizeof(MatScalar),&dk);CHKERRQ(ierr); 674 uik = dk + bs2; 675 work = uik + bs2; 676 ierr = PetscMalloc(bs*sizeof(PetscInt),&pivots);CHKERRQ(ierr); 677 678 ai = a->i; aj = a->j; aa = a->a; 679 680 /* for each row k */ 681 for (k = 0; k<mbs; k++){ 682 683 /*initialize k-th row with elements nonzero in row k of A */ 684 jmin = ai[k]; jmax = ai[k+1]; 685 ap = aa + jmin*bs2; 686 for (j = jmin; j < jmax; j++){ 687 vj = aj[j]; /* block col. index */ 688 rtmp_ptr = rtmp + vj*bs2; 689 for (i=0; i<bs2; i++) *rtmp_ptr++ = *ap++; 690 } 691 692 /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 693 ierr = PetscMemcpy(dk,rtmp+k*bs2,bs2*sizeof(MatScalar));CHKERRQ(ierr); 694 i = jl[k]; /* first row to be added to k_th row */ 695 696 while (i < k){ 697 nexti = jl[i]; /* next row to be added to k_th row */ 698 699 /* compute multiplier */ 700 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 701 702 /* uik = -inv(Di)*U_bar(i,k) */ 703 diag = ba + i*bs2; 704 u = ba + ili*bs2; 705 ierr = PetscMemzero(uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 706 Kernel_A_gets_A_minus_B_times_C(bs,uik,diag,u); 707 708 /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 709 Kernel_A_gets_A_plus_Btranspose_times_C(bs,dk,uik,u); 710 711 /* update -U(i,k) */ 712 ierr = PetscMemcpy(ba+ili*bs2,uik,bs2*sizeof(MatScalar));CHKERRQ(ierr); 713 714 /* add multiple of row i to k-th row ... */ 715 jmin = ili + 1; jmax = bi[i+1]; 716 if (jmin < jmax){ 717 for (j=jmin; j<jmax; j++) { 718 /* rtmp += -U(i,k)^T * U_bar(i,j) */ 719 rtmp_ptr = rtmp + bj[j]*bs2; 720 u = ba + j*bs2; 721 Kernel_A_gets_A_plus_Btranspose_times_C(bs,rtmp_ptr,uik,u); 722 } 723 724 /* ... add i to row list for next nonzero entry */ 725 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 726 j = bj[jmin]; 727 jl[i] = jl[j]; jl[j] = i; /* update jl */ 728 } 729 i = nexti; 730 } 731 732 /* save nonzero entries in k-th row of U ... */ 733 734 /* invert diagonal block */ 735 diag = ba+k*bs2; 736 ierr = PetscMemcpy(diag,dk,bs2*sizeof(MatScalar));CHKERRQ(ierr); 737 ierr = Kernel_A_gets_inverse_A(bs,diag,pivots,work);CHKERRQ(ierr); 738 739 jmin = bi[k]; jmax = bi[k+1]; 740 if (jmin < jmax) { 741 for (j=jmin; j<jmax; j++){ 742 vj = bj[j]; /* block col. index of U */ 743 u = ba + j*bs2; 744 rtmp_ptr = rtmp + vj*bs2; 745 for (k1=0; k1<bs2; k1++){ 746 *u++ = *rtmp_ptr; 747 *rtmp_ptr++ = 0.0; 748 } 749 } 750 751 /* ... add k to row list for first nonzero entry in k-th row */ 752 il[k] = jmin; 753 i = bj[jmin]; 754 jl[k] = jl[i]; jl[i] = k; 755 } 756 } 757 758 ierr = PetscFree(rtmp);CHKERRQ(ierr); 759 ierr = PetscFree(il);CHKERRQ(ierr); 760 ierr = PetscFree(dk);CHKERRQ(ierr); 761 ierr = PetscFree(pivots);CHKERRQ(ierr); 762 763 C->factor = FACTOR_CHOLESKY; 764 C->assembled = PETSC_TRUE; 765 C->preallocated = PETSC_TRUE; 766 PetscLogFlops(1.3333*bs*bs2*b->mbs); /* from inverting diagonal blocks */ 767 PetscFunctionReturn(0); 768 } 769 770 /* 771 Numeric U^T*D*U factorization for SBAIJ format. Modified from SNF of YSMP. 772 Version for blocks 2 by 2. 773 */ 774 #undef __FUNCT__ 775 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2" 776 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2(Mat A,MatFactorInfo *info,Mat *B) 777 { 778 Mat C = *B; 779 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 780 IS perm = b->row; 781 PetscErrorCode ierr; 782 PetscInt *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 783 PetscInt *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 784 MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 785 MatScalar *u,*diag,*rtmp,*rtmp_ptr; 786 787 PetscFunctionBegin; 788 /* initialization */ 789 /* il and jl record the first nonzero element in each row of the accessing 790 window U(0:k, k:mbs-1). 791 jl: list of rows to be added to uneliminated rows 792 i>= k: jl(i) is the first row to be added to row i 793 i< k: jl(i) is the row following row i in some list of rows 794 jl(i) = mbs indicates the end of a list 795 il(i): points to the first nonzero element in columns k,...,mbs-1 of 796 row i of U */ 797 ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 798 ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); 799 ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 800 jl = il + mbs; 801 for (i=0; i<mbs; i++) { 802 jl[i] = mbs; il[0] = 0; 803 } 804 ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr); 805 uik = dk + 4; 806 ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr); 807 808 /* check permutation */ 809 if (!a->permute){ 810 ai = a->i; aj = a->j; aa = a->a; 811 } else { 812 ai = a->inew; aj = a->jnew; 813 ierr = PetscMalloc(4*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr); 814 ierr = PetscMemcpy(aa,a->a,4*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr); 815 ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr); 816 ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr); 817 818 for (i=0; i<mbs; i++){ 819 jmin = ai[i]; jmax = ai[i+1]; 820 for (j=jmin; j<jmax; j++){ 821 while (a2anew[j] != j){ 822 k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k; 823 for (k1=0; k1<4; k1++){ 824 dk[k1] = aa[k*4+k1]; 825 aa[k*4+k1] = aa[j*4+k1]; 826 aa[j*4+k1] = dk[k1]; 827 } 828 } 829 /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */ 830 if (i > aj[j]){ 831 /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 832 ap = aa + j*4; /* ptr to the beginning of the block */ 833 dk[1] = ap[1]; /* swap ap[1] and ap[2] */ 834 ap[1] = ap[2]; 835 ap[2] = dk[1]; 836 } 837 } 838 } 839 ierr = PetscFree(a2anew);CHKERRQ(ierr); 840 } 841 842 /* for each row k */ 843 for (k = 0; k<mbs; k++){ 844 845 /*initialize k-th row with elements nonzero in row perm(k) of A */ 846 jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1]; 847 ap = aa + jmin*4; 848 for (j = jmin; j < jmax; j++){ 849 vj = perm_ptr[aj[j]]; /* block col. index */ 850 rtmp_ptr = rtmp + vj*4; 851 for (i=0; i<4; i++) *rtmp_ptr++ = *ap++; 852 } 853 854 /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 855 ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr); 856 i = jl[k]; /* first row to be added to k_th row */ 857 858 while (i < k){ 859 nexti = jl[i]; /* next row to be added to k_th row */ 860 861 /* compute multiplier */ 862 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 863 864 /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */ 865 diag = ba + i*4; 866 u = ba + ili*4; 867 uik[0] = -(diag[0]*u[0] + diag[2]*u[1]); 868 uik[1] = -(diag[1]*u[0] + diag[3]*u[1]); 869 uik[2] = -(diag[0]*u[2] + diag[2]*u[3]); 870 uik[3] = -(diag[1]*u[2] + diag[3]*u[3]); 871 872 /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */ 873 dk[0] += uik[0]*u[0] + uik[1]*u[1]; 874 dk[1] += uik[2]*u[0] + uik[3]*u[1]; 875 dk[2] += uik[0]*u[2] + uik[1]*u[3]; 876 dk[3] += uik[2]*u[2] + uik[3]*u[3]; 877 878 /* update -U(i,k): ba[ili] = uik */ 879 ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr); 880 881 /* add multiple of row i to k-th row ... */ 882 jmin = ili + 1; jmax = bi[i+1]; 883 if (jmin < jmax){ 884 for (j=jmin; j<jmax; j++) { 885 /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */ 886 rtmp_ptr = rtmp + bj[j]*4; 887 u = ba + j*4; 888 rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1]; 889 rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1]; 890 rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3]; 891 rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3]; 892 } 893 894 /* ... add i to row list for next nonzero entry */ 895 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 896 j = bj[jmin]; 897 jl[i] = jl[j]; jl[j] = i; /* update jl */ 898 } 899 i = nexti; 900 } 901 902 /* save nonzero entries in k-th row of U ... */ 903 904 /* invert diagonal block */ 905 diag = ba+k*4; 906 ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr); 907 ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr); 908 909 jmin = bi[k]; jmax = bi[k+1]; 910 if (jmin < jmax) { 911 for (j=jmin; j<jmax; j++){ 912 vj = bj[j]; /* block col. index of U */ 913 u = ba + j*4; 914 rtmp_ptr = rtmp + vj*4; 915 for (k1=0; k1<4; k1++){ 916 *u++ = *rtmp_ptr; 917 *rtmp_ptr++ = 0.0; 918 } 919 } 920 921 /* ... add k to row list for first nonzero entry in k-th row */ 922 il[k] = jmin; 923 i = bj[jmin]; 924 jl[k] = jl[i]; jl[i] = k; 925 } 926 } 927 928 ierr = PetscFree(rtmp);CHKERRQ(ierr); 929 ierr = PetscFree(il);CHKERRQ(ierr); 930 ierr = PetscFree(dk);CHKERRQ(ierr); 931 if (a->permute) { 932 ierr = PetscFree(aa);CHKERRQ(ierr); 933 } 934 ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr); 935 C->factor = FACTOR_CHOLESKY; 936 C->assembled = PETSC_TRUE; 937 C->preallocated = PETSC_TRUE; 938 PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */ 939 PetscFunctionReturn(0); 940 } 941 942 /* 943 Version for when blocks are 2 by 2 Using natural ordering 944 */ 945 #undef __FUNCT__ 946 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering" 947 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_2_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B) 948 { 949 Mat C = *B; 950 Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 951 PetscErrorCode ierr; 952 PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 953 PetscInt *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 954 MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 955 MatScalar *u,*diag,*rtmp,*rtmp_ptr; 956 957 PetscFunctionBegin; 958 /* initialization */ 959 /* il and jl record the first nonzero element in each row of the accessing 960 window U(0:k, k:mbs-1). 961 jl: list of rows to be added to uneliminated rows 962 i>= k: jl(i) is the first row to be added to row i 963 i< k: jl(i) is the row following row i in some list of rows 964 jl(i) = mbs indicates the end of a list 965 il(i): points to the first nonzero element in columns k,...,mbs-1 of 966 row i of U */ 967 ierr = PetscMalloc(4*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 968 ierr = PetscMemzero(rtmp,4*mbs*sizeof(MatScalar));CHKERRQ(ierr); 969 ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 970 jl = il + mbs; 971 for (i=0; i<mbs; i++) { 972 jl[i] = mbs; il[0] = 0; 973 } 974 ierr = PetscMalloc(8*sizeof(MatScalar),&dk);CHKERRQ(ierr); 975 uik = dk + 4; 976 977 ai = a->i; aj = a->j; aa = a->a; 978 979 /* for each row k */ 980 for (k = 0; k<mbs; k++){ 981 982 /*initialize k-th row with elements nonzero in row k of A */ 983 jmin = ai[k]; jmax = ai[k+1]; 984 ap = aa + jmin*4; 985 for (j = jmin; j < jmax; j++){ 986 vj = aj[j]; /* block col. index */ 987 rtmp_ptr = rtmp + vj*4; 988 for (i=0; i<4; i++) *rtmp_ptr++ = *ap++; 989 } 990 991 /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 992 ierr = PetscMemcpy(dk,rtmp+k*4,4*sizeof(MatScalar));CHKERRQ(ierr); 993 i = jl[k]; /* first row to be added to k_th row */ 994 995 while (i < k){ 996 nexti = jl[i]; /* next row to be added to k_th row */ 997 998 /* compute multiplier */ 999 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1000 1001 /* uik = -inv(Di)*U_bar(i,k): - ba[ili]*ba[i] */ 1002 diag = ba + i*4; 1003 u = ba + ili*4; 1004 uik[0] = -(diag[0]*u[0] + diag[2]*u[1]); 1005 uik[1] = -(diag[1]*u[0] + diag[3]*u[1]); 1006 uik[2] = -(diag[0]*u[2] + diag[2]*u[3]); 1007 uik[3] = -(diag[1]*u[2] + diag[3]*u[3]); 1008 1009 /* update D(k) += -U(i,k)^T * U_bar(i,k): dk += uik*ba[ili] */ 1010 dk[0] += uik[0]*u[0] + uik[1]*u[1]; 1011 dk[1] += uik[2]*u[0] + uik[3]*u[1]; 1012 dk[2] += uik[0]*u[2] + uik[1]*u[3]; 1013 dk[3] += uik[2]*u[2] + uik[3]*u[3]; 1014 1015 /* update -U(i,k): ba[ili] = uik */ 1016 ierr = PetscMemcpy(ba+ili*4,uik,4*sizeof(MatScalar));CHKERRQ(ierr); 1017 1018 /* add multiple of row i to k-th row ... */ 1019 jmin = ili + 1; jmax = bi[i+1]; 1020 if (jmin < jmax){ 1021 for (j=jmin; j<jmax; j++) { 1022 /* rtmp += -U(i,k)^T * U_bar(i,j): rtmp[bj[j]] += uik*ba[j]; */ 1023 rtmp_ptr = rtmp + bj[j]*4; 1024 u = ba + j*4; 1025 rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1]; 1026 rtmp_ptr[1] += uik[2]*u[0] + uik[3]*u[1]; 1027 rtmp_ptr[2] += uik[0]*u[2] + uik[1]*u[3]; 1028 rtmp_ptr[3] += uik[2]*u[2] + uik[3]*u[3]; 1029 } 1030 1031 /* ... add i to row list for next nonzero entry */ 1032 il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 1033 j = bj[jmin]; 1034 jl[i] = jl[j]; jl[j] = i; /* update jl */ 1035 } 1036 i = nexti; 1037 } 1038 1039 /* save nonzero entries in k-th row of U ... */ 1040 1041 /* invert diagonal block */ 1042 diag = ba+k*4; 1043 ierr = PetscMemcpy(diag,dk,4*sizeof(MatScalar));CHKERRQ(ierr); 1044 ierr = Kernel_A_gets_inverse_A_2(diag);CHKERRQ(ierr); 1045 1046 jmin = bi[k]; jmax = bi[k+1]; 1047 if (jmin < jmax) { 1048 for (j=jmin; j<jmax; j++){ 1049 vj = bj[j]; /* block col. index of U */ 1050 u = ba + j*4; 1051 rtmp_ptr = rtmp + vj*4; 1052 for (k1=0; k1<4; k1++){ 1053 *u++ = *rtmp_ptr; 1054 *rtmp_ptr++ = 0.0; 1055 } 1056 } 1057 1058 /* ... add k to row list for first nonzero entry in k-th row */ 1059 il[k] = jmin; 1060 i = bj[jmin]; 1061 jl[k] = jl[i]; jl[i] = k; 1062 } 1063 } 1064 1065 ierr = PetscFree(rtmp);CHKERRQ(ierr); 1066 ierr = PetscFree(il);CHKERRQ(ierr); 1067 ierr = PetscFree(dk);CHKERRQ(ierr); 1068 1069 C->factor = FACTOR_CHOLESKY; 1070 C->assembled = PETSC_TRUE; 1071 C->preallocated = PETSC_TRUE; 1072 PetscLogFlops(1.3333*8*b->mbs); /* from inverting diagonal blocks */ 1073 PetscFunctionReturn(0); 1074 } 1075 1076 /* 1077 Numeric U^T*D*U factorization for SBAIJ format. 1078 Version for blocks are 1 by 1. 1079 */ 1080 #undef __FUNCT__ 1081 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1" 1082 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1(Mat A,MatFactorInfo *info,Mat *B) 1083 { 1084 Mat C = *B; 1085 Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data; 1086 IS ip=b->row; 1087 PetscErrorCode ierr; 1088 PetscInt *rip,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j,*bcol; 1089 PetscInt *ai,*aj,*a2anew; 1090 PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 1091 MatScalar *rtmp,*ba=b->a,*bval,*aa,dk,uikdi; 1092 PetscReal zeropivot,rs,shiftnz; 1093 PetscTruth shiftpd; 1094 ChShift_Ctx sctx; 1095 PetscInt newshift; 1096 1097 PetscFunctionBegin; 1098 /* initialization */ 1099 shiftnz = info->shiftnz; 1100 shiftpd = info->shiftpd; 1101 zeropivot = info->zeropivot; 1102 1103 ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 1104 if (!a->permute){ 1105 ai = a->i; aj = a->j; aa = a->a; 1106 } else { 1107 ai = a->inew; aj = a->jnew; 1108 nz = ai[mbs]; 1109 ierr = PetscMalloc(nz*sizeof(MatScalar),&aa);CHKERRQ(ierr); 1110 a2anew = a->a2anew; 1111 bval = a->a; 1112 for (j=0; j<nz; j++){ 1113 aa[a2anew[j]] = *(bval++); 1114 } 1115 } 1116 1117 /* initialization */ 1118 /* il and jl record the first nonzero element in each row of the accessing 1119 window U(0:k, k:mbs-1). 1120 jl: list of rows to be added to uneliminated rows 1121 i>= k: jl(i) is the first row to be added to row i 1122 i< k: jl(i) is the row following row i in some list of rows 1123 jl(i) = mbs indicates the end of a list 1124 il(i): points to the first nonzero element in columns k,...,mbs-1 of 1125 row i of U */ 1126 nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 1127 ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 1128 jl = il + mbs; 1129 rtmp = (MatScalar*)(jl + mbs); 1130 1131 sctx.shift_amount = 0; 1132 sctx.nshift = 0; 1133 do { 1134 sctx.chshift = PETSC_FALSE; 1135 for (i=0; i<mbs; i++) { 1136 rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1137 } 1138 1139 for (k = 0; k<mbs; k++){ 1140 /*initialize k-th row by the perm[k]-th row of A */ 1141 jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 1142 bval = ba + bi[k]; 1143 for (j = jmin; j < jmax; j++){ 1144 col = rip[aj[j]]; 1145 rtmp[col] = aa[j]; 1146 *bval++ = 0.0; /* for in-place factorization */ 1147 } 1148 1149 /* shift the diagonal of the matrix */ 1150 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 1151 1152 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 1153 dk = rtmp[k]; 1154 i = jl[k]; /* first row to be added to k_th row */ 1155 1156 while (i < k){ 1157 nexti = jl[i]; /* next row to be added to k_th row */ 1158 1159 /* compute multiplier, update diag(k) and U(i,k) */ 1160 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1161 uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 1162 dk += uikdi*ba[ili]; 1163 ba[ili] = uikdi; /* -U(i,k) */ 1164 1165 /* add multiple of row i to k-th row */ 1166 jmin = ili + 1; jmax = bi[i+1]; 1167 if (jmin < jmax){ 1168 for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 1169 /* update il and jl for row i */ 1170 il[i] = jmin; 1171 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 1172 } 1173 i = nexti; 1174 } 1175 1176 /* shift the diagonals when zero pivot is detected */ 1177 /* compute rs=sum of abs(off-diagonal) */ 1178 rs = 0.0; 1179 jmin = bi[k]+1; 1180 nz = bi[k+1] - jmin; 1181 if (nz){ 1182 bcol = bj + jmin; 1183 while (nz--){ 1184 rs += PetscAbsScalar(rtmp[*bcol]); 1185 bcol++; 1186 } 1187 } 1188 1189 sctx.rs = rs; 1190 sctx.pv = dk; 1191 ierr = Mat_CholeskyCheckShift(info,&sctx,&newshift);CHKERRQ(ierr); 1192 if (newshift == 1){ 1193 break; /* sctx.shift_amount is updated */ 1194 } else if (newshift == -1){ 1195 SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs); 1196 } 1197 1198 /* copy data into U(k,:) */ 1199 ba[bi[k]] = 1.0/dk; /* U(k,k) */ 1200 jmin = bi[k]+1; jmax = bi[k+1]; 1201 if (jmin < jmax) { 1202 for (j=jmin; j<jmax; j++){ 1203 col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 1204 } 1205 /* add the k-th row into il and jl */ 1206 il[k] = jmin; 1207 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 1208 } 1209 } 1210 } while (sctx.chshift); 1211 ierr = PetscFree(il);CHKERRQ(ierr); 1212 if (a->permute){ierr = PetscFree(aa);CHKERRQ(ierr);} 1213 1214 ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 1215 C->factor = FACTOR_CHOLESKY; 1216 C->assembled = PETSC_TRUE; 1217 C->preallocated = PETSC_TRUE; 1218 PetscLogFlops(C->m); 1219 if (sctx.nshift){ 1220 if (shiftnz) { 1221 PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqSBAIJ_1: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount); 1222 } else if (shiftpd) { 1223 PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqSBAIJ_1: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount); 1224 } 1225 } 1226 PetscFunctionReturn(0); 1227 } 1228 1229 /* 1230 Version for when blocks are 1 by 1 Using natural ordering 1231 */ 1232 #undef __FUNCT__ 1233 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering" 1234 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering(Mat A,MatFactorInfo *info,Mat *B) 1235 { 1236 Mat C = *B; 1237 Mat_SeqSBAIJ *a=(Mat_SeqSBAIJ*)A->data,*b=(Mat_SeqSBAIJ *)C->data; 1238 PetscErrorCode ierr; 1239 PetscInt i,j,mbs = a->mbs; 1240 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 1241 PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz; 1242 MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 1243 PetscReal zeropivot,rs,shiftnz; 1244 PetscTruth shiftpd; 1245 ChShift_Ctx sctx; 1246 PetscInt newshift; 1247 1248 PetscFunctionBegin; 1249 /* initialization */ 1250 shiftnz = info->shiftnz; 1251 shiftpd = info->shiftpd; 1252 zeropivot = info->zeropivot; 1253 1254 /* il and jl record the first nonzero element in each row of the accessing 1255 window U(0:k, k:mbs-1). 1256 jl: list of rows to be added to uneliminated rows 1257 i>= k: jl(i) is the first row to be added to row i 1258 i< k: jl(i) is the row following row i in some list of rows 1259 jl(i) = mbs indicates the end of a list 1260 il(i): points to the first nonzero element in U(i,k:mbs-1) 1261 */ 1262 ierr = PetscMalloc(mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 1263 ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr); 1264 jl = il + mbs; 1265 1266 sctx.shift_amount = 0; 1267 sctx.nshift = 0; 1268 do { 1269 sctx.chshift = PETSC_FALSE; 1270 for (i=0; i<mbs; i++) { 1271 rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 1272 } 1273 1274 for (k = 0; k<mbs; k++){ 1275 /*initialize k-th row with elements nonzero in row perm(k) of A */ 1276 nz = ai[k+1] - ai[k]; 1277 acol = aj + ai[k]; 1278 aval = aa + ai[k]; 1279 bval = ba + bi[k]; 1280 while (nz -- ){ 1281 rtmp[*acol++] = *aval++; 1282 *bval++ = 0.0; /* for in-place factorization */ 1283 } 1284 1285 /* shift the diagonal of the matrix */ 1286 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 1287 1288 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 1289 dk = rtmp[k]; 1290 i = jl[k]; /* first row to be added to k_th row */ 1291 1292 while (i < k){ 1293 nexti = jl[i]; /* next row to be added to k_th row */ 1294 /* compute multiplier, update D(k) and U(i,k) */ 1295 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 1296 uikdi = - ba[ili]*ba[bi[i]]; 1297 dk += uikdi*ba[ili]; 1298 ba[ili] = uikdi; /* -U(i,k) */ 1299 1300 /* add multiple of row i to k-th row ... */ 1301 jmin = ili + 1; 1302 nz = bi[i+1] - jmin; 1303 if (nz > 0){ 1304 bcol = bj + jmin; 1305 bval = ba + jmin; 1306 while (nz --) rtmp[*bcol++] += uikdi*(*bval++); 1307 /* update il and jl for i-th row */ 1308 il[i] = jmin; 1309 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 1310 } 1311 i = nexti; 1312 } 1313 1314 /* shift the diagonals when zero pivot is detected */ 1315 /* compute rs=sum of abs(off-diagonal) */ 1316 rs = 0.0; 1317 jmin = bi[k]+1; 1318 nz = bi[k+1] - jmin; 1319 if (nz){ 1320 bcol = bj + jmin; 1321 while (nz--){ 1322 rs += PetscAbsScalar(rtmp[*bcol]); 1323 bcol++; 1324 } 1325 } 1326 1327 sctx.rs = rs; 1328 sctx.pv = dk; 1329 ierr = Mat_CholeskyCheckShift(info,&sctx,&newshift);CHKERRQ(ierr); 1330 if (newshift == 1){ 1331 break; /* sctx.shift_amount is updated */ 1332 } else if (newshift == -1){ 1333 SETERRQ4(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot row %D value %g tolerance %g * rs %g",k,PetscAbsScalar(dk),zeropivot,rs); 1334 } 1335 1336 /* copy data into U(k,:) */ 1337 ba[bi[k]] = 1.0/dk; 1338 jmin = bi[k]+1; 1339 nz = bi[k+1] - jmin; 1340 if (nz){ 1341 bcol = bj + jmin; 1342 bval = ba + jmin; 1343 while (nz--){ 1344 *bval++ = rtmp[*bcol]; 1345 rtmp[*bcol++] = 0.0; 1346 } 1347 /* add k-th row into il and jl */ 1348 il[k] = jmin; 1349 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 1350 } 1351 } /* end of for (k = 0; k<mbs; k++) */ 1352 } while (sctx.chshift); 1353 ierr = PetscFree(rtmp);CHKERRQ(ierr); 1354 ierr = PetscFree(il);CHKERRQ(ierr); 1355 1356 C->factor = FACTOR_CHOLESKY; 1357 C->assembled = PETSC_TRUE; 1358 C->preallocated = PETSC_TRUE; 1359 PetscLogFlops(C->m); 1360 if (sctx.nshift){ 1361 if (shiftnz) { 1362 PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering: number of shiftnz tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount); 1363 } else if (shiftpd) { 1364 PetscLogInfo(0,"MatCholeskyFactorNumeric_SeqSBAIJ_1_NaturalOrdering: number of shiftpd tries %D, shift_amount %g\n",sctx.nshift,sctx.shift_amount); 1365 } 1366 } 1367 PetscFunctionReturn(0); 1368 } 1369 1370 #undef __FUNCT__ 1371 #define __FUNCT__ "MatCholeskyFactor_SeqSBAIJ" 1372 PetscErrorCode MatCholeskyFactor_SeqSBAIJ(Mat A,IS perm,MatFactorInfo *info) 1373 { 1374 PetscErrorCode ierr; 1375 Mat C; 1376 1377 PetscFunctionBegin; 1378 ierr = MatCholeskyFactorSymbolic(A,perm,info,&C);CHKERRQ(ierr); 1379 ierr = MatCholeskyFactorNumeric(A,info,&C);CHKERRQ(ierr); 1380 ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 1381 PetscFunctionReturn(0); 1382 } 1383 1384 1385