1 #define PETSCMAT_DLL 2 3 /* 4 Factorization code for BAIJ format. 5 */ 6 #include "../src/mat/impls/baij/seq/baij.h" 7 #include "../src/mat/blockinvert.h" 8 9 /* ------------------------------------------------------------*/ 10 /* 11 Version for when blocks are 2 by 2 12 */ 13 #undef __FUNCT__ 14 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2" 15 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2(Mat B,Mat A,const MatFactorInfo *info) 16 { 17 Mat C = B; 18 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 19 IS isrow = b->row,isicol = b->icol; 20 PetscErrorCode ierr; 21 const PetscInt *r,*ic; 22 PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 23 PetscInt *ajtmpold,*ajtmp,nz,row; 24 PetscInt *diag_offset=b->diag,idx,*ai=a->i,*aj=a->j,*pj; 25 MatScalar *pv,*v,*rtmp,m1,m2,m3,m4,*pc,*w,*x,x1,x2,x3,x4; 26 MatScalar p1,p2,p3,p4; 27 MatScalar *ba = b->a,*aa = a->a; 28 PetscReal shift = info->shiftinblocks; 29 30 PetscFunctionBegin; 31 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 32 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 33 ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 34 35 for (i=0; i<n; i++) { 36 nz = bi[i+1] - bi[i]; 37 ajtmp = bj + bi[i]; 38 for (j=0; j<nz; j++) { 39 x = rtmp+4*ajtmp[j]; x[0] = x[1] = x[2] = x[3] = 0.0; 40 } 41 /* load in initial (unfactored row) */ 42 idx = r[i]; 43 nz = ai[idx+1] - ai[idx]; 44 ajtmpold = aj + ai[idx]; 45 v = aa + 4*ai[idx]; 46 for (j=0; j<nz; j++) { 47 x = rtmp+4*ic[ajtmpold[j]]; 48 x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 49 v += 4; 50 } 51 row = *ajtmp++; 52 while (row < i) { 53 pc = rtmp + 4*row; 54 p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 55 if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) { 56 pv = ba + 4*diag_offset[row]; 57 pj = bj + diag_offset[row] + 1; 58 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 59 pc[0] = m1 = p1*x1 + p3*x2; 60 pc[1] = m2 = p2*x1 + p4*x2; 61 pc[2] = m3 = p1*x3 + p3*x4; 62 pc[3] = m4 = p2*x3 + p4*x4; 63 nz = bi[row+1] - diag_offset[row] - 1; 64 pv += 4; 65 for (j=0; j<nz; j++) { 66 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 67 x = rtmp + 4*pj[j]; 68 x[0] -= m1*x1 + m3*x2; 69 x[1] -= m2*x1 + m4*x2; 70 x[2] -= m1*x3 + m3*x4; 71 x[3] -= m2*x3 + m4*x4; 72 pv += 4; 73 } 74 ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr); 75 } 76 row = *ajtmp++; 77 } 78 /* finished row so stick it into b->a */ 79 pv = ba + 4*bi[i]; 80 pj = bj + bi[i]; 81 nz = bi[i+1] - bi[i]; 82 for (j=0; j<nz; j++) { 83 x = rtmp+4*pj[j]; 84 pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 85 pv += 4; 86 } 87 /* invert diagonal block */ 88 w = ba + 4*diag_offset[i]; 89 ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr); 90 } 91 92 ierr = PetscFree(rtmp);CHKERRQ(ierr); 93 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 94 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 95 C->ops->solve = MatSolve_SeqBAIJ_2; 96 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2; 97 C->assembled = PETSC_TRUE; 98 ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 99 PetscFunctionReturn(0); 100 } 101 /* 102 Version for when blocks are 2 by 2 Using natural ordering 103 */ 104 #undef __FUNCT__ 105 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering" 106 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) 107 { 108 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 109 PetscErrorCode ierr; 110 PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 111 PetscInt *ajtmpold,*ajtmp,nz,row; 112 PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 113 MatScalar *pv,*v,*rtmp,*pc,*w,*x; 114 MatScalar p1,p2,p3,p4,m1,m2,m3,m4,x1,x2,x3,x4; 115 MatScalar *ba = b->a,*aa = a->a; 116 PetscReal shift = info->shiftinblocks; 117 118 PetscFunctionBegin; 119 ierr = PetscMalloc(4*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 120 for (i=0; i<n; i++) { 121 nz = bi[i+1] - bi[i]; 122 ajtmp = bj + bi[i]; 123 for (j=0; j<nz; j++) { 124 x = rtmp+4*ajtmp[j]; 125 x[0] = x[1] = x[2] = x[3] = 0.0; 126 } 127 /* load in initial (unfactored row) */ 128 nz = ai[i+1] - ai[i]; 129 ajtmpold = aj + ai[i]; 130 v = aa + 4*ai[i]; 131 for (j=0; j<nz; j++) { 132 x = rtmp+4*ajtmpold[j]; 133 x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 134 v += 4; 135 } 136 row = *ajtmp++; 137 while (row < i) { 138 pc = rtmp + 4*row; 139 p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 140 if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0) { 141 pv = ba + 4*diag_offset[row]; 142 pj = bj + diag_offset[row] + 1; 143 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 144 pc[0] = m1 = p1*x1 + p3*x2; 145 pc[1] = m2 = p2*x1 + p4*x2; 146 pc[2] = m3 = p1*x3 + p3*x4; 147 pc[3] = m4 = p2*x3 + p4*x4; 148 nz = bi[row+1] - diag_offset[row] - 1; 149 pv += 4; 150 for (j=0; j<nz; j++) { 151 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 152 x = rtmp + 4*pj[j]; 153 x[0] -= m1*x1 + m3*x2; 154 x[1] -= m2*x1 + m4*x2; 155 x[2] -= m1*x3 + m3*x4; 156 x[3] -= m2*x3 + m4*x4; 157 pv += 4; 158 } 159 ierr = PetscLogFlops(16.0*nz+12.0);CHKERRQ(ierr); 160 } 161 row = *ajtmp++; 162 } 163 /* finished row so stick it into b->a */ 164 pv = ba + 4*bi[i]; 165 pj = bj + bi[i]; 166 nz = bi[i+1] - bi[i]; 167 for (j=0; j<nz; j++) { 168 x = rtmp+4*pj[j]; 169 pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 170 /* 171 printf(" col %d:",pj[j]); 172 PetscInt j1; 173 for (j1=0; j1<4; j1++) printf(" %g,",*(pv+j1)); 174 printf("\n"); 175 */ 176 pv += 4; 177 } 178 /* invert diagonal block */ 179 w = ba + 4*diag_offset[i]; 180 /* 181 printf(" \n%d -th: diag: ",i); 182 for (j=0; j<4; j++){ 183 printf(" %g,",w[j]); 184 } 185 printf("\n----------------------------\n"); 186 */ 187 ierr = Kernel_A_gets_inverse_A_2(w,shift);CHKERRQ(ierr); 188 } 189 190 ierr = PetscFree(rtmp);CHKERRQ(ierr); 191 C->ops->solve = MatSolve_SeqBAIJ_2_NaturalOrdering; 192 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_2_NaturalOrdering; 193 C->assembled = PETSC_TRUE; 194 ierr = PetscLogFlops(1.3333*8*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 195 PetscFunctionReturn(0); 196 } 197 198 /* ----------------------------------------------------------- */ 199 /* 200 Version for when blocks are 1 by 1. 201 */ 202 #undef __FUNCT__ 203 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_1" 204 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_1(Mat C,Mat A,const MatFactorInfo *info) 205 { 206 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 207 IS isrow = b->row,isicol = b->icol; 208 PetscErrorCode ierr; 209 const PetscInt *r,*ic; 210 PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 211 PetscInt *ajtmpold,*ajtmp,nz,row,*ai = a->i,*aj = a->j; 212 PetscInt *diag_offset = b->diag,diag,*pj; 213 MatScalar *pv,*v,*rtmp,multiplier,*pc; 214 MatScalar *ba = b->a,*aa = a->a; 215 PetscTruth row_identity, col_identity; 216 217 PetscFunctionBegin; 218 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 219 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 220 ierr = PetscMalloc((n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 221 222 for (i=0; i<n; i++) { 223 nz = bi[i+1] - bi[i]; 224 ajtmp = bj + bi[i]; 225 for (j=0; j<nz; j++) rtmp[ajtmp[j]] = 0.0; 226 227 /* load in initial (unfactored row) */ 228 nz = ai[r[i]+1] - ai[r[i]]; 229 ajtmpold = aj + ai[r[i]]; 230 v = aa + ai[r[i]]; 231 for (j=0; j<nz; j++) rtmp[ic[ajtmpold[j]]] = v[j]; 232 233 row = *ajtmp++; 234 while (row < i) { 235 pc = rtmp + row; 236 if (*pc != 0.0) { 237 pv = ba + diag_offset[row]; 238 pj = bj + diag_offset[row] + 1; 239 multiplier = *pc * *pv++; 240 *pc = multiplier; 241 nz = bi[row+1] - diag_offset[row] - 1; 242 for (j=0; j<nz; j++) rtmp[pj[j]] -= multiplier * pv[j]; 243 ierr = PetscLogFlops(1.0+2.0*nz);CHKERRQ(ierr); 244 } 245 row = *ajtmp++; 246 } 247 /* finished row so stick it into b->a */ 248 pv = ba + bi[i]; 249 pj = bj + bi[i]; 250 nz = bi[i+1] - bi[i]; 251 for (j=0; j<nz; j++) {pv[j] = rtmp[pj[j]];} 252 diag = diag_offset[i] - bi[i]; 253 /* check pivot entry for current row */ 254 if (pv[diag] == 0.0) { 255 SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Zero pivot: row in original ordering %D in permuted ordering %D",r[i],i); 256 } 257 pv[diag] = 1.0/pv[diag]; 258 } 259 260 ierr = PetscFree(rtmp);CHKERRQ(ierr); 261 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 262 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 263 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 264 ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr); 265 if (row_identity && col_identity) { 266 C->ops->solve = MatSolve_SeqBAIJ_1_NaturalOrdering; 267 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1_NaturalOrdering; 268 } else { 269 C->ops->solve = MatSolve_SeqBAIJ_1; 270 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_1; 271 } 272 C->assembled = PETSC_TRUE; 273 ierr = PetscLogFlops(C->cmap->n);CHKERRQ(ierr); 274 PetscFunctionReturn(0); 275 } 276 277 EXTERN_C_BEGIN 278 #undef __FUNCT__ 279 #define __FUNCT__ "MatGetFactor_seqbaij_petsc" 280 PetscErrorCode MatGetFactor_seqbaij_petsc(Mat A,MatFactorType ftype,Mat *B) 281 { 282 PetscInt n = A->rmap->n; 283 PetscErrorCode ierr; 284 285 PetscFunctionBegin; 286 ierr = MatCreate(((PetscObject)A)->comm,B);CHKERRQ(ierr); 287 ierr = MatSetSizes(*B,n,n,n,n);CHKERRQ(ierr); 288 if (ftype == MAT_FACTOR_LU || ftype == MAT_FACTOR_ILU || ftype == MAT_FACTOR_ILUDT) { 289 ierr = MatSetType(*B,MATSEQBAIJ);CHKERRQ(ierr); 290 (*B)->ops->lufactorsymbolic = MatLUFactorSymbolic_SeqBAIJ; 291 (*B)->ops->ilufactorsymbolic = MatILUFactorSymbolic_SeqBAIJ; 292 (*B)->ops->iludtfactor = MatILUDTFactor_SeqBAIJ; 293 } else if (ftype == MAT_FACTOR_CHOLESKY || ftype == MAT_FACTOR_ICC) { 294 ierr = MatSetType(*B,MATSEQSBAIJ);CHKERRQ(ierr); 295 ierr = MatSeqSBAIJSetPreallocation(*B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 296 (*B)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqBAIJ; 297 (*B)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqBAIJ; 298 } else SETERRQ(PETSC_ERR_SUP,"Factor type not supported"); 299 (*B)->factor = ftype; 300 PetscFunctionReturn(0); 301 } 302 EXTERN_C_END 303 304 EXTERN_C_BEGIN 305 #undef __FUNCT__ 306 #define __FUNCT__ "MatGetFactorAvailable_seqaij_petsc" 307 PetscErrorCode MatGetFactorAvailable_seqbaij_petsc(Mat A,MatFactorType ftype,PetscTruth *flg) 308 { 309 PetscFunctionBegin; 310 *flg = PETSC_TRUE; 311 PetscFunctionReturn(0); 312 } 313 EXTERN_C_END 314 315 /* ----------------------------------------------------------- */ 316 #undef __FUNCT__ 317 #define __FUNCT__ "MatLUFactor_SeqBAIJ" 318 PetscErrorCode MatLUFactor_SeqBAIJ(Mat A,IS row,IS col,const MatFactorInfo *info) 319 { 320 PetscErrorCode ierr; 321 Mat C; 322 323 PetscFunctionBegin; 324 ierr = MatGetFactor(A,MAT_SOLVER_PETSC,MAT_FACTOR_LU,&C);CHKERRQ(ierr); 325 ierr = MatLUFactorSymbolic(C,A,row,col,info);CHKERRQ(ierr); 326 ierr = MatLUFactorNumeric(C,A,info);CHKERRQ(ierr); 327 A->ops->solve = C->ops->solve; 328 A->ops->solvetranspose = C->ops->solvetranspose; 329 ierr = MatHeaderCopy(A,C);CHKERRQ(ierr); 330 ierr = PetscLogObjectParent(A,((Mat_SeqBAIJ*)(A->data))->icol);CHKERRQ(ierr); 331 PetscFunctionReturn(0); 332 } 333 334 #include "../src/mat/impls/sbaij/seq/sbaij.h" 335 #undef __FUNCT__ 336 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N" 337 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N(Mat C,Mat A,const MatFactorInfo *info) 338 { 339 PetscErrorCode ierr; 340 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 341 Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 342 IS ip=b->row; 343 const PetscInt *rip; 344 PetscInt i,j,mbs=a->mbs,bs=A->rmap->bs,*bi=b->i,*bj=b->j,*bcol; 345 PetscInt *ai=a->i,*aj=a->j; 346 PetscInt k,jmin,jmax,*jl,*il,col,nexti,ili,nz; 347 MatScalar *rtmp,*ba=b->a,*bval,*aa=a->a,dk,uikdi; 348 PetscReal zeropivot,rs,shiftnz; 349 PetscReal shiftpd; 350 ChShift_Ctx sctx; 351 PetscInt newshift; 352 353 PetscFunctionBegin; 354 if (bs > 1) { 355 if (!a->sbaijMat){ 356 ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 357 } 358 ierr = (a->sbaijMat)->ops->choleskyfactornumeric(C,a->sbaijMat,info);CHKERRQ(ierr); 359 ierr = MatDestroy(a->sbaijMat);CHKERRQ(ierr); 360 a->sbaijMat = PETSC_NULL; 361 PetscFunctionReturn(0); 362 } 363 364 /* initialization */ 365 shiftnz = info->shiftnz; 366 shiftpd = info->shiftpd; 367 zeropivot = info->zeropivot; 368 369 ierr = ISGetIndices(ip,&rip);CHKERRQ(ierr); 370 nz = (2*mbs+1)*sizeof(PetscInt)+mbs*sizeof(MatScalar); 371 ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 372 jl = il + mbs; 373 rtmp = (MatScalar*)(jl + mbs); 374 375 sctx.shift_amount = 0; 376 sctx.nshift = 0; 377 do { 378 sctx.chshift = PETSC_FALSE; 379 for (i=0; i<mbs; i++) { 380 rtmp[i] = 0.0; jl[i] = mbs; il[0] = 0; 381 } 382 383 for (k = 0; k<mbs; k++){ 384 bval = ba + bi[k]; 385 /* initialize k-th row by the perm[k]-th row of A */ 386 jmin = ai[rip[k]]; jmax = ai[rip[k]+1]; 387 for (j = jmin; j < jmax; j++){ 388 col = rip[aj[j]]; 389 if (col >= k){ /* only take upper triangular entry */ 390 rtmp[col] = aa[j]; 391 *bval++ = 0.0; /* for in-place factorization */ 392 } 393 } 394 395 /* shift the diagonal of the matrix */ 396 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 397 398 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 399 dk = rtmp[k]; 400 i = jl[k]; /* first row to be added to k_th row */ 401 402 while (i < k){ 403 nexti = jl[i]; /* next row to be added to k_th row */ 404 405 /* compute multiplier, update diag(k) and U(i,k) */ 406 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 407 uikdi = - ba[ili]*ba[bi[i]]; /* diagonal(k) */ 408 dk += uikdi*ba[ili]; 409 ba[ili] = uikdi; /* -U(i,k) */ 410 411 /* add multiple of row i to k-th row */ 412 jmin = ili + 1; jmax = bi[i+1]; 413 if (jmin < jmax){ 414 for (j=jmin; j<jmax; j++) rtmp[bj[j]] += uikdi*ba[j]; 415 /* update il and jl for row i */ 416 il[i] = jmin; 417 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 418 } 419 i = nexti; 420 } 421 422 /* shift the diagonals when zero pivot is detected */ 423 /* compute rs=sum of abs(off-diagonal) */ 424 rs = 0.0; 425 jmin = bi[k]+1; 426 nz = bi[k+1] - jmin; 427 if (nz){ 428 bcol = bj + jmin; 429 while (nz--){ 430 rs += PetscAbsScalar(rtmp[*bcol]); 431 bcol++; 432 } 433 } 434 435 sctx.rs = rs; 436 sctx.pv = dk; 437 ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 438 if (newshift == 1) break; 439 440 /* copy data into U(k,:) */ 441 ba[bi[k]] = 1.0/dk; /* U(k,k) */ 442 jmin = bi[k]+1; jmax = bi[k+1]; 443 if (jmin < jmax) { 444 for (j=jmin; j<jmax; j++){ 445 col = bj[j]; ba[j] = rtmp[col]; rtmp[col] = 0.0; 446 } 447 /* add the k-th row into il and jl */ 448 il[k] = jmin; 449 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 450 } 451 } 452 } while (sctx.chshift); 453 ierr = PetscFree(il);CHKERRQ(ierr); 454 455 ierr = ISRestoreIndices(ip,&rip);CHKERRQ(ierr); 456 C->assembled = PETSC_TRUE; 457 C->preallocated = PETSC_TRUE; 458 ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 459 if (sctx.nshift){ 460 if (shiftpd) { 461 ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 462 } else if (shiftnz) { 463 ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 464 } 465 } 466 PetscFunctionReturn(0); 467 } 468 469 #undef __FUNCT__ 470 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering" 471 PetscErrorCode MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) 472 { 473 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 474 Mat_SeqSBAIJ *b=(Mat_SeqSBAIJ*)C->data; 475 PetscErrorCode ierr; 476 PetscInt i,j,am=a->mbs; 477 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 478 PetscInt k,jmin,*jl,*il,nexti,ili,*acol,*bcol,nz; 479 MatScalar *rtmp,*ba=b->a,*aa=a->a,dk,uikdi,*aval,*bval; 480 PetscReal zeropivot,rs,shiftnz; 481 PetscReal shiftpd; 482 ChShift_Ctx sctx; 483 PetscInt newshift; 484 485 PetscFunctionBegin; 486 /* initialization */ 487 shiftnz = info->shiftnz; 488 shiftpd = info->shiftpd; 489 zeropivot = info->zeropivot; 490 491 nz = (2*am+1)*sizeof(PetscInt)+am*sizeof(MatScalar); 492 ierr = PetscMalloc(nz,&il);CHKERRQ(ierr); 493 jl = il + am; 494 rtmp = (MatScalar*)(jl + am); 495 496 sctx.shift_amount = 0; 497 sctx.nshift = 0; 498 do { 499 sctx.chshift = PETSC_FALSE; 500 for (i=0; i<am; i++) { 501 rtmp[i] = 0.0; jl[i] = am; il[0] = 0; 502 } 503 504 for (k = 0; k<am; k++){ 505 /* initialize k-th row with elements nonzero in row perm(k) of A */ 506 nz = ai[k+1] - ai[k]; 507 acol = aj + ai[k]; 508 aval = aa + ai[k]; 509 bval = ba + bi[k]; 510 while (nz -- ){ 511 if (*acol < k) { /* skip lower triangular entries */ 512 acol++; aval++; 513 } else { 514 rtmp[*acol++] = *aval++; 515 *bval++ = 0.0; /* for in-place factorization */ 516 } 517 } 518 519 /* shift the diagonal of the matrix */ 520 if (sctx.nshift) rtmp[k] += sctx.shift_amount; 521 522 /* modify k-th row by adding in those rows i with U(i,k)!=0 */ 523 dk = rtmp[k]; 524 i = jl[k]; /* first row to be added to k_th row */ 525 526 while (i < k){ 527 nexti = jl[i]; /* next row to be added to k_th row */ 528 /* compute multiplier, update D(k) and U(i,k) */ 529 ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 530 uikdi = - ba[ili]*ba[bi[i]]; 531 dk += uikdi*ba[ili]; 532 ba[ili] = uikdi; /* -U(i,k) */ 533 534 /* add multiple of row i to k-th row ... */ 535 jmin = ili + 1; 536 nz = bi[i+1] - jmin; 537 if (nz > 0){ 538 bcol = bj + jmin; 539 bval = ba + jmin; 540 while (nz --) rtmp[*bcol++] += uikdi*(*bval++); 541 /* update il and jl for i-th row */ 542 il[i] = jmin; 543 j = bj[jmin]; jl[i] = jl[j]; jl[j] = i; 544 } 545 i = nexti; 546 } 547 548 /* shift the diagonals when zero pivot is detected */ 549 /* compute rs=sum of abs(off-diagonal) */ 550 rs = 0.0; 551 jmin = bi[k]+1; 552 nz = bi[k+1] - jmin; 553 if (nz){ 554 bcol = bj + jmin; 555 while (nz--){ 556 rs += PetscAbsScalar(rtmp[*bcol]); 557 bcol++; 558 } 559 } 560 561 sctx.rs = rs; 562 sctx.pv = dk; 563 ierr = MatCholeskyCheckShift_inline(info,sctx,k,newshift);CHKERRQ(ierr); 564 if (newshift == 1) break; /* sctx.shift_amount is updated */ 565 566 /* copy data into U(k,:) */ 567 ba[bi[k]] = 1.0/dk; 568 jmin = bi[k]+1; 569 nz = bi[k+1] - jmin; 570 if (nz){ 571 bcol = bj + jmin; 572 bval = ba + jmin; 573 while (nz--){ 574 *bval++ = rtmp[*bcol]; 575 rtmp[*bcol++] = 0.0; 576 } 577 /* add k-th row into il and jl */ 578 il[k] = jmin; 579 i = bj[jmin]; jl[k] = jl[i]; jl[i] = k; 580 } 581 } 582 } while (sctx.chshift); 583 ierr = PetscFree(il);CHKERRQ(ierr); 584 585 C->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 586 C->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 587 C->assembled = PETSC_TRUE; 588 C->preallocated = PETSC_TRUE; 589 ierr = PetscLogFlops(C->rmap->N);CHKERRQ(ierr); 590 if (sctx.nshift){ 591 if (shiftnz) { 592 ierr = PetscInfo2(A,"number of shiftnz tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 593 } else if (shiftpd) { 594 ierr = PetscInfo2(A,"number of shiftpd tries %D, shift_amount %G\n",sctx.nshift,sctx.shift_amount);CHKERRQ(ierr); 595 } 596 } 597 PetscFunctionReturn(0); 598 } 599 600 #include "petscbt.h" 601 #include "../src/mat/utils/freespace.h" 602 #undef __FUNCT__ 603 #define __FUNCT__ "MatICCFactorSymbolic_SeqBAIJ" 604 PetscErrorCode MatICCFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 605 { 606 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 607 Mat_SeqSBAIJ *b; 608 Mat B; 609 PetscErrorCode ierr; 610 PetscTruth perm_identity; 611 PetscInt reallocs=0,i,*ai=a->i,*aj=a->j,am=a->mbs,bs=A->rmap->bs,*ui; 612 const PetscInt *rip; 613 PetscInt jmin,jmax,nzk,k,j,*jl,prow,*il,nextprow; 614 PetscInt nlnk,*lnk,*lnk_lvl=PETSC_NULL,ncols,ncols_upper,*cols,*cols_lvl,*uj,**uj_ptr,**uj_lvl_ptr; 615 PetscReal fill=info->fill,levels=info->levels; 616 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 617 PetscFreeSpaceList free_space_lvl=PETSC_NULL,current_space_lvl=PETSC_NULL; 618 PetscBT lnkbt; 619 620 PetscFunctionBegin; 621 if (bs > 1){ 622 if (!a->sbaijMat){ 623 ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 624 } 625 (fact)->ops->iccfactorsymbolic = MatICCFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */ 626 ierr = MatICCFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr); 627 PetscFunctionReturn(0); 628 } 629 630 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 631 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 632 633 /* special case that simply copies fill pattern */ 634 if (!levels && perm_identity) { 635 ierr = MatMarkDiagonal_SeqBAIJ(A);CHKERRQ(ierr); 636 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 637 for (i=0; i<am; i++) { 638 ui[i] = ai[i+1] - a->diag[i]; /* ui: rowlengths - changes when !perm_identity */ 639 } 640 B = fact; 641 ierr = MatSeqSBAIJSetPreallocation(B,1,0,ui);CHKERRQ(ierr); 642 643 644 b = (Mat_SeqSBAIJ*)B->data; 645 uj = b->j; 646 for (i=0; i<am; i++) { 647 aj = a->j + a->diag[i]; 648 for (j=0; j<ui[i]; j++){ 649 *uj++ = *aj++; 650 } 651 b->ilen[i] = ui[i]; 652 } 653 ierr = PetscFree(ui);CHKERRQ(ierr); 654 B->factor = MAT_FACTOR_NONE; 655 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 656 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 657 B->factor = MAT_FACTOR_ICC; 658 659 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 660 PetscFunctionReturn(0); 661 } 662 663 /* initialization */ 664 ierr = PetscMalloc((am+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 665 ui[0] = 0; 666 ierr = PetscMalloc((2*am+1)*sizeof(PetscInt),&cols_lvl);CHKERRQ(ierr); 667 668 /* jl: linked list for storing indices of the pivot rows 669 il: il[i] points to the 1st nonzero entry of U(i,k:am-1) */ 670 ierr = PetscMalloc((2*am+1)*sizeof(PetscInt)+2*am*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 671 il = jl + am; 672 uj_ptr = (PetscInt**)(il + am); 673 uj_lvl_ptr = (PetscInt**)(uj_ptr + am); 674 for (i=0; i<am; i++){ 675 jl[i] = am; il[i] = 0; 676 } 677 678 /* create and initialize a linked list for storing column indices of the active row k */ 679 nlnk = am + 1; 680 ierr = PetscIncompleteLLCreate(am,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 681 682 /* initial FreeSpace size is fill*(ai[am]+1) */ 683 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space);CHKERRQ(ierr); 684 current_space = free_space; 685 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[am]+1)),&free_space_lvl);CHKERRQ(ierr); 686 current_space_lvl = free_space_lvl; 687 688 for (k=0; k<am; k++){ /* for each active row k */ 689 /* initialize lnk by the column indices of row rip[k] of A */ 690 nzk = 0; 691 ncols = ai[rip[k]+1] - ai[rip[k]]; 692 ncols_upper = 0; 693 cols = cols_lvl + am; 694 for (j=0; j<ncols; j++){ 695 i = rip[*(aj + ai[rip[k]] + j)]; 696 if (i >= k){ /* only take upper triangular entry */ 697 cols[ncols_upper] = i; 698 cols_lvl[ncols_upper] = -1; /* initialize level for nonzero entries */ 699 ncols_upper++; 700 } 701 } 702 ierr = PetscIncompleteLLAdd(ncols_upper,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 703 nzk += nlnk; 704 705 /* update lnk by computing fill-in for each pivot row to be merged in */ 706 prow = jl[k]; /* 1st pivot row */ 707 708 while (prow < k){ 709 nextprow = jl[prow]; 710 711 /* merge prow into k-th row */ 712 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:am-1) */ 713 jmax = ui[prow+1]; 714 ncols = jmax-jmin; 715 i = jmin - ui[prow]; 716 cols = uj_ptr[prow] + i; /* points to the 2nd nzero entry in U(prow,k:am-1) */ 717 for (j=0; j<ncols; j++) cols_lvl[j] = *(uj_lvl_ptr[prow] + i + j); 718 ierr = PetscIncompleteLLAddSorted(ncols,cols,levels,cols_lvl,am,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 719 nzk += nlnk; 720 721 /* update il and jl for prow */ 722 if (jmin < jmax){ 723 il[prow] = jmin; 724 j = *cols; jl[prow] = jl[j]; jl[j] = prow; 725 } 726 prow = nextprow; 727 } 728 729 /* if free space is not available, make more free space */ 730 if (current_space->local_remaining<nzk) { 731 i = am - k + 1; /* num of unfactored rows */ 732 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 733 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 734 ierr = PetscFreeSpaceGet(i,¤t_space_lvl);CHKERRQ(ierr); 735 reallocs++; 736 } 737 738 /* copy data into free_space and free_space_lvl, then initialize lnk */ 739 ierr = PetscIncompleteLLClean(am,am,nzk,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 740 741 /* add the k-th row into il and jl */ 742 if (nzk-1 > 0){ 743 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:am-1) */ 744 jl[k] = jl[i]; jl[i] = k; 745 il[k] = ui[k] + 1; 746 } 747 uj_ptr[k] = current_space->array; 748 uj_lvl_ptr[k] = current_space_lvl->array; 749 750 current_space->array += nzk; 751 current_space->local_used += nzk; 752 current_space->local_remaining -= nzk; 753 754 current_space_lvl->array += nzk; 755 current_space_lvl->local_used += nzk; 756 current_space_lvl->local_remaining -= nzk; 757 758 ui[k+1] = ui[k] + nzk; 759 } 760 761 #if defined(PETSC_USE_INFO) 762 if (ai[am] != 0) { 763 PetscReal af = ((PetscReal)(2*ui[am]-am))/((PetscReal)ai[am]); 764 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 765 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 766 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 767 } else { 768 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 769 } 770 #endif 771 772 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 773 ierr = PetscFree(jl);CHKERRQ(ierr); 774 ierr = PetscFree(cols_lvl);CHKERRQ(ierr); 775 776 /* destroy list of free space and other temporary array(s) */ 777 ierr = PetscMalloc((ui[am]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 778 ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 779 ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 780 ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 781 782 /* put together the new matrix in MATSEQSBAIJ format */ 783 B = fact; 784 ierr = MatSeqSBAIJSetPreallocation(B,1,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 785 786 b = (Mat_SeqSBAIJ*)B->data; 787 b->singlemalloc = PETSC_FALSE; 788 b->free_a = PETSC_TRUE; 789 b->free_ij = PETSC_TRUE; 790 ierr = PetscMalloc((ui[am]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 791 b->j = uj; 792 b->i = ui; 793 b->diag = 0; 794 b->ilen = 0; 795 b->imax = 0; 796 b->row = perm; 797 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 798 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 799 b->icol = perm; 800 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 801 ierr = PetscMalloc((am+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 802 ierr = PetscLogObjectMemory(B,(ui[am]-am)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 803 b->maxnz = b->nz = ui[am]; 804 805 B->info.factor_mallocs = reallocs; 806 B->info.fill_ratio_given = fill; 807 if (ai[am] != 0) { 808 B->info.fill_ratio_needed = ((PetscReal)ui[am])/((PetscReal)ai[am]); 809 } else { 810 B->info.fill_ratio_needed = 0.0; 811 } 812 if (perm_identity){ 813 B->ops->solve = MatSolve_SeqSBAIJ_1_NaturalOrdering; 814 B->ops->solvetranspose = MatSolve_SeqSBAIJ_1_NaturalOrdering; 815 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 816 } else { 817 (fact)->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 818 } 819 PetscFunctionReturn(0); 820 } 821 822 #undef __FUNCT__ 823 #define __FUNCT__ "MatCholeskyFactorSymbolic_SeqBAIJ" 824 PetscErrorCode MatCholeskyFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS perm,const MatFactorInfo *info) 825 { 826 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 827 Mat_SeqSBAIJ *b; 828 Mat B; 829 PetscErrorCode ierr; 830 PetscTruth perm_identity; 831 PetscReal fill = info->fill; 832 const PetscInt *rip; 833 PetscInt i,mbs=a->mbs,bs=A->rmap->bs,*ai=a->i,*aj=a->j,reallocs=0,prow; 834 PetscInt *jl,jmin,jmax,nzk,*ui,k,j,*il,nextprow; 835 PetscInt nlnk,*lnk,ncols,ncols_upper,*cols,*uj,**ui_ptr,*uj_ptr; 836 PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 837 PetscBT lnkbt; 838 839 PetscFunctionBegin; 840 if (bs > 1) { /* convert to seqsbaij */ 841 if (!a->sbaijMat){ 842 ierr = MatConvert(A,MATSEQSBAIJ,MAT_INITIAL_MATRIX,&a->sbaijMat);CHKERRQ(ierr); 843 } 844 (fact)->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_SeqSBAIJ; /* undue the change made in MatGetFactor_seqbaij_petsc */ 845 ierr = MatCholeskyFactorSymbolic(fact,a->sbaijMat,perm,info);CHKERRQ(ierr); 846 PetscFunctionReturn(0); 847 } 848 849 /* check whether perm is the identity mapping */ 850 ierr = ISIdentity(perm,&perm_identity);CHKERRQ(ierr); 851 if (!perm_identity) SETERRQ(PETSC_ERR_SUP,"Matrix reordering is not supported"); 852 ierr = ISGetIndices(perm,&rip);CHKERRQ(ierr); 853 854 /* initialization */ 855 ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&ui);CHKERRQ(ierr); 856 ui[0] = 0; 857 858 /* jl: linked list for storing indices of the pivot rows 859 il: il[i] points to the 1st nonzero entry of U(i,k:mbs-1) */ 860 ierr = PetscMalloc((3*mbs+1)*sizeof(PetscInt)+mbs*sizeof(PetscInt*),&jl);CHKERRQ(ierr); 861 il = jl + mbs; 862 cols = il + mbs; 863 ui_ptr = (PetscInt**)(cols + mbs); 864 for (i=0; i<mbs; i++){ 865 jl[i] = mbs; il[i] = 0; 866 } 867 868 /* create and initialize a linked list for storing column indices of the active row k */ 869 nlnk = mbs + 1; 870 ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 871 872 /* initial FreeSpace size is fill*(ai[mbs]+1) */ 873 ierr = PetscFreeSpaceGet((PetscInt)(fill*(ai[mbs]+1)),&free_space);CHKERRQ(ierr); 874 current_space = free_space; 875 876 for (k=0; k<mbs; k++){ /* for each active row k */ 877 /* initialize lnk by the column indices of row rip[k] of A */ 878 nzk = 0; 879 ncols = ai[rip[k]+1] - ai[rip[k]]; 880 ncols_upper = 0; 881 for (j=0; j<ncols; j++){ 882 i = rip[*(aj + ai[rip[k]] + j)]; 883 if (i >= k){ /* only take upper triangular entry */ 884 cols[ncols_upper] = i; 885 ncols_upper++; 886 } 887 } 888 ierr = PetscLLAdd(ncols_upper,cols,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 889 nzk += nlnk; 890 891 /* update lnk by computing fill-in for each pivot row to be merged in */ 892 prow = jl[k]; /* 1st pivot row */ 893 894 while (prow < k){ 895 nextprow = jl[prow]; 896 /* merge prow into k-th row */ 897 jmin = il[prow] + 1; /* index of the 2nd nzero entry in U(prow,k:mbs-1) */ 898 jmax = ui[prow+1]; 899 ncols = jmax-jmin; 900 uj_ptr = ui_ptr[prow] + jmin - ui[prow]; /* points to the 2nd nzero entry in U(prow,k:mbs-1) */ 901 ierr = PetscLLAddSorted(ncols,uj_ptr,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 902 nzk += nlnk; 903 904 /* update il and jl for prow */ 905 if (jmin < jmax){ 906 il[prow] = jmin; 907 j = *uj_ptr; jl[prow] = jl[j]; jl[j] = prow; 908 } 909 prow = nextprow; 910 } 911 912 /* if free space is not available, make more free space */ 913 if (current_space->local_remaining<nzk) { 914 i = mbs - k + 1; /* num of unfactored rows */ 915 i = PetscMin(i*nzk, i*(i-1)); /* i*nzk, i*(i-1): estimated and max additional space needed */ 916 ierr = PetscFreeSpaceGet(i,¤t_space);CHKERRQ(ierr); 917 reallocs++; 918 } 919 920 /* copy data into free space, then initialize lnk */ 921 ierr = PetscLLClean(mbs,mbs,nzk,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 922 923 /* add the k-th row into il and jl */ 924 if (nzk-1 > 0){ 925 i = current_space->array[1]; /* col value of the first nonzero element in U(k, k+1:mbs-1) */ 926 jl[k] = jl[i]; jl[i] = k; 927 il[k] = ui[k] + 1; 928 } 929 ui_ptr[k] = current_space->array; 930 current_space->array += nzk; 931 current_space->local_used += nzk; 932 current_space->local_remaining -= nzk; 933 934 ui[k+1] = ui[k] + nzk; 935 } 936 937 #if defined(PETSC_USE_INFO) 938 if (ai[mbs] != 0) { 939 PetscReal af = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 940 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %G needed %G\n",reallocs,fill,af);CHKERRQ(ierr); 941 ierr = PetscInfo1(A,"Run with -pc_factor_fill %G or use \n",af);CHKERRQ(ierr); 942 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%G) for best performance.\n",af);CHKERRQ(ierr); 943 } else { 944 ierr = PetscInfo(A,"Empty matrix.\n");CHKERRQ(ierr); 945 } 946 #endif 947 948 ierr = ISRestoreIndices(perm,&rip);CHKERRQ(ierr); 949 ierr = PetscFree(jl);CHKERRQ(ierr); 950 951 /* destroy list of free space and other temporary array(s) */ 952 ierr = PetscMalloc((ui[mbs]+1)*sizeof(PetscInt),&uj);CHKERRQ(ierr); 953 ierr = PetscFreeSpaceContiguous(&free_space,uj);CHKERRQ(ierr); 954 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 955 956 /* put together the new matrix in MATSEQSBAIJ format */ 957 B = fact; 958 ierr = MatSeqSBAIJSetPreallocation(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 959 960 b = (Mat_SeqSBAIJ*)B->data; 961 b->singlemalloc = PETSC_FALSE; 962 b->free_a = PETSC_TRUE; 963 b->free_ij = PETSC_TRUE; 964 ierr = PetscMalloc((ui[mbs]+1)*sizeof(MatScalar),&b->a);CHKERRQ(ierr); 965 b->j = uj; 966 b->i = ui; 967 b->diag = 0; 968 b->ilen = 0; 969 b->imax = 0; 970 b->row = perm; 971 b->pivotinblocks = PETSC_FALSE; /* need to get from MatFactorInfo */ 972 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 973 b->icol = perm; 974 ierr = PetscObjectReference((PetscObject)perm);CHKERRQ(ierr); 975 ierr = PetscMalloc((mbs+1)*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 976 ierr = PetscLogObjectMemory(B,(ui[mbs]-mbs)*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 977 b->maxnz = b->nz = ui[mbs]; 978 979 B->info.factor_mallocs = reallocs; 980 B->info.fill_ratio_given = fill; 981 if (ai[mbs] != 0) { 982 B->info.fill_ratio_needed = ((PetscReal)ui[mbs])/((PetscReal)ai[mbs]); 983 } else { 984 B->info.fill_ratio_needed = 0.0; 985 } 986 if (perm_identity){ 987 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N_NaturalOrdering; 988 } else { 989 B->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_SeqBAIJ_N; 990 } 991 PetscFunctionReturn(0); 992 } 993 994 /* --------------------------------------------------------- */ 995 #undef __FUNCT__ 996 #define __FUNCT__ "MatSolve_SeqBAIJ_N_NaturalOrdering_newdatastruct" 997 PetscErrorCode MatSolve_SeqBAIJ_N_NaturalOrdering_newdatastruct(Mat A,Vec bb,Vec xx) 998 { 999 Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 1000 PetscErrorCode ierr; 1001 const PetscInt *ai=a->i,*aj=a->j,*vi; 1002 PetscInt i,n=a->mbs; 1003 PetscInt nz,bs=A->rmap->bs,bs2=a->bs2,*adiag=a->diag; 1004 MatScalar *aa=a->a,*v; 1005 PetscScalar *x,*b,*s,*t,*ls; 1006 1007 PetscFunctionBegin; 1008 /* printf("MatSolve_SeqBAIJ_NaturalOrdering_iludt..bs %d\n",bs); */ 1009 ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1010 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1011 t = a->solve_work; 1012 1013 /* forward solve the lower triangular */ 1014 ierr = PetscMemcpy(t,b,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy 1st block of b to t */ 1015 1016 for (i=1; i<n; i++) { 1017 v = aa + bs2*ai[i]; 1018 vi = aj + ai[i]; 1019 nz = ai[i+1] - ai[i]; 1020 s = t + bs*i; 1021 ierr = PetscMemcpy(s,b+bs*i,bs*sizeof(PetscScalar));CHKERRQ(ierr); /* copy i_th block of b to t */ 1022 while (nz--) { 1023 Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*(*vi++)); 1024 v += bs2; 1025 } 1026 } 1027 1028 /* backward solve the upper triangular */ 1029 ls = a->solve_work + A->cmap->n; 1030 for (i=n-1; i>=0; i--){ 1031 v = aa + bs2*ai[2*n-i]; 1032 vi = aj + ai[2*n-i]; 1033 nz = ai[2*n-i +1] - ai[2*n-i]-1; 1034 ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 1035 while (nz--) { 1036 Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*(*vi++)); 1037 v += bs2; 1038 } 1039 Kernel_w_gets_A_times_v(bs,ls,aa+bs2*adiag[i],t+i*bs); /* *inv(diagonal[i]) */ 1040 ierr = PetscMemcpy(x+i*bs,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 1041 } 1042 1043 ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 1044 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1045 ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr); 1046 PetscFunctionReturn(0); 1047 } 1048 1049 #undef __FUNCT__ 1050 #define __FUNCT__ "MatSolve_SeqBAIJ_N_newdatastruct" 1051 PetscErrorCode MatSolve_SeqBAIJ_N_newdatastruct(Mat A,Vec bb,Vec xx) 1052 { 1053 Mat_SeqBAIJ *a=(Mat_SeqBAIJ *)A->data; 1054 IS iscol=a->col,isrow=a->row; 1055 PetscErrorCode ierr; 1056 const PetscInt *r,*c,*rout,*cout,*ai=a->i,*aj=a->j,*vi; 1057 PetscInt i,n=a->mbs; 1058 PetscInt nz,bs=A->rmap->bs,bs2=a->bs2,*adiag=a->diag; 1059 MatScalar *aa=a->a,*v; 1060 PetscScalar *x,*b,*s,*t,*ls; 1061 1062 PetscFunctionBegin; 1063 ierr = VecGetArray(bb,&b);CHKERRQ(ierr); 1064 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 1065 t = a->solve_work; 1066 1067 ierr = ISGetIndices(isrow,&rout);CHKERRQ(ierr); r = rout; 1068 ierr = ISGetIndices(iscol,&cout);CHKERRQ(ierr); c = cout + (n-1); 1069 1070 /* forward solve the lower triangular */ 1071 ierr = PetscMemcpy(t,b+bs*(*r++),bs*sizeof(PetscScalar));CHKERRQ(ierr); 1072 for (i=1; i<n; i++) { 1073 v = aa + bs2*ai[i]; 1074 vi = aj + ai[i]; 1075 nz = ai[i+1] - ai[i]; 1076 s = t + bs*i; 1077 ierr = PetscMemcpy(s,b+bs*(*r++),bs*sizeof(PetscScalar));CHKERRQ(ierr); 1078 while (nz--) { 1079 Kernel_v_gets_v_minus_A_times_w(bs,s,v,t+bs*(*vi++)); 1080 v += bs2; 1081 } 1082 } 1083 1084 /* backward solve the upper triangular */ 1085 ls = a->solve_work + A->cmap->n; 1086 for (i=n-1; i>=0; i--){ 1087 v = aa + bs2*(adiag[i+1] + 1); 1088 vi = aj + a->diag[i+1] + 1; 1089 nz = adiag[i] - adiag[i+1] - 1; 1090 ierr = PetscMemcpy(ls,t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 1091 while (nz--) { 1092 Kernel_v_gets_v_minus_A_times_w(bs,ls,v,t+bs*(*vi++)); 1093 v += bs2; 1094 } 1095 Kernel_w_gets_A_times_v(bs,ls,aa+bs2*adiag[i],t+i*bs); /* *inv(diagonal[i]) */ 1096 ierr = PetscMemcpy(x + bs*(*c--),t+i*bs,bs*sizeof(PetscScalar));CHKERRQ(ierr); 1097 } 1098 ierr = ISRestoreIndices(isrow,&rout);CHKERRQ(ierr); 1099 ierr = ISRestoreIndices(iscol,&cout);CHKERRQ(ierr); 1100 ierr = VecRestoreArray(bb,&b);CHKERRQ(ierr); 1101 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 1102 ierr = PetscLogFlops(2.0*(a->bs2)*(a->nz) - A->rmap->bs*A->cmap->n);CHKERRQ(ierr); 1103 PetscFunctionReturn(0); 1104 } 1105 1106 #undef __FUNCT__ 1107 #define __FUNCT__ "BlockAbs_privat" 1108 PetscErrorCode BlockAbs_private(PetscInt nbs,PetscInt bs2,PetscScalar *blockarray,PetscReal *absarray) 1109 { 1110 PetscErrorCode ierr; 1111 PetscInt i,j; 1112 PetscFunctionBegin; 1113 ierr = PetscMemzero(absarray,(nbs+1)*sizeof(PetscReal));CHKERRQ(ierr); 1114 for (i=0; i<nbs; i++){ 1115 for (j=0; j<bs2; j++){ 1116 if (absarray[i] < PetscAbsScalar(blockarray[i*nbs+j])) absarray[i] = PetscAbsScalar(blockarray[i*nbs+j]); 1117 } 1118 } 1119 PetscFunctionReturn(0); 1120 } 1121 1122 #undef __FUNCT__ 1123 #define __FUNCT__ "MatILUDTFactor_SeqBAIJ" 1124 PetscErrorCode MatILUDTFactor_SeqBAIJ(Mat A,IS isrow,IS iscol,const MatFactorInfo *info,Mat *fact) 1125 { 1126 Mat B = *fact; 1127 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b; 1128 IS isicol; 1129 PetscErrorCode ierr; 1130 const PetscInt *r,*ic; 1131 PetscInt i,mbs=a->mbs,bs=A->rmap->bs,bs2=a->bs2,*ai=a->i,*aj=a->j,*ajtmp,*adiag; 1132 PetscInt *bi,*bj,*bdiag; 1133 1134 PetscInt row,nzi,nzi_bl,nzi_bu,*im,dtcount,nzi_al,nzi_au; 1135 PetscInt nlnk,*lnk; 1136 PetscBT lnkbt; 1137 PetscTruth row_identity,icol_identity,both_identity; 1138 MatScalar *aatmp,*pv,*batmp,*ba,*rtmp,*pc,*multiplier,*vtmp; 1139 const PetscInt *ics; 1140 PetscInt j,nz,*pj,*bjtmp,k,ncut,*jtmp; 1141 1142 PetscReal dt=info->dt; /* shift=info->shiftinblocks; */ 1143 PetscInt nnz_max; 1144 PetscTruth missing; 1145 PetscReal *vtmp_abs; 1146 MatScalar *v_work; 1147 PetscInt *v_pivots; 1148 1149 PetscFunctionBegin; 1150 /* ------- symbolic factorization, can be reused ---------*/ 1151 ierr = MatMissingDiagonal(A,&missing,&i);CHKERRQ(ierr); 1152 if (missing) SETERRQ1(PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",i); 1153 adiag=a->diag; 1154 1155 ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 1156 1157 /* bdiag is location of diagonal in factor */ 1158 ierr = PetscMalloc((mbs+1)*sizeof(PetscInt),&bdiag);CHKERRQ(ierr); 1159 1160 /* allocate row pointers bi */ 1161 ierr = PetscMalloc((2*mbs+2)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 1162 1163 /* allocate bj and ba; max num of nonzero entries is (ai[n]+2*n*dtcount+2) */ 1164 dtcount = (PetscInt)info->dtcount; 1165 if (dtcount > mbs-1) dtcount = mbs-1; 1166 nnz_max = ai[mbs]+2*mbs*dtcount +2; 1167 /* printf("MatILUDTFactor_SeqBAIJ, bs %d, ai[mbs] %d, nnz_max %d, dtcount %d\n",bs,ai[mbs],nnz_max,dtcount); */ 1168 ierr = PetscMalloc(nnz_max*sizeof(PetscInt),&bj);CHKERRQ(ierr); 1169 nnz_max = nnz_max*bs2; 1170 ierr = PetscMalloc(nnz_max*sizeof(MatScalar),&ba);CHKERRQ(ierr); 1171 1172 /* put together the new matrix */ 1173 ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(B,bs,MAT_SKIP_ALLOCATION,PETSC_NULL);CHKERRQ(ierr); 1174 ierr = PetscLogObjectParent(B,isicol);CHKERRQ(ierr); 1175 b = (Mat_SeqBAIJ*)(B)->data; 1176 b->free_a = PETSC_TRUE; 1177 b->free_ij = PETSC_TRUE; 1178 b->singlemalloc = PETSC_FALSE; 1179 b->a = ba; 1180 b->j = bj; 1181 b->i = bi; 1182 b->diag = bdiag; 1183 b->ilen = 0; 1184 b->imax = 0; 1185 b->row = isrow; 1186 b->col = iscol; 1187 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 1188 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 1189 b->icol = isicol; 1190 ierr = PetscMalloc((bs*(mbs+1))*sizeof(PetscScalar),&b->solve_work);CHKERRQ(ierr); 1191 1192 ierr = PetscLogObjectMemory(B,nnz_max*(sizeof(PetscInt)+sizeof(MatScalar)));CHKERRQ(ierr); 1193 b->maxnz = nnz_max/bs2; 1194 1195 (B)->factor = MAT_FACTOR_ILUDT; 1196 (B)->info.factor_mallocs = 0; 1197 (B)->info.fill_ratio_given = ((PetscReal)nnz_max)/((PetscReal)(ai[mbs]*bs2)); 1198 CHKMEMQ; 1199 /* ------- end of symbolic factorization ---------*/ 1200 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 1201 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 1202 ics = ic; 1203 1204 /* linked list for storing column indices of the active row */ 1205 nlnk = mbs + 1; 1206 ierr = PetscLLCreate(mbs,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1207 1208 /* im: used by PetscLLAddSortedLU(); jtmp: working array for column indices of active row */ 1209 ierr = PetscMalloc((2*mbs+1)*sizeof(PetscInt),&im);CHKERRQ(ierr); 1210 jtmp = im + mbs; 1211 /* rtmp, vtmp: working arrays for sparse and contiguous row entries of active row */ 1212 ierr = PetscMalloc((2*mbs*bs2+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 1213 vtmp = rtmp + bs2*mbs; 1214 ierr = PetscMalloc((mbs+1)*sizeof(PetscReal),&vtmp_abs);CHKERRQ(ierr); 1215 1216 ierr = PetscMalloc(bs*sizeof(PetscInt) + (bs+bs2)*sizeof(MatScalar),&v_work);CHKERRQ(ierr); 1217 multiplier = v_work + bs; 1218 v_pivots = (PetscInt*)(multiplier + bs2); 1219 1220 bi[0] = 0; 1221 bdiag[0] = (nnz_max/bs2)-1; /* location of diagonal in factor B */ 1222 bi[2*mbs+1] = bdiag[0]+1; /* endof bj and ba array */ 1223 for (i=0; i<mbs; i++) { 1224 /* copy initial fill into linked list */ 1225 nzi = 0; /* nonzeros for active row i */ 1226 nzi = ai[r[i]+1] - ai[r[i]]; 1227 if (!nzi) SETERRQ2(PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 1228 nzi_al = adiag[r[i]] - ai[r[i]]; 1229 nzi_au = ai[r[i]+1] - adiag[r[i]] -1; 1230 /* printf("row %d, nzi_al/au %d %d\n",i,nzi_al,nzi_au); */ 1231 1232 /* load in initial unfactored row */ 1233 ajtmp = aj + ai[r[i]]; 1234 ierr = PetscLLAddPerm(nzi,ajtmp,ic,mbs,nlnk,lnk,lnkbt);CHKERRQ(ierr); 1235 ierr = PetscMemzero(rtmp,(mbs*bs2+1)*sizeof(PetscScalar));CHKERRQ(ierr); 1236 aatmp = a->a + bs2*ai[r[i]]; 1237 for (j=0; j<nzi; j++) { 1238 ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],aatmp+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1239 } 1240 1241 /* add pivot rows into linked list */ 1242 row = lnk[mbs]; 1243 while (row < i) { 1244 nzi_bl = bi[row+1] - bi[row] + 1; 1245 bjtmp = bj + bdiag[row+1]+1; /* points to 1st column next to the diagonal in U */ 1246 ierr = PetscLLAddSortedLU(bjtmp,row,nlnk,lnk,lnkbt,i,nzi_bl,im);CHKERRQ(ierr); 1247 nzi += nlnk; 1248 row = lnk[row]; 1249 } 1250 1251 /* copy data from lnk into jtmp, then initialize lnk */ 1252 ierr = PetscLLClean(mbs,mbs,nzi,lnk,jtmp,lnkbt);CHKERRQ(ierr); 1253 1254 /* numerical factorization */ 1255 bjtmp = jtmp; 1256 row = *bjtmp++; /* 1st pivot row */ 1257 1258 while (row < i) { 1259 pc = rtmp + bs2*row; 1260 pv = ba + bs2*bdiag[row]; /* inv(diag) of the pivot row */ 1261 Kernel_A_gets_A_times_B(bs,pc,pv,multiplier); /* pc= multiplier = pc*inv(diag[row]) */ 1262 ierr = BlockAbs_private(1,bs2,pc,vtmp_abs);CHKERRQ(ierr); 1263 if (vtmp_abs[0] > dt){ /* apply tolerance dropping rule */ 1264 pj = bj + bdiag[row+1] + 1; /* point to 1st entry of U(row,:) */ 1265 pv = ba + bs2*(bdiag[row+1] + 1); 1266 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:), excluding diagonal */ 1267 for (j=0; j<nz; j++){ 1268 Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); 1269 } 1270 /* ierr = PetscLogFlops(bslog*(nz+1.0)-bs);CHKERRQ(ierr); */ 1271 } 1272 row = *bjtmp++; 1273 } 1274 1275 /* copy sparse rtmp into contiguous vtmp; separate L and U part */ 1276 nzi_bl = 0; j = 0; 1277 while (jtmp[j] < i){ /* L-part. Note: jtmp is sorted */ 1278 ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1279 nzi_bl++; j++; 1280 } 1281 nzi_bu = nzi - nzi_bl -1; 1282 /* printf("nzi %d, nzi_bl %d, nzi_bu %d\n",nzi,nzi_bl,nzi_bu); */ 1283 1284 while (j < nzi){ /* U-part */ 1285 ierr = PetscMemcpy(vtmp+bs2*j,rtmp+bs2*jtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1286 /* 1287 printf(" col %d: ",jtmp[j]); 1288 for (j1=0; j1<bs2; j1++) printf(" %g",*(vtmp+bs2*j+j1)); 1289 printf(" \n"); 1290 */ 1291 j++; 1292 } 1293 1294 ierr = BlockAbs_private(nzi,bs2,vtmp,vtmp_abs);CHKERRQ(ierr); 1295 /* 1296 printf(" row %d, nzi %d, vtmp_abs\n",i,nzi); 1297 for (j1=0; j1<nzi; j1++) printf(" (%d %g),",jtmp[j1],vtmp_abs[j1]); 1298 printf(" \n"); 1299 */ 1300 bjtmp = bj + bi[i]; 1301 batmp = ba + bs2*bi[i]; 1302 /* apply level dropping rule to L part */ 1303 ncut = nzi_al + dtcount; 1304 if (ncut < nzi_bl){ 1305 ierr = PetscSortSplitReal(ncut,nzi_bl,vtmp_abs,jtmp);CHKERRQ(ierr); 1306 ierr = PetscSortIntWithScalarArray(ncut,jtmp,vtmp);CHKERRQ(ierr); 1307 } else { 1308 ncut = nzi_bl; 1309 } 1310 for (j=0; j<ncut; j++){ 1311 bjtmp[j] = jtmp[j]; 1312 ierr = PetscMemcpy(batmp+bs2*j,rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1313 /* 1314 printf(" col %d: ",bjtmp[j]); 1315 for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*j+j1)); 1316 printf("\n"); 1317 */ 1318 } 1319 bi[i+1] = bi[i] + ncut; 1320 nzi = ncut + 1; 1321 1322 /* apply level dropping rule to U part */ 1323 ncut = nzi_au + dtcount; 1324 if (ncut < nzi_bu){ 1325 ierr = PetscSortSplitReal(ncut,nzi_bu,vtmp_abs+nzi_bl+1,jtmp+nzi_bl+1);CHKERRQ(ierr); 1326 ierr = PetscSortIntWithScalarArray(ncut,jtmp+nzi_bl+1,vtmp+nzi_bl+1);CHKERRQ(ierr); 1327 } else { 1328 ncut = nzi_bu; 1329 } 1330 nzi += ncut; 1331 1332 /* mark bdiagonal */ 1333 bdiag[i+1] = bdiag[i] - (ncut + 1); 1334 bi[2*mbs - i] = bi[2*mbs - i +1] - (ncut + 1); 1335 1336 bjtmp = bj + bdiag[i]; 1337 batmp = ba + bs2*bdiag[i]; 1338 ierr = PetscMemcpy(batmp,rtmp+bs2*i,bs2*sizeof(MatScalar));CHKERRQ(ierr); 1339 *bjtmp = i; 1340 /* 1341 printf(" diag %d: ",*bjtmp); 1342 for (j=0; j<bs2; j++){ 1343 printf(" %g,",batmp[j]); 1344 } 1345 printf("\n"); 1346 */ 1347 bjtmp = bj + bdiag[i+1]+1; 1348 batmp = ba + (bdiag[i+1]+1)*bs2; 1349 1350 for (k=0; k<ncut; k++){ 1351 bjtmp[k] = jtmp[nzi_bl+1+k]; 1352 ierr = PetscMemcpy(batmp+bs2*k,rtmp+bs2*bjtmp[k],bs2*sizeof(MatScalar));CHKERRQ(ierr); 1353 /* 1354 printf(" col %d:",bjtmp[k]); 1355 for (j1=0; j1<bs2; j1++) printf(" %g,",*(batmp+bs2*k+j1)); 1356 printf("\n"); 1357 */ 1358 } 1359 1360 im[i] = nzi; /* used by PetscLLAddSortedLU() */ 1361 1362 /* invert diagonal block for simplier triangular solves - add shift??? */ 1363 batmp = ba + bs2*bdiag[i]; 1364 ierr = Kernel_A_gets_inverse_A(bs,batmp,v_pivots,v_work);CHKERRQ(ierr); 1365 } /* for (i=0; i<mbs; i++) */ 1366 1367 /* printf("end of L %d, beginning of U %d\n",bi[mbs],bdiag[mbs]); */ 1368 if (bi[mbs] >= bdiag[mbs]) SETERRQ2(PETSC_ERR_ARG_SIZ,"end of L array %d cannot >= the beginning of U array %d",bi[mbs],bdiag[mbs]); 1369 1370 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 1371 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 1372 1373 ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1374 1375 ierr = PetscFree(im);CHKERRQ(ierr); 1376 ierr = PetscFree(rtmp);CHKERRQ(ierr); 1377 1378 ierr = PetscLogFlops(bs2*B->cmap->n);CHKERRQ(ierr); 1379 b->maxnz = b->nz = bi[mbs] + bdiag[0] - bdiag[mbs]; 1380 1381 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 1382 ierr = ISIdentity(isicol,&icol_identity);CHKERRQ(ierr); 1383 both_identity = (PetscTruth) (row_identity && icol_identity); 1384 if (row_identity && icol_identity) { 1385 B->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering_newdatastruct; 1386 } else { 1387 B->ops->solve = MatSolve_SeqBAIJ_N_newdatastruct; 1388 } 1389 1390 B->ops->solveadd = 0; 1391 B->ops->solvetranspose = 0; 1392 B->ops->solvetransposeadd = 0; 1393 B->ops->matsolve = 0; 1394 B->assembled = PETSC_TRUE; 1395 B->preallocated = PETSC_TRUE; 1396 PetscFunctionReturn(0); 1397 } 1398