Actual source code: baijsolvtrannat1.c
1: #include <../src/mat/impls/baij/seq/baij.h>
3: PetscErrorCode MatSolveTranspose_SeqBAIJ_1_NaturalOrdering(Mat A,Vec bb,Vec xx)
4: {
5: Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data;
6: PetscErrorCode ierr;
7: const PetscInt *adiag = a->diag,*ai = a->i,*aj = a->j,*vi;
8: PetscInt i,n = a->mbs,j;
9: PetscInt nz;
10: PetscScalar *x,*tmp,s1;
11: const MatScalar *aa = a->a,*v;
12: const PetscScalar *b;
15: VecGetArrayRead(bb,&b);
16: VecGetArray(xx,&x);
17: tmp = a->solve_work;
19: /* copy the b into temp work space according to permutation */
20: for (i=0; i<n; i++) tmp[i] = b[i];
22: /* forward solve the U^T */
23: for (i=0; i<n; i++) {
24: v = aa + adiag[i+1] + 1;
25: vi = aj + adiag[i+1] + 1;
26: nz = adiag[i] - adiag[i+1] - 1;
27: s1 = tmp[i];
28: s1 *= v[nz]; /* multiply by inverse of diagonal entry */
29: for (j=0; j<nz; j++) tmp[vi[j]] -= s1*v[j];
30: tmp[i] = s1;
31: }
33: /* backward solve the L^T */
34: for (i=n-1; i>=0; i--) {
35: v = aa + ai[i];
36: vi = aj + ai[i];
37: nz = ai[i+1] - ai[i];
38: s1 = tmp[i];
39: for (j=0; j<nz; j++) tmp[vi[j]] -= s1*v[j];
40: }
42: /* copy tmp into x according to permutation */
43: for (i=0; i<n; i++) x[i] = tmp[i];
45: VecRestoreArrayRead(bb,&b);
46: VecRestoreArray(xx,&x);
48: PetscLogFlops(2.0*a->nz-A->cmap->n);
49: return(0);
50: }
52: PetscErrorCode MatSolveTranspose_SeqBAIJ_1_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)
53: {
54: Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data;
55: PetscErrorCode ierr;
56: PetscInt i,nz;
57: const PetscInt *diag = a->diag,n=a->mbs,*vi,*ai=a->i,*aj=a->j;
58: const MatScalar *aa =a->a,*v;
59: PetscScalar s1,*x;
62: VecCopy(bb,xx);
63: VecGetArray(xx,&x);
65: /* forward solve the U^T */
66: for (i=0; i<n; i++) {
68: v = aa + diag[i];
69: /* multiply by the inverse of the block diagonal */
70: s1 = (*v++)*x[i];
71: vi = aj + diag[i] + 1;
72: nz = ai[i+1] - diag[i] - 1;
73: while (nz--) {
74: x[*vi++] -= (*v++)*s1;
75: }
76: x[i] = s1;
77: }
78: /* backward solve the L^T */
79: for (i=n-1; i>=0; i--) {
80: v = aa + diag[i] - 1;
81: vi = aj + diag[i] - 1;
82: nz = diag[i] - ai[i];
83: s1 = x[i];
84: while (nz--) {
85: x[*vi--] -= (*v--)*s1;
86: }
87: }
88: VecRestoreArray(xx,&x);
89: PetscLogFlops(2.0*(a->nz) - A->cmap->n);
90: return(0);
91: }