1 #include <../src/mat/impls/baij/seq/baij.h> 2 #include <petsc/private/kernels/blockinvert.h> 3 4 PetscErrorCode MatSolveTranspose_SeqBAIJ_1(Mat A, Vec bb, Vec xx) { 5 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 6 IS iscol = a->col, isrow = a->row; 7 const PetscInt *rout, *cout, *r, *c, *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; 13 14 PetscFunctionBegin; 15 PetscCall(VecGetArrayRead(bb, &b)); 16 PetscCall(VecGetArray(xx, &x)); 17 tmp = a->solve_work; 18 19 PetscCall(ISGetIndices(isrow, &rout)); 20 r = rout; 21 PetscCall(ISGetIndices(iscol, &cout)); 22 c = cout; 23 24 /* copy the b into temp work space according to permutation */ 25 for (i = 0; i < n; i++) tmp[i] = b[c[i]]; 26 27 /* forward solve the U^T */ 28 for (i = 0; i < n; i++) { 29 v = aa + adiag[i + 1] + 1; 30 vi = aj + adiag[i + 1] + 1; 31 nz = adiag[i] - adiag[i + 1] - 1; 32 s1 = tmp[i]; 33 s1 *= v[nz]; /* multiply by inverse of diagonal entry */ 34 for (j = 0; j < nz; j++) tmp[vi[j]] -= s1 * v[j]; 35 tmp[i] = s1; 36 } 37 38 /* backward solve the L^T */ 39 for (i = n - 1; i >= 0; i--) { 40 v = aa + ai[i]; 41 vi = aj + ai[i]; 42 nz = ai[i + 1] - ai[i]; 43 s1 = tmp[i]; 44 for (j = 0; j < nz; j++) tmp[vi[j]] -= s1 * v[j]; 45 } 46 47 /* copy tmp into x according to permutation */ 48 for (i = 0; i < n; i++) x[r[i]] = tmp[i]; 49 50 PetscCall(ISRestoreIndices(isrow, &rout)); 51 PetscCall(ISRestoreIndices(iscol, &cout)); 52 PetscCall(VecRestoreArrayRead(bb, &b)); 53 PetscCall(VecRestoreArray(xx, &x)); 54 55 PetscCall(PetscLogFlops(2.0 * a->nz - A->cmap->n)); 56 PetscFunctionReturn(0); 57 } 58 59 PetscErrorCode MatSolveTranspose_SeqBAIJ_1_inplace(Mat A, Vec bb, Vec xx) { 60 Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 61 IS iscol = a->col, isrow = a->row; 62 const PetscInt *r, *c, *rout, *cout; 63 const PetscInt *diag = a->diag, n = a->mbs, *vi, *ai = a->i, *aj = a->j; 64 PetscInt i, nz; 65 const MatScalar *aa = a->a, *v; 66 PetscScalar s1, *x, *t; 67 const PetscScalar *b; 68 69 PetscFunctionBegin; 70 PetscCall(VecGetArrayRead(bb, &b)); 71 PetscCall(VecGetArray(xx, &x)); 72 t = a->solve_work; 73 74 PetscCall(ISGetIndices(isrow, &rout)); 75 r = rout; 76 PetscCall(ISGetIndices(iscol, &cout)); 77 c = cout; 78 79 /* copy the b into temp work space according to permutation */ 80 for (i = 0; i < n; i++) t[i] = b[c[i]]; 81 82 /* forward solve the U^T */ 83 for (i = 0; i < n; i++) { 84 v = aa + diag[i]; 85 /* multiply by the inverse of the block diagonal */ 86 s1 = (*v++) * t[i]; 87 vi = aj + diag[i] + 1; 88 nz = ai[i + 1] - diag[i] - 1; 89 while (nz--) { t[*vi++] -= (*v++) * s1; } 90 t[i] = s1; 91 } 92 /* backward solve the L^T */ 93 for (i = n - 1; i >= 0; i--) { 94 v = aa + diag[i] - 1; 95 vi = aj + diag[i] - 1; 96 nz = diag[i] - ai[i]; 97 s1 = t[i]; 98 while (nz--) { t[*vi--] -= (*v--) * s1; } 99 } 100 101 /* copy t into x according to permutation */ 102 for (i = 0; i < n; i++) x[r[i]] = t[i]; 103 104 PetscCall(ISRestoreIndices(isrow, &rout)); 105 PetscCall(ISRestoreIndices(iscol, &cout)); 106 PetscCall(VecRestoreArrayRead(bb, &b)); 107 PetscCall(VecRestoreArray(xx, &x)); 108 PetscCall(PetscLogFlops(2.0 * (a->nz) - A->cmap->n)); 109 PetscFunctionReturn(0); 110 } 111