xref: /petsc/src/mat/impls/baij/seq/baijsolvtrannat5.c (revision 503c0ea9b45bcfbcebbb1ea5341243bbc69f0bea)
1 #include <../src/mat/impls/baij/seq/baij.h>
2 
3 PetscErrorCode MatSolveTranspose_SeqBAIJ_5_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)
4 {
5   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ*)A->data;
6   const PetscInt  *diag=a->diag,n=a->mbs,*vi,*ai=a->i,*aj=a->j;
7   PetscInt        i,nz,idx,idt,oidx;
8   const MatScalar *aa=a->a,*v;
9   PetscScalar     s1,s2,s3,s4,s5,x1,x2,x3,x4,x5,*x;
10 
11   PetscFunctionBegin;
12   PetscCall(VecCopy(bb,xx));
13   PetscCall(VecGetArray(xx,&x));
14 
15   /* forward solve the U^T */
16   idx = 0;
17   for (i=0; i<n; i++) {
18 
19     v = aa + 25*diag[i];
20     /* multiply by the inverse of the block diagonal */
21     x1 = x[idx];   x2 = x[1+idx]; x3    = x[2+idx]; x4 = x[3+idx]; x5 = x[4+idx];
22     s1 = v[0]*x1  +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
23     s2 = v[5]*x1  +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
24     s3 = v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
25     s4 = v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
26     s5 = v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
27     v += 25;
28 
29     vi = aj + diag[i] + 1;
30     nz = ai[i+1] - diag[i] - 1;
31     while (nz--) {
32       oidx       = 5*(*vi++);
33       x[oidx]   -= v[0]*s1  +  v[1]*s2 +  v[2]*s3 +  v[3]*s4 +  v[4]*s5;
34       x[oidx+1] -= v[5]*s1  +  v[6]*s2 +  v[7]*s3 +  v[8]*s4 +  v[9]*s5;
35       x[oidx+2] -= v[10]*s1 + v[11]*s2 + v[12]*s3 + v[13]*s4 + v[14]*s5;
36       x[oidx+3] -= v[15]*s1 + v[16]*s2 + v[17]*s3 + v[18]*s4 + v[19]*s5;
37       x[oidx+4] -= v[20]*s1 + v[21]*s2 + v[22]*s3 + v[23]*s4 + v[24]*s5;
38       v         += 25;
39     }
40     x[idx] = s1;x[1+idx] = s2; x[2+idx] = s3;x[3+idx] = s4; x[4+idx] = s5;
41     idx   += 5;
42   }
43   /* backward solve the L^T */
44   for (i=n-1; i>=0; i--) {
45     v   = aa + 25*diag[i] - 25;
46     vi  = aj + diag[i] - 1;
47     nz  = diag[i] - ai[i];
48     idt = 5*i;
49     s1  = x[idt];  s2 = x[1+idt]; s3 = x[2+idt];s4 = x[3+idt]; s5 = x[4+idt];
50     while (nz--) {
51       idx       = 5*(*vi--);
52       x[idx]   -=  v[0]*s1 +  v[1]*s2 +  v[2]*s3 +  v[3]*s4 +  v[4]*s5;
53       x[idx+1] -=  v[5]*s1 +  v[6]*s2 +  v[7]*s3 +  v[8]*s4 +  v[9]*s5;
54       x[idx+2] -= v[10]*s1 + v[11]*s2 + v[12]*s3 + v[13]*s4 + v[14]*s5;
55       x[idx+3] -= v[15]*s1 + v[16]*s2 + v[17]*s3 + v[18]*s4 + v[19]*s5;
56       x[idx+4] -= v[20]*s1 + v[21]*s2 + v[22]*s3 + v[23]*s4 + v[24]*s5;
57       v        -= 25;
58     }
59   }
60   PetscCall(VecRestoreArray(xx,&x));
61   PetscCall(PetscLogFlops(2.0*25*(a->nz) - 5.0*A->cmap->n));
62   PetscFunctionReturn(0);
63 }
64 
65 PetscErrorCode MatSolveTranspose_SeqBAIJ_5_NaturalOrdering(Mat A,Vec bb,Vec xx)
66 {
67   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ*)A->data;
68   const PetscInt  n=a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag;
69   PetscInt        nz,idx,idt,j,i,oidx;
70   const PetscInt  bs =A->rmap->bs,bs2=a->bs2;
71   const MatScalar *aa=a->a,*v;
72   PetscScalar     s1,s2,s3,s4,s5,x1,x2,x3,x4,x5,*x;
73 
74   PetscFunctionBegin;
75   PetscCall(VecCopy(bb,xx));
76   PetscCall(VecGetArray(xx,&x));
77 
78   /* forward solve the U^T */
79   idx = 0;
80   for (i=0; i<n; i++) {
81     v = aa + bs2*diag[i];
82     /* multiply by the inverse of the block diagonal */
83     x1 = x[idx];   x2 = x[1+idx];  x3 = x[2+idx];  x4 = x[3+idx];
84     x5 = x[4+idx];
85     s1 =  v[0]*x1   +  v[1]*x2   + v[2]*x3   + v[3]*x4   + v[4]*x5;
86     s2 =  v[5]*x1   +  v[6]*x2   + v[7]*x3   + v[8]*x4   + v[9]*x5;
87     s3 =  v[10]*x1  +  v[11]*x2  + v[12]*x3  + v[13]*x4  + v[14]*x5;
88     s4 =  v[15]*x1  +  v[16]*x2  + v[17]*x3  + v[18]*x4  + v[19]*x5;
89     s5 =  v[20]*x1  +  v[21]*x2  + v[22]*x3  + v[23]*x4   + v[24]*x5;
90     v -= bs2;
91 
92     vi = aj + diag[i] - 1;
93     nz = diag[i] - diag[i+1] - 1;
94     for (j=0; j>-nz; j--) {
95       oidx       = bs*vi[j];
96       x[oidx]   -=  v[0]*s1   +  v[1]*s2   + v[2]*s3   + v[3]*s4   + v[4]*s5;
97       x[oidx+1] -=  v[5]*s1   +  v[6]*s2   + v[7]*s3   + v[8]*s4   + v[9]*s5;
98       x[oidx+2] -=  v[10]*s1  +  v[11]*s2  + v[12]*s3  + v[13]*s4  + v[14]*s5;
99       x[oidx+3] -=  v[15]*s1  +  v[16]*s2  + v[17]*s3  + v[18]*s4  + v[19]*s5;
100       x[oidx+4] -=  v[20]*s1  +  v[21]*s2  + v[22]*s3  + v[23]*s4   + v[24]*s5;
101       v         -= bs2;
102     }
103     x[idx] = s1;x[1+idx] = s2;  x[2+idx] = s3;  x[3+idx] = s4; x[4+idx] = s5;
104     idx   += bs;
105   }
106   /* backward solve the L^T */
107   for (i=n-1; i>=0; i--) {
108     v   = aa + bs2*ai[i];
109     vi  = aj + ai[i];
110     nz  = ai[i+1] - ai[i];
111     idt = bs*i;
112     s1  = x[idt];  s2 = x[1+idt];  s3 = x[2+idt];  s4 = x[3+idt];  s5 = x[4+idt];
113     for (j=0; j<nz; j++) {
114       idx       = bs*vi[j];
115       x[idx]   -=  v[0]*s1   +  v[1]*s2   + v[2]*s3   + v[3]*s4   + v[4]*s5;
116       x[idx+1] -=  v[5]*s1   +  v[6]*s2   + v[7]*s3   + v[8]*s4   + v[9]*s5;
117       x[idx+2] -=  v[10]*s1  +  v[11]*s2  + v[12]*s3  + v[13]*s4  + v[14]*s5;
118       x[idx+3] -=  v[15]*s1  +  v[16]*s2  + v[17]*s3  + v[18]*s4  + v[19]*s5;
119       x[idx+4] -=  v[20]*s1  +  v[21]*s2  + v[22]*s3  + v[23]*s4   + v[24]*s5;
120       v        += bs2;
121     }
122   }
123   PetscCall(VecRestoreArray(xx,&x));
124   PetscCall(PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n));
125   PetscFunctionReturn(0);
126 }
127