Actual source code: baijsolvtran5.c
petsc-3.13.6 2020-09-29
1: #include <../src/mat/impls/baij/seq/baij.h>
2: #include <petsc/private/kernels/blockinvert.h>
4: PetscErrorCode MatSolveTranspose_SeqBAIJ_5_inplace(Mat A,Vec bb,Vec xx)
5: {
6: Mat_SeqBAIJ *a =(Mat_SeqBAIJ*)A->data;
7: IS iscol=a->col,isrow=a->row;
8: PetscErrorCode ierr;
9: const PetscInt *r,*c,*rout,*cout;
10: const PetscInt *diag=a->diag,n=a->mbs,*vi,*ai=a->i,*aj=a->j;
11: PetscInt i,nz,idx,idt,ii,ic,ir,oidx;
12: const MatScalar *aa=a->a,*v;
13: PetscScalar s1,s2,s3,s4,s5,x1,x2,x3,x4,x5,*x,*t;
14: const PetscScalar *b;
17: VecGetArrayRead(bb,&b);
18: VecGetArray(xx,&x);
19: t = a->solve_work;
21: ISGetIndices(isrow,&rout); r = rout;
22: ISGetIndices(iscol,&cout); c = cout;
24: /* copy the b into temp work space according to permutation */
25: ii = 0;
26: for (i=0; i<n; i++) {
27: ic = 5*c[i];
28: t[ii] = b[ic];
29: t[ii+1] = b[ic+1];
30: t[ii+2] = b[ic+2];
31: t[ii+3] = b[ic+3];
32: t[ii+4] = b[ic+4];
33: ii += 5;
34: }
36: /* forward solve the U^T */
37: idx = 0;
38: for (i=0; i<n; i++) {
40: v = aa + 25*diag[i];
41: /* multiply by the inverse of the block diagonal */
42: x1 = t[idx]; x2 = t[1+idx]; x3 = t[2+idx]; x4 = t[3+idx]; x5 = t[4+idx];
43: s1 = v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5;
44: s2 = v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5;
45: s3 = v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
46: s4 = v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
47: s5 = v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
48: v += 25;
50: vi = aj + diag[i] + 1;
51: nz = ai[i+1] - diag[i] - 1;
52: while (nz--) {
53: oidx = 5*(*vi++);
54: t[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5;
55: t[oidx+1] -= v[5]*s1 + v[6]*s2 + v[7]*s3 + v[8]*s4 + v[9]*s5;
56: t[oidx+2] -= v[10]*s1 + v[11]*s2 + v[12]*s3 + v[13]*s4 + v[14]*s5;
57: t[oidx+3] -= v[15]*s1 + v[16]*s2 + v[17]*s3 + v[18]*s4 + v[19]*s5;
58: t[oidx+4] -= v[20]*s1 + v[21]*s2 + v[22]*s3 + v[23]*s4 + v[24]*s5;
59: v += 25;
60: }
61: t[idx] = s1;t[1+idx] = s2; t[2+idx] = s3;t[3+idx] = s4; t[4+idx] = s5;
62: idx += 5;
63: }
64: /* backward solve the L^T */
65: for (i=n-1; i>=0; i--) {
66: v = aa + 25*diag[i] - 25;
67: vi = aj + diag[i] - 1;
68: nz = diag[i] - ai[i];
69: idt = 5*i;
70: s1 = t[idt]; s2 = t[1+idt]; s3 = t[2+idt];s4 = t[3+idt]; s5 = t[4+idt];
71: while (nz--) {
72: idx = 5*(*vi--);
73: t[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5;
74: t[idx+1] -= v[5]*s1 + v[6]*s2 + v[7]*s3 + v[8]*s4 + v[9]*s5;
75: t[idx+2] -= v[10]*s1 + v[11]*s2 + v[12]*s3 + v[13]*s4 + v[14]*s5;
76: t[idx+3] -= v[15]*s1 + v[16]*s2 + v[17]*s3 + v[18]*s4 + v[19]*s5;
77: t[idx+4] -= v[20]*s1 + v[21]*s2 + v[22]*s3 + v[23]*s4 + v[24]*s5;
78: v -= 25;
79: }
80: }
82: /* copy t into x according to permutation */
83: ii = 0;
84: for (i=0; i<n; i++) {
85: ir = 5*r[i];
86: x[ir] = t[ii];
87: x[ir+1] = t[ii+1];
88: x[ir+2] = t[ii+2];
89: x[ir+3] = t[ii+3];
90: x[ir+4] = t[ii+4];
91: ii += 5;
92: }
94: ISRestoreIndices(isrow,&rout);
95: ISRestoreIndices(iscol,&cout);
96: VecRestoreArrayRead(bb,&b);
97: VecRestoreArray(xx,&x);
98: PetscLogFlops(2.0*25*(a->nz) - 5.0*A->cmap->n);
99: return(0);
100: }
102: PetscErrorCode MatSolveTranspose_SeqBAIJ_5(Mat A,Vec bb,Vec xx)
103: {
104: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data;
105: PetscErrorCode ierr;
106: IS iscol=a->col,isrow=a->row;
107: const PetscInt n =a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag;
108: const PetscInt *r,*c,*rout,*cout;
109: PetscInt nz,idx,idt,j,i,oidx,ii,ic,ir;
110: const PetscInt bs =A->rmap->bs,bs2=a->bs2;
111: const MatScalar *aa=a->a,*v;
112: PetscScalar s1,s2,s3,s4,s5,x1,x2,x3,x4,x5,*x,*t;
113: const PetscScalar *b;
116: VecGetArrayRead(bb,&b);
117: VecGetArray(xx,&x);
118: t = a->solve_work;
120: ISGetIndices(isrow,&rout); r = rout;
121: ISGetIndices(iscol,&cout); c = cout;
123: /* copy b into temp work space according to permutation */
124: for (i=0; i<n; i++) {
125: ii = bs*i; ic = bs*c[i];
126: t[ii] = b[ic]; t[ii+1] = b[ic+1]; t[ii+2] = b[ic+2]; t[ii+3] = b[ic+3];
127: t[ii+4] = b[ic+4];
128: }
130: /* forward solve the U^T */
131: idx = 0;
132: for (i=0; i<n; i++) {
133: v = aa + bs2*diag[i];
134: /* multiply by the inverse of the block diagonal */
135: x1 = t[idx]; x2 = t[1+idx]; x3 = t[2+idx]; x4 = t[3+idx]; x5 = t[4+idx];
136: s1 = v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5;
137: s2 = v[5]*x1 + v[6]*x2 + v[7]*x3 + v[8]*x4 + v[9]*x5;
138: s3 = v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
139: s4 = v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
140: s5 = v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
141: v -= bs2;
143: vi = aj + diag[i] - 1;
144: nz = diag[i] - diag[i+1] - 1;
145: for (j=0; j>-nz; j--) {
146: oidx = bs*vi[j];
147: t[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5;
148: t[oidx+1] -= v[5]*s1 + v[6]*s2 + v[7]*s3 + v[8]*s4 + v[9]*s5;
149: t[oidx+2] -= v[10]*s1 + v[11]*s2 + v[12]*s3 + v[13]*s4 + v[14]*s5;
150: t[oidx+3] -= v[15]*s1 + v[16]*s2 + v[17]*s3 + v[18]*s4 + v[19]*s5;
151: t[oidx+4] -= v[20]*s1 + v[21]*s2 + v[22]*s3 + v[23]*s4 + v[24]*s5;
152: v -= bs2;
153: }
154: t[idx] = s1;t[1+idx] = s2; t[2+idx] = s3; t[3+idx] = s4; t[4+idx] =s5;
155: idx += bs;
156: }
157: /* backward solve the L^T */
158: for (i=n-1; i>=0; i--) {
159: v = aa + bs2*ai[i];
160: vi = aj + ai[i];
161: nz = ai[i+1] - ai[i];
162: idt = bs*i;
163: s1 = t[idt]; s2 = t[1+idt]; s3 = t[2+idt]; s4 = t[3+idt]; s5 = t[4+idt];
164: for (j=0; j<nz; j++) {
165: idx = bs*vi[j];
166: t[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5;
167: t[idx+1] -= v[5]*s1 + v[6]*s2 + v[7]*s3 + v[8]*s4 + v[9]*s5;
168: t[idx+2] -= v[10]*s1 + v[11]*s2 + v[12]*s3 + v[13]*s4 + v[14]*s5;
169: t[idx+3] -= v[15]*s1 + v[16]*s2 + v[17]*s3 + v[18]*s4 + v[19]*s5;
170: t[idx+4] -= v[20]*s1 + v[21]*s2 + v[22]*s3 + v[23]*s4 + v[24]*s5;
171: v += bs2;
172: }
173: }
175: /* copy t into x according to permutation */
176: for (i=0; i<n; i++) {
177: ii = bs*i; ir = bs*r[i];
178: x[ir] = t[ii]; x[ir+1] = t[ii+1]; x[ir+2] = t[ii+2]; x[ir+3] = t[ii+3];
179: x[ir+4] = t[ii+4];
180: }
182: ISRestoreIndices(isrow,&rout);
183: ISRestoreIndices(iscol,&cout);
184: VecRestoreArrayRead(bb,&b);
185: VecRestoreArray(xx,&x);
186: PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);
187: return(0);
188: }