Actual source code: baijfact5.c


  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:       Version for when blocks are 7 by 7
  9: */
 10: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_7_inplace(Mat C,Mat A,const MatFactorInfo *info)
 11: {
 12:   Mat_SeqBAIJ    *a    = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
 13:   IS             isrow = b->row,isicol = b->icol;
 14:   const PetscInt *r,*ic,*bi = b->i,*bj = b->j,*ajtmp,*diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj,*ajtmpold;
 15:   PetscInt       i,j,n = a->mbs,nz,row,idx;
 16:   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
 17:   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
 18:   MatScalar      p5,p6,p7,p8,p9,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,x16;
 19:   MatScalar      x17,x18,x19,x20,x21,x22,x23,x24,x25,p10,p11,p12,p13,p14;
 20:   MatScalar      p15,p16,p17,p18,p19,p20,p21,p22,p23,p24,p25,m10,m11,m12;
 21:   MatScalar      m13,m14,m15,m16,m17,m18,m19,m20,m21,m22,m23,m24,m25;
 22:   MatScalar      p26,p27,p28,p29,p30,p31,p32,p33,p34,p35,p36;
 23:   MatScalar      p37,p38,p39,p40,p41,p42,p43,p44,p45,p46,p47,p48,p49;
 24:   MatScalar      x26,x27,x28,x29,x30,x31,x32,x33,x34,x35,x36;
 25:   MatScalar      x37,x38,x39,x40,x41,x42,x43,x44,x45,x46,x47,x48,x49;
 26:   MatScalar      m26,m27,m28,m29,m30,m31,m32,m33,m34,m35,m36;
 27:   MatScalar      m37,m38,m39,m40,m41,m42,m43,m44,m45,m46,m47,m48,m49;
 28:   MatScalar      *ba   = b->a,*aa = a->a;
 29:   PetscReal      shift = info->shiftamount;
 30:   PetscBool      allowzeropivot,zeropivotdetected;

 32:   allowzeropivot = PetscNot(A->erroriffailure);
 33:   ISGetIndices(isrow,&r);
 34:   ISGetIndices(isicol,&ic);
 35:   PetscMalloc1(49*(n+1),&rtmp);

 37:   for (i=0; i<n; i++) {
 38:     nz    = bi[i+1] - bi[i];
 39:     ajtmp = bj + bi[i];
 40:     for  (j=0; j<nz; j++) {
 41:       x     = rtmp+49*ajtmp[j];
 42:       x[0]  = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
 43:       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
 44:       x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = x[25] = 0.0;
 45:       x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0;
 46:       x[34] = x[35] = x[36] = x[37] = x[38] = x[39] = x[40] = x[41] = 0.0;
 47:       x[42] = x[43] = x[44] = x[45] = x[46] = x[47] = x[48] = 0.0;
 48:     }
 49:     /* load in initial (unfactored row) */
 50:     idx      = r[i];
 51:     nz       = ai[idx+1] - ai[idx];
 52:     ajtmpold = aj + ai[idx];
 53:     v        = aa + 49*ai[idx];
 54:     for (j=0; j<nz; j++) {
 55:       x     = rtmp+49*ic[ajtmpold[j]];
 56:       x[0]  =  v[0];  x[1] =  v[1];  x[2] =  v[2];  x[3] =  v[3];
 57:       x[4]  =  v[4];  x[5] =  v[5];  x[6] =  v[6];  x[7] =  v[7];
 58:       x[8]  =  v[8];  x[9] =  v[9];  x[10] = v[10]; x[11] = v[11];
 59:       x[12] = v[12]; x[13] = v[13]; x[14] = v[14]; x[15] = v[15];
 60:       x[16] = v[16]; x[17] = v[17]; x[18] = v[18]; x[19] = v[19];
 61:       x[20] = v[20]; x[21] = v[21]; x[22] = v[22]; x[23] = v[23];
 62:       x[24] = v[24]; x[25] = v[25]; x[26] = v[26]; x[27] = v[27];
 63:       x[28] = v[28]; x[29] = v[29]; x[30] = v[30]; x[31] = v[31];
 64:       x[32] = v[32]; x[33] = v[33]; x[34] = v[34]; x[35] = v[35];
 65:       x[36] = v[36]; x[37] = v[37]; x[38] = v[38]; x[39] = v[39];
 66:       x[40] = v[40]; x[41] = v[41]; x[42] = v[42]; x[43] = v[43];
 67:       x[44] = v[44]; x[45] = v[45]; x[46] = v[46]; x[47] = v[47];
 68:       x[48] = v[48];
 69:       v    += 49;
 70:     }
 71:     row = *ajtmp++;
 72:     while (row < i) {
 73:       pc  =  rtmp + 49*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];
 76:       p9  = pc[8];  p10 = pc[9];  p11 = pc[10]; p12 = pc[11];
 77:       p13 = pc[12]; p14 = pc[13]; p15 = pc[14]; p16 = pc[15];
 78:       p17 = pc[16]; p18 = pc[17]; p19 = pc[18]; p20 = pc[19];
 79:       p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23];
 80:       p25 = pc[24]; p26 = pc[25]; p27 = pc[26]; p28 = pc[27];
 81:       p29 = pc[28]; p30 = pc[29]; p31 = pc[30]; p32 = pc[31];
 82:       p33 = pc[32]; p34 = pc[33]; p35 = pc[34]; p36 = pc[35];
 83:       p37 = pc[36]; p38 = pc[37]; p39 = pc[38]; p40 = pc[39];
 84:       p41 = pc[40]; p42 = pc[41]; p43 = pc[42]; p44 = pc[43];
 85:       p45 = pc[44]; p46 = pc[45]; p47 = pc[46]; p48 = pc[47];
 86:       p49 = pc[48];
 87:       if (p1  != 0.0 || p2  != 0.0 || p3  != 0.0 || p4  != 0.0 ||
 88:           p5  != 0.0 || p6  != 0.0 || p7  != 0.0 || p8  != 0.0 ||
 89:           p9  != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 ||
 90:           p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 ||
 91:           p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 ||
 92:           p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 ||
 93:           p25 != 0.0 || p26 != 0.0 || p27 != 0.0 || p28 != 0.0 ||
 94:           p29 != 0.0 || p30 != 0.0 || p31 != 0.0 || p32 != 0.0 ||
 95:           p33 != 0.0 || p34 != 0.0 || p35 != 0.0 || p36 != 0.0 ||
 96:           p37 != 0.0 || p38 != 0.0 || p39 != 0.0 || p40 != 0.0 ||
 97:           p41 != 0.0 || p42 != 0.0 || p43 != 0.0 || p44 != 0.0 ||
 98:           p45 != 0.0 || p46 != 0.0 || p47 != 0.0 || p48 != 0.0 ||
 99:           p49 != 0.0) {
100:         pv    = ba + 49*diag_offset[row];
101:         pj    = bj + diag_offset[row] + 1;
102:         x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
103:         x5    = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];
104:         x9    = pv[8];  x10 = pv[9];  x11 = pv[10]; x12 = pv[11];
105:         x13   = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
106:         x17   = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
107:         x21   = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
108:         x25   = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
109:         x29   = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
110:         x33   = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
111:         x37   = pv[36]; x38 = pv[37]; x39 = pv[38]; x40 = pv[39];
112:         x41   = pv[40]; x42 = pv[41]; x43 = pv[42]; x44 = pv[43];
113:         x45   = pv[44]; x46 = pv[45]; x47 = pv[46]; x48 = pv[47];
114:         x49   = pv[48];
115:         pc[0] = m1  = p1*x1  + p8*x2   + p15*x3  + p22*x4  + p29*x5  + p36*x6 + p43*x7;
116:         pc[1] = m2  = p2*x1  + p9*x2   + p16*x3  + p23*x4  + p30*x5  + p37*x6 + p44*x7;
117:         pc[2] = m3  = p3*x1  + p10*x2  + p17*x3  + p24*x4  + p31*x5  + p38*x6 + p45*x7;
118:         pc[3] = m4  = p4*x1  + p11*x2  + p18*x3  + p25*x4  + p32*x5  + p39*x6 + p46*x7;
119:         pc[4] = m5  = p5*x1  + p12*x2  + p19*x3  + p26*x4  + p33*x5  + p40*x6 + p47*x7;
120:         pc[5] = m6  = p6*x1  + p13*x2  + p20*x3  + p27*x4  + p34*x5  + p41*x6 + p48*x7;
121:         pc[6] = m7  = p7*x1  + p14*x2  + p21*x3  + p28*x4  + p35*x5  + p42*x6 + p49*x7;

123:         pc[7]  = m8  = p1*x8  + p8*x9   + p15*x10 + p22*x11 + p29*x12 + p36*x13 + p43*x14;
124:         pc[8]  = m9  = p2*x8  + p9*x9   + p16*x10 + p23*x11 + p30*x12 + p37*x13 + p44*x14;
125:         pc[9]  = m10 = p3*x8  + p10*x9  + p17*x10 + p24*x11 + p31*x12 + p38*x13 + p45*x14;
126:         pc[10] = m11 = p4*x8  + p11*x9  + p18*x10 + p25*x11 + p32*x12 + p39*x13 + p46*x14;
127:         pc[11] = m12 = p5*x8  + p12*x9  + p19*x10 + p26*x11 + p33*x12 + p40*x13 + p47*x14;
128:         pc[12] = m13 = p6*x8  + p13*x9  + p20*x10 + p27*x11 + p34*x12 + p41*x13 + p48*x14;
129:         pc[13] = m14 = p7*x8  + p14*x9  + p21*x10 + p28*x11 + p35*x12 + p42*x13 + p49*x14;

131:         pc[14] = m15 = p1*x15 + p8*x16  + p15*x17 + p22*x18 + p29*x19 + p36*x20 + p43*x21;
132:         pc[15] = m16 = p2*x15 + p9*x16  + p16*x17 + p23*x18 + p30*x19 + p37*x20 + p44*x21;
133:         pc[16] = m17 = p3*x15 + p10*x16 + p17*x17 + p24*x18 + p31*x19 + p38*x20 + p45*x21;
134:         pc[17] = m18 = p4*x15 + p11*x16 + p18*x17 + p25*x18 + p32*x19 + p39*x20 + p46*x21;
135:         pc[18] = m19 = p5*x15 + p12*x16 + p19*x17 + p26*x18 + p33*x19 + p40*x20 + p47*x21;
136:         pc[19] = m20 = p6*x15 + p13*x16 + p20*x17 + p27*x18 + p34*x19 + p41*x20 + p48*x21;
137:         pc[20] = m21 = p7*x15 + p14*x16 + p21*x17 + p28*x18 + p35*x19 + p42*x20 + p49*x21;

139:         pc[21] = m22 = p1*x22 + p8*x23  + p15*x24 + p22*x25 + p29*x26 + p36*x27 + p43*x28;
140:         pc[22] = m23 = p2*x22 + p9*x23  + p16*x24 + p23*x25 + p30*x26 + p37*x27 + p44*x28;
141:         pc[23] = m24 = p3*x22 + p10*x23 + p17*x24 + p24*x25 + p31*x26 + p38*x27 + p45*x28;
142:         pc[24] = m25 = p4*x22 + p11*x23 + p18*x24 + p25*x25 + p32*x26 + p39*x27 + p46*x28;
143:         pc[25] = m26 = p5*x22 + p12*x23 + p19*x24 + p26*x25 + p33*x26 + p40*x27 + p47*x28;
144:         pc[26] = m27 = p6*x22 + p13*x23 + p20*x24 + p27*x25 + p34*x26 + p41*x27 + p48*x28;
145:         pc[27] = m28 = p7*x22 + p14*x23 + p21*x24 + p28*x25 + p35*x26 + p42*x27 + p49*x28;

147:         pc[28] = m29 = p1*x29 + p8*x30  + p15*x31 + p22*x32 + p29*x33 + p36*x34 + p43*x35;
148:         pc[29] = m30 = p2*x29 + p9*x30  + p16*x31 + p23*x32 + p30*x33 + p37*x34 + p44*x35;
149:         pc[30] = m31 = p3*x29 + p10*x30 + p17*x31 + p24*x32 + p31*x33 + p38*x34 + p45*x35;
150:         pc[31] = m32 = p4*x29 + p11*x30 + p18*x31 + p25*x32 + p32*x33 + p39*x34 + p46*x35;
151:         pc[32] = m33 = p5*x29 + p12*x30 + p19*x31 + p26*x32 + p33*x33 + p40*x34 + p47*x35;
152:         pc[33] = m34 = p6*x29 + p13*x30 + p20*x31 + p27*x32 + p34*x33 + p41*x34 + p48*x35;
153:         pc[34] = m35 = p7*x29 + p14*x30 + p21*x31 + p28*x32 + p35*x33 + p42*x34 + p49*x35;

155:         pc[35] = m36 = p1*x36 + p8*x37  + p15*x38 + p22*x39 + p29*x40 + p36*x41 + p43*x42;
156:         pc[36] = m37 = p2*x36 + p9*x37  + p16*x38 + p23*x39 + p30*x40 + p37*x41 + p44*x42;
157:         pc[37] = m38 = p3*x36 + p10*x37 + p17*x38 + p24*x39 + p31*x40 + p38*x41 + p45*x42;
158:         pc[38] = m39 = p4*x36 + p11*x37 + p18*x38 + p25*x39 + p32*x40 + p39*x41 + p46*x42;
159:         pc[39] = m40 = p5*x36 + p12*x37 + p19*x38 + p26*x39 + p33*x40 + p40*x41 + p47*x42;
160:         pc[40] = m41 = p6*x36 + p13*x37 + p20*x38 + p27*x39 + p34*x40 + p41*x41 + p48*x42;
161:         pc[41] = m42 = p7*x36 + p14*x37 + p21*x38 + p28*x39 + p35*x40 + p42*x41 + p49*x42;

163:         pc[42] = m43 = p1*x43 + p8*x44  + p15*x45 + p22*x46 + p29*x47 + p36*x48 + p43*x49;
164:         pc[43] = m44 = p2*x43 + p9*x44  + p16*x45 + p23*x46 + p30*x47 + p37*x48 + p44*x49;
165:         pc[44] = m45 = p3*x43 + p10*x44 + p17*x45 + p24*x46 + p31*x47 + p38*x48 + p45*x49;
166:         pc[45] = m46 = p4*x43 + p11*x44 + p18*x45 + p25*x46 + p32*x47 + p39*x48 + p46*x49;
167:         pc[46] = m47 = p5*x43 + p12*x44 + p19*x45 + p26*x46 + p33*x47 + p40*x48 + p47*x49;
168:         pc[47] = m48 = p6*x43 + p13*x44 + p20*x45 + p27*x46 + p34*x47 + p41*x48 + p48*x49;
169:         pc[48] = m49 = p7*x43 + p14*x44 + p21*x45 + p28*x46 + p35*x47 + p42*x48 + p49*x49;

171:         nz  = bi[row+1] - diag_offset[row] - 1;
172:         pv += 49;
173:         for (j=0; j<nz; j++) {
174:           x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
175:           x5    = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];
176:           x9    = pv[8];  x10 = pv[9];  x11 = pv[10]; x12 = pv[11];
177:           x13   = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
178:           x17   = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
179:           x21   = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
180:           x25   = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
181:           x29   = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
182:           x33   = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
183:           x37   = pv[36]; x38 = pv[37]; x39 = pv[38]; x40 = pv[39];
184:           x41   = pv[40]; x42 = pv[41]; x43 = pv[42]; x44 = pv[43];
185:           x45   = pv[44]; x46 = pv[45]; x47 = pv[46]; x48 = pv[47];
186:           x49   = pv[48];
187:           x     = rtmp + 49*pj[j];
188:           x[0] -= m1*x1  + m8*x2   + m15*x3  + m22*x4  + m29*x5  + m36*x6 + m43*x7;
189:           x[1] -= m2*x1  + m9*x2   + m16*x3  + m23*x4  + m30*x5  + m37*x6 + m44*x7;
190:           x[2] -= m3*x1  + m10*x2  + m17*x3  + m24*x4  + m31*x5  + m38*x6 + m45*x7;
191:           x[3] -= m4*x1  + m11*x2  + m18*x3  + m25*x4  + m32*x5  + m39*x6 + m46*x7;
192:           x[4] -= m5*x1  + m12*x2  + m19*x3  + m26*x4  + m33*x5  + m40*x6 + m47*x7;
193:           x[5] -= m6*x1  + m13*x2  + m20*x3  + m27*x4  + m34*x5  + m41*x6 + m48*x7;
194:           x[6] -= m7*x1  + m14*x2  + m21*x3  + m28*x4  + m35*x5  + m42*x6 + m49*x7;

196:           x[7]  -= m1*x8  + m8*x9   + m15*x10 + m22*x11 + m29*x12 + m36*x13 + m43*x14;
197:           x[8]  -= m2*x8  + m9*x9   + m16*x10 + m23*x11 + m30*x12 + m37*x13 + m44*x14;
198:           x[9]  -= m3*x8  + m10*x9  + m17*x10 + m24*x11 + m31*x12 + m38*x13 + m45*x14;
199:           x[10] -= m4*x8  + m11*x9  + m18*x10 + m25*x11 + m32*x12 + m39*x13 + m46*x14;
200:           x[11] -= m5*x8  + m12*x9  + m19*x10 + m26*x11 + m33*x12 + m40*x13 + m47*x14;
201:           x[12] -= m6*x8  + m13*x9  + m20*x10 + m27*x11 + m34*x12 + m41*x13 + m48*x14;
202:           x[13] -= m7*x8  + m14*x9  + m21*x10 + m28*x11 + m35*x12 + m42*x13 + m49*x14;

204:           x[14] -= m1*x15 + m8*x16  + m15*x17 + m22*x18 + m29*x19 + m36*x20 + m43*x21;
205:           x[15] -= m2*x15 + m9*x16  + m16*x17 + m23*x18 + m30*x19 + m37*x20 + m44*x21;
206:           x[16] -= m3*x15 + m10*x16 + m17*x17 + m24*x18 + m31*x19 + m38*x20 + m45*x21;
207:           x[17] -= m4*x15 + m11*x16 + m18*x17 + m25*x18 + m32*x19 + m39*x20 + m46*x21;
208:           x[18] -= m5*x15 + m12*x16 + m19*x17 + m26*x18 + m33*x19 + m40*x20 + m47*x21;
209:           x[19] -= m6*x15 + m13*x16 + m20*x17 + m27*x18 + m34*x19 + m41*x20 + m48*x21;
210:           x[20] -= m7*x15 + m14*x16 + m21*x17 + m28*x18 + m35*x19 + m42*x20 + m49*x21;

212:           x[21] -= m1*x22 + m8*x23  + m15*x24 + m22*x25 + m29*x26 + m36*x27 + m43*x28;
213:           x[22] -= m2*x22 + m9*x23  + m16*x24 + m23*x25 + m30*x26 + m37*x27 + m44*x28;
214:           x[23] -= m3*x22 + m10*x23 + m17*x24 + m24*x25 + m31*x26 + m38*x27 + m45*x28;
215:           x[24] -= m4*x22 + m11*x23 + m18*x24 + m25*x25 + m32*x26 + m39*x27 + m46*x28;
216:           x[25] -= m5*x22 + m12*x23 + m19*x24 + m26*x25 + m33*x26 + m40*x27 + m47*x28;
217:           x[26] -= m6*x22 + m13*x23 + m20*x24 + m27*x25 + m34*x26 + m41*x27 + m48*x28;
218:           x[27] -= m7*x22 + m14*x23 + m21*x24 + m28*x25 + m35*x26 + m42*x27 + m49*x28;

220:           x[28] -= m1*x29 + m8*x30  + m15*x31 + m22*x32 + m29*x33 + m36*x34 + m43*x35;
221:           x[29] -= m2*x29 + m9*x30  + m16*x31 + m23*x32 + m30*x33 + m37*x34 + m44*x35;
222:           x[30] -= m3*x29 + m10*x30 + m17*x31 + m24*x32 + m31*x33 + m38*x34 + m45*x35;
223:           x[31] -= m4*x29 + m11*x30 + m18*x31 + m25*x32 + m32*x33 + m39*x34 + m46*x35;
224:           x[32] -= m5*x29 + m12*x30 + m19*x31 + m26*x32 + m33*x33 + m40*x34 + m47*x35;
225:           x[33] -= m6*x29 + m13*x30 + m20*x31 + m27*x32 + m34*x33 + m41*x34 + m48*x35;
226:           x[34] -= m7*x29 + m14*x30 + m21*x31 + m28*x32 + m35*x33 + m42*x34 + m49*x35;

228:           x[35] -= m1*x36 + m8*x37  + m15*x38 + m22*x39 + m29*x40 + m36*x41 + m43*x42;
229:           x[36] -= m2*x36 + m9*x37  + m16*x38 + m23*x39 + m30*x40 + m37*x41 + m44*x42;
230:           x[37] -= m3*x36 + m10*x37 + m17*x38 + m24*x39 + m31*x40 + m38*x41 + m45*x42;
231:           x[38] -= m4*x36 + m11*x37 + m18*x38 + m25*x39 + m32*x40 + m39*x41 + m46*x42;
232:           x[39] -= m5*x36 + m12*x37 + m19*x38 + m26*x39 + m33*x40 + m40*x41 + m47*x42;
233:           x[40] -= m6*x36 + m13*x37 + m20*x38 + m27*x39 + m34*x40 + m41*x41 + m48*x42;
234:           x[41] -= m7*x36 + m14*x37 + m21*x38 + m28*x39 + m35*x40 + m42*x41 + m49*x42;

236:           x[42] -= m1*x43 + m8*x44  + m15*x45 + m22*x46 + m29*x47 + m36*x48 + m43*x49;
237:           x[43] -= m2*x43 + m9*x44  + m16*x45 + m23*x46 + m30*x47 + m37*x48 + m44*x49;
238:           x[44] -= m3*x43 + m10*x44 + m17*x45 + m24*x46 + m31*x47 + m38*x48 + m45*x49;
239:           x[45] -= m4*x43 + m11*x44 + m18*x45 + m25*x46 + m32*x47 + m39*x48 + m46*x49;
240:           x[46] -= m5*x43 + m12*x44 + m19*x45 + m26*x46 + m33*x47 + m40*x48 + m47*x49;
241:           x[47] -= m6*x43 + m13*x44 + m20*x45 + m27*x46 + m34*x47 + m41*x48 + m48*x49;
242:           x[48] -= m7*x43 + m14*x44 + m21*x45 + m28*x46 + m35*x47 + m42*x48 + m49*x49;
243:           pv    += 49;
244:         }
245:         PetscLogFlops(686.0*nz+637.0);
246:       }
247:       row = *ajtmp++;
248:     }
249:     /* finished row so stick it into b->a */
250:     pv = ba + 49*bi[i];
251:     pj = bj + bi[i];
252:     nz = bi[i+1] - bi[i];
253:     for (j=0; j<nz; j++) {
254:       x      = rtmp+49*pj[j];
255:       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
256:       pv[4]  = x[4];  pv[5]  = x[5];  pv[6]  = x[6];  pv[7]  = x[7];
257:       pv[8]  = x[8];  pv[9]  = x[9];  pv[10] = x[10]; pv[11] = x[11];
258:       pv[12] = x[12]; pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15];
259:       pv[16] = x[16]; pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19];
260:       pv[20] = x[20]; pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23];
261:       pv[24] = x[24]; pv[25] = x[25]; pv[26] = x[26]; pv[27] = x[27];
262:       pv[28] = x[28]; pv[29] = x[29]; pv[30] = x[30]; pv[31] = x[31];
263:       pv[32] = x[32]; pv[33] = x[33]; pv[34] = x[34]; pv[35] = x[35];
264:       pv[36] = x[36]; pv[37] = x[37]; pv[38] = x[38]; pv[39] = x[39];
265:       pv[40] = x[40]; pv[41] = x[41]; pv[42] = x[42]; pv[43] = x[43];
266:       pv[44] = x[44]; pv[45] = x[45]; pv[46] = x[46]; pv[47] = x[47];
267:       pv[48] = x[48];
268:       pv    += 49;
269:     }
270:     /* invert diagonal block */
271:     w    = ba + 49*diag_offset[i];
272:     PetscKernel_A_gets_inverse_A_7(w,shift,allowzeropivot,&zeropivotdetected);
273:     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
274:   }

276:   PetscFree(rtmp);
277:   ISRestoreIndices(isicol,&ic);
278:   ISRestoreIndices(isrow,&r);

280:   C->ops->solve          = MatSolve_SeqBAIJ_7_inplace;
281:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_7_inplace;
282:   C->assembled           = PETSC_TRUE;

284:   PetscLogFlops(1.333333333333*7*7*7*b->mbs); /* from inverting diagonal blocks */
285:   return 0;
286: }

288: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_7(Mat B,Mat A,const MatFactorInfo *info)
289: {
290:   Mat            C     =B;
291:   Mat_SeqBAIJ    *a    =(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
292:   IS             isrow = b->row,isicol = b->icol;
293:   const PetscInt *r,*ic;
294:   PetscInt       i,j,k,nz,nzL,row;
295:   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
296:   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
297:   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
298:   PetscInt       flg;
299:   PetscReal      shift = info->shiftamount;
300:   PetscBool      allowzeropivot,zeropivotdetected;

302:   allowzeropivot = PetscNot(A->erroriffailure);
303:   ISGetIndices(isrow,&r);
304:   ISGetIndices(isicol,&ic);

306:   /* generate work space needed by the factorization */
307:   PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);
308:   PetscArrayzero(rtmp,bs2*n);

310:   for (i=0; i<n; i++) {
311:     /* zero rtmp */
312:     /* L part */
313:     nz    = bi[i+1] - bi[i];
314:     bjtmp = bj + bi[i];
315:     for  (j=0; j<nz; j++) {
316:       PetscArrayzero(rtmp+bs2*bjtmp[j],bs2);
317:     }

319:     /* U part */
320:     nz    = bdiag[i] - bdiag[i+1];
321:     bjtmp = bj + bdiag[i+1]+1;
322:     for  (j=0; j<nz; j++) {
323:       PetscArrayzero(rtmp+bs2*bjtmp[j],bs2);
324:     }

326:     /* load in initial (unfactored row) */
327:     nz    = ai[r[i]+1] - ai[r[i]];
328:     ajtmp = aj + ai[r[i]];
329:     v     = aa + bs2*ai[r[i]];
330:     for (j=0; j<nz; j++) {
331:       PetscArraycpy(rtmp+bs2*ic[ajtmp[j]],v+bs2*j,bs2);
332:     }

334:     /* elimination */
335:     bjtmp = bj + bi[i];
336:     nzL   = bi[i+1] - bi[i];
337:     for (k=0; k < nzL; k++) {
338:       row = bjtmp[k];
339:       pc  = rtmp + bs2*row;
340:       for (flg=0,j=0; j<bs2; j++) {
341:         if (pc[j]!=0.0) {
342:           flg = 1;
343:           break;
344:         }
345:       }
346:       if (flg) {
347:         pv = b->a + bs2*bdiag[row];
348:         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
349:         PetscKernel_A_gets_A_times_B_7(pc,pv,mwork);

351:         pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
352:         pv = b->a + bs2*(bdiag[row+1]+1);
353:         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
354:         for (j=0; j<nz; j++) {
355:           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
356:           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
357:           v    = rtmp + bs2*pj[j];
358:           PetscKernel_A_gets_A_minus_B_times_C_7(v,pc,pv);
359:           pv  += bs2;
360:         }
361:         PetscLogFlops(686.0*nz+637); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
362:       }
363:     }

365:     /* finished row so stick it into b->a */
366:     /* L part */
367:     pv = b->a + bs2*bi[i];
368:     pj = b->j + bi[i];
369:     nz = bi[i+1] - bi[i];
370:     for (j=0; j<nz; j++) {
371:       PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2);
372:     }

374:     /* Mark diagonal and invert diagonal for simpler triangular solves */
375:     pv   = b->a + bs2*bdiag[i];
376:     pj   = b->j + bdiag[i];
377:     PetscArraycpy(pv,rtmp+bs2*pj[0],bs2);
378:     PetscKernel_A_gets_inverse_A_7(pv,shift,allowzeropivot,&zeropivotdetected);
379:     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;

381:     /* U part */
382:     pv = b->a + bs2*(bdiag[i+1]+1);
383:     pj = b->j + bdiag[i+1]+1;
384:     nz = bdiag[i] - bdiag[i+1] - 1;
385:     for (j=0; j<nz; j++) {
386:       PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2);
387:     }
388:   }

390:   PetscFree2(rtmp,mwork);
391:   ISRestoreIndices(isicol,&ic);
392:   ISRestoreIndices(isrow,&r);

394:   C->ops->solve          = MatSolve_SeqBAIJ_7;
395:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_7;
396:   C->assembled           = PETSC_TRUE;

398:   PetscLogFlops(1.333333333333*7*7*7*n); /* from inverting diagonal blocks */
399:   return 0;
400: }

402: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering_inplace(Mat C,Mat A,const MatFactorInfo *info)
403: {
404:   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ*)C->data;
405:   PetscInt       i,j,n = a->mbs,*bi = b->i,*bj = b->j;
406:   PetscInt       *ajtmpold,*ajtmp,nz,row;
407:   PetscInt       *diag_offset = b->diag,*ai=a->i,*aj=a->j,*pj;
408:   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
409:   MatScalar      x1,x2,x3,x4,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15;
410:   MatScalar      x16,x17,x18,x19,x20,x21,x22,x23,x24,x25;
411:   MatScalar      p1,p2,p3,p4,p5,p6,p7,p8,p9,p10,p11,p12,p13,p14,p15;
412:   MatScalar      p16,p17,p18,p19,p20,p21,p22,p23,p24,p25;
413:   MatScalar      m1,m2,m3,m4,m5,m6,m7,m8,m9,m10,m11,m12,m13,m14,m15;
414:   MatScalar      m16,m17,m18,m19,m20,m21,m22,m23,m24,m25;
415:   MatScalar      p26,p27,p28,p29,p30,p31,p32,p33,p34,p35,p36;
416:   MatScalar      p37,p38,p39,p40,p41,p42,p43,p44,p45,p46,p47,p48,p49;
417:   MatScalar      x26,x27,x28,x29,x30,x31,x32,x33,x34,x35,x36;
418:   MatScalar      x37,x38,x39,x40,x41,x42,x43,x44,x45,x46,x47,x48,x49;
419:   MatScalar      m26,m27,m28,m29,m30,m31,m32,m33,m34,m35,m36;
420:   MatScalar      m37,m38,m39,m40,m41,m42,m43,m44,m45,m46,m47,m48,m49;
421:   MatScalar      *ba   = b->a,*aa = a->a;
422:   PetscReal      shift = info->shiftamount;
423:   PetscBool      allowzeropivot,zeropivotdetected;

425:   allowzeropivot = PetscNot(A->erroriffailure);
426:   PetscMalloc1(49*(n+1),&rtmp);
427:   for (i=0; i<n; i++) {
428:     nz    = bi[i+1] - bi[i];
429:     ajtmp = bj + bi[i];
430:     for  (j=0; j<nz; j++) {
431:       x     = rtmp+49*ajtmp[j];
432:       x[0]  = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
433:       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
434:       x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = x[25] = 0.0;
435:       x[26] = x[27] = x[28] = x[29] = x[30] = x[31] = x[32] = x[33] = 0.0;
436:       x[34] = x[35] = x[36] = x[37] = x[38] = x[39] = x[40] = x[41] = 0.0;
437:       x[42] = x[43] = x[44] = x[45] = x[46] = x[47] = x[48] = 0.0;
438:     }
439:     /* load in initial (unfactored row) */
440:     nz       = ai[i+1] - ai[i];
441:     ajtmpold = aj + ai[i];
442:     v        = aa + 49*ai[i];
443:     for (j=0; j<nz; j++) {
444:       x     = rtmp+49*ajtmpold[j];
445:       x[0]  =  v[0];  x[1] =  v[1];  x[2] =  v[2];  x[3] =  v[3];
446:       x[4]  =  v[4];  x[5] =  v[5];  x[6] =  v[6];  x[7] =  v[7];
447:       x[8]  =  v[8];  x[9] =  v[9];  x[10] = v[10]; x[11] = v[11];
448:       x[12] = v[12]; x[13] = v[13]; x[14] = v[14]; x[15] = v[15];
449:       x[16] = v[16]; x[17] = v[17]; x[18] = v[18]; x[19] = v[19];
450:       x[20] = v[20]; x[21] = v[21]; x[22] = v[22]; x[23] = v[23];
451:       x[24] = v[24]; x[25] = v[25]; x[26] = v[26]; x[27] = v[27];
452:       x[28] = v[28]; x[29] = v[29]; x[30] = v[30]; x[31] = v[31];
453:       x[32] = v[32]; x[33] = v[33]; x[34] = v[34]; x[35] = v[35];
454:       x[36] = v[36]; x[37] = v[37]; x[38] = v[38]; x[39] = v[39];
455:       x[40] = v[40]; x[41] = v[41]; x[42] = v[42]; x[43] = v[43];
456:       x[44] = v[44]; x[45] = v[45]; x[46] = v[46]; x[47] = v[47];
457:       x[48] = v[48];
458:       v    += 49;
459:     }
460:     row = *ajtmp++;
461:     while (row < i) {
462:       pc  = rtmp + 49*row;
463:       p1  = pc[0];  p2  = pc[1];  p3  = pc[2];  p4  = pc[3];
464:       p5  = pc[4];  p6  = pc[5];  p7  = pc[6];  p8  = pc[7];
465:       p9  = pc[8];  p10 = pc[9];  p11 = pc[10]; p12 = pc[11];
466:       p13 = pc[12]; p14 = pc[13]; p15 = pc[14]; p16 = pc[15];
467:       p17 = pc[16]; p18 = pc[17]; p19 = pc[18]; p20 = pc[19];
468:       p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23];
469:       p25 = pc[24]; p26 = pc[25]; p27 = pc[26]; p28 = pc[27];
470:       p29 = pc[28]; p30 = pc[29]; p31 = pc[30]; p32 = pc[31];
471:       p33 = pc[32]; p34 = pc[33]; p35 = pc[34]; p36 = pc[35];
472:       p37 = pc[36]; p38 = pc[37]; p39 = pc[38]; p40 = pc[39];
473:       p41 = pc[40]; p42 = pc[41]; p43 = pc[42]; p44 = pc[43];
474:       p45 = pc[44]; p46 = pc[45]; p47 = pc[46]; p48 = pc[47];
475:       p49 = pc[48];
476:       if (p1  != 0.0 || p2  != 0.0 || p3  != 0.0 || p4  != 0.0 ||
477:           p5  != 0.0 || p6  != 0.0 || p7  != 0.0 || p8  != 0.0 ||
478:           p9  != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 ||
479:           p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 ||
480:           p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 ||
481:           p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 ||
482:           p25 != 0.0 || p26 != 0.0 || p27 != 0.0 || p28 != 0.0 ||
483:           p29 != 0.0 || p30 != 0.0 || p31 != 0.0 || p32 != 0.0 ||
484:           p33 != 0.0 || p34 != 0.0 || p35 != 0.0 || p36 != 0.0 ||
485:           p37 != 0.0 || p38 != 0.0 || p39 != 0.0 || p40 != 0.0 ||
486:           p41 != 0.0 || p42 != 0.0 || p43 != 0.0 || p44 != 0.0 ||
487:           p45 != 0.0 || p46 != 0.0 || p47 != 0.0 || p48 != 0.0 ||
488:           p49 != 0.0) {
489:         pv    = ba + 49*diag_offset[row];
490:         pj    = bj + diag_offset[row] + 1;
491:         x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
492:         x5    = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];
493:         x9    = pv[8];  x10 = pv[9];  x11 = pv[10]; x12 = pv[11];
494:         x13   = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
495:         x17   = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
496:         x21   = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
497:         x25   = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
498:         x29   = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
499:         x33   = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
500:         x37   = pv[36]; x38 = pv[37]; x39 = pv[38]; x40 = pv[39];
501:         x41   = pv[40]; x42 = pv[41]; x43 = pv[42]; x44 = pv[43];
502:         x45   = pv[44]; x46 = pv[45]; x47 = pv[46]; x48 = pv[47];
503:         x49   = pv[48];
504:         pc[0] = m1  = p1*x1  + p8*x2   + p15*x3  + p22*x4  + p29*x5  + p36*x6 + p43*x7;
505:         pc[1] = m2  = p2*x1  + p9*x2   + p16*x3  + p23*x4  + p30*x5  + p37*x6 + p44*x7;
506:         pc[2] = m3  = p3*x1  + p10*x2  + p17*x3  + p24*x4  + p31*x5  + p38*x6 + p45*x7;
507:         pc[3] = m4  = p4*x1  + p11*x2  + p18*x3  + p25*x4  + p32*x5  + p39*x6 + p46*x7;
508:         pc[4] = m5  = p5*x1  + p12*x2  + p19*x3  + p26*x4  + p33*x5  + p40*x6 + p47*x7;
509:         pc[5] = m6  = p6*x1  + p13*x2  + p20*x3  + p27*x4  + p34*x5  + p41*x6 + p48*x7;
510:         pc[6] = m7  = p7*x1  + p14*x2  + p21*x3  + p28*x4  + p35*x5  + p42*x6 + p49*x7;

512:         pc[7]  = m8  = p1*x8  + p8*x9   + p15*x10 + p22*x11 + p29*x12 + p36*x13 + p43*x14;
513:         pc[8]  = m9  = p2*x8  + p9*x9   + p16*x10 + p23*x11 + p30*x12 + p37*x13 + p44*x14;
514:         pc[9]  = m10 = p3*x8  + p10*x9  + p17*x10 + p24*x11 + p31*x12 + p38*x13 + p45*x14;
515:         pc[10] = m11 = p4*x8  + p11*x9  + p18*x10 + p25*x11 + p32*x12 + p39*x13 + p46*x14;
516:         pc[11] = m12 = p5*x8  + p12*x9  + p19*x10 + p26*x11 + p33*x12 + p40*x13 + p47*x14;
517:         pc[12] = m13 = p6*x8  + p13*x9  + p20*x10 + p27*x11 + p34*x12 + p41*x13 + p48*x14;
518:         pc[13] = m14 = p7*x8  + p14*x9  + p21*x10 + p28*x11 + p35*x12 + p42*x13 + p49*x14;

520:         pc[14] = m15 = p1*x15 + p8*x16  + p15*x17 + p22*x18 + p29*x19 + p36*x20 + p43*x21;
521:         pc[15] = m16 = p2*x15 + p9*x16  + p16*x17 + p23*x18 + p30*x19 + p37*x20 + p44*x21;
522:         pc[16] = m17 = p3*x15 + p10*x16 + p17*x17 + p24*x18 + p31*x19 + p38*x20 + p45*x21;
523:         pc[17] = m18 = p4*x15 + p11*x16 + p18*x17 + p25*x18 + p32*x19 + p39*x20 + p46*x21;
524:         pc[18] = m19 = p5*x15 + p12*x16 + p19*x17 + p26*x18 + p33*x19 + p40*x20 + p47*x21;
525:         pc[19] = m20 = p6*x15 + p13*x16 + p20*x17 + p27*x18 + p34*x19 + p41*x20 + p48*x21;
526:         pc[20] = m21 = p7*x15 + p14*x16 + p21*x17 + p28*x18 + p35*x19 + p42*x20 + p49*x21;

528:         pc[21] = m22 = p1*x22 + p8*x23  + p15*x24 + p22*x25 + p29*x26 + p36*x27 + p43*x28;
529:         pc[22] = m23 = p2*x22 + p9*x23  + p16*x24 + p23*x25 + p30*x26 + p37*x27 + p44*x28;
530:         pc[23] = m24 = p3*x22 + p10*x23 + p17*x24 + p24*x25 + p31*x26 + p38*x27 + p45*x28;
531:         pc[24] = m25 = p4*x22 + p11*x23 + p18*x24 + p25*x25 + p32*x26 + p39*x27 + p46*x28;
532:         pc[25] = m26 = p5*x22 + p12*x23 + p19*x24 + p26*x25 + p33*x26 + p40*x27 + p47*x28;
533:         pc[26] = m27 = p6*x22 + p13*x23 + p20*x24 + p27*x25 + p34*x26 + p41*x27 + p48*x28;
534:         pc[27] = m28 = p7*x22 + p14*x23 + p21*x24 + p28*x25 + p35*x26 + p42*x27 + p49*x28;

536:         pc[28] = m29 = p1*x29 + p8*x30  + p15*x31 + p22*x32 + p29*x33 + p36*x34 + p43*x35;
537:         pc[29] = m30 = p2*x29 + p9*x30  + p16*x31 + p23*x32 + p30*x33 + p37*x34 + p44*x35;
538:         pc[30] = m31 = p3*x29 + p10*x30 + p17*x31 + p24*x32 + p31*x33 + p38*x34 + p45*x35;
539:         pc[31] = m32 = p4*x29 + p11*x30 + p18*x31 + p25*x32 + p32*x33 + p39*x34 + p46*x35;
540:         pc[32] = m33 = p5*x29 + p12*x30 + p19*x31 + p26*x32 + p33*x33 + p40*x34 + p47*x35;
541:         pc[33] = m34 = p6*x29 + p13*x30 + p20*x31 + p27*x32 + p34*x33 + p41*x34 + p48*x35;
542:         pc[34] = m35 = p7*x29 + p14*x30 + p21*x31 + p28*x32 + p35*x33 + p42*x34 + p49*x35;

544:         pc[35] = m36 = p1*x36 + p8*x37  + p15*x38 + p22*x39 + p29*x40 + p36*x41 + p43*x42;
545:         pc[36] = m37 = p2*x36 + p9*x37  + p16*x38 + p23*x39 + p30*x40 + p37*x41 + p44*x42;
546:         pc[37] = m38 = p3*x36 + p10*x37 + p17*x38 + p24*x39 + p31*x40 + p38*x41 + p45*x42;
547:         pc[38] = m39 = p4*x36 + p11*x37 + p18*x38 + p25*x39 + p32*x40 + p39*x41 + p46*x42;
548:         pc[39] = m40 = p5*x36 + p12*x37 + p19*x38 + p26*x39 + p33*x40 + p40*x41 + p47*x42;
549:         pc[40] = m41 = p6*x36 + p13*x37 + p20*x38 + p27*x39 + p34*x40 + p41*x41 + p48*x42;
550:         pc[41] = m42 = p7*x36 + p14*x37 + p21*x38 + p28*x39 + p35*x40 + p42*x41 + p49*x42;

552:         pc[42] = m43 = p1*x43 + p8*x44  + p15*x45 + p22*x46 + p29*x47 + p36*x48 + p43*x49;
553:         pc[43] = m44 = p2*x43 + p9*x44  + p16*x45 + p23*x46 + p30*x47 + p37*x48 + p44*x49;
554:         pc[44] = m45 = p3*x43 + p10*x44 + p17*x45 + p24*x46 + p31*x47 + p38*x48 + p45*x49;
555:         pc[45] = m46 = p4*x43 + p11*x44 + p18*x45 + p25*x46 + p32*x47 + p39*x48 + p46*x49;
556:         pc[46] = m47 = p5*x43 + p12*x44 + p19*x45 + p26*x46 + p33*x47 + p40*x48 + p47*x49;
557:         pc[47] = m48 = p6*x43 + p13*x44 + p20*x45 + p27*x46 + p34*x47 + p41*x48 + p48*x49;
558:         pc[48] = m49 = p7*x43 + p14*x44 + p21*x45 + p28*x46 + p35*x47 + p42*x48 + p49*x49;

560:         nz  = bi[row+1] - diag_offset[row] - 1;
561:         pv += 49;
562:         for (j=0; j<nz; j++) {
563:           x1    = pv[0];  x2  = pv[1];  x3  = pv[2];  x4  = pv[3];
564:           x5    = pv[4];  x6  = pv[5];  x7  = pv[6];  x8  = pv[7];
565:           x9    = pv[8];  x10 = pv[9];  x11 = pv[10]; x12 = pv[11];
566:           x13   = pv[12]; x14 = pv[13]; x15 = pv[14]; x16 = pv[15];
567:           x17   = pv[16]; x18 = pv[17]; x19 = pv[18]; x20 = pv[19];
568:           x21   = pv[20]; x22 = pv[21]; x23 = pv[22]; x24 = pv[23];
569:           x25   = pv[24]; x26 = pv[25]; x27 = pv[26]; x28 = pv[27];
570:           x29   = pv[28]; x30 = pv[29]; x31 = pv[30]; x32 = pv[31];
571:           x33   = pv[32]; x34 = pv[33]; x35 = pv[34]; x36 = pv[35];
572:           x37   = pv[36]; x38 = pv[37]; x39 = pv[38]; x40 = pv[39];
573:           x41   = pv[40]; x42 = pv[41]; x43 = pv[42]; x44 = pv[43];
574:           x45   = pv[44]; x46 = pv[45]; x47 = pv[46]; x48 = pv[47];
575:           x49   = pv[48];
576:           x     = rtmp + 49*pj[j];
577:           x[0] -= m1*x1  + m8*x2   + m15*x3  + m22*x4  + m29*x5  + m36*x6 + m43*x7;
578:           x[1] -= m2*x1  + m9*x2   + m16*x3  + m23*x4  + m30*x5  + m37*x6 + m44*x7;
579:           x[2] -= m3*x1  + m10*x2  + m17*x3  + m24*x4  + m31*x5  + m38*x6 + m45*x7;
580:           x[3] -= m4*x1  + m11*x2  + m18*x3  + m25*x4  + m32*x5  + m39*x6 + m46*x7;
581:           x[4] -= m5*x1  + m12*x2  + m19*x3  + m26*x4  + m33*x5  + m40*x6 + m47*x7;
582:           x[5] -= m6*x1  + m13*x2  + m20*x3  + m27*x4  + m34*x5  + m41*x6 + m48*x7;
583:           x[6] -= m7*x1  + m14*x2  + m21*x3  + m28*x4  + m35*x5  + m42*x6 + m49*x7;

585:           x[7]  -= m1*x8  + m8*x9   + m15*x10 + m22*x11 + m29*x12 + m36*x13 + m43*x14;
586:           x[8]  -= m2*x8  + m9*x9   + m16*x10 + m23*x11 + m30*x12 + m37*x13 + m44*x14;
587:           x[9]  -= m3*x8  + m10*x9  + m17*x10 + m24*x11 + m31*x12 + m38*x13 + m45*x14;
588:           x[10] -= m4*x8  + m11*x9  + m18*x10 + m25*x11 + m32*x12 + m39*x13 + m46*x14;
589:           x[11] -= m5*x8  + m12*x9  + m19*x10 + m26*x11 + m33*x12 + m40*x13 + m47*x14;
590:           x[12] -= m6*x8  + m13*x9  + m20*x10 + m27*x11 + m34*x12 + m41*x13 + m48*x14;
591:           x[13] -= m7*x8  + m14*x9  + m21*x10 + m28*x11 + m35*x12 + m42*x13 + m49*x14;

593:           x[14] -= m1*x15 + m8*x16  + m15*x17 + m22*x18 + m29*x19 + m36*x20 + m43*x21;
594:           x[15] -= m2*x15 + m9*x16  + m16*x17 + m23*x18 + m30*x19 + m37*x20 + m44*x21;
595:           x[16] -= m3*x15 + m10*x16 + m17*x17 + m24*x18 + m31*x19 + m38*x20 + m45*x21;
596:           x[17] -= m4*x15 + m11*x16 + m18*x17 + m25*x18 + m32*x19 + m39*x20 + m46*x21;
597:           x[18] -= m5*x15 + m12*x16 + m19*x17 + m26*x18 + m33*x19 + m40*x20 + m47*x21;
598:           x[19] -= m6*x15 + m13*x16 + m20*x17 + m27*x18 + m34*x19 + m41*x20 + m48*x21;
599:           x[20] -= m7*x15 + m14*x16 + m21*x17 + m28*x18 + m35*x19 + m42*x20 + m49*x21;

601:           x[21] -= m1*x22 + m8*x23  + m15*x24 + m22*x25 + m29*x26 + m36*x27 + m43*x28;
602:           x[22] -= m2*x22 + m9*x23  + m16*x24 + m23*x25 + m30*x26 + m37*x27 + m44*x28;
603:           x[23] -= m3*x22 + m10*x23 + m17*x24 + m24*x25 + m31*x26 + m38*x27 + m45*x28;
604:           x[24] -= m4*x22 + m11*x23 + m18*x24 + m25*x25 + m32*x26 + m39*x27 + m46*x28;
605:           x[25] -= m5*x22 + m12*x23 + m19*x24 + m26*x25 + m33*x26 + m40*x27 + m47*x28;
606:           x[26] -= m6*x22 + m13*x23 + m20*x24 + m27*x25 + m34*x26 + m41*x27 + m48*x28;
607:           x[27] -= m7*x22 + m14*x23 + m21*x24 + m28*x25 + m35*x26 + m42*x27 + m49*x28;

609:           x[28] -= m1*x29 + m8*x30  + m15*x31 + m22*x32 + m29*x33 + m36*x34 + m43*x35;
610:           x[29] -= m2*x29 + m9*x30  + m16*x31 + m23*x32 + m30*x33 + m37*x34 + m44*x35;
611:           x[30] -= m3*x29 + m10*x30 + m17*x31 + m24*x32 + m31*x33 + m38*x34 + m45*x35;
612:           x[31] -= m4*x29 + m11*x30 + m18*x31 + m25*x32 + m32*x33 + m39*x34 + m46*x35;
613:           x[32] -= m5*x29 + m12*x30 + m19*x31 + m26*x32 + m33*x33 + m40*x34 + m47*x35;
614:           x[33] -= m6*x29 + m13*x30 + m20*x31 + m27*x32 + m34*x33 + m41*x34 + m48*x35;
615:           x[34] -= m7*x29 + m14*x30 + m21*x31 + m28*x32 + m35*x33 + m42*x34 + m49*x35;

617:           x[35] -= m1*x36 + m8*x37  + m15*x38 + m22*x39 + m29*x40 + m36*x41 + m43*x42;
618:           x[36] -= m2*x36 + m9*x37  + m16*x38 + m23*x39 + m30*x40 + m37*x41 + m44*x42;
619:           x[37] -= m3*x36 + m10*x37 + m17*x38 + m24*x39 + m31*x40 + m38*x41 + m45*x42;
620:           x[38] -= m4*x36 + m11*x37 + m18*x38 + m25*x39 + m32*x40 + m39*x41 + m46*x42;
621:           x[39] -= m5*x36 + m12*x37 + m19*x38 + m26*x39 + m33*x40 + m40*x41 + m47*x42;
622:           x[40] -= m6*x36 + m13*x37 + m20*x38 + m27*x39 + m34*x40 + m41*x41 + m48*x42;
623:           x[41] -= m7*x36 + m14*x37 + m21*x38 + m28*x39 + m35*x40 + m42*x41 + m49*x42;

625:           x[42] -= m1*x43 + m8*x44  + m15*x45 + m22*x46 + m29*x47 + m36*x48 + m43*x49;
626:           x[43] -= m2*x43 + m9*x44  + m16*x45 + m23*x46 + m30*x47 + m37*x48 + m44*x49;
627:           x[44] -= m3*x43 + m10*x44 + m17*x45 + m24*x46 + m31*x47 + m38*x48 + m45*x49;
628:           x[45] -= m4*x43 + m11*x44 + m18*x45 + m25*x46 + m32*x47 + m39*x48 + m46*x49;
629:           x[46] -= m5*x43 + m12*x44 + m19*x45 + m26*x46 + m33*x47 + m40*x48 + m47*x49;
630:           x[47] -= m6*x43 + m13*x44 + m20*x45 + m27*x46 + m34*x47 + m41*x48 + m48*x49;
631:           x[48] -= m7*x43 + m14*x44 + m21*x45 + m28*x46 + m35*x47 + m42*x48 + m49*x49;
632:           pv    += 49;
633:         }
634:         PetscLogFlops(686.0*nz+637.0);
635:       }
636:       row = *ajtmp++;
637:     }
638:     /* finished row so stick it into b->a */
639:     pv = ba + 49*bi[i];
640:     pj = bj + bi[i];
641:     nz = bi[i+1] - bi[i];
642:     for (j=0; j<nz; j++) {
643:       x      = rtmp+49*pj[j];
644:       pv[0]  = x[0];  pv[1]  = x[1];  pv[2]  = x[2];  pv[3]  = x[3];
645:       pv[4]  = x[4];  pv[5]  = x[5];  pv[6]  = x[6];  pv[7]  = x[7];
646:       pv[8]  = x[8];  pv[9]  = x[9];  pv[10] = x[10]; pv[11] = x[11];
647:       pv[12] = x[12]; pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15];
648:       pv[16] = x[16]; pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19];
649:       pv[20] = x[20]; pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23];
650:       pv[24] = x[24]; pv[25] = x[25]; pv[26] = x[26]; pv[27] = x[27];
651:       pv[28] = x[28]; pv[29] = x[29]; pv[30] = x[30]; pv[31] = x[31];
652:       pv[32] = x[32]; pv[33] = x[33]; pv[34] = x[34]; pv[35] = x[35];
653:       pv[36] = x[36]; pv[37] = x[37]; pv[38] = x[38]; pv[39] = x[39];
654:       pv[40] = x[40]; pv[41] = x[41]; pv[42] = x[42]; pv[43] = x[43];
655:       pv[44] = x[44]; pv[45] = x[45]; pv[46] = x[46]; pv[47] = x[47];
656:       pv[48] = x[48];
657:       pv    += 49;
658:     }
659:     /* invert diagonal block */
660:     w    = ba + 49*diag_offset[i];
661:     PetscKernel_A_gets_inverse_A_7(w,shift,allowzeropivot,&zeropivotdetected);
662:     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
663:   }

665:   PetscFree(rtmp);

667:   C->ops->solve          = MatSolve_SeqBAIJ_7_NaturalOrdering_inplace;
668:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_7_NaturalOrdering_inplace;
669:   C->assembled           = PETSC_TRUE;

671:   PetscLogFlops(1.333333333333*7*7*7*b->mbs); /* from inverting diagonal blocks */
672:   return 0;
673: }

675: PetscErrorCode MatLUFactorNumeric_SeqBAIJ_7_NaturalOrdering(Mat B,Mat A,const MatFactorInfo *info)
676: {
677:   Mat            C =B;
678:   Mat_SeqBAIJ    *a=(Mat_SeqBAIJ*)A->data,*b=(Mat_SeqBAIJ*)C->data;
679:   PetscInt       i,j,k,nz,nzL,row;
680:   const PetscInt n=a->mbs,*ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j;
681:   const PetscInt *ajtmp,*bjtmp,*bdiag=b->diag,*pj,bs2=a->bs2;
682:   MatScalar      *rtmp,*pc,*mwork,*v,*pv,*aa=a->a;
683:   PetscInt       flg;
684:   PetscReal      shift = info->shiftamount;
685:   PetscBool      allowzeropivot,zeropivotdetected;

687:   allowzeropivot = PetscNot(A->erroriffailure);

689:   /* generate work space needed by the factorization */
690:   PetscMalloc2(bs2*n,&rtmp,bs2,&mwork);
691:   PetscArrayzero(rtmp,bs2*n);

693:   for (i=0; i<n; i++) {
694:     /* zero rtmp */
695:     /* L part */
696:     nz    = bi[i+1] - bi[i];
697:     bjtmp = bj + bi[i];
698:     for  (j=0; j<nz; j++) {
699:       PetscArrayzero(rtmp+bs2*bjtmp[j],bs2);
700:     }

702:     /* U part */
703:     nz    = bdiag[i] - bdiag[i+1];
704:     bjtmp = bj + bdiag[i+1]+1;
705:     for  (j=0; j<nz; j++) {
706:       PetscArrayzero(rtmp+bs2*bjtmp[j],bs2);
707:     }

709:     /* load in initial (unfactored row) */
710:     nz    = ai[i+1] - ai[i];
711:     ajtmp = aj + ai[i];
712:     v     = aa + bs2*ai[i];
713:     for (j=0; j<nz; j++) {
714:       PetscArraycpy(rtmp+bs2*ajtmp[j],v+bs2*j,bs2);
715:     }

717:     /* elimination */
718:     bjtmp = bj + bi[i];
719:     nzL   = bi[i+1] - bi[i];
720:     for (k=0; k < nzL; k++) {
721:       row = bjtmp[k];
722:       pc  = rtmp + bs2*row;
723:       for (flg=0,j=0; j<bs2; j++) {
724:         if (pc[j]!=0.0) {
725:           flg = 1;
726:           break;
727:         }
728:       }
729:       if (flg) {
730:         pv = b->a + bs2*bdiag[row];
731:         /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */
732:         PetscKernel_A_gets_A_times_B_7(pc,pv,mwork);

734:         pj = b->j + bdiag[row+1]+1; /* beginning of U(row,:) */
735:         pv = b->a + bs2*(bdiag[row+1]+1);
736:         nz = bdiag[row] - bdiag[row+1] - 1; /* num of entries inU(row,:), excluding diag */
737:         for (j=0; j<nz; j++) {
738:           /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */
739:           /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */
740:           v    = rtmp + bs2*pj[j];
741:           PetscKernel_A_gets_A_minus_B_times_C_7(v,pc,pv);
742:           pv  += bs2;
743:         }
744:         PetscLogFlops(686.0*nz+637); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */
745:       }
746:     }

748:     /* finished row so stick it into b->a */
749:     /* L part */
750:     pv = b->a + bs2*bi[i];
751:     pj = b->j + bi[i];
752:     nz = bi[i+1] - bi[i];
753:     for (j=0; j<nz; j++) {
754:       PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2);
755:     }

757:     /* Mark diagonal and invert diagonal for simpler triangular solves */
758:     pv   = b->a + bs2*bdiag[i];
759:     pj   = b->j + bdiag[i];
760:     PetscArraycpy(pv,rtmp+bs2*pj[0],bs2);
761:     PetscKernel_A_gets_inverse_A_7(pv,shift,allowzeropivot,&zeropivotdetected);
762:     if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;

764:     /* U part */
765:     pv = b->a + bs2*(bdiag[i+1]+1);
766:     pj = b->j + bdiag[i+1]+1;
767:     nz = bdiag[i] - bdiag[i+1] - 1;
768:     for (j=0; j<nz; j++) {
769:       PetscArraycpy(pv+bs2*j,rtmp+bs2*pj[j],bs2);
770:     }
771:   }
772:   PetscFree2(rtmp,mwork);

774:   C->ops->solve          = MatSolve_SeqBAIJ_7_NaturalOrdering;
775:   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_7_NaturalOrdering;
776:   C->assembled           = PETSC_TRUE;

778:   PetscLogFlops(1.333333333333*7*7*7*n); /* from inverting diagonal blocks */
779:   return 0;
780: }