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 Version for when blocks are 3 by 3 11 */ 12 #undef __FUNCT__ 13 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3" 14 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3(Mat C,Mat A,const MatFactorInfo *info) 15 { 16 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 17 IS isrow = b->row,isicol = b->icol; 18 PetscErrorCode ierr; 19 const PetscInt *r,*ic; 20 PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 21 PetscInt *ajtmpold,*ajtmp,nz,row,*ai=a->i,*aj=a->j; 22 PetscInt *diag_offset = b->diag,idx,*pj; 23 MatScalar *pv,*v,*rtmp,*pc,*w,*x; 24 MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; 25 MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9; 26 MatScalar *ba = b->a,*aa = a->a; 27 PetscReal shift = info->shiftinblocks; 28 29 PetscFunctionBegin; 30 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 31 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 32 ierr = PetscMalloc(9*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 33 34 for (i=0; i<n; i++) { 35 nz = bi[i+1] - bi[i]; 36 ajtmp = bj + bi[i]; 37 for (j=0; j<nz; j++) { 38 x = rtmp + 9*ajtmp[j]; 39 x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = 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 + 9*ai[idx]; 46 for (j=0; j<nz; j++) { 47 x = rtmp + 9*ic[ajtmpold[j]]; 48 x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 49 x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8]; 50 v += 9; 51 } 52 row = *ajtmp++; 53 while (row < i) { 54 pc = rtmp + 9*row; 55 p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 56 p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8]; 57 if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || 58 p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) { 59 pv = ba + 9*diag_offset[row]; 60 pj = bj + diag_offset[row] + 1; 61 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 62 x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 63 pc[0] = m1 = p1*x1 + p4*x2 + p7*x3; 64 pc[1] = m2 = p2*x1 + p5*x2 + p8*x3; 65 pc[2] = m3 = p3*x1 + p6*x2 + p9*x3; 66 67 pc[3] = m4 = p1*x4 + p4*x5 + p7*x6; 68 pc[4] = m5 = p2*x4 + p5*x5 + p8*x6; 69 pc[5] = m6 = p3*x4 + p6*x5 + p9*x6; 70 71 pc[6] = m7 = p1*x7 + p4*x8 + p7*x9; 72 pc[7] = m8 = p2*x7 + p5*x8 + p8*x9; 73 pc[8] = m9 = p3*x7 + p6*x8 + p9*x9; 74 nz = bi[row+1] - diag_offset[row] - 1; 75 pv += 9; 76 for (j=0; j<nz; j++) { 77 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 78 x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 79 x = rtmp + 9*pj[j]; 80 x[0] -= m1*x1 + m4*x2 + m7*x3; 81 x[1] -= m2*x1 + m5*x2 + m8*x3; 82 x[2] -= m3*x1 + m6*x2 + m9*x3; 83 84 x[3] -= m1*x4 + m4*x5 + m7*x6; 85 x[4] -= m2*x4 + m5*x5 + m8*x6; 86 x[5] -= m3*x4 + m6*x5 + m9*x6; 87 88 x[6] -= m1*x7 + m4*x8 + m7*x9; 89 x[7] -= m2*x7 + m5*x8 + m8*x9; 90 x[8] -= m3*x7 + m6*x8 + m9*x9; 91 pv += 9; 92 } 93 ierr = PetscLogFlops(54.0*nz+36.0);CHKERRQ(ierr); 94 } 95 row = *ajtmp++; 96 } 97 /* finished row so stick it into b->a */ 98 pv = ba + 9*bi[i]; 99 pj = bj + bi[i]; 100 nz = bi[i+1] - bi[i]; 101 for (j=0; j<nz; j++) { 102 x = rtmp + 9*pj[j]; 103 pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 104 pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8]; 105 pv += 9; 106 } 107 /* invert diagonal block */ 108 w = ba + 9*diag_offset[i]; 109 ierr = Kernel_A_gets_inverse_A_3(w,shift);CHKERRQ(ierr); 110 } 111 112 ierr = PetscFree(rtmp);CHKERRQ(ierr); 113 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 114 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 115 C->ops->solve = MatSolve_SeqBAIJ_3; 116 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3; 117 C->assembled = PETSC_TRUE; 118 ierr = PetscLogFlops(1.3333*27*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 119 PetscFunctionReturn(0); 120 } 121 122 /* MatLUFactorNumeric_SeqBAIJ_3_newdatastruct - 123 copied from MatLUFactorNumeric_SeqBAIJ_N_newdatastruct() and manually re-implemented 124 Kernel_A_gets_A_times_B() 125 Kernel_A_gets_A_minus_B_times_C() 126 Kernel_A_gets_inverse_A() 127 */ 128 #undef __FUNCT__ 129 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_newdatastruct" 130 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_newdatastruct(Mat B,Mat A,const MatFactorInfo *info) 131 { 132 Mat C=B; 133 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 134 IS isrow = b->row,isicol = b->icol; 135 PetscErrorCode ierr; 136 const PetscInt *r,*ic,*ics; 137 PetscInt i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 138 PetscInt *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj; 139 MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 140 PetscInt bs2 = a->bs2,flg; 141 PetscReal shift = info->shiftinblocks; 142 143 PetscFunctionBegin; 144 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 145 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 146 147 /* generate work space needed by the factorization */ 148 ierr = PetscMalloc((bs2*n+bs2+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 149 mwork = rtmp + bs2*n; 150 ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 151 ics = ic; 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 = bi[2*n-i+1] - bi[2*n-i]; 164 bjtmp = bj + bi[2*n-i]; 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++) { if (pc[j]!=0.0) { flg = 1; break; }} 184 if (flg) { 185 pv = b->a + bs2*bdiag[row]; 186 /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 187 ierr = Kernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr); 188 189 pj = b->j + bi[2*n-row]; /* begining of U(row,:) */ 190 pv = b->a + bs2*bi[2*n-row]; 191 nz = bi[2*n-row+1] - bi[2*n-row] - 1; /* num of entries inU(row,:), excluding diag */ 192 for (j=0; j<nz; j++) { 193 /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 194 /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 195 v = rtmp + bs2*pj[j]; 196 ierr = Kernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr); 197 pv += bs2; 198 } 199 ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 200 } 201 } 202 203 /* finished row so stick it into b->a */ 204 /* L part */ 205 pv = b->a + bs2*bi[i] ; 206 pj = b->j + bi[i] ; 207 nz = bi[i+1] - bi[i]; 208 for (j=0; j<nz; j++) { 209 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 210 } 211 212 /* Mark diagonal and invert diagonal for simplier triangular solves */ 213 pv = b->a + bs2*bdiag[i]; 214 pj = b->j + bdiag[i]; 215 ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 216 /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 217 ierr = Kernel_A_gets_inverse_A_3(pv,shift);CHKERRQ(ierr); 218 219 /* U part */ 220 pv = b->a + bs2*bi[2*n-i]; 221 pj = b->j + bi[2*n-i]; 222 nz = bi[2*n-i+1] - bi[2*n-i] - 1; 223 for (j=0; j<nz; j++){ 224 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 225 } 226 } 227 228 ierr = PetscFree(rtmp);CHKERRQ(ierr); 229 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 230 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 231 232 C->assembled = PETSC_TRUE; 233 ierr = PetscLogFlops(1.3333*bs2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 234 PetscFunctionReturn(0); 235 } 236 237 #undef __FUNCT__ 238 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_newdatastruct_v2" 239 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_newdatastruct_v2(Mat B,Mat A,const MatFactorInfo *info) 240 { 241 Mat C=B; 242 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 243 IS isrow = b->row,isicol = b->icol; 244 PetscErrorCode ierr; 245 const PetscInt *r,*ic,*ics; 246 PetscInt i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 247 PetscInt *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj; 248 MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 249 PetscInt bs2 = a->bs2,flg; 250 PetscReal shift = info->shiftinblocks; 251 252 PetscFunctionBegin; 253 ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr); 254 ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr); 255 256 /* generate work space needed by the factorization */ 257 ierr = PetscMalloc((bs2*n+bs2+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 258 mwork = rtmp + bs2*n; 259 ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 260 ics = ic; 261 262 for (i=0; i<n; i++){ 263 /* zero rtmp */ 264 /* L part */ 265 nz = bi[i+1] - bi[i]; 266 bjtmp = bj + bi[i]; 267 for (j=0; j<nz; j++){ 268 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 269 } 270 271 /* U part */ 272 nz = bdiag[i] - bdiag[i+1]; 273 bjtmp = bj + bdiag[i+1]+1; 274 for (j=0; j<nz; j++){ 275 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 276 } 277 278 /* load in initial (unfactored row) */ 279 nz = ai[r[i]+1] - ai[r[i]]; 280 ajtmp = aj + ai[r[i]]; 281 v = aa + bs2*ai[r[i]]; 282 for (j=0; j<nz; j++) { 283 ierr = PetscMemcpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 284 } 285 286 /* elimination */ 287 bjtmp = bj + bi[i]; 288 nzL = bi[i+1] - bi[i]; 289 for(k = 0;k < nzL;k++){ 290 row = bjtmp[k]; 291 pc = rtmp + bs2*row; 292 for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }} 293 if (flg) { 294 pv = b->a + bs2*bdiag[row]; 295 /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 296 ierr = Kernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr); 297 298 pj = b->j + bdiag[row+1] + 1; /* beginning of U(row,:) */ 299 pv = b->a + bs2*(bdiag[row+1]+1); 300 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */ 301 for (j=0; j<nz; j++) { 302 /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 303 /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 304 v = rtmp + bs2*pj[j]; 305 ierr = Kernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr); 306 pv += bs2; 307 } 308 ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 309 } 310 } 311 312 /* finished row so stick it into b->a */ 313 /* L part */ 314 pv = b->a + bs2*bi[i] ; 315 pj = b->j + bi[i] ; 316 nz = bi[i+1] - bi[i]; 317 for (j=0; j<nz; j++) { 318 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 319 } 320 321 /* Mark diagonal and invert diagonal for simplier triangular solves */ 322 pv = b->a + bs2*bdiag[i]; 323 pj = b->j + bdiag[i]; 324 ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 325 /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 326 ierr = Kernel_A_gets_inverse_A_3(pv,shift);CHKERRQ(ierr); 327 328 /* U part */ 329 pj = b->j + bdiag[i+1] + 1; 330 pv = b->a + bs2*(bdiag[i+1]+1); 331 nz = bdiag[i] - bdiag[i+1] - 1; 332 for (j=0; j<nz; j++){ 333 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 334 } 335 } 336 337 ierr = PetscFree(rtmp);CHKERRQ(ierr); 338 ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr); 339 ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr); 340 341 C->assembled = PETSC_TRUE; 342 ierr = PetscLogFlops(1.3333*bs2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 343 PetscFunctionReturn(0); 344 } 345 346 #undef __FUNCT__ 347 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering" 348 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) 349 { 350 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data; 351 PetscErrorCode ierr; 352 PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j; 353 PetscInt *ajtmpold,*ajtmp,nz,row; 354 PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 355 MatScalar *pv,*v,*rtmp,*pc,*w,*x; 356 MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; 357 MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9; 358 MatScalar *ba = b->a,*aa = a->a; 359 PetscReal shift = info->shiftinblocks; 360 361 PetscFunctionBegin; 362 ierr = PetscMalloc(9*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 363 364 for (i=0; i<n; i++) { 365 nz = bi[i+1] - bi[i]; 366 ajtmp = bj + bi[i]; 367 for (j=0; j<nz; j++) { 368 x = rtmp+9*ajtmp[j]; 369 x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = 0.0; 370 } 371 /* load in initial (unfactored row) */ 372 nz = ai[i+1] - ai[i]; 373 ajtmpold = aj + ai[i]; 374 v = aa + 9*ai[i]; 375 for (j=0; j<nz; j++) { 376 x = rtmp+9*ajtmpold[j]; 377 x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 378 x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8]; 379 v += 9; 380 } 381 row = *ajtmp++; 382 while (row < i) { 383 pc = rtmp + 9*row; 384 p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 385 p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8]; 386 if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || 387 p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0) { 388 pv = ba + 9*diag_offset[row]; 389 pj = bj + diag_offset[row] + 1; 390 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 391 x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 392 pc[0] = m1 = p1*x1 + p4*x2 + p7*x3; 393 pc[1] = m2 = p2*x1 + p5*x2 + p8*x3; 394 pc[2] = m3 = p3*x1 + p6*x2 + p9*x3; 395 396 pc[3] = m4 = p1*x4 + p4*x5 + p7*x6; 397 pc[4] = m5 = p2*x4 + p5*x5 + p8*x6; 398 pc[5] = m6 = p3*x4 + p6*x5 + p9*x6; 399 400 pc[6] = m7 = p1*x7 + p4*x8 + p7*x9; 401 pc[7] = m8 = p2*x7 + p5*x8 + p8*x9; 402 pc[8] = m9 = p3*x7 + p6*x8 + p9*x9; 403 404 nz = bi[row+1] - diag_offset[row] - 1; 405 pv += 9; 406 for (j=0; j<nz; j++) { 407 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 408 x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 409 x = rtmp + 9*pj[j]; 410 x[0] -= m1*x1 + m4*x2 + m7*x3; 411 x[1] -= m2*x1 + m5*x2 + m8*x3; 412 x[2] -= m3*x1 + m6*x2 + m9*x3; 413 414 x[3] -= m1*x4 + m4*x5 + m7*x6; 415 x[4] -= m2*x4 + m5*x5 + m8*x6; 416 x[5] -= m3*x4 + m6*x5 + m9*x6; 417 418 x[6] -= m1*x7 + m4*x8 + m7*x9; 419 x[7] -= m2*x7 + m5*x8 + m8*x9; 420 x[8] -= m3*x7 + m6*x8 + m9*x9; 421 pv += 9; 422 } 423 ierr = PetscLogFlops(54.0*nz+36.0);CHKERRQ(ierr); 424 } 425 row = *ajtmp++; 426 } 427 /* finished row so stick it into b->a */ 428 pv = ba + 9*bi[i]; 429 pj = bj + bi[i]; 430 nz = bi[i+1] - bi[i]; 431 for (j=0; j<nz; j++) { 432 x = rtmp+9*pj[j]; 433 pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 434 pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8]; 435 pv += 9; 436 } 437 /* invert diagonal block */ 438 w = ba + 9*diag_offset[i]; 439 ierr = Kernel_A_gets_inverse_A_3(w,shift);CHKERRQ(ierr); 440 } 441 442 ierr = PetscFree(rtmp);CHKERRQ(ierr); 443 C->ops->solve = MatSolve_SeqBAIJ_3_NaturalOrdering; 444 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_3_NaturalOrdering; 445 C->assembled = PETSC_TRUE; 446 ierr = PetscLogFlops(1.3333*27*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 447 PetscFunctionReturn(0); 448 } 449 450 /* 451 MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_newdatastruct - 452 copied from MatLUFactorNumeric_SeqBAIJ_2_NaturalOrdering_newdatastruct() 453 */ 454 #undef __FUNCT__ 455 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_newdatastruct" 456 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_newdatastruct(Mat B,Mat A,const MatFactorInfo *info) 457 { 458 Mat C=B; 459 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 460 PetscErrorCode ierr; 461 PetscInt i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 462 PetscInt *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj; 463 MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 464 PetscInt bs2 = a->bs2,flg; 465 PetscReal shift = info->shiftinblocks; 466 467 PetscFunctionBegin; 468 /* generate work space needed by the factorization */ 469 ierr = PetscMalloc((bs2*n+bs2+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 470 mwork = rtmp + bs2*n; 471 ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 472 473 for (i=0; i<n; i++){ 474 /* zero rtmp */ 475 /* L part */ 476 nz = bi[i+1] - bi[i]; 477 bjtmp = bj + bi[i]; 478 for (j=0; j<nz; j++){ 479 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 480 } 481 482 /* U part */ 483 nz = bi[2*n-i+1] - bi[2*n-i]; 484 bjtmp = bj + bi[2*n-i]; 485 for (j=0; j<nz; j++){ 486 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 487 } 488 489 /* load in initial (unfactored row) */ 490 nz = ai[i+1] - ai[i]; 491 ajtmp = aj + ai[i]; 492 v = aa + bs2*ai[i]; 493 for (j=0; j<nz; j++) { 494 ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 495 } 496 497 /* elimination */ 498 bjtmp = bj + bi[i]; 499 nzL = bi[i+1] - bi[i]; 500 for(k=0;k<nzL;k++){ 501 row = bjtmp[k]; 502 pc = rtmp + bs2*row; 503 for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }} 504 if (flg) { 505 pv = b->a + bs2*bdiag[row]; 506 /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 507 ierr = Kernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr); 508 509 pj = b->j + bi[2*n-row]; /* begining of U(row,:) */ 510 pv = b->a + bs2*bi[2*n-row]; 511 nz = bi[2*n-row+1] - bi[2*n-row] - 1; /* num of entries inU(row,:), excluding diag */ 512 for (j=0; j<nz; j++) { 513 /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 514 /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 515 v = rtmp + bs2*pj[j]; 516 ierr = Kernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr); 517 pv += bs2; 518 } 519 ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 520 } 521 } 522 523 /* finished row so stick it into b->a */ 524 /* L part */ 525 pv = b->a + bs2*bi[i] ; 526 pj = b->j + bi[i] ; 527 nz = bi[i+1] - bi[i]; 528 for (j=0; j<nz; j++) { 529 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 530 } 531 532 /* Mark diagonal and invert diagonal for simplier triangular solves */ 533 pv = b->a + bs2*bdiag[i]; 534 pj = b->j + bdiag[i]; 535 ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 536 /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 537 ierr = Kernel_A_gets_inverse_A_3(pv,shift);CHKERRQ(ierr); 538 539 /* U part */ 540 pv = b->a + bs2*bi[2*n-i]; 541 pj = b->j + bi[2*n-i]; 542 nz = bi[2*n-i+1] - bi[2*n-i] - 1; 543 for (j=0; j<nz; j++){ 544 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 545 } 546 } 547 548 ierr = PetscFree(rtmp);CHKERRQ(ierr); 549 C->assembled = PETSC_TRUE; 550 ierr = PetscLogFlops(1.3333*bs2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 551 PetscFunctionReturn(0); 552 } 553 554 #undef __FUNCT__ 555 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_newdatastruct_v2" 556 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_3_NaturalOrdering_newdatastruct_v2(Mat B,Mat A,const MatFactorInfo *info) 557 { 558 Mat C=B; 559 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ *)C->data; 560 PetscErrorCode ierr; 561 PetscInt i,j,k,n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 562 PetscInt *ajtmp,*bjtmp,nz,nzL,row,*bdiag=b->diag,*pj; 563 MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a; 564 PetscInt bs2 = a->bs2,flg; 565 PetscReal shift = info->shiftinblocks; 566 567 PetscFunctionBegin; 568 /* generate work space needed by the factorization */ 569 ierr = PetscMalloc((bs2*n+bs2+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 570 mwork = rtmp + bs2*n; 571 ierr = PetscMemzero(rtmp,bs2*n*sizeof(MatScalar));CHKERRQ(ierr); 572 573 for (i=0; i<n; i++){ 574 /* zero rtmp */ 575 /* L part */ 576 nz = bi[i+1] - bi[i]; 577 bjtmp = bj + bi[i]; 578 for (j=0; j<nz; j++){ 579 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 580 } 581 582 /* U part */ 583 nz = bdiag[i] - bdiag[i+1]; 584 bjtmp = bj + bdiag[i+1] + 1; 585 for (j=0; j<nz; j++){ 586 ierr = PetscMemzero(rtmp+bs2*bjtmp[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 587 } 588 589 /* load in initial (unfactored row) */ 590 nz = ai[i+1] - ai[i]; 591 ajtmp = aj + ai[i]; 592 v = aa + bs2*ai[i]; 593 for (j=0; j<nz; j++) { 594 ierr = PetscMemcpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2*sizeof(MatScalar));CHKERRQ(ierr); 595 } 596 597 /* elimination */ 598 bjtmp = bj + bi[i]; 599 nzL = bi[i+1] - bi[i]; 600 for(k=0;k<nzL;k++){ 601 row = bjtmp[k]; 602 pc = rtmp + bs2*row; 603 for (flg=0,j=0; j<bs2; j++) { if (pc[j]!=0.0) { flg = 1; break; }} 604 if (flg) { 605 pv = b->a + bs2*bdiag[row]; 606 /* Kernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 607 ierr = Kernel_A_gets_A_times_B_3(pc,pv,mwork);CHKERRQ(ierr); 608 609 pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */ 610 pv = b->a + bs2*(bdiag[row+1]+1); 611 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries in U(row,:) excluding diag */ 612 for (j=0; j<nz; j++) { 613 /* Kernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 614 /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 615 v = rtmp + bs2*pj[j]; 616 ierr = Kernel_A_gets_A_minus_B_times_C_3(v,pc,pv);CHKERRQ(ierr); 617 pv += bs2; 618 } 619 ierr = PetscLogFlops(54*nz+45);CHKERRQ(ierr); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 620 } 621 } 622 623 /* finished row so stick it into b->a */ 624 /* L part */ 625 pv = b->a + bs2*bi[i] ; 626 pj = b->j + bi[i] ; 627 nz = bi[i+1] - bi[i]; 628 for (j=0; j<nz; j++) { 629 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 630 } 631 632 /* Mark diagonal and invert diagonal for simplier triangular solves */ 633 pv = b->a + bs2*bdiag[i]; 634 pj = b->j + bdiag[i]; 635 ierr = PetscMemcpy(pv,rtmp+bs2*pj[0],bs2*sizeof(MatScalar));CHKERRQ(ierr); 636 /* ierr = Kernel_A_gets_inverse_A(bs,pv,v_pivots,v_work);CHKERRQ(ierr); */ 637 ierr = Kernel_A_gets_inverse_A_3(pv,shift);CHKERRQ(ierr); 638 639 /* U part */ 640 pv = b->a + bs2*(bdiag[i+1]+1); 641 pj = b->j + bdiag[i+1]+1; 642 nz = bdiag[i] - bdiag[i+1] - 1; 643 for (j=0; j<nz; j++){ 644 ierr = PetscMemcpy(pv+bs2*j,rtmp+bs2*pj[j],bs2*sizeof(MatScalar));CHKERRQ(ierr); 645 } 646 } 647 648 ierr = PetscFree(rtmp);CHKERRQ(ierr); 649 C->assembled = PETSC_TRUE; 650 ierr = PetscLogFlops(1.3333*bs2*n);CHKERRQ(ierr); /* from inverting diagonal blocks */ 651 PetscFunctionReturn(0); 652 } 653 654