Actual source code: baijsolvtrannat7.c

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

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

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

 15:   /* forward solve the U^T */
 16:   idx = 0;
 17:   for (i = 0; i < n; i++) {
 18:     v = aa + 49 * diag[i];
 19:     /* multiply by the inverse of the block diagonal */
 20:     x1 = x[idx];
 21:     x2 = x[1 + idx];
 22:     x3 = x[2 + idx];
 23:     x4 = x[3 + idx];
 24:     x5 = x[4 + idx];
 25:     x6 = x[5 + idx];
 26:     x7 = x[6 + idx];
 27:     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5 + v[5] * x6 + v[6] * x7;
 28:     s2 = v[7] * x1 + v[8] * x2 + v[9] * x3 + v[10] * x4 + v[11] * x5 + v[12] * x6 + v[13] * x7;
 29:     s3 = v[14] * x1 + v[15] * x2 + v[16] * x3 + v[17] * x4 + v[18] * x5 + v[19] * x6 + v[20] * x7;
 30:     s4 = v[21] * x1 + v[22] * x2 + v[23] * x3 + v[24] * x4 + v[25] * x5 + v[26] * x6 + v[27] * x7;
 31:     s5 = v[28] * x1 + v[29] * x2 + v[30] * x3 + v[31] * x4 + v[32] * x5 + v[33] * x6 + v[34] * x7;
 32:     s6 = v[35] * x1 + v[36] * x2 + v[37] * x3 + v[38] * x4 + v[39] * x5 + v[40] * x6 + v[41] * x7;
 33:     s7 = v[42] * x1 + v[43] * x2 + v[44] * x3 + v[45] * x4 + v[46] * x5 + v[47] * x6 + v[48] * x7;
 34:     v += 49;

 36:     vi = aj + diag[i] + 1;
 37:     nz = ai[i + 1] - diag[i] - 1;
 38:     while (nz--) {
 39:       oidx = 7 * (*vi++);
 40:       x[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6 + v[6] * s7;
 41:       x[oidx + 1] -= v[7] * s1 + v[8] * s2 + v[9] * s3 + v[10] * s4 + v[11] * s5 + v[12] * s6 + v[13] * s7;
 42:       x[oidx + 2] -= v[14] * s1 + v[15] * s2 + v[16] * s3 + v[17] * s4 + v[18] * s5 + v[19] * s6 + v[20] * s7;
 43:       x[oidx + 3] -= v[21] * s1 + v[22] * s2 + v[23] * s3 + v[24] * s4 + v[25] * s5 + v[26] * s6 + v[27] * s7;
 44:       x[oidx + 4] -= v[28] * s1 + v[29] * s2 + v[30] * s3 + v[31] * s4 + v[32] * s5 + v[33] * s6 + v[34] * s7;
 45:       x[oidx + 5] -= v[35] * s1 + v[36] * s2 + v[37] * s3 + v[38] * s4 + v[39] * s5 + v[40] * s6 + v[41] * s7;
 46:       x[oidx + 6] -= v[42] * s1 + v[43] * s2 + v[44] * s3 + v[45] * s4 + v[46] * s5 + v[47] * s6 + v[48] * s7;
 47:       v += 49;
 48:     }
 49:     x[idx]     = s1;
 50:     x[1 + idx] = s2;
 51:     x[2 + idx] = s3;
 52:     x[3 + idx] = s4;
 53:     x[4 + idx] = s5;
 54:     x[5 + idx] = s6;
 55:     x[6 + idx] = s7;
 56:     idx += 7;
 57:   }
 58:   /* backward solve the L^T */
 59:   for (i = n - 1; i >= 0; i--) {
 60:     v   = aa + 49 * diag[i] - 49;
 61:     vi  = aj + diag[i] - 1;
 62:     nz  = diag[i] - ai[i];
 63:     idt = 7 * i;
 64:     s1  = x[idt];
 65:     s2  = x[1 + idt];
 66:     s3  = x[2 + idt];
 67:     s4  = x[3 + idt];
 68:     s5  = x[4 + idt];
 69:     s6  = x[5 + idt];
 70:     s7  = x[6 + idt];
 71:     while (nz--) {
 72:       idx = 7 * (*vi--);
 73:       x[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6 + v[6] * s7;
 74:       x[idx + 1] -= v[7] * s1 + v[8] * s2 + v[9] * s3 + v[10] * s4 + v[11] * s5 + v[12] * s6 + v[13] * s7;
 75:       x[idx + 2] -= v[14] * s1 + v[15] * s2 + v[16] * s3 + v[17] * s4 + v[18] * s5 + v[19] * s6 + v[20] * s7;
 76:       x[idx + 3] -= v[21] * s1 + v[22] * s2 + v[23] * s3 + v[24] * s4 + v[25] * s5 + v[26] * s6 + v[27] * s7;
 77:       x[idx + 4] -= v[28] * s1 + v[29] * s2 + v[30] * s3 + v[31] * s4 + v[32] * s5 + v[33] * s6 + v[34] * s7;
 78:       x[idx + 5] -= v[35] * s1 + v[36] * s2 + v[37] * s3 + v[38] * s4 + v[39] * s5 + v[40] * s6 + v[41] * s7;
 79:       x[idx + 6] -= v[42] * s1 + v[43] * s2 + v[44] * s3 + v[45] * s4 + v[46] * s5 + v[47] * s6 + v[48] * s7;
 80:       v -= 49;
 81:     }
 82:   }
 83:   PetscCall(VecRestoreArray(xx, &x));
 84:   PetscCall(PetscLogFlops(2.0 * 49 * (a->nz) - 7.0 * A->cmap->n));
 85:   PetscFunctionReturn(PETSC_SUCCESS);
 86: }
 87: PetscErrorCode MatSolveTranspose_SeqBAIJ_7_NaturalOrdering(Mat A, Vec bb, Vec xx)
 88: {
 89:   Mat_SeqBAIJ     *a = (Mat_SeqBAIJ *)A->data;
 90:   const PetscInt   n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag;
 91:   PetscInt         nz, idx, idt, j, i, oidx;
 92:   const PetscInt   bs = A->rmap->bs, bs2 = a->bs2;
 93:   const MatScalar *aa = a->a, *v;
 94:   PetscScalar      s1, s2, s3, s4, s5, s6, s7, x1, x2, x3, x4, x5, x6, x7, *x;

 96:   PetscFunctionBegin;
 97:   PetscCall(VecCopy(bb, xx));
 98:   PetscCall(VecGetArray(xx, &x));

100:   /* forward solve the U^T */
101:   idx = 0;
102:   for (i = 0; i < n; i++) {
103:     v = aa + bs2 * diag[i];
104:     /* multiply by the inverse of the block diagonal */
105:     x1 = x[idx];
106:     x2 = x[1 + idx];
107:     x3 = x[2 + idx];
108:     x4 = x[3 + idx];
109:     x5 = x[4 + idx];
110:     x6 = x[5 + idx];
111:     x7 = x[6 + idx];
112:     s1 = v[0] * x1 + v[1] * x2 + v[2] * x3 + v[3] * x4 + v[4] * x5 + v[5] * x6 + v[6] * x7;
113:     s2 = v[7] * x1 + v[8] * x2 + v[9] * x3 + v[10] * x4 + v[11] * x5 + v[12] * x6 + v[13] * x7;
114:     s3 = v[14] * x1 + v[15] * x2 + v[16] * x3 + v[17] * x4 + v[18] * x5 + v[19] * x6 + v[20] * x7;
115:     s4 = v[21] * x1 + v[22] * x2 + v[23] * x3 + v[24] * x4 + v[25] * x5 + v[26] * x6 + v[27] * x7;
116:     s5 = v[28] * x1 + v[29] * x2 + v[30] * x3 + v[31] * x4 + v[32] * x5 + v[33] * x6 + v[34] * x7;
117:     s6 = v[35] * x1 + v[36] * x2 + v[37] * x3 + v[38] * x4 + v[39] * x5 + v[40] * x6 + v[41] * x7;
118:     s7 = v[42] * x1 + v[43] * x2 + v[44] * x3 + v[45] * x4 + v[46] * x5 + v[47] * x6 + v[48] * x7;
119:     v -= bs2;
120:     vi = aj + diag[i] - 1;
121:     nz = diag[i] - diag[i + 1] - 1;
122:     for (j = 0; j > -nz; j--) {
123:       oidx = bs * vi[j];
124:       x[oidx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6 + v[6] * s7;
125:       x[oidx + 1] -= v[7] * s1 + v[8] * s2 + v[9] * s3 + v[10] * s4 + v[11] * s5 + v[12] * s6 + v[13] * s7;
126:       x[oidx + 2] -= v[14] * s1 + v[15] * s2 + v[16] * s3 + v[17] * s4 + v[18] * s5 + v[19] * s6 + v[20] * s7;
127:       x[oidx + 3] -= v[21] * s1 + v[22] * s2 + v[23] * s3 + v[24] * s4 + v[25] * s5 + v[26] * s6 + v[27] * s7;
128:       x[oidx + 4] -= v[28] * s1 + v[29] * s2 + v[30] * s3 + v[31] * s4 + v[32] * s5 + v[33] * s6 + v[34] * s7;
129:       x[oidx + 5] -= v[35] * s1 + v[36] * s2 + v[37] * s3 + v[38] * s4 + v[39] * s5 + v[40] * s6 + v[41] * s7;
130:       x[oidx + 6] -= v[42] * s1 + v[43] * s2 + v[44] * s3 + v[45] * s4 + v[46] * s5 + v[47] * s6 + v[48] * s7;
131:       v -= bs2;
132:     }
133:     x[idx]     = s1;
134:     x[1 + idx] = s2;
135:     x[2 + idx] = s3;
136:     x[3 + idx] = s4;
137:     x[4 + idx] = s5;
138:     x[5 + idx] = s6;
139:     x[6 + idx] = s7;
140:     idx += bs;
141:   }
142:   /* backward solve the L^T */
143:   for (i = n - 1; i >= 0; i--) {
144:     v   = aa + bs2 * ai[i];
145:     vi  = aj + ai[i];
146:     nz  = ai[i + 1] - ai[i];
147:     idt = bs * i;
148:     s1  = x[idt];
149:     s2  = x[1 + idt];
150:     s3  = x[2 + idt];
151:     s4  = x[3 + idt];
152:     s5  = x[4 + idt];
153:     s6  = x[5 + idt];
154:     s7  = x[6 + idt];
155:     for (j = 0; j < nz; j++) {
156:       idx = bs * vi[j];
157:       x[idx] -= v[0] * s1 + v[1] * s2 + v[2] * s3 + v[3] * s4 + v[4] * s5 + v[5] * s6 + v[6] * s7;
158:       x[idx + 1] -= v[7] * s1 + v[8] * s2 + v[9] * s3 + v[10] * s4 + v[11] * s5 + v[12] * s6 + v[13] * s7;
159:       x[idx + 2] -= v[14] * s1 + v[15] * s2 + v[16] * s3 + v[17] * s4 + v[18] * s5 + v[19] * s6 + v[20] * s7;
160:       x[idx + 3] -= v[21] * s1 + v[22] * s2 + v[23] * s3 + v[24] * s4 + v[25] * s5 + v[26] * s6 + v[27] * s7;
161:       x[idx + 4] -= v[28] * s1 + v[29] * s2 + v[30] * s3 + v[31] * s4 + v[32] * s5 + v[33] * s6 + v[34] * s7;
162:       x[idx + 5] -= v[35] * s1 + v[36] * s2 + v[37] * s3 + v[38] * s4 + v[39] * s5 + v[40] * s6 + v[41] * s7;
163:       x[idx + 6] -= v[42] * s1 + v[43] * s2 + v[44] * s3 + v[45] * s4 + v[46] * s5 + v[47] * s6 + v[48] * s7;
164:       v += bs2;
165:     }
166:   }
167:   PetscCall(VecRestoreArray(xx, &x));
168:   PetscCall(PetscLogFlops(2.0 * bs2 * (a->nz) - bs * A->cmap->n));
169:   PetscFunctionReturn(PETSC_SUCCESS);
170: }