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