Actual source code: baijsolvnat6.c
1: #include <../src/mat/impls/baij/seq/baij.h>
2: #include <petsc/private/kernels/blockinvert.h>
4: PetscErrorCode MatSolve_SeqBAIJ_6_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)
5: {
6: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
7: PetscInt i,nz,idx,idt,jdx;
8: PetscErrorCode ierr;
9: const PetscInt *diag = a->diag,*vi,n=a->mbs,*ai=a->i,*aj=a->j;
10: const MatScalar *aa =a->a,*v;
11: PetscScalar *x,s1,s2,s3,s4,s5,s6,x1,x2,x3,x4,x5,x6;
12: const PetscScalar *b;
15: VecGetArrayRead(bb,&b);
16: VecGetArray(xx,&x);
17: /* forward solve the lower triangular */
18: idx = 0;
19: x[0] = b[idx]; x[1] = b[1+idx]; x[2] = b[2+idx];
20: x[3] = b[3+idx]; x[4] = b[4+idx]; x[5] = b[5+idx];
21: for (i=1; i<n; i++) {
22: v = aa + 36*ai[i];
23: vi = aj + ai[i];
24: nz = diag[i] - ai[i];
25: idx = 6*i;
26: s1 = b[idx]; s2 = b[1+idx]; s3 = b[2+idx];
27: s4 = b[3+idx]; s5 = b[4+idx]; s6 = b[5+idx];
28: while (nz--) {
29: jdx = 6*(*vi++);
30: x1 = x[jdx]; x2 = x[1+jdx]; x3 = x[2+jdx];
31: x4 = x[3+jdx]; x5 = x[4+jdx]; x6 = x[5+jdx];
32: s1 -= v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6;
33: s2 -= v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6;
34: s3 -= v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6;
35: s4 -= v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6;
36: s5 -= v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6;
37: s6 -= v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6;
38: v += 36;
39: }
40: x[idx] = s1;
41: x[1+idx] = s2;
42: x[2+idx] = s3;
43: x[3+idx] = s4;
44: x[4+idx] = s5;
45: x[5+idx] = s6;
46: }
47: /* backward solve the upper triangular */
48: for (i=n-1; i>=0; i--) {
49: v = aa + 36*diag[i] + 36;
50: vi = aj + diag[i] + 1;
51: nz = ai[i+1] - diag[i] - 1;
52: idt = 6*i;
53: s1 = x[idt]; s2 = x[1+idt];
54: s3 = x[2+idt]; s4 = x[3+idt];
55: s5 = x[4+idt]; s6 = x[5+idt];
56: while (nz--) {
57: idx = 6*(*vi++);
58: x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx];
59: x4 = x[3+idx]; x5 = x[4+idx]; x6 = x[5+idx];
60: s1 -= v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6;
61: s2 -= v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6;
62: s3 -= v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6;
63: s4 -= v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6;
64: s5 -= v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6;
65: s6 -= v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6;
66: v += 36;
67: }
68: v = aa + 36*diag[i];
69: x[idt] = v[0]*s1 + v[6]*s2 + v[12]*s3 + v[18]*s4 + v[24]*s5 + v[30]*s6;
70: x[1+idt] = v[1]*s1 + v[7]*s2 + v[13]*s3 + v[19]*s4 + v[25]*s5 + v[31]*s6;
71: x[2+idt] = v[2]*s1 + v[8]*s2 + v[14]*s3 + v[20]*s4 + v[26]*s5 + v[32]*s6;
72: x[3+idt] = v[3]*s1 + v[9]*s2 + v[15]*s3 + v[21]*s4 + v[27]*s5 + v[33]*s6;
73: x[4+idt] = v[4]*s1 + v[10]*s2 + v[16]*s3 + v[22]*s4 + v[28]*s5 + v[34]*s6;
74: x[5+idt] = v[5]*s1 + v[11]*s2 + v[17]*s3 + v[23]*s4 + v[29]*s5 + v[35]*s6;
75: }
77: VecRestoreArrayRead(bb,&b);
78: VecRestoreArray(xx,&x);
79: PetscLogFlops(2.0*36*(a->nz) - 6.0*A->cmap->n);
80: return(0);
81: }
83: PetscErrorCode MatSolve_SeqBAIJ_6_NaturalOrdering(Mat A,Vec bb,Vec xx)
84: {
85: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
86: const PetscInt n =a->mbs,*vi,*ai=a->i,*aj=a->j,*adiag=a->diag;
87: PetscErrorCode ierr;
88: PetscInt i,k,nz,idx,jdx,idt;
89: const PetscInt bs = A->rmap->bs,bs2 = a->bs2;
90: const MatScalar *aa=a->a,*v;
91: PetscScalar *x;
92: const PetscScalar *b;
93: PetscScalar s1,s2,s3,s4,s5,s6,x1,x2,x3,x4,x5,x6;
96: VecGetArrayRead(bb,&b);
97: VecGetArray(xx,&x);
98: /* forward solve the lower triangular */
99: idx = 0;
100: x[0] = b[idx]; x[1] = b[1+idx];x[2] = b[2+idx];x[3] = b[3+idx];
101: x[4] = b[4+idx];x[5] = b[5+idx];
102: for (i=1; i<n; i++) {
103: v = aa + bs2*ai[i];
104: vi = aj + ai[i];
105: nz = ai[i+1] - ai[i];
106: idx = bs*i;
107: s1 = b[idx];s2 = b[1+idx];s3 = b[2+idx];s4 = b[3+idx];
108: s5 = b[4+idx];s6 = b[5+idx];
109: for (k=0; k<nz; k++) {
110: jdx = bs*vi[k];
111: x1 = x[jdx];x2 = x[1+jdx]; x3 =x[2+jdx];x4 =x[3+jdx];
112: x5 = x[4+jdx]; x6 = x[5+jdx];
113: s1 -= v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6;
114: s2 -= v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6;
115: s3 -= v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6;
116: s4 -= v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6;
117: s5 -= v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6;
118: s6 -= v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6;
119: v += bs2;
120: }
122: x[idx] = s1;
123: x[1+idx] = s2;
124: x[2+idx] = s3;
125: x[3+idx] = s4;
126: x[4+idx] = s5;
127: x[5+idx] = s6;
128: }
130: /* backward solve the upper triangular */
131: for (i=n-1; i>=0; i--) {
132: v = aa + bs2*(adiag[i+1]+1);
133: vi = aj + adiag[i+1]+1;
134: nz = adiag[i] - adiag[i+1]-1;
135: idt = bs*i;
136: s1 = x[idt]; s2 = x[1+idt];s3 = x[2+idt];s4 = x[3+idt];
137: s5 = x[4+idt];s6 = x[5+idt];
138: for (k=0; k<nz; k++) {
139: idx = bs*vi[k];
140: x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx];x4 = x[3+idx];
141: x5 = x[4+idx];x6 = x[5+idx];
142: s1 -= v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6;
143: s2 -= v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6;
144: s3 -= v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6;
145: s4 -= v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6;
146: s5 -= v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6;
147: s6 -= v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6;
148: v += bs2;
149: }
150: /* x = inv_diagonal*x */
151: x[idt] = v[0]*s1 + v[6]*s2 + v[12]*s3 + v[18]*s4 + v[24]*s5 + v[30]*s6;
152: x[1+idt] = v[1]*s1 + v[7]*s2 + v[13]*s3 + v[19]*s4 + v[25]*s5 + v[31]*s6;
153: x[2+idt] = v[2]*s1 + v[8]*s2 + v[14]*s3 + v[20]*s4 + v[26]*s5 + v[32]*s6;
154: x[3+idt] = v[3]*s1 + v[9]*s2 + v[15]*s3 + v[21]*s4 + v[27]*s5 + v[33]*s6;
155: x[4+idt] = v[4]*s1 + v[10]*s2 + v[16]*s3 + v[22]*s4 + v[28]*s5 + v[34]*s6;
156: x[5+idt] = v[5]*s1 + v[11]*s2 + v[17]*s3 + v[23]*s4 + v[29]*s5 + v[35]*s6;
157: }
159: VecRestoreArrayRead(bb,&b);
160: VecRestoreArray(xx,&x);
161: PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);
162: return(0);
163: }