xref: /petsc/src/mat/impls/baij/seq/baijsolvnat11.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 /* Block operations are done by accessing one column at at time */
5*2c733ed4SBarry Smith /* Default MatSolve for block size 11 */
6*2c733ed4SBarry Smith 
7*2c733ed4SBarry Smith PetscErrorCode MatSolve_SeqBAIJ_11_NaturalOrdering(Mat A,Vec bb,Vec xx)
8*2c733ed4SBarry Smith {
9*2c733ed4SBarry Smith   Mat_SeqBAIJ       *a=(Mat_SeqBAIJ*)A->data;
10*2c733ed4SBarry Smith   PetscErrorCode    ierr;
11*2c733ed4SBarry Smith   const PetscInt    n=a->mbs,*ai=a->i,*aj=a->j,*adiag=a->diag,*vi,bs=A->rmap->bs,bs2=a->bs2;
12*2c733ed4SBarry Smith   PetscInt          i,k,nz,idx,idt,m;
13*2c733ed4SBarry Smith   const MatScalar   *aa=a->a,*v;
14*2c733ed4SBarry Smith   PetscScalar       s[11];
15*2c733ed4SBarry Smith   PetscScalar       *x,xv;
16*2c733ed4SBarry Smith   const PetscScalar *b;
17*2c733ed4SBarry Smith 
18*2c733ed4SBarry Smith   PetscFunctionBegin;
19*2c733ed4SBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
20*2c733ed4SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
21*2c733ed4SBarry Smith 
22*2c733ed4SBarry Smith   /* forward solve the lower triangular */
23*2c733ed4SBarry Smith   for (i=0; i<n; i++) {
24*2c733ed4SBarry Smith     v         = aa + bs2*ai[i];
25*2c733ed4SBarry Smith     vi        = aj + ai[i];
26*2c733ed4SBarry Smith     nz        = ai[i+1] - ai[i];
27*2c733ed4SBarry Smith     idt       = bs*i;
28*2c733ed4SBarry Smith     x[idt]    = b[idt];    x[1+idt]  = b[1+idt];  x[2+idt]  = b[2+idt];  x[3+idt]  = b[3+idt];  x[4+idt]  = b[4+idt];
29*2c733ed4SBarry Smith     x[5+idt]  = b[5+idt];  x[6+idt]  = b[6+idt];  x[7+idt]  = b[7+idt];  x[8+idt]  = b[8+idt];  x[9+idt] = b[9+idt];
30*2c733ed4SBarry Smith     x[10+idt] = b[10+idt];
31*2c733ed4SBarry Smith     for (m=0; m<nz; m++) {
32*2c733ed4SBarry Smith       idx = bs*vi[m];
33*2c733ed4SBarry Smith       for (k=0; k<11; k++) {
34*2c733ed4SBarry Smith         xv         = x[k + idx];
35*2c733ed4SBarry Smith         x[idt]    -= v[0]*xv;
36*2c733ed4SBarry Smith         x[1+idt]  -= v[1]*xv;
37*2c733ed4SBarry Smith         x[2+idt]  -= v[2]*xv;
38*2c733ed4SBarry Smith         x[3+idt]  -= v[3]*xv;
39*2c733ed4SBarry Smith         x[4+idt]  -= v[4]*xv;
40*2c733ed4SBarry Smith         x[5+idt]  -= v[5]*xv;
41*2c733ed4SBarry Smith         x[6+idt]  -= v[6]*xv;
42*2c733ed4SBarry Smith         x[7+idt]  -= v[7]*xv;
43*2c733ed4SBarry Smith         x[8+idt]  -= v[8]*xv;
44*2c733ed4SBarry Smith         x[9+idt]  -= v[9]*xv;
45*2c733ed4SBarry Smith         x[10+idt] -= v[10]*xv;
46*2c733ed4SBarry Smith         v         += 11;
47*2c733ed4SBarry Smith       }
48*2c733ed4SBarry Smith     }
49*2c733ed4SBarry Smith   }
50*2c733ed4SBarry Smith   /* backward solve the upper triangular */
51*2c733ed4SBarry Smith   for (i=n-1; i>=0; i--) {
52*2c733ed4SBarry Smith     v     = aa + bs2*(adiag[i+1]+1);
53*2c733ed4SBarry Smith     vi    = aj + adiag[i+1]+1;
54*2c733ed4SBarry Smith     nz    = adiag[i] - adiag[i+1] - 1;
55*2c733ed4SBarry Smith     idt   = bs*i;
56*2c733ed4SBarry Smith     s[0]  = x[idt];    s[1]  = x[1+idt];  s[2]  = x[2+idt];  s[3]  = x[3+idt];  s[4]  = x[4+idt];
57*2c733ed4SBarry Smith     s[5]  = x[5+idt];  s[6]  = x[6+idt];  s[7]  = x[7+idt];  s[8]  = x[8+idt];  s[9]  = x[9+idt];
58*2c733ed4SBarry Smith     s[10] = x[10+idt];
59*2c733ed4SBarry Smith 
60*2c733ed4SBarry Smith     for (m=0; m<nz; m++) {
61*2c733ed4SBarry Smith       idx = bs*vi[m];
62*2c733ed4SBarry Smith       for (k=0; k<11; k++) {
63*2c733ed4SBarry Smith         xv     = x[k + idx];
64*2c733ed4SBarry Smith         s[0]  -= v[0]*xv;
65*2c733ed4SBarry Smith         s[1]  -= v[1]*xv;
66*2c733ed4SBarry Smith         s[2]  -= v[2]*xv;
67*2c733ed4SBarry Smith         s[3]  -= v[3]*xv;
68*2c733ed4SBarry Smith         s[4]  -= v[4]*xv;
69*2c733ed4SBarry Smith         s[5]  -= v[5]*xv;
70*2c733ed4SBarry Smith         s[6]  -= v[6]*xv;
71*2c733ed4SBarry Smith         s[7]  -= v[7]*xv;
72*2c733ed4SBarry Smith         s[8]  -= v[8]*xv;
73*2c733ed4SBarry Smith         s[9]  -= v[9]*xv;
74*2c733ed4SBarry Smith         s[10] -= v[10]*xv;
75*2c733ed4SBarry Smith         v     += 11;
76*2c733ed4SBarry Smith       }
77*2c733ed4SBarry Smith     }
78*2c733ed4SBarry Smith     ierr = PetscMemzero(x+idt,bs*sizeof(MatScalar));CHKERRQ(ierr);
79*2c733ed4SBarry Smith     for (k=0; k<11; k++) {
80*2c733ed4SBarry Smith       x[idt]    += v[0]*s[k];
81*2c733ed4SBarry Smith       x[1+idt]  += v[1]*s[k];
82*2c733ed4SBarry Smith       x[2+idt]  += v[2]*s[k];
83*2c733ed4SBarry Smith       x[3+idt]  += v[3]*s[k];
84*2c733ed4SBarry Smith       x[4+idt]  += v[4]*s[k];
85*2c733ed4SBarry Smith       x[5+idt]  += v[5]*s[k];
86*2c733ed4SBarry Smith       x[6+idt]  += v[6]*s[k];
87*2c733ed4SBarry Smith       x[7+idt]  += v[7]*s[k];
88*2c733ed4SBarry Smith       x[8+idt]  += v[8]*s[k];
89*2c733ed4SBarry Smith       x[9+idt]  += v[9]*s[k];
90*2c733ed4SBarry Smith       x[10+idt] += v[10]*s[k];
91*2c733ed4SBarry Smith       v         += 11;
92*2c733ed4SBarry Smith     }
93*2c733ed4SBarry Smith   }
94*2c733ed4SBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
95*2c733ed4SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
96*2c733ed4SBarry Smith   ierr = PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);CHKERRQ(ierr);
97*2c733ed4SBarry Smith   PetscFunctionReturn(0);
98*2c733ed4SBarry Smith }
99