xref: /petsc/src/mat/impls/baij/seq/baijsolvnat6.c (revision 503c0ea9b45bcfbcebbb1ea5341243bbc69f0bea)
1 #include <../src/mat/impls/baij/seq/baij.h>
2 #include <petsc/private/kernels/blockinvert.h>
3 
4 PetscErrorCode MatSolve_SeqBAIJ_6_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)
5 {
6   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
7   PetscInt          i,nz,idx,idt,jdx;
8   const PetscInt    *diag = a->diag,*vi,n=a->mbs,*ai=a->i,*aj=a->j;
9   const MatScalar   *aa   =a->a,*v;
10   PetscScalar       *x,s1,s2,s3,s4,s5,s6,x1,x2,x3,x4,x5,x6;
11   const PetscScalar *b;
12 
13   PetscFunctionBegin;
14   PetscCall(VecGetArrayRead(bb,&b));
15   PetscCall(VecGetArray(xx,&x));
16   /* forward solve the lower triangular */
17   idx  = 0;
18   x[0] = b[idx];   x[1] = b[1+idx]; x[2] = b[2+idx];
19   x[3] = b[3+idx]; x[4] = b[4+idx]; x[5] = b[5+idx];
20   for (i=1; i<n; i++) {
21     v   =  aa + 36*ai[i];
22     vi  =  aj + ai[i];
23     nz  =  diag[i] - ai[i];
24     idx =  6*i;
25     s1  =  b[idx];   s2 = b[1+idx]; s3 = b[2+idx];
26     s4  =  b[3+idx]; s5 = b[4+idx]; s6 = b[5+idx];
27     while (nz--) {
28       jdx = 6*(*vi++);
29       x1  = x[jdx];   x2 = x[1+jdx]; x3 = x[2+jdx];
30       x4  = x[3+jdx]; x5 = x[4+jdx]; x6 = x[5+jdx];
31       s1 -= v[0]*x1 + v[6]*x2  + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6;
32       s2 -= v[1]*x1 + v[7]*x2  + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6;
33       s3 -= v[2]*x1 + v[8]*x2  + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6;
34       s4 -= v[3]*x1 + v[9]*x2  + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6;
35       s5 -= v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6;
36       s6 -= v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6;
37       v  += 36;
38     }
39     x[idx]   = s1;
40     x[1+idx] = s2;
41     x[2+idx] = s3;
42     x[3+idx] = s4;
43     x[4+idx] = s5;
44     x[5+idx] = s6;
45   }
46   /* backward solve the upper triangular */
47   for (i=n-1; i>=0; i--) {
48     v   = aa + 36*diag[i] + 36;
49     vi  = aj + diag[i] + 1;
50     nz  = ai[i+1] - diag[i] - 1;
51     idt = 6*i;
52     s1  = x[idt];   s2 = x[1+idt];
53     s3  = x[2+idt]; s4 = x[3+idt];
54     s5  = x[4+idt]; s6 = x[5+idt];
55     while (nz--) {
56       idx = 6*(*vi++);
57       x1  = x[idx];   x2 = x[1+idx]; x3 = x[2+idx];
58       x4  = x[3+idx]; x5 = x[4+idx]; x6 = x[5+idx];
59       s1 -= v[0]*x1 + v[6]*x2  + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6;
60       s2 -= v[1]*x1 + v[7]*x2  + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6;
61       s3 -= v[2]*x1 + v[8]*x2  + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6;
62       s4 -= v[3]*x1 + v[9]*x2  + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6;
63       s5 -= v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6;
64       s6 -= v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6;
65       v  += 36;
66     }
67     v        = aa + 36*diag[i];
68     x[idt]   = v[0]*s1 + v[6]*s2  + v[12]*s3 + v[18]*s4 + v[24]*s5 + v[30]*s6;
69     x[1+idt] = v[1]*s1 + v[7]*s2  + v[13]*s3 + v[19]*s4 + v[25]*s5 + v[31]*s6;
70     x[2+idt] = v[2]*s1 + v[8]*s2  + v[14]*s3 + v[20]*s4 + v[26]*s5 + v[32]*s6;
71     x[3+idt] = v[3]*s1 + v[9]*s2  + v[15]*s3 + v[21]*s4 + v[27]*s5 + v[33]*s6;
72     x[4+idt] = v[4]*s1 + v[10]*s2 + v[16]*s3 + v[22]*s4 + v[28]*s5 + v[34]*s6;
73     x[5+idt] = v[5]*s1 + v[11]*s2 + v[17]*s3 + v[23]*s4 + v[29]*s5 + v[35]*s6;
74   }
75 
76   PetscCall(VecRestoreArrayRead(bb,&b));
77   PetscCall(VecRestoreArray(xx,&x));
78   PetscCall(PetscLogFlops(2.0*36*(a->nz) - 6.0*A->cmap->n));
79   PetscFunctionReturn(0);
80 }
81 
82 PetscErrorCode MatSolve_SeqBAIJ_6_NaturalOrdering(Mat A,Vec bb,Vec xx)
83 {
84   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
85   const PetscInt    n  =a->mbs,*vi,*ai=a->i,*aj=a->j,*adiag=a->diag;
86   PetscInt          i,k,nz,idx,jdx,idt;
87   const PetscInt    bs = A->rmap->bs,bs2 = a->bs2;
88   const MatScalar   *aa=a->a,*v;
89   PetscScalar       *x;
90   const PetscScalar *b;
91   PetscScalar       s1,s2,s3,s4,s5,s6,x1,x2,x3,x4,x5,x6;
92 
93   PetscFunctionBegin;
94   PetscCall(VecGetArrayRead(bb,&b));
95   PetscCall(VecGetArray(xx,&x));
96   /* forward solve the lower triangular */
97   idx  = 0;
98   x[0] = b[idx]; x[1] = b[1+idx];x[2] = b[2+idx];x[3] = b[3+idx];
99   x[4] = b[4+idx];x[5] = b[5+idx];
100   for (i=1; i<n; i++) {
101     v   = aa + bs2*ai[i];
102     vi  = aj + ai[i];
103     nz  = ai[i+1] - ai[i];
104     idx = bs*i;
105     s1  = b[idx];s2 = b[1+idx];s3 = b[2+idx];s4 = b[3+idx];
106     s5  = b[4+idx];s6 = b[5+idx];
107     for (k=0; k<nz; k++) {
108       jdx = bs*vi[k];
109       x1  = x[jdx];x2 = x[1+jdx]; x3 =x[2+jdx];x4 =x[3+jdx];
110       x5  = x[4+jdx]; x6 = x[5+jdx];
111       s1 -= v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4  + v[24]*x5 + v[30]*x6;
112       s2 -= v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4  + v[25]*x5 + v[31]*x6;
113       s3 -= v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4  + v[26]*x5 + v[32]*x6;
114       s4 -= v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4  + v[27]*x5 + v[33]*x6;
115       s5 -= v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4  + v[28]*x5 + v[34]*x6;
116       s6 -= v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4  + v[29]*x5 + v[35]*x6;
117       v  +=  bs2;
118     }
119 
120     x[idx]   = s1;
121     x[1+idx] = s2;
122     x[2+idx] = s3;
123     x[3+idx] = s4;
124     x[4+idx] = s5;
125     x[5+idx] = s6;
126   }
127 
128   /* backward solve the upper triangular */
129   for (i=n-1; i>=0; i--) {
130     v   = aa + bs2*(adiag[i+1]+1);
131     vi  = aj + adiag[i+1]+1;
132     nz  = adiag[i] - adiag[i+1]-1;
133     idt = bs*i;
134     s1  = x[idt];  s2 = x[1+idt];s3 = x[2+idt];s4 = x[3+idt];
135     s5  = x[4+idt];s6 = x[5+idt];
136     for (k=0; k<nz; k++) {
137       idx = bs*vi[k];
138       x1  = x[idx];   x2 = x[1+idx]; x3 = x[2+idx];x4 = x[3+idx];
139       x5  = x[4+idx];x6 = x[5+idx];
140       s1 -= v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4  + v[24]*x5 + v[30]*x6;
141       s2 -= v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4  + v[25]*x5 + v[31]*x6;
142       s3 -= v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4  + v[26]*x5 + v[32]*x6;
143       s4 -= v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4  + v[27]*x5 + v[33]*x6;
144       s5 -= v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4  + v[28]*x5 + v[34]*x6;
145       s6 -= v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4  + v[29]*x5 + v[35]*x6;
146       v  +=  bs2;
147     }
148     /* x = inv_diagonal*x */
149     x[idt]   = v[0]*s1 + v[6]*s2 + v[12]*s3 + v[18]*s4 + v[24]*s5 + v[30]*s6;
150     x[1+idt] = v[1]*s1 + v[7]*s2 + v[13]*s3 + v[19]*s4 + v[25]*s5 + v[31]*s6;
151     x[2+idt] = v[2]*s1 + v[8]*s2 + v[14]*s3 + v[20]*s4 + v[26]*s5 + v[32]*s6;
152     x[3+idt] = v[3]*s1 + v[9]*s2 + v[15]*s3 + v[21]*s4 + v[27]*s5 + v[33]*s6;
153     x[4+idt] = v[4]*s1 + v[10]*s2 + v[16]*s3 + v[22]*s4 + v[28]*s5 + v[34]*s6;
154     x[5+idt] = v[5]*s1 + v[11]*s2 + v[17]*s3 + v[23]*s4 + v[29]*s5 + v[35]*s6;
155   }
156 
157   PetscCall(VecRestoreArrayRead(bb,&b));
158   PetscCall(VecRestoreArray(xx,&x));
159   PetscCall(PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n));
160   PetscFunctionReturn(0);
161 }
162