Actual source code: baijsolvtrannat3.c

  1: #include <../src/mat/impls/baij/seq/baij.h>

  3: PetscErrorCode MatSolveTranspose_SeqBAIJ_3_NaturalOrdering_inplace(Mat A, Vec bb, Vec xx)
  4: {
  5:   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *)A->data;
  6:   const PetscInt   n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag;
  7:   PetscInt         i, nz, idx, idt, oidx;
  8:   const MatScalar *aa = a->a, *v;
  9:   PetscScalar      s1, s2, s3, x1, x2, x3, *x;

 11:   VecCopy(bb, xx);
 12:   VecGetArray(xx, &x);

 14:   /* forward solve the U^T */
 15:   idx = 0;
 16:   for (i = 0; i < n; i++) {
 17:     v = aa + 9 * diag[i];
 18:     /* multiply by the inverse of the block diagonal */
 19:     x1 = x[idx];
 20:     x2 = x[1 + idx];
 21:     x3 = x[2 + idx];
 22:     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3;
 23:     s2 = v[3] * x1 + v[4] * x2 + v[5] * x3;
 24:     s3 = v[6] * x1 + v[7] * x2 + v[8] * x3;
 25:     v += 9;

 27:     vi = aj + diag[i] + 1;
 28:     nz = ai[i + 1] - diag[i] - 1;
 29:     while (nz--) {
 30:       oidx = 3 * (*vi++);
 31:       x[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3;
 32:       x[oidx + 1] -= v[3] * s1 + v[4] * s2 + v[5] * s3;
 33:       x[oidx + 2] -= v[6] * s1 + v[7] * s2 + v[8] * s3;
 34:       v += 9;
 35:     }
 36:     x[idx]     = s1;
 37:     x[1 + idx] = s2;
 38:     x[2 + idx] = s3;
 39:     idx += 3;
 40:   }
 41:   /* backward solve the L^T */
 42:   for (i = n - 1; i >= 0; i--) {
 43:     v   = aa + 9 * diag[i] - 9;
 44:     vi  = aj + diag[i] - 1;
 45:     nz  = diag[i] - ai[i];
 46:     idt = 3 * i;
 47:     s1  = x[idt];
 48:     s2  = x[1 + idt];
 49:     s3  = x[2 + idt];
 50:     while (nz--) {
 51:       idx = 3 * (*vi--);
 52:       x[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3;
 53:       x[idx + 1] -= v[3] * s1 + v[4] * s2 + v[5] * s3;
 54:       x[idx + 2] -= v[6] * s1 + v[7] * s2 + v[8] * s3;
 55:       v -= 9;
 56:     }
 57:   }
 58:   VecRestoreArray(xx, &x);
 59:   PetscLogFlops(2.0 * 9.0 * (a->nz) - 3.0 * A->cmap->n);
 60:   return 0;
 61: }

 63: PetscErrorCode MatSolveTranspose_SeqBAIJ_3_NaturalOrdering(Mat A, Vec bb, Vec xx)
 64: {
 65:   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *)A->data;
 66:   const PetscInt   n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag;
 67:   PetscInt         nz, idx, idt, j, i, oidx;
 68:   const PetscInt   bs = A->rmap->bs, bs2 = a->bs2;
 69:   const MatScalar *aa = a->a, *v;
 70:   PetscScalar      s1, s2, s3, x1, x2, x3, *x;

 72:   VecCopy(bb, xx);
 73:   VecGetArray(xx, &x);

 75:   /* forward solve the U^T */
 76:   idx = 0;
 77:   for (i = 0; i < n; i++) {
 78:     v = aa + bs2 * diag[i];
 79:     /* multiply by the inverse of the block diagonal */
 80:     x1 = x[idx];
 81:     x2 = x[1 + idx];
 82:     x3 = x[2 + idx];
 83:     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3;
 84:     s2 = v[3] * x1 + v[4] * x2 + v[5] * x3;
 85:     s3 = v[6] * x1 + v[7] * x2 + v[8] * x3;
 86:     v -= bs2;

 88:     vi = aj + diag[i] - 1;
 89:     nz = diag[i] - diag[i + 1] - 1;
 90:     for (j = 0; j > -nz; j--) {
 91:       oidx = bs * vi[j];
 92:       x[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3;
 93:       x[oidx + 1] -= v[3] * s1 + v[4] * s2 + v[5] * s3;
 94:       x[oidx + 2] -= v[6] * s1 + v[7] * s2 + v[8] * s3;
 95:       v -= bs2;
 96:     }
 97:     x[idx]     = s1;
 98:     x[1 + idx] = s2;
 99:     x[2 + idx] = s3;
100:     idx += bs;
101:   }
102:   /* backward solve the L^T */
103:   for (i = n - 1; i >= 0; i--) {
104:     v   = aa + bs2 * ai[i];
105:     vi  = aj + ai[i];
106:     nz  = ai[i + 1] - ai[i];
107:     idt = bs * i;
108:     s1  = x[idt];
109:     s2  = x[1 + idt];
110:     s3  = x[2 + idt];
111:     for (j = 0; j < nz; j++) {
112:       idx = bs * vi[j];
113:       x[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3;
114:       x[idx + 1] -= v[3] * s1 + v[4] * s2 + v[5] * s3;
115:       x[idx + 2] -= v[6] * s1 + v[7] * s2 + v[8] * s3;
116:       v += bs2;
117:     }
118:   }
119:   VecRestoreArray(xx, &x);
120:   PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n);
121:   return 0;
122: }