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