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