xref: /petsc/src/mat/impls/baij/seq/baijsolvtrannat3.c (revision 2c733ed45a445b0336611d812c866924feeaee2c)
1*2c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h>
2*2c733ed4SBarry Smith 
3*2c733ed4SBarry Smith PetscErrorCode MatSolveTranspose_SeqBAIJ_3_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)
4*2c733ed4SBarry Smith {
5*2c733ed4SBarry Smith   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ*)A->data;
6*2c733ed4SBarry Smith   PetscErrorCode  ierr;
7*2c733ed4SBarry Smith   const PetscInt  n=a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag;
8*2c733ed4SBarry Smith   PetscInt        i,nz,idx,idt,oidx;
9*2c733ed4SBarry Smith   const MatScalar *aa=a->a,*v;
10*2c733ed4SBarry Smith   PetscScalar     s1,s2,s3,x1,x2,x3,*x;
11*2c733ed4SBarry Smith 
12*2c733ed4SBarry Smith   PetscFunctionBegin;
13*2c733ed4SBarry Smith   ierr = VecCopy(bb,xx);CHKERRQ(ierr);
14*2c733ed4SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
15*2c733ed4SBarry Smith 
16*2c733ed4SBarry Smith   /* forward solve the U^T */
17*2c733ed4SBarry Smith   idx = 0;
18*2c733ed4SBarry Smith   for (i=0; i<n; i++) {
19*2c733ed4SBarry Smith 
20*2c733ed4SBarry Smith     v = aa + 9*diag[i];
21*2c733ed4SBarry Smith     /* multiply by the inverse of the block diagonal */
22*2c733ed4SBarry Smith     x1 = x[idx];   x2 = x[1+idx]; x3    = x[2+idx];
23*2c733ed4SBarry Smith     s1 = v[0]*x1  +  v[1]*x2 +  v[2]*x3;
24*2c733ed4SBarry Smith     s2 = v[3]*x1  +  v[4]*x2 +  v[5]*x3;
25*2c733ed4SBarry Smith     s3 = v[6]*x1  +  v[7]*x2 + v[8]*x3;
26*2c733ed4SBarry Smith     v += 9;
27*2c733ed4SBarry Smith 
28*2c733ed4SBarry Smith     vi = aj + diag[i] + 1;
29*2c733ed4SBarry Smith     nz = ai[i+1] - diag[i] - 1;
30*2c733ed4SBarry Smith     while (nz--) {
31*2c733ed4SBarry Smith       oidx       = 3*(*vi++);
32*2c733ed4SBarry Smith       x[oidx]   -= v[0]*s1  +  v[1]*s2 +  v[2]*s3;
33*2c733ed4SBarry Smith       x[oidx+1] -= v[3]*s1  +  v[4]*s2 +  v[5]*s3;
34*2c733ed4SBarry Smith       x[oidx+2] -= v[6]*s1 + v[7]*s2 + v[8]*s3;
35*2c733ed4SBarry Smith       v         += 9;
36*2c733ed4SBarry Smith     }
37*2c733ed4SBarry Smith     x[idx] = s1;x[1+idx] = s2; x[2+idx] = s3;
38*2c733ed4SBarry Smith     idx   += 3;
39*2c733ed4SBarry Smith   }
40*2c733ed4SBarry Smith   /* backward solve the L^T */
41*2c733ed4SBarry Smith   for (i=n-1; i>=0; i--) {
42*2c733ed4SBarry Smith     v   = aa + 9*diag[i] - 9;
43*2c733ed4SBarry Smith     vi  = aj + diag[i] - 1;
44*2c733ed4SBarry Smith     nz  = diag[i] - ai[i];
45*2c733ed4SBarry Smith     idt = 3*i;
46*2c733ed4SBarry Smith     s1  = x[idt];  s2 = x[1+idt]; s3 = x[2+idt];
47*2c733ed4SBarry Smith     while (nz--) {
48*2c733ed4SBarry Smith       idx       = 3*(*vi--);
49*2c733ed4SBarry Smith       x[idx]   -=  v[0]*s1 +  v[1]*s2 +  v[2]*s3;
50*2c733ed4SBarry Smith       x[idx+1] -=  v[3]*s1 +  v[4]*s2 +  v[5]*s3;
51*2c733ed4SBarry Smith       x[idx+2] -= v[6]*s1 + v[7]*s2 + v[8]*s3;
52*2c733ed4SBarry Smith       v        -= 9;
53*2c733ed4SBarry Smith     }
54*2c733ed4SBarry Smith   }
55*2c733ed4SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
56*2c733ed4SBarry Smith   ierr = PetscLogFlops(2.0*9.0*(a->nz) - 3.0*A->cmap->n);CHKERRQ(ierr);
57*2c733ed4SBarry Smith   PetscFunctionReturn(0);
58*2c733ed4SBarry Smith }
59*2c733ed4SBarry Smith 
60*2c733ed4SBarry Smith PetscErrorCode MatSolveTranspose_SeqBAIJ_3_NaturalOrdering(Mat A,Vec bb,Vec xx)
61*2c733ed4SBarry Smith {
62*2c733ed4SBarry Smith   Mat_SeqBAIJ     *a=(Mat_SeqBAIJ*)A->data;
63*2c733ed4SBarry Smith   PetscErrorCode  ierr;
64*2c733ed4SBarry Smith   const PetscInt  n=a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag;
65*2c733ed4SBarry Smith   PetscInt        nz,idx,idt,j,i,oidx;
66*2c733ed4SBarry Smith   const PetscInt  bs =A->rmap->bs,bs2=a->bs2;
67*2c733ed4SBarry Smith   const MatScalar *aa=a->a,*v;
68*2c733ed4SBarry Smith   PetscScalar     s1,s2,s3,x1,x2,x3,*x;
69*2c733ed4SBarry Smith 
70*2c733ed4SBarry Smith   PetscFunctionBegin;
71*2c733ed4SBarry Smith   ierr = VecCopy(bb,xx);CHKERRQ(ierr);
72*2c733ed4SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
73*2c733ed4SBarry Smith 
74*2c733ed4SBarry Smith   /* forward solve the U^T */
75*2c733ed4SBarry Smith   idx = 0;
76*2c733ed4SBarry Smith   for (i=0; i<n; i++) {
77*2c733ed4SBarry Smith     v = aa + bs2*diag[i];
78*2c733ed4SBarry Smith     /* multiply by the inverse of the block diagonal */
79*2c733ed4SBarry Smith     x1 = x[idx];   x2 = x[1+idx];  x3 = x[2+idx];
80*2c733ed4SBarry Smith     s1 = v[0]*x1  +  v[1]*x2  + v[2]*x3;
81*2c733ed4SBarry Smith     s2 = v[3]*x1  +  v[4]*x2  + v[5]*x3;
82*2c733ed4SBarry Smith     s3 = v[6]*x1  +  v[7]*x2  + v[8]*x3;
83*2c733ed4SBarry Smith     v -= bs2;
84*2c733ed4SBarry Smith 
85*2c733ed4SBarry Smith     vi = aj + diag[i] - 1;
86*2c733ed4SBarry Smith     nz = diag[i] - diag[i+1] - 1;
87*2c733ed4SBarry Smith     for (j=0; j>-nz; j--) {
88*2c733ed4SBarry Smith       oidx       = bs*vi[j];
89*2c733ed4SBarry Smith       x[oidx]   -= v[0]*s1  +  v[1]*s2  + v[2]*s3;
90*2c733ed4SBarry Smith       x[oidx+1] -= v[3]*s1  +  v[4]*s2  + v[5]*s3;
91*2c733ed4SBarry Smith       x[oidx+2] -= v[6]*s1  +  v[7]*s2  + v[8]*s3;
92*2c733ed4SBarry Smith       v         -= bs2;
93*2c733ed4SBarry Smith     }
94*2c733ed4SBarry Smith     x[idx] = s1;x[1+idx] = s2;  x[2+idx] = s3;
95*2c733ed4SBarry Smith     idx   += bs;
96*2c733ed4SBarry Smith   }
97*2c733ed4SBarry Smith   /* backward solve the L^T */
98*2c733ed4SBarry Smith   for (i=n-1; i>=0; i--) {
99*2c733ed4SBarry Smith     v   = aa + bs2*ai[i];
100*2c733ed4SBarry Smith     vi  = aj + ai[i];
101*2c733ed4SBarry Smith     nz  = ai[i+1] - ai[i];
102*2c733ed4SBarry Smith     idt = bs*i;
103*2c733ed4SBarry Smith     s1  = x[idt];  s2 = x[1+idt];  s3 = x[2+idt];
104*2c733ed4SBarry Smith     for (j=0; j<nz; j++) {
105*2c733ed4SBarry Smith       idx       = bs*vi[j];
106*2c733ed4SBarry Smith       x[idx]   -= v[0]*s1  +  v[1]*s2  + v[2]*s3;
107*2c733ed4SBarry Smith       x[idx+1] -= v[3]*s1  +  v[4]*s2  + v[5]*s3;
108*2c733ed4SBarry Smith       x[idx+2] -= v[6]*s1  +  v[7]*s2  + v[8]*s3;
109*2c733ed4SBarry Smith       v        += bs2;
110*2c733ed4SBarry Smith     }
111*2c733ed4SBarry Smith   }
112*2c733ed4SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
113*2c733ed4SBarry Smith   ierr = PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);CHKERRQ(ierr);
114*2c733ed4SBarry Smith   PetscFunctionReturn(0);
115*2c733ed4SBarry Smith }
116