xref: /petsc/src/mat/impls/baij/seq/baijfact9.c (revision a4efd8ea4bb0600a5aef92b2c23cbfa314827193)
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/inline/ilu.h"
8 
9 /* ------------------------------------------------------------*/
10 /*
11       Version for when blocks are 5 by 5
12 */
13 #undef __FUNCT__
14 #define __FUNCT__ "MatLUFactorNumeric_SeqBAIJ_5"
15 PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5(Mat A,MatFactorInfo *info,Mat *B)
16 {
17   Mat            C = *B;
18   Mat_SeqBAIJ    *a = (Mat_SeqBAIJ*)A->data,*b = (Mat_SeqBAIJ *)C->data;
19   IS             isrow = b->row,isicol = b->icol;
20   PetscErrorCode ierr;
21   PetscInt       *r,*ic,i,j,n = a->mbs,*bi = b->i,*bj = b->j;
22   PetscInt       *ajtmpold,*ajtmp,nz,row;
23   PetscInt       *diag_offset = b->diag,idx,*ai=a->i,*aj=a->j,*pj;
24   MatScalar      *pv,*v,*rtmp,*pc,*w,*x;
25   MatScalar      p1,p2,p3,p4,m1,m2,m3,m4,m5,m6,m7,m8,m9,x1,x2,x3,x4;
26   MatScalar      p5,p6,p7,p8,p9,x5,x6,x7,x8,x9,x10,x11,x12,x13,x14,x15,x16;
27   MatScalar      x17,x18,x19,x20,x21,x22,x23,x24,x25,p10,p11,p12,p13,p14;
28   MatScalar      p15,p16,p17,p18,p19,p20,p21,p22,p23,p24,p25,m10,m11,m12;
29   MatScalar      m13,m14,m15,m16,m17,m18,m19,m20,m21,m22,m23,m24,m25;
30   MatScalar      *ba = b->a,*aa = a->a;
31   PetscReal      shift = info->shiftinblocks;
32 
33   PetscFunctionBegin;
34   ierr = ISGetIndices(isrow,&r);CHKERRQ(ierr);
35   ierr = ISGetIndices(isicol,&ic);CHKERRQ(ierr);
36   ierr = PetscMalloc(25*(n+1)*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
37 
38   for (i=0; i<n; i++) {
39     nz    = bi[i+1] - bi[i];
40     ajtmp = bj + bi[i];
41     for  (j=0; j<nz; j++) {
42       x = rtmp+25*ajtmp[j];
43       x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0;
44       x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0;
45       x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = 0.0;
46     }
47     /* load in initial (unfactored row) */
48     idx      = r[i];
49     nz       = ai[idx+1] - ai[idx];
50     ajtmpold = aj + ai[idx];
51     v        = aa + 25*ai[idx];
52     for (j=0; j<nz; j++) {
53       x    = rtmp+25*ic[ajtmpold[j]];
54       x[0] = v[0]; x[1] = v[1]; x[2] = v[2]; x[3] = v[3];
55       x[4] = v[4]; x[5] = v[5]; x[6] = v[6]; x[7] = v[7]; x[8] = v[8];
56       x[9] = v[9]; x[10] = v[10]; x[11] = v[11]; x[12] = v[12]; x[13] = v[13];
57       x[14] = v[14]; x[15] = v[15]; x[16] = v[16]; x[17] = v[17];
58       x[18] = v[18]; x[19] = v[19]; x[20] = v[20]; x[21] = v[21];
59       x[22] = v[22]; x[23] = v[23]; x[24] = v[24];
60       v    += 25;
61     }
62     row = *ajtmp++;
63     while (row < i) {
64       pc = rtmp + 25*row;
65       p1 = pc[0]; p2 = pc[1]; p3 = pc[2]; p4 = pc[3];
66       p5 = pc[4]; p6 = pc[5]; p7 = pc[6]; p8 = pc[7]; p9 = pc[8];
67       p10 = pc[9]; p11 = pc[10]; p12 = pc[11]; p13 = pc[12]; p14 = pc[13];
68       p15 = pc[14]; p16 = pc[15]; p17 = pc[16]; p18 = pc[17]; p19 = pc[18];
69       p20 = pc[19]; p21 = pc[20]; p22 = pc[21]; p23 = pc[22]; p24 = pc[23];
70       p25 = pc[24];
71       if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 ||
72           p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 ||
73           p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0
74           || p16 != 0.0 || p17 != 0.0 || p18 != 0.0 || p19 != 0.0 ||
75           p20 != 0.0 || p21 != 0.0 || p22 != 0.0 || p23 != 0.0 ||
76           p24 != 0.0 || p25 != 0.0) {
77         pv = ba + 25*diag_offset[row];
78         pj = bj + diag_offset[row] + 1;
79         x1 = pv[0]; x2 = pv[1]; x3 = pv[2]; x4 = pv[3];
80         x5 = pv[4]; x6 = pv[5]; x7 = pv[6]; x8 = pv[7]; x9 = pv[8];
81         x10 = pv[9]; x11 = pv[10]; x12 = pv[11]; x13 = pv[12]; x14 = pv[13];
82         x15 = pv[14]; x16 = pv[15]; x17 = pv[16]; x18 = pv[17];
83         x19 = pv[18]; x20 = pv[19]; x21 = pv[20]; x22 = pv[21];
84         x23 = pv[22]; x24 = pv[23]; x25 = pv[24];
85         pc[0] = m1 = p1*x1 + p6*x2  + p11*x3 + p16*x4 + p21*x5;
86         pc[1] = m2 = p2*x1 + p7*x2  + p12*x3 + p17*x4 + p22*x5;
87         pc[2] = m3 = p3*x1 + p8*x2  + p13*x3 + p18*x4 + p23*x5;
88         pc[3] = m4 = p4*x1 + p9*x2  + p14*x3 + p19*x4 + p24*x5;
89         pc[4] = m5 = p5*x1 + p10*x2 + p15*x3 + p20*x4 + p25*x5;
90 
91         pc[5] = m6 = p1*x6 + p6*x7  + p11*x8 + p16*x9 + p21*x10;
92         pc[6] = m7 = p2*x6 + p7*x7  + p12*x8 + p17*x9 + p22*x10;
93         pc[7] = m8 = p3*x6 + p8*x7  + p13*x8 + p18*x9 + p23*x10;
94         pc[8] = m9 = p4*x6 + p9*x7  + p14*x8 + p19*x9 + p24*x10;
95         pc[9] = m10 = p5*x6 + p10*x7 + p15*x8 + p20*x9 + p25*x10;
96 
97         pc[10] = m11 = p1*x11 + p6*x12  + p11*x13 + p16*x14 + p21*x15;
98         pc[11] = m12 = p2*x11 + p7*x12  + p12*x13 + p17*x14 + p22*x15;
99         pc[12] = m13 = p3*x11 + p8*x12  + p13*x13 + p18*x14 + p23*x15;
100         pc[13] = m14 = p4*x11 + p9*x12  + p14*x13 + p19*x14 + p24*x15;
101         pc[14] = m15 = p5*x11 + p10*x12 + p15*x13 + p20*x14 + p25*x15;
102 
103         pc[15] = m16 = p1*x16 + p6*x17  + p11*x18 + p16*x19 + p21*x20;
104         pc[16] = m17 = p2*x16 + p7*x17  + p12*x18 + p17*x19 + p22*x20;
105         pc[17] = m18 = p3*x16 + p8*x17  + p13*x18 + p18*x19 + p23*x20;
106         pc[18] = m19 = p4*x16 + p9*x17  + p14*x18 + p19*x19 + p24*x20;
107         pc[19] = m20 = p5*x16 + p10*x17 + p15*x18 + p20*x19 + p25*x20;
108 
109         pc[20] = m21 = p1*x21 + p6*x22  + p11*x23 + p16*x24 + p21*x25;
110         pc[21] = m22 = p2*x21 + p7*x22  + p12*x23 + p17*x24 + p22*x25;
111         pc[22] = m23 = p3*x21 + p8*x22  + p13*x23 + p18*x24 + p23*x25;
112         pc[23] = m24 = p4*x21 + p9*x22  + p14*x23 + p19*x24 + p24*x25;
113         pc[24] = m25 = p5*x21 + p10*x22 + p15*x23 + p20*x24 + p25*x25;
114 
115         nz = bi[row+1] - diag_offset[row] - 1;
116         pv += 25;
117         for (j=0; j<nz; j++) {
118           x1   = pv[0];  x2 = pv[1];   x3  = pv[2];  x4  = pv[3];
119           x5   = pv[4];  x6 = pv[5];   x7  = pv[6];  x8  = pv[7]; x9 = pv[8];
120           x10  = pv[9];  x11 = pv[10]; x12 = pv[11]; x13 = pv[12];
121           x14  = pv[13]; x15 = pv[14]; x16 = pv[15]; x17 = pv[16];
122           x18  = pv[17]; x19 = pv[18]; x20 = pv[19]; x21 = pv[20];
123           x22  = pv[21]; x23 = pv[22]; x24 = pv[23]; x25 = pv[24];
124           x    = rtmp + 25*pj[j];
125           x[0] -= m1*x1 + m6*x2  + m11*x3 + m16*x4 + m21*x5;
126           x[1] -= m2*x1 + m7*x2  + m12*x3 + m17*x4 + m22*x5;
127           x[2] -= m3*x1 + m8*x2  + m13*x3 + m18*x4 + m23*x5;
128           x[3] -= m4*x1 + m9*x2  + m14*x3 + m19*x4 + m24*x5;
129           x[4] -= m5*x1 + m10*x2 + m15*x3 + m20*x4 + m25*x5;
130 
131           x[5] -= m1*x6 + m6*x7  + m11*x8 + m16*x9 + m21*x10;
132           x[6] -= m2*x6 + m7*x7  + m12*x8 + m17*x9 + m22*x10;
133           x[7] -= m3*x6 + m8*x7  + m13*x8 + m18*x9 + m23*x10;
134           x[8] -= m4*x6 + m9*x7  + m14*x8 + m19*x9 + m24*x10;
135           x[9] -= m5*x6 + m10*x7 + m15*x8 + m20*x9 + m25*x10;
136 
137           x[10] -= m1*x11 + m6*x12  + m11*x13 + m16*x14 + m21*x15;
138           x[11] -= m2*x11 + m7*x12  + m12*x13 + m17*x14 + m22*x15;
139           x[12] -= m3*x11 + m8*x12  + m13*x13 + m18*x14 + m23*x15;
140           x[13] -= m4*x11 + m9*x12  + m14*x13 + m19*x14 + m24*x15;
141           x[14] -= m5*x11 + m10*x12 + m15*x13 + m20*x14 + m25*x15;
142 
143           x[15] -= m1*x16 + m6*x17  + m11*x18 + m16*x19 + m21*x20;
144           x[16] -= m2*x16 + m7*x17  + m12*x18 + m17*x19 + m22*x20;
145           x[17] -= m3*x16 + m8*x17  + m13*x18 + m18*x19 + m23*x20;
146           x[18] -= m4*x16 + m9*x17  + m14*x18 + m19*x19 + m24*x20;
147           x[19] -= m5*x16 + m10*x17 + m15*x18 + m20*x19 + m25*x20;
148 
149           x[20] -= m1*x21 + m6*x22  + m11*x23 + m16*x24 + m21*x25;
150           x[21] -= m2*x21 + m7*x22  + m12*x23 + m17*x24 + m22*x25;
151           x[22] -= m3*x21 + m8*x22  + m13*x23 + m18*x24 + m23*x25;
152           x[23] -= m4*x21 + m9*x22  + m14*x23 + m19*x24 + m24*x25;
153           x[24] -= m5*x21 + m10*x22 + m15*x23 + m20*x24 + m25*x25;
154 
155           pv   += 25;
156         }
157         ierr = PetscLogFlops(250*nz+225);CHKERRQ(ierr);
158       }
159       row = *ajtmp++;
160     }
161     /* finished row so stick it into b->a */
162     pv = ba + 25*bi[i];
163     pj = bj + bi[i];
164     nz = bi[i+1] - bi[i];
165     for (j=0; j<nz; j++) {
166       x     = rtmp+25*pj[j];
167       pv[0] = x[0]; pv[1] = x[1]; pv[2] = x[2]; pv[3] = x[3];
168       pv[4] = x[4]; pv[5] = x[5]; pv[6] = x[6]; pv[7] = x[7]; pv[8] = x[8];
169       pv[9] = x[9]; pv[10] = x[10]; pv[11] = x[11]; pv[12] = x[12];
170       pv[13] = x[13]; pv[14] = x[14]; pv[15] = x[15]; pv[16] = x[16];
171       pv[17] = x[17]; pv[18] = x[18]; pv[19] = x[19]; pv[20] = x[20];
172       pv[21] = x[21]; pv[22] = x[22]; pv[23] = x[23]; pv[24] = x[24];
173       pv   += 25;
174     }
175     /* invert diagonal block */
176     w = ba + 25*diag_offset[i];
177     ierr = Kernel_A_gets_inverse_A_5(w,shift);CHKERRQ(ierr);
178   }
179 
180   ierr = PetscFree(rtmp);CHKERRQ(ierr);
181   ierr = ISRestoreIndices(isicol,&ic);CHKERRQ(ierr);
182   ierr = ISRestoreIndices(isrow,&r);CHKERRQ(ierr);
183   C->ops->solve          = MatSolve_SeqBAIJ_5;
184   C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5;
185   C->assembled = PETSC_TRUE;
186   ierr = PetscLogFlops(1.3333*125*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
187   PetscFunctionReturn(0);
188 }
189