1 2 /* 3 Factorization code for BAIJ format. 4 */ 5 6 #include <../src/mat/impls/baij/seq/baij.h> 7 #include <petsc/private/kernels/blockinvert.h> 8 #include <petscbt.h> 9 #include <../src/mat/utils/freespace.h> 10 11 /* ----------------------------------------------------------------*/ 12 extern PetscErrorCode MatDuplicateNoCreate_SeqBAIJ(Mat,Mat,MatDuplicateOption,PetscBool); 13 14 #undef __FUNCT__ 15 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering" 16 /* 17 This is not much faster than MatLUFactorNumeric_SeqBAIJ_N() but the solve is faster at least sometimes 18 */ 19 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_15_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info) 20 { 21 Mat C =B; 22 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data; 23 PetscErrorCode ierr; 24 PetscInt i,j,k,ipvt[15]; 25 const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j,*ajtmp,*bjtmp,*bdiag=b->diag,*pj; 26 PetscInt nz,nzL,row; 27 MatScalar *rtmp,*pc,*mwork,*pv,*vv,work[225]; 28 const MatScalar *v,*aa=a->a; 29 PetscInt bs2 = a->bs2,bs=A->rmap->bs,flg; 30 PetscInt sol_ver; 31 PetscBool zeropivotdetected; 32 33 PetscFunctionBegin; 34 ierr = PetscOptionsGetInt(NULL,((PetscObject)A)->prefix,"-sol_ver",&sol_ver,NULL);CHKERRQ(ierr); 35 36 /* generate work space needed by the factorization */ 37 ierr = PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);CHKERRQ(ierr); 38 ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 39 40 for (i=0; i<n; i++) { 41 /* zero rtmp */ 42 /* L part */ 43 nz = bi[i+1] - bi[i]; 44 bjtmp = bj + bi[i]; 45 for (j=0; j<nz; j++) { 46 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 47 } 48 49 /* U part */ 50 nz = bdiag[i] - bdiag[i+1]; 51 bjtmp = bj + bdiag[i+1]+1; 52 for (j=0; j<nz; j++) { 53 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 54 } 55 56 /* load in initial (unfactored row) */ 57 nz = ai[i+1] - ai[i]; 58 ajtmp = aj + ai[i]; 59 v = aa + bs2*ai[i]; 60 for (j=0; j<nz; j++) { 61 ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 62 } 63 64 /* elimination */ 65 bjtmp = bj + bi[i]; 66 nzL = bi[i+1] - bi[i]; 67 for (k=0; k < nzL; k++) { 68 row = bjtmp[k]; 69 pc = rtmp + bs2*row; 70 for (flg=0,j=0; j<bs2; j++) { 71 if (pc[j]!=0.0) { 72 flg = 1; 73 break; 74 } 75 } 76 if (flg) { 77 pv = b->a + bs2*bdiag[row]; 78 PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); 79 /*ierr = PetscKernel_A_gets_A_times_B_15(pc,pv,mwork);CHKERRQ(ierr);*/ 80 pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */ 81 pv = b->a + bs2*(bdiag[row+1]+1); 82 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */ 83 for (j=0; j<nz; j++) { 84 vv = rtmp + bs2*pj[j]; 85 PetscKernel_A_gets_A_minus_B_times_C(bs,vv,pc,pv); 86 /* ierr = PetscKernel_A_gets_A_minus_B_times_C_15(vv,pc,pv);CHKERRQ(ierr); */ 87 pv += bs2; 88 } 89 ierr = PetscLogFlops(2*bs2*bs*(nz+1)-bs2);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 90 } 91 } 92 93 /* finished row so stick it into b->a */ 94 /* L part */ 95 pv = b->a + bs2*bi[i]; 96 pj = b->j + bi[i]; 97 nz = bi[i+1] - bi[i]; 98 for (j=0; j<nz; j++) { 99 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 100 } 101 102 /* Mark diagonal and invert diagonal for simplier triangular solves */ 103 pv = b->a + bs2*bdiag[i]; 104 pj = b->j + bdiag[i]; 105 ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 106 /* PetscKernel_A_gets_inverse_A(bs,pv,pivots,work); */ 107 ierr = PetscKernel_A_gets_inverse_A_15(pv,ipvt,work,info->shiftamount,!A->erroriffailure,&zeropivotdetected);CHKERRQ(ierr); 108 if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 109 110 /* U part */ 111 pv = b->a + bs2*(bdiag[i+1]+1); 112 pj = b->j + bdiag[i+1]+1; 113 nz = bdiag[i] - bdiag[i+1] - 1; 114 for (j=0; j<nz; j++) { 115 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 116 } 117 } 118 119 ierr = PetscFree2(rtmp,mwork);CHKERRQ(ierr); 120 121 C->ops->solve = MatSolve_SeqBAIJ_15_NaturalOrdering_ver1; 122 C->ops->solvetranspose = MatSolve_SeqBAIJ_N_NaturalOrdering; 123 C->assembled = PETSC_TRUE; 124 125 ierr = PetscLogFlops(1.333333333333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 126 PetscFunctionReturn(0); 127 } 128 129 #undef __FUNCT__ 130 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_N" 131 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_N(Mat B,Mat A,const MatFactorInfo *info) 132 { 133 Mat C =B; 134 Mat_SeqBAIJ *a =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data; 135 IS isrow = b->row,isicol = b->icol; 136 PetscErrorCode ierr; 137 const PetscInt *r,*ic; 138 PetscInt i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 139 PetscInt *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj; 140 MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 141 PetscInt bs=A->rmap->bs,bs2 = a->bs2,*v_pivots,flg; 142 MatScalar *v_work; 143 PetscBool col_identity,row_identity,both_identity; 144 145 PetscFunctionBegin; 146 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 147 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 148 149 ierr = PetscMalloc1(bs2*n,&rtmp);CHKERRQ(ierr); 150 ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 151 152 /* generate work space needed by dense LU factorization */ 153 ierr = PetscMalloc3(bs,&v_work,bs2,&mwork,bs,&v_pivots);CHKERRQ(ierr); 154 155 for (i=0; i<n; i++) { 156 /* zero rtmp */ 157 /* L part */ 158 nz = bi[i+1] - bi[i]; 159 bjtmp = bj + bi[i]; 160 for (j=0; j<nz; j++) { 161 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 162 } 163 164 /* U part */ 165 nz = bdiag[i] - bdiag[i+1]; 166 bjtmp = bj + bdiag[i+1]+1; 167 for (j=0; j<nz; j++) { 168 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 169 } 170 171 /* load in initial (unfactored row) */ 172 nz = ai[r[i]+1] - ai[r[i]]; 173 ajtmp = aj + ai[r[i]]; 174 v = aa + bs2*ai[r[i]]; 175 for (j=0; j<nz; j++) { 176 ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 177 } 178 179 /* elimination */ 180 bjtmp = bj + bi[i]; 181 nzL = bi[i+1] - bi[i]; 182 for (k=0; k < nzL; k++) { 183 row = bjtmp[k]; 184 pc = rtmp + bs2*row; 185 for (flg=0,j=0; j<bs2; j++) { 186 if (pc[j]!=0.0) { 187 flg = 1; 188 break; 189 } 190 } 191 if (flg) { 192 pv = b->a + bs2*bdiag[row]; 193 PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); /* *pc = *pc * (*pv); */ 194 pj = b->j + bdiag[row+1]+1; /* begining of U(row,:) */ 195 pv = b->a + bs2*(bdiag[row+1]+1); 196 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */ 197 for (j=0; j<nz; j++) { 198 PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); 199 } 200 ierr = PetscLogFlops(2*bs2*bs*(nz+1)-bs2);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 201 } 202 } 203 204 /* finished row so stick it into b->a */ 205 /* L part */ 206 pv = b->a + bs2*bi[i]; 207 pj = b->j + bi[i]; 208 nz = bi[i+1] - bi[i]; 209 for (j=0; j<nz; j++) { 210 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 211 } 212 213 /* Mark diagonal and invert diagonal for simplier triangular solves */ 214 pv = b->a + bs2*bdiag[i]; 215 pj = b->j + bdiag[i]; 216 /* if (*pj != i)SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"row %d != *pj %d",i,*pj); */ 217 ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 218 ierr = PetscKernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); 219 220 /* U part */ 221 pv = b->a + bs2*(bdiag[i+1]+1); 222 pj = b->j + bdiag[i+1]+1; 223 nz = bdiag[i] - bdiag[i+1] - 1; 224 for (j=0; j<nz; j++) { 225 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 226 } 227 } 228 229 ierr = PetscFree(rtmp);CHKERRQ(ierr); 230 ierr = PetscFree3(v_work,mwork,v_pivots);CHKERRQ(ierr); 231 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 232 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 233 234 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 235 ierr = ISIdentity(isicol,&col_identity);CHKERRQ(ierr); 236 237 both_identity = (PetscBool) (row_identity && col_identity); 238 if (both_identity) { 239 C->ops->solve = MatSolve_SeqBAIJ_N_NaturalOrdering; 240 } else { 241 C->ops->solve = MatSolve_SeqBAIJ_N; 242 } 243 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_N; 244 245 C->assembled = PETSC_TRUE; 246 247 ierr = PetscLogFlops(1.333333333333*bs*bs2*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 248 PetscFunctionReturn(0); 249 } 250 251 /* 252 ilu(0) with natural ordering under new data structure. 253 See MatILUFactorSymbolic_SeqAIJ_ilu0() for detailed description 254 because this code is almost identical to MatILUFactorSymbolic_SeqAIJ_ilu0_inplace(). 255 */ 256 257 #undef __FUNCT__ 258 #define __FUNCT__ "MatILUFactorSymbolic_SeqBAIJ_ilu0" 259 PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_ilu0(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 260 { 261 262 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 263 PetscErrorCode ierr; 264 PetscInt n=a->mbs,*ai=a->i,*aj,*adiag=a->diag,bs2 = a->bs2; 265 PetscInt i,j,nz,*bi,*bj,*bdiag,bi_temp; 266 267 PetscFunctionBegin; 268 ierr = MatDuplicateNoCreate_SeqBAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_FALSE);CHKERRQ(ierr); 269 b = (Mat_SeqBAIJ*)(fact)->data; 270 271 /* allocate matrix arrays for new data structure */ 272 ierr = PetscMalloc3(bs2*ai[n]+1,&b->a,ai[n]+1,&b->j,n+1,&b->i);CHKERRQ(ierr); 273 ierr = PetscLogObjectMemory((PetscObject)fact,ai[n]*(bs2*sizeof(PetscScalar)+sizeof(PetscInt))+(n+1)*sizeof(PetscInt));CHKERRQ(ierr); 274 275 b->singlemalloc = PETSC_TRUE; 276 b->free_a = PETSC_TRUE; 277 b->free_ij = PETSC_TRUE; 278 fact->preallocated = PETSC_TRUE; 279 fact->assembled = PETSC_TRUE; 280 if (!b->diag) { 281 ierr = PetscMalloc1(n+1,&b->diag);CHKERRQ(ierr); 282 ierr = PetscLogObjectMemory((PetscObject)fact,(n+1)*sizeof(PetscInt));CHKERRQ(ierr); 283 } 284 bdiag = b->diag; 285 286 if (n > 0) { 287 ierr = PetscMemzero(b->a,bs2*ai[n]*sizeof(MatScalar));CHKERRQ(ierr); 288 } 289 290 /* set bi and bj with new data structure */ 291 bi = b->i; 292 bj = b->j; 293 294 /* L part */ 295 bi[0] = 0; 296 for (i=0; i<n; i++) { 297 nz = adiag[i] - ai[i]; 298 bi[i+1] = bi[i] + nz; 299 aj = a->j + ai[i]; 300 for (j=0; j<nz; j++) { 301 *bj = aj[j]; bj++; 302 } 303 } 304 305 /* U part */ 306 bi_temp = bi[n]; 307 bdiag[n] = bi[n]-1; 308 for (i=n-1; i>=0; i--) { 309 nz = ai[i+1] - adiag[i] - 1; 310 bi_temp = bi_temp + nz + 1; 311 aj = a->j + adiag[i] + 1; 312 for (j=0; j<nz; j++) { 313 *bj = aj[j]; bj++; 314 } 315 /* diag[i] */ 316 *bj = i; bj++; 317 bdiag[i] = bi_temp - 1; 318 } 319 PetscFunctionReturn(0); 320 } 321 322 #undef __FUNCT__ 323 #define __FUNCT__ "MatILUFactorSymbolic_SeqBAIJ" 324 PetscErrorCode MatILUFactorSymbolic_SeqBAIJ(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 325 { 326 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 327 IS isicol; 328 PetscErrorCode ierr; 329 const PetscInt *r,*ic; 330 PetscInt n=a->mbs,*ai=a->i,*aj=a->j,d; 331 PetscInt *bi,*cols,nnz,*cols_lvl; 332 PetscInt *bdiag,prow,fm,nzbd,reallocs=0,dcount=0; 333 PetscInt i,levels,diagonal_fill; 334 PetscBool col_identity,row_identity,both_identity; 335 PetscReal f; 336 PetscInt nlnk,*lnk,*lnk_lvl=NULL; 337 PetscBT lnkbt; 338 PetscInt nzi,*bj,**bj_ptr,**bjlvl_ptr; 339 PetscFreeSpaceList free_space =NULL,current_space=NULL; 340 PetscFreeSpaceList free_space_lvl=NULL,current_space_lvl=NULL; 341 PetscBool missing; 342 PetscInt bs=A->rmap->bs,bs2=a->bs2; 343 344 PetscFunctionBegin; 345 if (A->rmap->n != A->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Must be square matrix, rows %D columns %D",A->rmap->n,A->cmap->n); 346 if (bs>1) { /* check shifttype */ 347 if (info->shifttype == MAT_SHIFT_NONZERO || info->shifttype == MAT_SHIFT_POSITIVE_DEFINITE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Only MAT_SHIFT_NONE and MAT_SHIFT_INBLOCKS are supported for BAIJ matrix"); 348 } 349 350 ierr = MatMissingDiagonal(A,&missing,&d);CHKERRQ(ierr); 351 if (missing) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix is missing diagonal entry %D",d); 352 353 f = info->fill; 354 levels = (PetscInt)info->levels; 355 diagonal_fill = (PetscInt)info->diagonal_fill; 356 357 ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 358 359 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 360 ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 361 362 both_identity = (PetscBool) (row_identity && col_identity); 363 364 if (!levels && both_identity) { 365 /* special case: ilu(0) with natural ordering */ 366 ierr = MatILUFactorSymbolic_SeqBAIJ_ilu0(fact,A,isrow,iscol,info);CHKERRQ(ierr); 367 ierr = MatSeqBAIJSetNumericFactorization(fact,both_identity);CHKERRQ(ierr); 368 369 fact->factortype = MAT_FACTOR_ILU; 370 (fact)->info.factor_mallocs = 0; 371 (fact)->info.fill_ratio_given = info->fill; 372 (fact)->info.fill_ratio_needed = 1.0; 373 374 b = (Mat_SeqBAIJ*)(fact)->data; 375 b->row = isrow; 376 b->col = iscol; 377 b->icol = isicol; 378 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 379 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 380 b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 381 382 ierr = PetscMalloc1((n+1)*bs,&b->solve_work);CHKERRQ(ierr); 383 PetscFunctionReturn(0); 384 } 385 386 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 387 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 388 389 /* get new row pointers */ 390 ierr = PetscMalloc1(n+1,&bi);CHKERRQ(ierr); 391 bi[0] = 0; 392 /* bdiag is location of diagonal in factor */ 393 ierr = PetscMalloc1(n+1,&bdiag);CHKERRQ(ierr); 394 bdiag[0] = 0; 395 396 ierr = PetscMalloc2(n,&bj_ptr,n,&bjlvl_ptr);CHKERRQ(ierr); 397 398 /* create a linked list for storing column indices of the active row */ 399 nlnk = n + 1; 400 ierr = PetscIncompleteLLCreate(n,n,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 401 402 /* initial FreeSpace size is f*(ai[n]+1) */ 403 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space);CHKERRQ(ierr); 404 current_space = free_space; 405 ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(f,ai[n]+1),&free_space_lvl);CHKERRQ(ierr); 406 current_space_lvl = free_space_lvl; 407 408 for (i=0; i<n; i++) { 409 nzi = 0; 410 /* copy current row into linked list */ 411 nnz = ai[r[i]+1] - ai[r[i]]; 412 if (!nnz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[i],i); 413 cols = aj + ai[r[i]]; 414 lnk[i] = -1; /* marker to indicate if diagonal exists */ 415 ierr = PetscIncompleteLLInit(nnz,cols,n,ic,nlnk,lnk,lnk_lvl,lnkbt);CHKERRQ(ierr); 416 nzi += nlnk; 417 418 /* make sure diagonal entry is included */ 419 if (diagonal_fill && lnk[i] == -1) { 420 fm = n; 421 while (lnk[fm] < i) fm = lnk[fm]; 422 lnk[i] = lnk[fm]; /* insert diagonal into linked list */ 423 lnk[fm] = i; 424 lnk_lvl[i] = 0; 425 nzi++; dcount++; 426 } 427 428 /* add pivot rows into the active row */ 429 nzbd = 0; 430 prow = lnk[n]; 431 while (prow < i) { 432 nnz = bdiag[prow]; 433 cols = bj_ptr[prow] + nnz + 1; 434 cols_lvl = bjlvl_ptr[prow] + nnz + 1; 435 nnz = bi[prow+1] - bi[prow] - nnz - 1; 436 437 ierr = PetscILULLAddSorted(nnz,cols,levels,cols_lvl,prow,nlnk,lnk,lnk_lvl,lnkbt,prow);CHKERRQ(ierr); 438 nzi += nlnk; 439 prow = lnk[prow]; 440 nzbd++; 441 } 442 bdiag[i] = nzbd; 443 bi[i+1] = bi[i] + nzi; 444 445 /* if free space is not available, make more free space */ 446 if (current_space->local_remaining<nzi) { 447 nnz = PetscIntMultTruncate(2,PetscIntMultTruncate(nzi,(n - i))); /* estimated and max additional space needed */ 448 ierr = PetscFreeSpaceGet(nnz,¤t_space);CHKERRQ(ierr); 449 ierr = PetscFreeSpaceGet(nnz,¤t_space_lvl);CHKERRQ(ierr); 450 reallocs++; 451 } 452 453 /* copy data into free_space and free_space_lvl, then initialize lnk */ 454 ierr = PetscIncompleteLLClean(n,n,nzi,lnk,lnk_lvl,current_space->array,current_space_lvl->array,lnkbt);CHKERRQ(ierr); 455 456 bj_ptr[i] = current_space->array; 457 bjlvl_ptr[i] = current_space_lvl->array; 458 459 /* make sure the active row i has diagonal entry */ 460 if (*(bj_ptr[i]+bdiag[i]) != i) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\ntry running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",i); 461 462 current_space->array += nzi; 463 current_space->local_used += nzi; 464 current_space->local_remaining -= nzi; 465 466 current_space_lvl->array += nzi; 467 current_space_lvl->local_used += nzi; 468 current_space_lvl->local_remaining -= nzi; 469 } 470 471 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 472 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 473 474 /* copy free_space into bj and free free_space; set bi, bj, bdiag in new datastructure; */ 475 ierr = PetscMalloc1(bi[n]+1,&bj);CHKERRQ(ierr); 476 ierr = PetscFreeSpaceContiguous_LU(&free_space,bj,n,bi,bdiag);CHKERRQ(ierr); 477 478 ierr = PetscIncompleteLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 479 ierr = PetscFreeSpaceDestroy(free_space_lvl);CHKERRQ(ierr); 480 ierr = PetscFree2(bj_ptr,bjlvl_ptr);CHKERRQ(ierr); 481 482 #if defined(PETSC_USE_INFO) 483 { 484 PetscReal af = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]); 485 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocs,(double)f,(double)af);CHKERRQ(ierr); 486 ierr = PetscInfo1(A,"Run with -[sub_]pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr); 487 ierr = PetscInfo1(A,"PCFactorSetFill([sub]pc,%g);\n",(double)af);CHKERRQ(ierr); 488 ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 489 if (diagonal_fill) { 490 ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals\n",dcount);CHKERRQ(ierr); 491 } 492 } 493 #endif 494 495 /* put together the new matrix */ 496 ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 497 ierr = PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol);CHKERRQ(ierr); 498 499 b = (Mat_SeqBAIJ*)(fact)->data; 500 b->free_a = PETSC_TRUE; 501 b->free_ij = PETSC_TRUE; 502 b->singlemalloc = PETSC_FALSE; 503 504 ierr = PetscMalloc1(bs2*(bdiag[0]+1),&b->a);CHKERRQ(ierr); 505 506 b->j = bj; 507 b->i = bi; 508 b->diag = bdiag; 509 b->free_diag = PETSC_TRUE; 510 b->ilen = 0; 511 b->imax = 0; 512 b->row = isrow; 513 b->col = iscol; 514 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 515 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 516 b->icol = isicol; 517 518 ierr = PetscMalloc1(bs*n+bs,&b->solve_work);CHKERRQ(ierr); 519 /* In b structure: Free imax, ilen, old a, old j. 520 Allocate bdiag, solve_work, new a, new j */ 521 ierr = PetscLogObjectMemory((PetscObject)fact,(bdiag[0]+1) * (sizeof(PetscInt)+bs2*sizeof(PetscScalar)));CHKERRQ(ierr); 522 b->maxnz = b->nz = bdiag[0]+1; 523 524 fact->info.factor_mallocs = reallocs; 525 fact->info.fill_ratio_given = f; 526 fact->info.fill_ratio_needed = ((PetscReal)(bdiag[0]+1))/((PetscReal)ai[n]); 527 528 ierr = MatSeqBAIJSetNumericFactorization(fact,both_identity);CHKERRQ(ierr); 529 PetscFunctionReturn(0); 530 } 531 532 /* 533 This code is virtually identical to MatILUFactorSymbolic_SeqAIJ 534 except that the data structure of Mat_SeqAIJ is slightly different. 535 Not a good example of code reuse. 536 */ 537 #undef __FUNCT__ 538 #define __FUNCT__ "MatILUFactorSymbolic_SeqBAIJ_inplace" 539 PetscErrorCode MatILUFactorSymbolic_SeqBAIJ_inplace(Mat fact,Mat A,IS isrow,IS iscol,const MatFactorInfo *info) 540 { 541 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b; 542 IS isicol; 543 PetscErrorCode ierr; 544 const PetscInt *r,*ic,*ai = a->i,*aj = a->j,*xi; 545 PetscInt prow,n = a->mbs,*ainew,*ajnew,jmax,*fill,nz,*im,*ajfill,*flev,*xitmp; 546 PetscInt *dloc,idx,row,m,fm,nzf,nzi,reallocate = 0,dcount = 0; 547 PetscInt incrlev,nnz,i,bs = A->rmap->bs,bs2 = a->bs2,levels,diagonal_fill,dd; 548 PetscBool col_identity,row_identity,both_identity,flg; 549 PetscReal f; 550 551 PetscFunctionBegin; 552 ierr = MatMissingDiagonal_SeqBAIJ(A,&flg,&dd);CHKERRQ(ierr); 553 if (flg) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Matrix A is missing diagonal entry in row %D",dd); 554 555 f = info->fill; 556 levels = (PetscInt)info->levels; 557 diagonal_fill = (PetscInt)info->diagonal_fill; 558 559 ierr = ISInvertPermutation(iscol,PETSC_DECIDE,&isicol);CHKERRQ(ierr); 560 561 ierr = ISIdentity(isrow,&row_identity);CHKERRQ(ierr); 562 ierr = ISIdentity(iscol,&col_identity);CHKERRQ(ierr); 563 both_identity = (PetscBool) (row_identity && col_identity); 564 565 if (!levels && both_identity) { /* special case copy the nonzero structure */ 566 ierr = MatDuplicateNoCreate_SeqBAIJ(fact,A,MAT_DO_NOT_COPY_VALUES,PETSC_TRUE);CHKERRQ(ierr); 567 ierr = MatSeqBAIJSetNumericFactorization_inplace(fact,both_identity);CHKERRQ(ierr); 568 569 fact->factortype = MAT_FACTOR_ILU; 570 b = (Mat_SeqBAIJ*)fact->data; 571 b->row = isrow; 572 b->col = iscol; 573 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 574 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 575 b->icol = isicol; 576 b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 577 578 ierr = PetscMalloc1((n+1)*bs,&b->solve_work);CHKERRQ(ierr); 579 PetscFunctionReturn(0); 580 } 581 582 /* general case perform the symbolic factorization */ 583 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 584 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 585 586 /* get new row pointers */ 587 ierr = PetscMalloc1(n+1,&ainew);CHKERRQ(ierr); 588 ainew[0] = 0; 589 /* don't know how many column pointers are needed so estimate */ 590 jmax = (PetscInt)(f*ai[n] + 1); 591 ierr = PetscMalloc1(jmax,&ajnew);CHKERRQ(ierr); 592 /* ajfill is level of fill for each fill entry */ 593 ierr = PetscMalloc1(jmax,&ajfill);CHKERRQ(ierr); 594 /* fill is a linked list of nonzeros in active row */ 595 ierr = PetscMalloc1(n+1,&fill);CHKERRQ(ierr); 596 /* im is level for each filled value */ 597 ierr = PetscMalloc1(n+1,&im);CHKERRQ(ierr); 598 /* dloc is location of diagonal in factor */ 599 ierr = PetscMalloc1(n+1,&dloc);CHKERRQ(ierr); 600 dloc[0] = 0; 601 for (prow=0; prow<n; prow++) { 602 603 /* copy prow into linked list */ 604 nzf = nz = ai[r[prow]+1] - ai[r[prow]]; 605 if (!nz) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Empty row in matrix: row in original ordering %D in permuted ordering %D",r[prow],prow); 606 xi = aj + ai[r[prow]]; 607 fill[n] = n; 608 fill[prow] = -1; /* marker for diagonal entry */ 609 while (nz--) { 610 fm = n; 611 idx = ic[*xi++]; 612 do { 613 m = fm; 614 fm = fill[m]; 615 } while (fm < idx); 616 fill[m] = idx; 617 fill[idx] = fm; 618 im[idx] = 0; 619 } 620 621 /* make sure diagonal entry is included */ 622 if (diagonal_fill && fill[prow] == -1) { 623 fm = n; 624 while (fill[fm] < prow) fm = fill[fm]; 625 fill[prow] = fill[fm]; /* insert diagonal into linked list */ 626 fill[fm] = prow; 627 im[prow] = 0; 628 nzf++; 629 dcount++; 630 } 631 632 nzi = 0; 633 row = fill[n]; 634 while (row < prow) { 635 incrlev = im[row] + 1; 636 nz = dloc[row]; 637 xi = ajnew + ainew[row] + nz + 1; 638 flev = ajfill + ainew[row] + nz + 1; 639 nnz = ainew[row+1] - ainew[row] - nz - 1; 640 fm = row; 641 while (nnz-- > 0) { 642 idx = *xi++; 643 if (*flev + incrlev > levels) { 644 flev++; 645 continue; 646 } 647 do { 648 m = fm; 649 fm = fill[m]; 650 } while (fm < idx); 651 if (fm != idx) { 652 im[idx] = *flev + incrlev; 653 fill[m] = idx; 654 fill[idx] = fm; 655 fm = idx; 656 nzf++; 657 } else if (im[idx] > *flev + incrlev) im[idx] = *flev+incrlev; 658 flev++; 659 } 660 row = fill[row]; 661 nzi++; 662 } 663 /* copy new filled row into permanent storage */ 664 ainew[prow+1] = ainew[prow] + nzf; 665 if (ainew[prow+1] > jmax) { 666 667 /* estimate how much additional space we will need */ 668 /* use the strategy suggested by David Hysom <hysom@perch-t.icase.edu> */ 669 /* just double the memory each time */ 670 PetscInt maxadd = jmax; 671 /* maxadd = (int)(((f*ai[n]+1)*(n-prow+5))/n); */ 672 if (maxadd < nzf) maxadd = (n-prow)*(nzf+1); 673 jmax += maxadd; 674 675 /* allocate a longer ajnew and ajfill */ 676 ierr = PetscMalloc1(jmax,&xitmp);CHKERRQ(ierr); 677 ierr = PetscMemcpy(xitmp,ajnew,ainew[prow]*sizeof(PetscInt));CHKERRQ(ierr); 678 ierr = PetscFree(ajnew);CHKERRQ(ierr); 679 ajnew = xitmp; 680 ierr = PetscMalloc1(jmax,&xitmp);CHKERRQ(ierr); 681 ierr = PetscMemcpy(xitmp,ajfill,ainew[prow]*sizeof(PetscInt));CHKERRQ(ierr); 682 ierr = PetscFree(ajfill);CHKERRQ(ierr); 683 ajfill = xitmp; 684 reallocate++; /* count how many reallocations are needed */ 685 } 686 xitmp = ajnew + ainew[prow]; 687 flev = ajfill + ainew[prow]; 688 dloc[prow] = nzi; 689 fm = fill[n]; 690 while (nzf--) { 691 *xitmp++ = fm; 692 *flev++ = im[fm]; 693 fm = fill[fm]; 694 } 695 /* make sure row has diagonal entry */ 696 if (ajnew[ainew[prow]+dloc[prow]] != prow) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_MAT_LU_ZRPVT,"Row %D has missing diagonal in factored matrix\n\ 697 try running with -pc_factor_nonzeros_along_diagonal or -pc_factor_diagonal_fill",prow); 698 } 699 ierr = PetscFree(ajfill);CHKERRQ(ierr); 700 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 701 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 702 ierr = PetscFree(fill);CHKERRQ(ierr); 703 ierr = PetscFree(im);CHKERRQ(ierr); 704 705 #if defined(PETSC_USE_INFO) 706 { 707 PetscReal af = ((PetscReal)ainew[n])/((PetscReal)ai[n]); 708 ierr = PetscInfo3(A,"Reallocs %D Fill ratio:given %g needed %g\n",reallocate,(double)f,(double)af);CHKERRQ(ierr); 709 ierr = PetscInfo1(A,"Run with -pc_factor_fill %g or use \n",(double)af);CHKERRQ(ierr); 710 ierr = PetscInfo1(A,"PCFactorSetFill(pc,%g);\n",(double)af);CHKERRQ(ierr); 711 ierr = PetscInfo(A,"for best performance.\n");CHKERRQ(ierr); 712 if (diagonal_fill) { 713 ierr = PetscInfo1(A,"Detected and replaced %D missing diagonals\n",dcount);CHKERRQ(ierr); 714 } 715 } 716 #endif 717 718 /* put together the new matrix */ 719 ierr = MatSeqBAIJSetPreallocation_SeqBAIJ(fact,bs,MAT_SKIP_ALLOCATION,NULL);CHKERRQ(ierr); 720 ierr = PetscLogObjectParent((PetscObject)fact,(PetscObject)isicol);CHKERRQ(ierr); 721 b = (Mat_SeqBAIJ*)fact->data; 722 723 b->free_a = PETSC_TRUE; 724 b->free_ij = PETSC_TRUE; 725 b->singlemalloc = PETSC_FALSE; 726 727 ierr = PetscMalloc1(bs2*ainew[n],&b->a);CHKERRQ(ierr); 728 729 b->j = ajnew; 730 b->i = ainew; 731 for (i=0; i<n; i++) dloc[i] += ainew[i]; 732 b->diag = dloc; 733 b->free_diag = PETSC_TRUE; 734 b->ilen = 0; 735 b->imax = 0; 736 b->row = isrow; 737 b->col = iscol; 738 b->pivotinblocks = (info->pivotinblocks) ? PETSC_TRUE : PETSC_FALSE; 739 740 ierr = PetscObjectReference((PetscObject)isrow);CHKERRQ(ierr); 741 ierr = PetscObjectReference((PetscObject)iscol);CHKERRQ(ierr); 742 b->icol = isicol; 743 ierr = PetscMalloc1(bs*n+bs,&b->solve_work);CHKERRQ(ierr); 744 /* In b structure: Free imax, ilen, old a, old j. 745 Allocate dloc, solve_work, new a, new j */ 746 ierr = PetscLogObjectMemory((PetscObject)fact,(ainew[n]-n)*(sizeof(PetscInt))+bs2*ainew[n]*sizeof(PetscScalar));CHKERRQ(ierr); 747 b->maxnz = b->nz = ainew[n]; 748 749 fact->info.factor_mallocs = reallocate; 750 fact->info.fill_ratio_given = f; 751 fact->info.fill_ratio_needed = ((PetscReal)ainew[n])/((PetscReal)ai[prow]); 752 753 ierr = MatSeqBAIJSetNumericFactorization_inplace(fact,both_identity);CHKERRQ(ierr); 754 PetscFunctionReturn(0); 755 } 756 757 #undef __FUNCT__ 758 #define __FUNCT__ "MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE" 759 PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE(Mat A) 760 { 761 /* Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; */ 762 /* int i,*AJ=a->j,nz=a->nz; */ 763 764 PetscFunctionBegin; 765 /* Undo Column scaling */ 766 /* while (nz--) { */ 767 /* AJ[i] = AJ[i]/4; */ 768 /* } */ 769 /* This should really invoke a push/pop logic, but we don't have that yet. */ 770 A->ops->setunfactored = NULL; 771 PetscFunctionReturn(0); 772 } 773 774 #undef __FUNCT__ 775 #define __FUNCT__ "MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj" 776 PetscErrorCode MatSetUnfactored_SeqBAIJ_4_NaturalOrdering_SSE_usj(Mat A) 777 { 778 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 779 PetscInt *AJ=a->j,nz=a->nz; 780 unsigned short *aj=(unsigned short*)AJ; 781 782 PetscFunctionBegin; 783 /* Is this really necessary? */ 784 while (nz--) { 785 AJ[nz] = (int)((unsigned int)aj[nz]); /* First extend, then convert to signed. */ 786 } 787 A->ops->setunfactored = NULL; 788 PetscFunctionReturn(0); 789 } 790 791 792