Actual source code: baijsolvtran6.c

  1: #include <../src/mat/impls/baij/seq/baij.h>
  2: #include <petsc/private/kernels/blockinvert.h>

  4: PetscErrorCode MatSolveTranspose_SeqBAIJ_6_inplace(Mat A,Vec bb,Vec xx)
  5: {
  6:   Mat_SeqBAIJ       *a   =(Mat_SeqBAIJ*)A->data;
  7:   IS                iscol=a->col,isrow=a->row;
  8:   const PetscInt    *r,*c,*rout,*cout;
  9:   const PetscInt    *diag=a->diag,n=a->mbs,*vi,*ai=a->i,*aj=a->j;
 10:   PetscInt          i,nz,idx,idt,ii,ic,ir,oidx;
 11:   const MatScalar   *aa=a->a,*v;
 12:   PetscScalar       s1,s2,s3,s4,s5,s6,x1,x2,x3,x4,x5,x6,*x,*t;
 13:   const PetscScalar *b;

 15:   VecGetArrayRead(bb,&b);
 16:   VecGetArray(xx,&x);
 17:   t    = a->solve_work;

 19:   ISGetIndices(isrow,&rout); r = rout;
 20:   ISGetIndices(iscol,&cout); c = cout;

 22:   /* copy the b into temp work space according to permutation */
 23:   ii = 0;
 24:   for (i=0; i<n; i++) {
 25:     ic      = 6*c[i];
 26:     t[ii]   = b[ic];
 27:     t[ii+1] = b[ic+1];
 28:     t[ii+2] = b[ic+2];
 29:     t[ii+3] = b[ic+3];
 30:     t[ii+4] = b[ic+4];
 31:     t[ii+5] = b[ic+5];
 32:     ii     += 6;
 33:   }

 35:   /* forward solve the U^T */
 36:   idx = 0;
 37:   for (i=0; i<n; i++) {

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

 51:     vi = aj + diag[i] + 1;
 52:     nz = ai[i+1] - diag[i] - 1;
 53:     while (nz--) {
 54:       oidx       = 6*(*vi++);
 55:       t[oidx]   -= v[0]*s1  +  v[1]*s2 +  v[2]*s3 +  v[3]*s4 +  v[4]*s5 +  v[5]*s6;
 56:       t[oidx+1] -= v[6]*s1  +  v[7]*s2 +  v[8]*s3 +  v[9]*s4 + v[10]*s5 + v[11]*s6;
 57:       t[oidx+2] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4 + v[16]*s5 + v[17]*s6;
 58:       t[oidx+3] -= v[18]*s1 + v[19]*s2 + v[20]*s3 + v[21]*s4 + v[22]*s5 + v[23]*s6;
 59:       t[oidx+4] -= v[24]*s1 + v[25]*s2 + v[26]*s3 + v[27]*s4 + v[28]*s5 + v[29]*s6;
 60:       t[oidx+5] -= v[30]*s1 + v[31]*s2 + v[32]*s3 + v[33]*s4 + v[34]*s5 + v[35]*s6;
 61:       v         += 36;
 62:     }
 63:     t[idx]   = s1;t[1+idx] = s2; t[2+idx] = s3;t[3+idx] = s4; t[4+idx] = s5;
 64:     t[5+idx] = s6;
 65:     idx     += 6;
 66:   }
 67:   /* backward solve the L^T */
 68:   for (i=n-1; i>=0; i--) {
 69:     v   = aa + 36*diag[i] - 36;
 70:     vi  = aj + diag[i] - 1;
 71:     nz  = diag[i] - ai[i];
 72:     idt = 6*i;
 73:     s1  = t[idt];  s2 = t[1+idt]; s3 = t[2+idt];s4 = t[3+idt]; s5 = t[4+idt];
 74:     s6  = t[5+idt];
 75:     while (nz--) {
 76:       idx       = 6*(*vi--);
 77:       t[idx]   -=  v[0]*s1 +  v[1]*s2 +  v[2]*s3 +  v[3]*s4 +  v[4]*s5 +  v[5]*s6;
 78:       t[idx+1] -=  v[6]*s1 +  v[7]*s2 +  v[8]*s3 +  v[9]*s4 + v[10]*s5 + v[11]*s6;
 79:       t[idx+2] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4 + v[16]*s5 + v[17]*s6;
 80:       t[idx+3] -= v[18]*s1 + v[19]*s2 + v[20]*s3 + v[21]*s4 + v[22]*s5 + v[23]*s6;
 81:       t[idx+4] -= v[24]*s1 + v[25]*s2 + v[26]*s3 + v[27]*s4 + v[28]*s5 + v[29]*s6;
 82:       t[idx+5] -= v[30]*s1 + v[31]*s2 + v[32]*s3 + v[33]*s4 + v[34]*s5 + v[35]*s6;
 83:       v        -= 36;
 84:     }
 85:   }

 87:   /* copy t into x according to permutation */
 88:   ii = 0;
 89:   for (i=0; i<n; i++) {
 90:     ir      = 6*r[i];
 91:     x[ir]   = t[ii];
 92:     x[ir+1] = t[ii+1];
 93:     x[ir+2] = t[ii+2];
 94:     x[ir+3] = t[ii+3];
 95:     x[ir+4] = t[ii+4];
 96:     x[ir+5] = t[ii+5];
 97:     ii     += 6;
 98:   }

100:   ISRestoreIndices(isrow,&rout);
101:   ISRestoreIndices(iscol,&cout);
102:   VecRestoreArrayRead(bb,&b);
103:   VecRestoreArray(xx,&x);
104:   PetscLogFlops(2.0*36*(a->nz) - 6.0*A->cmap->n);
105:   return 0;
106: }

108: PetscErrorCode MatSolveTranspose_SeqBAIJ_6(Mat A,Vec bb,Vec xx)
109: {
110:   Mat_SeqBAIJ       *a=(Mat_SeqBAIJ*)A->data;
111:   IS                iscol=a->col,isrow=a->row;
112:   const PetscInt    n    =a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag;
113:   const PetscInt    *r,*c,*rout,*cout;
114:   PetscInt          nz,idx,idt,j,i,oidx,ii,ic,ir;
115:   const PetscInt    bs =A->rmap->bs,bs2=a->bs2;
116:   const MatScalar   *aa=a->a,*v;
117:   PetscScalar       s1,s2,s3,s4,s5,s6,x1,x2,x3,x4,x5,x6,*x,*t;
118:   const PetscScalar *b;

120:   VecGetArrayRead(bb,&b);
121:   VecGetArray(xx,&x);
122:   t    = a->solve_work;

124:   ISGetIndices(isrow,&rout); r = rout;
125:   ISGetIndices(iscol,&cout); c = cout;

127:   /* copy b into temp work space according to permutation */
128:   for (i=0; i<n; i++) {
129:     ii      = bs*i; ic = bs*c[i];
130:     t[ii]   = b[ic]; t[ii+1] = b[ic+1]; t[ii+2] = b[ic+2]; t[ii+3] = b[ic+3];
131:     t[ii+4] = b[ic+4];  t[ii+5] = b[ic+5];
132:   }

134:   /* forward solve the U^T */
135:   idx = 0;
136:   for (i=0; i<n; i++) {
137:     v = aa + bs2*diag[i];
138:     /* multiply by the inverse of the block diagonal */
139:     x1 = t[idx];   x2 = t[1+idx]; x3    = t[2+idx]; x4 = t[3+idx]; x5 = t[4+idx];
140:     x6 = t[5+idx];
141:     s1 = v[0]*x1  +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5 +  v[5]*x6;
142:     s2 = v[6]*x1  +  v[7]*x2 +  v[8]*x3 +  v[9]*x4 + v[10]*x5 + v[11]*x6;
143:     s3 = v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4 + v[16]*x5 + v[17]*x6;
144:     s4 = v[18]*x1 + v[19]*x2 + v[20]*x3 + v[21]*x4 + v[22]*x5 + v[23]*x6;
145:     s5 = v[24]*x1 + v[25]*x2 + v[26]*x3 + v[27]*x4 + v[28]*x5 + v[29]*x6;
146:     s6 = v[30]*x1 + v[31]*x2 + v[32]*x3 + v[33]*x4 + v[34]*x5 + v[35]*x6;
147:     v -= bs2;

149:     vi = aj + diag[i] - 1;
150:     nz = diag[i] - diag[i+1] - 1;
151:     for (j=0; j>-nz; j--) {
152:       oidx       = bs*vi[j];
153:       t[oidx]   -= v[0]*s1  +  v[1]*s2 +  v[2]*s3 +  v[3]*s4 +  v[4]*s5 +  v[5]*s6;
154:       t[oidx+1] -= v[6]*s1  +  v[7]*s2 +  v[8]*s3 +  v[9]*s4 + v[10]*s5 + v[11]*s6;
155:       t[oidx+2] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4 + v[16]*s5 + v[17]*s6;
156:       t[oidx+3] -= v[18]*s1 + v[19]*s2 + v[20]*s3 + v[21]*s4 + v[22]*s5 + v[23]*s6;
157:       t[oidx+4] -= v[24]*s1 + v[25]*s2 + v[26]*s3 + v[27]*s4 + v[28]*s5 + v[29]*s6;
158:       t[oidx+5] -= v[30]*s1 + v[31]*s2 + v[32]*s3 + v[33]*s4 + v[34]*s5 + v[35]*s6;
159:       v         -= bs2;
160:     }
161:     t[idx]   = s1;t[1+idx] = s2;  t[2+idx] = s3;  t[3+idx] = s4; t[4+idx] =s5;
162:     t[5+idx] = s6;
163:     idx     += bs;
164:   }
165:   /* backward solve the L^T */
166:   for (i=n-1; i>=0; i--) {
167:     v   = aa + bs2*ai[i];
168:     vi  = aj + ai[i];
169:     nz  = ai[i+1] - ai[i];
170:     idt = bs*i;
171:     s1  = t[idt];  s2 = t[1+idt];  s3 = t[2+idt];  s4 = t[3+idt]; s5 = t[4+idt];
172:     s6  = t[5+idt];
173:     for (j=0; j<nz; j++) {
174:       idx       = bs*vi[j];
175:       t[idx]   -= v[0]*s1  +  v[1]*s2 +  v[2]*s3 +  v[3]*s4 +  v[4]*s5 +  v[5]*s6;
176:       t[idx+1] -= v[6]*s1  +  v[7]*s2 +  v[8]*s3 +  v[9]*s4 + v[10]*s5 + v[11]*s6;
177:       t[idx+2] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4 + v[16]*s5 + v[17]*s6;
178:       t[idx+3] -= v[18]*s1 + v[19]*s2 + v[20]*s3 + v[21]*s4 + v[22]*s5 + v[23]*s6;
179:       t[idx+4] -= v[24]*s1 + v[25]*s2 + v[26]*s3 + v[27]*s4 + v[28]*s5 + v[29]*s6;
180:       t[idx+5] -= v[30]*s1 + v[31]*s2 + v[32]*s3 + v[33]*s4 + v[34]*s5 + v[35]*s6;
181:       v        += bs2;
182:     }
183:   }

185:   /* copy t into x according to permutation */
186:   for (i=0; i<n; i++) {
187:     ii      = bs*i;  ir = bs*r[i];
188:     x[ir]   = t[ii];  x[ir+1] = t[ii+1]; x[ir+2] = t[ii+2];  x[ir+3] = t[ii+3];
189:     x[ir+4] = t[ii+4];  x[ir+5] = t[ii+5];
190:   }

192:   ISRestoreIndices(isrow,&rout);
193:   ISRestoreIndices(iscol,&cout);
194:   VecRestoreArrayRead(bb,&b);
195:   VecRestoreArray(xx,&x);
196:   PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);
197:   return 0;
198: }