1 2 /* 3 Factorization code for BAIJ format. 4 */ 5 #include <../src/mat/impls/baij/seq/baij.h> 6 #include <petsc/private/kernels/blockinvert.h> 7 8 /* ------------------------------------------------------------*/ 9 /* 10 Version for when blocks are 5 by 5 11 */ 12 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_inplace(Mat C,Mat A,const MatFactorInfo *info) 13 { 14 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data; 15 IS isrow = b->row,isicol = b->icol; 16 const PetscInt *r,*ic,*bi = b->i,*bj = b->j,*ajtmpold,*ajtmp; 17 PetscInt i,j,n = a->mbs,nz,row,idx,ipvt[5]; 18 const PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 19 MatScalar *w,*pv,*rtmp,*x,*pc; 20 const MatScalar *v,*aa = a->a; 21 MatScalar p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4; 22 MatScalar p5,p6,p7,p8,p9,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,x16; 23 MatScalar x17,x18,x19,x20,x21,x22,x23,x24,x25,p10,p11,p12,p13,p14; 24 MatScalar p15,p16,p17,p18,p19,p20,p21,p22,p23,p24,p25,m10,m11,m12; 25 MatScalar m13,m14,m15,m16,m17,m18,m19,m20,m21,m22,m23,m24,m25; 26 MatScalar *ba = b->a,work[25]; 27 PetscReal shift = info->shiftamount; 28 PetscBool allowzeropivot,zeropivotdetected; 29 30 PetscFunctionBegin; 31 allowzeropivot = PetscNot(A->erroriffailure); 32 PetscCall(ISGetIndices(isrow,&r)); 33 PetscCall(ISGetIndices(isicol,&ic)); 34 PetscCall(PetscMalloc1(25*(n+1),&rtmp)); 35 36 #define PETSC_USE_MEMZERO 1 37 #define PETSC_USE_MEMCPY 1 38 39 for (i=0; i<n; i++) { 40 nz = bi[i+1] - bi[i]; 41 ajtmp = bj + bi[i]; 42 for (j=0; j<nz; j++) { 43 #if defined(PETSC_USE_MEMZERO) 44 PetscCall(PetscArrayzero(rtmp+25*ajtmp[j],25)); 45 #else 46 x = rtmp+25*ajtmp[j]; 47 x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0; 48 x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0; 49 x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = 0.0; 50 #endif 51 } 52 /* load in initial (unfactored row) */ 53 idx = r[i]; 54 nz = ai[idx+1] - ai[idx]; 55 ajtmpold = aj + ai[idx]; 56 v = aa + 25*ai[idx]; 57 for (j=0; j<nz; j++) { 58 #if defined(PETSC_USE_MEMCPY) 59 PetscCall(PetscArraycpy(rtmp+25*ic[ajtmpold[j]],v,25)); 60 #else 61 x = rtmp+25*ic[ajtmpold[j]]; 62 x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 63 x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8]; 64 x[9] = v[9]; x[10] = v[10]; x[11] = v[11]; x[12] = v[12]; x[13] = v[13]; 65 x[14] = v[14]; x[15] = v[15]; x[16] = v[16]; x[17] = v[17]; 66 x[18] = v[18]; x[19] = v[19]; x[20] = v[20]; x[21] = v[21]; 67 x[22] = v[22]; x[23] = v[23]; x[24] = v[24]; 68 #endif 69 v += 25; 70 } 71 row = *ajtmp++; 72 while (row < i) { 73 pc = rtmp + 25*row; 74 p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 75 p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8]; 76 p10 = pc[9]; p11 = pc[10]; p12 = pc[11]; p13 = pc[12]; p14 = pc[13]; 77 p15 = pc[14]; p16 = pc[15]; p17 = pc[16]; p18 = pc[17]; p19 = pc[18]; 78 p20 = pc[19]; p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23]; 79 p25 = pc[24]; 80 if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || 81 p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 || 82 p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0 83 || p16 != 0.0 || p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || 84 p20 != 0.0 || p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || 85 p24 != 0.0 || p25 != 0.0) { 86 pv = ba + 25*diag_offset[row]; 87 pj = bj + diag_offset[row] + 1; 88 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 89 x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 90 x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; x14 = pv[13]; 91 x15 = pv[14]; x16 = pv[15]; x17 = pv[16]; x18 = pv[17]; 92 x19 = pv[18]; x20 = pv[19]; x21 = pv[20]; x22 = pv[21]; 93 x23 = pv[22]; x24 = pv[23]; x25 = pv[24]; 94 pc[0] = m1 = p1*x1 + p6*x2 + p11*x3 + p16*x4 + p21*x5; 95 pc[1] = m2 = p2*x1 + p7*x2 + p12*x3 + p17*x4 + p22*x5; 96 pc[2] = m3 = p3*x1 + p8*x2 + p13*x3 + p18*x4 + p23*x5; 97 pc[3] = m4 = p4*x1 + p9*x2 + p14*x3 + p19*x4 + p24*x5; 98 pc[4] = m5 = p5*x1 + p10*x2 + p15*x3 + p20*x4 + p25*x5; 99 100 pc[5] = m6 = p1*x6 + p6*x7 + p11*x8 + p16*x9 + p21*x10; 101 pc[6] = m7 = p2*x6 + p7*x7 + p12*x8 + p17*x9 + p22*x10; 102 pc[7] = m8 = p3*x6 + p8*x7 + p13*x8 + p18*x9 + p23*x10; 103 pc[8] = m9 = p4*x6 + p9*x7 + p14*x8 + p19*x9 + p24*x10; 104 pc[9] = m10 = p5*x6 + p10*x7 + p15*x8 + p20*x9 + p25*x10; 105 106 pc[10] = m11 = p1*x11 + p6*x12 + p11*x13 + p16*x14 + p21*x15; 107 pc[11] = m12 = p2*x11 + p7*x12 + p12*x13 + p17*x14 + p22*x15; 108 pc[12] = m13 = p3*x11 + p8*x12 + p13*x13 + p18*x14 + p23*x15; 109 pc[13] = m14 = p4*x11 + p9*x12 + p14*x13 + p19*x14 + p24*x15; 110 pc[14] = m15 = p5*x11 + p10*x12 + p15*x13 + p20*x14 + p25*x15; 111 112 pc[15] = m16 = p1*x16 + p6*x17 + p11*x18 + p16*x19 + p21*x20; 113 pc[16] = m17 = p2*x16 + p7*x17 + p12*x18 + p17*x19 + p22*x20; 114 pc[17] = m18 = p3*x16 + p8*x17 + p13*x18 + p18*x19 + p23*x20; 115 pc[18] = m19 = p4*x16 + p9*x17 + p14*x18 + p19*x19 + p24*x20; 116 pc[19] = m20 = p5*x16 + p10*x17 + p15*x18 + p20*x19 + p25*x20; 117 118 pc[20] = m21 = p1*x21 + p6*x22 + p11*x23 + p16*x24 + p21*x25; 119 pc[21] = m22 = p2*x21 + p7*x22 + p12*x23 + p17*x24 + p22*x25; 120 pc[22] = m23 = p3*x21 + p8*x22 + p13*x23 + p18*x24 + p23*x25; 121 pc[23] = m24 = p4*x21 + p9*x22 + p14*x23 + p19*x24 + p24*x25; 122 pc[24] = m25 = p5*x21 + p10*x22 + p15*x23 + p20*x24 + p25*x25; 123 124 nz = bi[row+1] - diag_offset[row] - 1; 125 pv += 25; 126 for (j=0; j<nz; j++) { 127 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 128 x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 129 x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; 130 x14 = pv[13]; x15 = pv[14]; x16 = pv[15]; x17 = pv[16]; 131 x18 = pv[17]; x19 = pv[18]; x20 = pv[19]; x21 = pv[20]; 132 x22 = pv[21]; x23 = pv[22]; x24 = pv[23]; x25 = pv[24]; 133 x = rtmp + 25*pj[j]; 134 x[0] -= m1*x1 + m6*x2 + m11*x3 + m16*x4 + m21*x5; 135 x[1] -= m2*x1 + m7*x2 + m12*x3 + m17*x4 + m22*x5; 136 x[2] -= m3*x1 + m8*x2 + m13*x3 + m18*x4 + m23*x5; 137 x[3] -= m4*x1 + m9*x2 + m14*x3 + m19*x4 + m24*x5; 138 x[4] -= m5*x1 + m10*x2 + m15*x3 + m20*x4 + m25*x5; 139 140 x[5] -= m1*x6 + m6*x7 + m11*x8 + m16*x9 + m21*x10; 141 x[6] -= m2*x6 + m7*x7 + m12*x8 + m17*x9 + m22*x10; 142 x[7] -= m3*x6 + m8*x7 + m13*x8 + m18*x9 + m23*x10; 143 x[8] -= m4*x6 + m9*x7 + m14*x8 + m19*x9 + m24*x10; 144 x[9] -= m5*x6 + m10*x7 + m15*x8 + m20*x9 + m25*x10; 145 146 x[10] -= m1*x11 + m6*x12 + m11*x13 + m16*x14 + m21*x15; 147 x[11] -= m2*x11 + m7*x12 + m12*x13 + m17*x14 + m22*x15; 148 x[12] -= m3*x11 + m8*x12 + m13*x13 + m18*x14 + m23*x15; 149 x[13] -= m4*x11 + m9*x12 + m14*x13 + m19*x14 + m24*x15; 150 x[14] -= m5*x11 + m10*x12 + m15*x13 + m20*x14 + m25*x15; 151 152 x[15] -= m1*x16 + m6*x17 + m11*x18 + m16*x19 + m21*x20; 153 x[16] -= m2*x16 + m7*x17 + m12*x18 + m17*x19 + m22*x20; 154 x[17] -= m3*x16 + m8*x17 + m13*x18 + m18*x19 + m23*x20; 155 x[18] -= m4*x16 + m9*x17 + m14*x18 + m19*x19 + m24*x20; 156 x[19] -= m5*x16 + m10*x17 + m15*x18 + m20*x19 + m25*x20; 157 158 x[20] -= m1*x21 + m6*x22 + m11*x23 + m16*x24 + m21*x25; 159 x[21] -= m2*x21 + m7*x22 + m12*x23 + m17*x24 + m22*x25; 160 x[22] -= m3*x21 + m8*x22 + m13*x23 + m18*x24 + m23*x25; 161 x[23] -= m4*x21 + m9*x22 + m14*x23 + m19*x24 + m24*x25; 162 x[24] -= m5*x21 + m10*x22 + m15*x23 + m20*x24 + m25*x25; 163 164 pv += 25; 165 } 166 PetscCall(PetscLogFlops(250.0*nz+225.0)); 167 } 168 row = *ajtmp++; 169 } 170 /* finished row so stick it into b->a */ 171 pv = ba + 25*bi[i]; 172 pj = bj + bi[i]; 173 nz = bi[i+1] - bi[i]; 174 for (j=0; j<nz; j++) { 175 #if defined(PETSC_USE_MEMCPY) 176 PetscCall(PetscArraycpy(pv,rtmp+25*pj[j],25)); 177 #else 178 x = rtmp+25*pj[j]; 179 pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 180 pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8]; 181 pv[9] = x[9]; pv[10] = x[10]; pv[11] = x[11]; pv[12] = x[12]; 182 pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15]; pv[16] = x[16]; 183 pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19]; pv[20] = x[20]; 184 pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23]; pv[24] = x[24]; 185 #endif 186 pv += 25; 187 } 188 /* invert diagonal block */ 189 w = ba + 25*diag_offset[i]; 190 PetscCall(PetscKernel_A_gets_inverse_A_5(w,ipvt,work,shift,allowzeropivot,&zeropivotdetected)); 191 if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 192 } 193 194 PetscCall(PetscFree(rtmp)); 195 PetscCall(ISRestoreIndices(isicol,&ic)); 196 PetscCall(ISRestoreIndices(isrow,&r)); 197 198 C->ops->solve = MatSolve_SeqBAIJ_5_inplace; 199 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_inplace; 200 C->assembled = PETSC_TRUE; 201 202 PetscCall(PetscLogFlops(1.333333333333*5*5*5*b->mbs)); /* from inverting diagonal blocks */ 203 PetscFunctionReturn(0); 204 } 205 206 /* MatLUFactorNumeric_SeqBAIJ_5 - 207 copied from MatLUFactorNumeric_SeqBAIJ_N_inplace() and manually re-implemented 208 PetscKernel_A_gets_A_times_B() 209 PetscKernel_A_gets_A_minus_B_times_C() 210 PetscKernel_A_gets_inverse_A() 211 */ 212 213 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5(Mat B,Mat A,const MatFactorInfo *info) 214 { 215 Mat C =B; 216 Mat_SeqBAIJ *a =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data; 217 IS isrow = b->row,isicol = b->icol; 218 const PetscInt *r,*ic; 219 PetscInt i,j,k,nz,nzL,row; 220 const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 221 const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2; 222 MatScalar *rtmp,*pc,*mwork,*v,*pv,*aa=a->a,work[25]; 223 PetscInt flg,ipvt[5]; 224 PetscReal shift = info->shiftamount; 225 PetscBool allowzeropivot,zeropivotdetected; 226 227 PetscFunctionBegin; 228 allowzeropivot = PetscNot(A->erroriffailure); 229 PetscCall(ISGetIndices(isrow,&r)); 230 PetscCall(ISGetIndices(isicol,&ic)); 231 232 /* generate work space needed by the factorization */ 233 PetscCall(PetscMalloc2(bs2*n,&rtmp,bs2,&mwork)); 234 PetscCall(PetscArrayzero(rtmp,bs2*n)); 235 236 for (i=0; i<n; i++) { 237 /* zero rtmp */ 238 /* L part */ 239 nz = bi[i+1] - bi[i]; 240 bjtmp = bj + bi[i]; 241 for (j=0; j<nz; j++) { 242 PetscCall(PetscArrayzero(rtmp+bs2*bjtmp[j],bs2)); 243 } 244 245 /* U part */ 246 nz = bdiag[i] - bdiag[i+1]; 247 bjtmp = bj + bdiag[i+1]+1; 248 for (j=0; j<nz; j++) { 249 PetscCall(PetscArrayzero(rtmp+bs2*bjtmp[j],bs2)); 250 } 251 252 /* load in initial (unfactored row) */ 253 nz = ai[r[i]+1] - ai[r[i]]; 254 ajtmp = aj + ai[r[i]]; 255 v = aa + bs2*ai[r[i]]; 256 for (j=0; j<nz; j++) { 257 PetscCall(PetscArraycpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2)); 258 } 259 260 /* elimination */ 261 bjtmp = bj + bi[i]; 262 nzL = bi[i+1] - bi[i]; 263 for (k=0; k < nzL; k++) { 264 row = bjtmp[k]; 265 pc = rtmp + bs2*row; 266 for (flg=0,j=0; j<bs2; j++) { 267 if (pc[j]!=0.0) { 268 flg = 1; 269 break; 270 } 271 } 272 if (flg) { 273 pv = b->a + bs2*bdiag[row]; 274 /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 275 PetscCall(PetscKernel_A_gets_A_times_B_5(pc,pv,mwork)); 276 277 pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */ 278 pv = b->a + bs2*(bdiag[row+1]+1); 279 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */ 280 for (j=0; j<nz; j++) { 281 /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 282 /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 283 v = rtmp + bs2*pj[j]; 284 PetscCall(PetscKernel_A_gets_A_minus_B_times_C_5(v,pc,pv)); 285 pv += bs2; 286 } 287 PetscCall(PetscLogFlops(250.0*nz+225)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 288 } 289 } 290 291 /* finished row so stick it into b->a */ 292 /* L part */ 293 pv = b->a + bs2*bi[i]; 294 pj = b->j + bi[i]; 295 nz = bi[i+1] - bi[i]; 296 for (j=0; j<nz; j++) { 297 PetscCall(PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2)); 298 } 299 300 /* Mark diagonal and invert diagonal for simpler triangular solves */ 301 pv = b->a + bs2*bdiag[i]; 302 pj = b->j + bdiag[i]; 303 PetscCall(PetscArraycpy(pv,rtmp+bs2*pj[0],bs2)); 304 PetscCall(PetscKernel_A_gets_inverse_A_5(pv,ipvt,work,shift,allowzeropivot,&zeropivotdetected)); 305 if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 306 307 /* U part */ 308 pv = b->a + bs2*(bdiag[i+1]+1); 309 pj = b->j + bdiag[i+1]+1; 310 nz = bdiag[i] - bdiag[i+1] - 1; 311 for (j=0; j<nz; j++) { 312 PetscCall(PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2)); 313 } 314 } 315 316 PetscCall(PetscFree2(rtmp,mwork)); 317 PetscCall(ISRestoreIndices(isicol,&ic)); 318 PetscCall(ISRestoreIndices(isrow,&r)); 319 320 C->ops->solve = MatSolve_SeqBAIJ_5; 321 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5; 322 C->assembled = PETSC_TRUE; 323 324 PetscCall(PetscLogFlops(1.333333333333*5*5*5*n)); /* from inverting diagonal blocks */ 325 PetscFunctionReturn(0); 326 } 327 328 /* 329 Version for when blocks are 5 by 5 Using natural ordering 330 */ 331 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info) 332 { 333 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data; 334 PetscInt i,j,n = a->mbs,*bi = b->i,*bj = b->j,ipvt[5]; 335 PetscInt *ajtmpold,*ajtmp,nz,row; 336 PetscInt *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj; 337 MatScalar *pv,*v,*rtmp,*pc,*w,*x; 338 MatScalar x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15; 339 MatScalar x16,x17,x18,x19,x20,x21,x22,x23,x24,x25; 340 MatScalar p1,p2,p3,p4,p5,p6,p7,p8,p9,p10,p11,p12,p13,p14,p15; 341 MatScalar p16,p17,p18,p19,p20,p21,p22,p23,p24,p25; 342 MatScalar m1,m2,m3,m4,m5,m6,m7,m8,m9,m10,m11,m12,m13,m14,m15; 343 MatScalar m16,m17,m18,m19,m20,m21,m22,m23,m24,m25; 344 MatScalar *ba = b->a,*aa = a->a,work[25]; 345 PetscReal shift = info->shiftamount; 346 PetscBool allowzeropivot,zeropivotdetected; 347 348 PetscFunctionBegin; 349 allowzeropivot = PetscNot(A->erroriffailure); 350 PetscCall(PetscMalloc1(25*(n+1),&rtmp)); 351 for (i=0; i<n; i++) { 352 nz = bi[i+1] - bi[i]; 353 ajtmp = bj + bi[i]; 354 for (j=0; j<nz; j++) { 355 x = rtmp+25*ajtmp[j]; 356 x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0; 357 x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = 0.0; 358 x[16] = x[17] = x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = 0.0; 359 } 360 /* load in initial (unfactored row) */ 361 nz = ai[i+1] - ai[i]; 362 ajtmpold = aj + ai[i]; 363 v = aa + 25*ai[i]; 364 for (j=0; j<nz; j++) { 365 x = rtmp+25*ajtmpold[j]; 366 x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3]; 367 x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8]; 368 x[9] = v[9]; x[10] = v[10]; x[11] = v[11]; x[12] = v[12]; x[13] = v[13]; 369 x[14] = v[14]; x[15] = v[15]; x[16] = v[16]; x[17] = v[17]; x[18] = v[18]; 370 x[19] = v[19]; x[20] = v[20]; x[21] = v[21]; x[22] = v[22]; x[23] = v[23]; 371 x[24] = v[24]; 372 v += 25; 373 } 374 row = *ajtmp++; 375 while (row < i) { 376 pc = rtmp + 25*row; 377 p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3]; 378 p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8]; 379 p10 = pc[9]; p11 = pc[10]; p12 = pc[11]; p13 = pc[12]; p14 = pc[13]; 380 p15 = pc[14]; p16 = pc[15]; p17 = pc[16]; p18 = pc[17]; 381 p19 = pc[18]; p20 = pc[19]; p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; 382 p24 = pc[23]; p25 = pc[24]; 383 if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || 384 p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 || 385 p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0 386 || p16 != 0.0 || p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 387 || p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 || p25 != 0.0) { 388 pv = ba + 25*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 x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; x14 = pv[13]; 393 x15 = pv[14]; x16 = pv[15]; x17 = pv[16]; x18 = pv[17]; x19 = pv[18]; 394 x20 = pv[19]; x21 = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23]; 395 x25 = pv[24]; 396 pc[0] = m1 = p1*x1 + p6*x2 + p11*x3 + p16*x4 + p21*x5; 397 pc[1] = m2 = p2*x1 + p7*x2 + p12*x3 + p17*x4 + p22*x5; 398 pc[2] = m3 = p3*x1 + p8*x2 + p13*x3 + p18*x4 + p23*x5; 399 pc[3] = m4 = p4*x1 + p9*x2 + p14*x3 + p19*x4 + p24*x5; 400 pc[4] = m5 = p5*x1 + p10*x2 + p15*x3 + p20*x4 + p25*x5; 401 402 pc[5] = m6 = p1*x6 + p6*x7 + p11*x8 + p16*x9 + p21*x10; 403 pc[6] = m7 = p2*x6 + p7*x7 + p12*x8 + p17*x9 + p22*x10; 404 pc[7] = m8 = p3*x6 + p8*x7 + p13*x8 + p18*x9 + p23*x10; 405 pc[8] = m9 = p4*x6 + p9*x7 + p14*x8 + p19*x9 + p24*x10; 406 pc[9] = m10 = p5*x6 + p10*x7 + p15*x8 + p20*x9 + p25*x10; 407 408 pc[10] = m11 = p1*x11 + p6*x12 + p11*x13 + p16*x14 + p21*x15; 409 pc[11] = m12 = p2*x11 + p7*x12 + p12*x13 + p17*x14 + p22*x15; 410 pc[12] = m13 = p3*x11 + p8*x12 + p13*x13 + p18*x14 + p23*x15; 411 pc[13] = m14 = p4*x11 + p9*x12 + p14*x13 + p19*x14 + p24*x15; 412 pc[14] = m15 = p5*x11 + p10*x12 + p15*x13 + p20*x14 + p25*x15; 413 414 pc[15] = m16 = p1*x16 + p6*x17 + p11*x18 + p16*x19 + p21*x20; 415 pc[16] = m17 = p2*x16 + p7*x17 + p12*x18 + p17*x19 + p22*x20; 416 pc[17] = m18 = p3*x16 + p8*x17 + p13*x18 + p18*x19 + p23*x20; 417 pc[18] = m19 = p4*x16 + p9*x17 + p14*x18 + p19*x19 + p24*x20; 418 pc[19] = m20 = p5*x16 + p10*x17 + p15*x18 + p20*x19 + p25*x20; 419 420 pc[20] = m21 = p1*x21 + p6*x22 + p11*x23 + p16*x24 + p21*x25; 421 pc[21] = m22 = p2*x21 + p7*x22 + p12*x23 + p17*x24 + p22*x25; 422 pc[22] = m23 = p3*x21 + p8*x22 + p13*x23 + p18*x24 + p23*x25; 423 pc[23] = m24 = p4*x21 + p9*x22 + p14*x23 + p19*x24 + p24*x25; 424 pc[24] = m25 = p5*x21 + p10*x22 + p15*x23 + p20*x24 + p25*x25; 425 426 nz = bi[row+1] - diag_offset[row] - 1; 427 pv += 25; 428 for (j=0; j<nz; j++) { 429 x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3]; 430 x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8]; 431 x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; 432 x14 = pv[13]; x15 = pv[14]; x16 = pv[15]; x17 = pv[16]; x18 = pv[17]; 433 x19 = pv[18]; x20 = pv[19]; x21 = pv[20]; x22 = pv[21]; x23 = pv[22]; 434 x24 = pv[23]; x25 = pv[24]; 435 x = rtmp + 25*pj[j]; 436 x[0] -= m1*x1 + m6*x2 + m11*x3 + m16*x4 + m21*x5; 437 x[1] -= m2*x1 + m7*x2 + m12*x3 + m17*x4 + m22*x5; 438 x[2] -= m3*x1 + m8*x2 + m13*x3 + m18*x4 + m23*x5; 439 x[3] -= m4*x1 + m9*x2 + m14*x3 + m19*x4 + m24*x5; 440 x[4] -= m5*x1 + m10*x2 + m15*x3 + m20*x4 + m25*x5; 441 442 x[5] -= m1*x6 + m6*x7 + m11*x8 + m16*x9 + m21*x10; 443 x[6] -= m2*x6 + m7*x7 + m12*x8 + m17*x9 + m22*x10; 444 x[7] -= m3*x6 + m8*x7 + m13*x8 + m18*x9 + m23*x10; 445 x[8] -= m4*x6 + m9*x7 + m14*x8 + m19*x9 + m24*x10; 446 x[9] -= m5*x6 + m10*x7 + m15*x8 + m20*x9 + m25*x10; 447 448 x[10] -= m1*x11 + m6*x12 + m11*x13 + m16*x14 + m21*x15; 449 x[11] -= m2*x11 + m7*x12 + m12*x13 + m17*x14 + m22*x15; 450 x[12] -= m3*x11 + m8*x12 + m13*x13 + m18*x14 + m23*x15; 451 x[13] -= m4*x11 + m9*x12 + m14*x13 + m19*x14 + m24*x15; 452 x[14] -= m5*x11 + m10*x12 + m15*x13 + m20*x14 + m25*x15; 453 454 x[15] -= m1*x16 + m6*x17 + m11*x18 + m16*x19 + m21*x20; 455 x[16] -= m2*x16 + m7*x17 + m12*x18 + m17*x19 + m22*x20; 456 x[17] -= m3*x16 + m8*x17 + m13*x18 + m18*x19 + m23*x20; 457 x[18] -= m4*x16 + m9*x17 + m14*x18 + m19*x19 + m24*x20; 458 x[19] -= m5*x16 + m10*x17 + m15*x18 + m20*x19 + m25*x20; 459 460 x[20] -= m1*x21 + m6*x22 + m11*x23 + m16*x24 + m21*x25; 461 x[21] -= m2*x21 + m7*x22 + m12*x23 + m17*x24 + m22*x25; 462 x[22] -= m3*x21 + m8*x22 + m13*x23 + m18*x24 + m23*x25; 463 x[23] -= m4*x21 + m9*x22 + m14*x23 + m19*x24 + m24*x25; 464 x[24] -= m5*x21 + m10*x22 + m15*x23 + m20*x24 + m25*x25; 465 pv += 25; 466 } 467 PetscCall(PetscLogFlops(250.0*nz+225.0)); 468 } 469 row = *ajtmp++; 470 } 471 /* finished row so stick it into b->a */ 472 pv = ba + 25*bi[i]; 473 pj = bj + bi[i]; 474 nz = bi[i+1] - bi[i]; 475 for (j=0; j<nz; j++) { 476 x = rtmp+25*pj[j]; 477 pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3]; 478 pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8]; 479 pv[9] = x[9]; pv[10] = x[10]; pv[11] = x[11]; pv[12] = x[12]; 480 pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15]; pv[16] = x[16]; pv[17] = x[17]; 481 pv[18] = x[18]; pv[19] = x[19]; pv[20] = x[20]; pv[21] = x[21]; pv[22] = x[22]; 482 pv[23] = x[23]; pv[24] = x[24]; 483 pv += 25; 484 } 485 /* invert diagonal block */ 486 w = ba + 25*diag_offset[i]; 487 PetscCall(PetscKernel_A_gets_inverse_A_5(w,ipvt,work,shift,allowzeropivot,&zeropivotdetected)); 488 if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 489 } 490 491 PetscCall(PetscFree(rtmp)); 492 493 C->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering_inplace; 494 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering_inplace; 495 C->assembled = PETSC_TRUE; 496 497 PetscCall(PetscLogFlops(1.333333333333*5*5*5*b->mbs)); /* from inverting diagonal blocks */ 498 PetscFunctionReturn(0); 499 } 500 501 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info) 502 { 503 Mat C =B; 504 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data; 505 PetscInt i,j,k,nz,nzL,row; 506 const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 507 const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2; 508 MatScalar *rtmp,*pc,*mwork,*v,*vv,*pv,*aa=a->a,work[25]; 509 PetscInt flg,ipvt[5]; 510 PetscReal shift = info->shiftamount; 511 PetscBool allowzeropivot,zeropivotdetected; 512 513 PetscFunctionBegin; 514 allowzeropivot = PetscNot(A->erroriffailure); 515 516 /* generate work space needed by the factorization */ 517 PetscCall(PetscMalloc2(bs2*n,&rtmp,bs2,&mwork)); 518 PetscCall(PetscArrayzero(rtmp,bs2*n)); 519 520 for (i=0; i<n; i++) { 521 /* zero rtmp */ 522 /* L part */ 523 nz = bi[i+1] - bi[i]; 524 bjtmp = bj + bi[i]; 525 for (j=0; j<nz; j++) { 526 PetscCall(PetscArrayzero(rtmp+bs2*bjtmp[j],bs2)); 527 } 528 529 /* U part */ 530 nz = bdiag[i] - bdiag[i+1]; 531 bjtmp = bj + bdiag[i+1]+1; 532 for (j=0; j<nz; j++) { 533 PetscCall(PetscArrayzero(rtmp+bs2*bjtmp[j],bs2)); 534 } 535 536 /* load in initial (unfactored row) */ 537 nz = ai[i+1] - ai[i]; 538 ajtmp = aj + ai[i]; 539 v = aa + bs2*ai[i]; 540 for (j=0; j<nz; j++) { 541 PetscCall(PetscArraycpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2)); 542 } 543 544 /* elimination */ 545 bjtmp = bj + bi[i]; 546 nzL = bi[i+1] - bi[i]; 547 for (k=0; k < nzL; k++) { 548 row = bjtmp[k]; 549 pc = rtmp + bs2*row; 550 for (flg=0,j=0; j<bs2; j++) { 551 if (pc[j]!=0.0) { 552 flg = 1; 553 break; 554 } 555 } 556 if (flg) { 557 pv = b->a + bs2*bdiag[row]; 558 /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 559 PetscCall(PetscKernel_A_gets_A_times_B_5(pc,pv,mwork)); 560 561 pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */ 562 pv = b->a + bs2*(bdiag[row+1]+1); 563 nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */ 564 for (j=0; j<nz; j++) { 565 /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 566 /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 567 vv = rtmp + bs2*pj[j]; 568 PetscCall(PetscKernel_A_gets_A_minus_B_times_C_5(vv,pc,pv)); 569 pv += bs2; 570 } 571 PetscCall(PetscLogFlops(250.0*nz+225)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 572 } 573 } 574 575 /* finished row so stick it into b->a */ 576 /* L part */ 577 pv = b->a + bs2*bi[i]; 578 pj = b->j + bi[i]; 579 nz = bi[i+1] - bi[i]; 580 for (j=0; j<nz; j++) { 581 PetscCall(PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2)); 582 } 583 584 /* Mark diagonal and invert diagonal for simpler triangular solves */ 585 pv = b->a + bs2*bdiag[i]; 586 pj = b->j + bdiag[i]; 587 PetscCall(PetscArraycpy(pv,rtmp+bs2*pj[0],bs2)); 588 PetscCall(PetscKernel_A_gets_inverse_A_5(pv,ipvt,work,shift,allowzeropivot,&zeropivotdetected)); 589 if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 590 591 /* U part */ 592 pv = b->a + bs2*(bdiag[i+1]+1); 593 pj = b->j + bdiag[i+1]+1; 594 nz = bdiag[i] - bdiag[i+1] - 1; 595 for (j=0; j<nz; j++) { 596 PetscCall(PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2)); 597 } 598 } 599 PetscCall(PetscFree2(rtmp,mwork)); 600 601 C->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering; 602 C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering; 603 C->assembled = PETSC_TRUE; 604 605 PetscCall(PetscLogFlops(1.333333333333*5*5*5*n)); /* from inverting diagonal blocks */ 606 PetscFunctionReturn(0); 607 } 608