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