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