Actual source code: baijsolvtran4.c
petsc-3.14.6 2021-03-30
1: #include <../src/mat/impls/baij/seq/baij.h>
2: #include <petsc/private/kernels/blockinvert.h>
4: PetscErrorCode MatSolveTranspose_SeqBAIJ_4_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,x1,x2,x3,x4,*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 = 4*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: ii += 4;
33: }
35: /* forward solve the U^T */
36: idx = 0;
37: for (i=0; i<n; i++) {
39: v = aa + 16*diag[i];
40: /* multiply by the inverse of the block diagonal */
41: x1 = t[idx]; x2 = t[1+idx]; x3 = t[2+idx]; x4 = t[3+idx];
42: s1 = v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4;
43: s2 = v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4;
44: s3 = v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4;
45: s4 = v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
46: v += 16;
48: vi = aj + diag[i] + 1;
49: nz = ai[i+1] - diag[i] - 1;
50: while (nz--) {
51: oidx = 4*(*vi++);
52: t[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4;
53: t[oidx+1] -= v[4]*s1 + v[5]*s2 + v[6]*s3 + v[7]*s4;
54: t[oidx+2] -= v[8]*s1 + v[9]*s2 + v[10]*s3 + v[11]*s4;
55: t[oidx+3] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4;
56: v += 16;
57: }
58: t[idx] = s1;t[1+idx] = s2; t[2+idx] = s3;t[3+idx] = s4;
59: idx += 4;
60: }
61: /* backward solve the L^T */
62: for (i=n-1; i>=0; i--) {
63: v = aa + 16*diag[i] - 16;
64: vi = aj + diag[i] - 1;
65: nz = diag[i] - ai[i];
66: idt = 4*i;
67: s1 = t[idt]; s2 = t[1+idt]; s3 = t[2+idt];s4 = t[3+idt];
68: while (nz--) {
69: idx = 4*(*vi--);
70: t[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4;
71: t[idx+1] -= v[4]*s1 + v[5]*s2 + v[6]*s3 + v[7]*s4;
72: t[idx+2] -= v[8]*s1 + v[9]*s2 + v[10]*s3 + v[11]*s4;
73: t[idx+3] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4;
74: v -= 16;
75: }
76: }
78: /* copy t into x according to permutation */
79: ii = 0;
80: for (i=0; i<n; i++) {
81: ir = 4*r[i];
82: x[ir] = t[ii];
83: x[ir+1] = t[ii+1];
84: x[ir+2] = t[ii+2];
85: x[ir+3] = t[ii+3];
86: ii += 4;
87: }
89: ISRestoreIndices(isrow,&rout);
90: ISRestoreIndices(iscol,&cout);
91: VecRestoreArrayRead(bb,&b);
92: VecRestoreArray(xx,&x);
93: PetscLogFlops(2.0*16*(a->nz) - 4.0*A->cmap->n);
94: return(0);
95: }
97: PetscErrorCode MatSolveTranspose_SeqBAIJ_4(Mat A,Vec bb,Vec xx)
98: {
99: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data;
100: PetscErrorCode ierr;
101: IS iscol=a->col,isrow=a->row;
102: const PetscInt n =a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag;
103: const PetscInt *r,*c,*rout,*cout;
104: PetscInt nz,idx,idt,j,i,oidx,ii,ic,ir;
105: const PetscInt bs =A->rmap->bs,bs2=a->bs2;
106: const MatScalar *aa=a->a,*v;
107: PetscScalar s1,s2,s3,s4,x1,x2,x3,x4,*x,*t;
108: const PetscScalar *b;
111: VecGetArrayRead(bb,&b);
112: VecGetArray(xx,&x);
113: t = a->solve_work;
115: ISGetIndices(isrow,&rout); r = rout;
116: ISGetIndices(iscol,&cout); c = cout;
118: /* copy b into temp work space according to permutation */
119: for (i=0; i<n; i++) {
120: ii = bs*i; ic = bs*c[i];
121: t[ii] = b[ic]; t[ii+1] = b[ic+1]; t[ii+2] = b[ic+2]; t[ii+3] = b[ic+3];
122: }
124: /* forward solve the U^T */
125: idx = 0;
126: for (i=0; i<n; i++) {
127: v = aa + bs2*diag[i];
128: /* multiply by the inverse of the block diagonal */
129: x1 = t[idx]; x2 = t[1+idx]; x3 = t[2+idx]; x4 = t[3+idx];
130: s1 = v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4;
131: s2 = v[4]*x1 + v[5]*x2 + v[6]*x3 + v[7]*x4;
132: s3 = v[8]*x1 + v[9]*x2 + v[10]*x3 + v[11]*x4;
133: s4 = v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4;
134: v -= bs2;
136: vi = aj + diag[i] - 1;
137: nz = diag[i] - diag[i+1] - 1;
138: for (j=0; j>-nz; j--) {
139: oidx = bs*vi[j];
140: t[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4;
141: t[oidx+1] -= v[4]*s1 + v[5]*s2 + v[6]*s3 + v[7]*s4;
142: t[oidx+2] -= v[8]*s1 + v[9]*s2 + v[10]*s3 + v[11]*s4;
143: t[oidx+3] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4;
144: v -= bs2;
145: }
146: t[idx] = s1;t[1+idx] = s2; t[2+idx] = s3; t[3+idx] = s4;
147: idx += bs;
148: }
149: /* backward solve the L^T */
150: for (i=n-1; i>=0; i--) {
151: v = aa + bs2*ai[i];
152: vi = aj + ai[i];
153: nz = ai[i+1] - ai[i];
154: idt = bs*i;
155: s1 = t[idt]; s2 = t[1+idt]; s3 = t[2+idt]; s4 = t[3+idt];
156: for (j=0; j<nz; j++) {
157: idx = bs*vi[j];
158: t[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4;
159: t[idx+1] -= v[4]*s1 + v[5]*s2 + v[6]*s3 + v[7]*s4;
160: t[idx+2] -= v[8]*s1 + v[9]*s2 + v[10]*s3 + v[11]*s4;
161: t[idx+3] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4;
162: v += bs2;
163: }
164: }
166: /* copy t into x according to permutation */
167: for (i=0; i<n; i++) {
168: ii = bs*i; ir = bs*r[i];
169: x[ir] = t[ii]; x[ir+1] = t[ii+1]; x[ir+2] = t[ii+2]; x[ir+3] = t[ii+3];
170: }
172: ISRestoreIndices(isrow,&rout);
173: ISRestoreIndices(iscol,&cout);
174: VecRestoreArrayRead(bb,&b);
175: VecRestoreArray(xx,&x);
176: PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);
177: return(0);
178: }