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