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