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