Actual source code: baijsolvtran5.c

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

  4: PetscErrorCode MatSolveTranspose_SeqBAIJ_5_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,x1,x2,x3,x4,x5,*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      = 5*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:     ii     += 5;
 32:   }

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

 38:     v = aa + 25*diag[i];
 39:     /* multiply by the inverse of the block diagonal */
 40:     x1 = t[idx];   x2 = t[1+idx]; x3    = t[2+idx]; x4 = t[3+idx]; x5 = t[4+idx];
 41:     s1 = v[0]*x1  +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
 42:     s2 = v[5]*x1  +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
 43:     s3 = v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
 44:     s4 = v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
 45:     s5 = v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
 46:     v += 25;

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

 80:   /* copy t into x according to permutation */
 81:   ii = 0;
 82:   for (i=0; i<n; i++) {
 83:     ir      = 5*r[i];
 84:     x[ir]   = t[ii];
 85:     x[ir+1] = t[ii+1];
 86:     x[ir+2] = t[ii+2];
 87:     x[ir+3] = t[ii+3];
 88:     x[ir+4] = t[ii+4];
 89:     ii     += 5;
 90:   }

 92:   ISRestoreIndices(isrow,&rout);
 93:   ISRestoreIndices(iscol,&cout);
 94:   VecRestoreArrayRead(bb,&b);
 95:   VecRestoreArray(xx,&x);
 96:   PetscLogFlops(2.0*25*(a->nz) - 5.0*A->cmap->n);
 97:   return 0;
 98: }

100: PetscErrorCode MatSolveTranspose_SeqBAIJ_5(Mat A,Vec bb,Vec xx)
101: {
102:   Mat_SeqBAIJ       *a=(Mat_SeqBAIJ*)A->data;
103:   IS                iscol=a->col,isrow=a->row;
104:   const PetscInt    n    =a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag;
105:   const PetscInt    *r,*c,*rout,*cout;
106:   PetscInt          nz,idx,idt,j,i,oidx,ii,ic,ir;
107:   const PetscInt    bs =A->rmap->bs,bs2=a->bs2;
108:   const MatScalar   *aa=a->a,*v;
109:   PetscScalar       s1,s2,s3,s4,s5,x1,x2,x3,x4,x5,*x,*t;
110:   const PetscScalar *b;

112:   VecGetArrayRead(bb,&b);
113:   VecGetArray(xx,&x);
114:   t    = a->solve_work;

116:   ISGetIndices(isrow,&rout); r = rout;
117:   ISGetIndices(iscol,&cout); c = cout;

119:   /* copy b into temp work space according to permutation */
120:   for (i=0; i<n; i++) {
121:     ii      = bs*i; ic = bs*c[i];
122:     t[ii]   = b[ic]; t[ii+1] = b[ic+1]; t[ii+2] = b[ic+2]; t[ii+3] = b[ic+3];
123:     t[ii+4] = b[ic+4];
124:   }

126:   /* forward solve the U^T */
127:   idx = 0;
128:   for (i=0; i<n; i++) {
129:     v = aa + bs2*diag[i];
130:     /* multiply by the inverse of the block diagonal */
131:     x1 = t[idx];   x2 = t[1+idx]; x3    = t[2+idx]; x4 = t[3+idx]; x5 = t[4+idx];
132:     s1 = v[0]*x1  +  v[1]*x2 +  v[2]*x3 +  v[3]*x4 +  v[4]*x5;
133:     s2 = v[5]*x1  +  v[6]*x2 +  v[7]*x3 +  v[8]*x4 +  v[9]*x5;
134:     s3 = v[10]*x1 + v[11]*x2 + v[12]*x3 + v[13]*x4 + v[14]*x5;
135:     s4 = v[15]*x1 + v[16]*x2 + v[17]*x3 + v[18]*x4 + v[19]*x5;
136:     s5 = v[20]*x1 + v[21]*x2 + v[22]*x3 + v[23]*x4 + v[24]*x5;
137:     v -= bs2;

139:     vi = aj + diag[i] - 1;
140:     nz = diag[i] - diag[i+1] - 1;
141:     for (j=0; j>-nz; j--) {
142:       oidx       = bs*vi[j];
143:       t[oidx]   -= v[0]*s1  +  v[1]*s2 +  v[2]*s3 +  v[3]*s4 +  v[4]*s5;
144:       t[oidx+1] -= v[5]*s1  +  v[6]*s2 +  v[7]*s3 +  v[8]*s4 +  v[9]*s5;
145:       t[oidx+2] -= v[10]*s1 + v[11]*s2 + v[12]*s3 + v[13]*s4 + v[14]*s5;
146:       t[oidx+3] -= v[15]*s1 + v[16]*s2 + v[17]*s3 + v[18]*s4 + v[19]*s5;
147:       t[oidx+4] -= v[20]*s1 + v[21]*s2 + v[22]*s3 + v[23]*s4 + v[24]*s5;
148:       v         -= bs2;
149:     }
150:     t[idx] = s1;t[1+idx] = s2;  t[2+idx] = s3;  t[3+idx] = s4; t[4+idx] =s5;
151:     idx   += bs;
152:   }
153:   /* backward solve the L^T */
154:   for (i=n-1; i>=0; i--) {
155:     v   = aa + bs2*ai[i];
156:     vi  = aj + ai[i];
157:     nz  = ai[i+1] - ai[i];
158:     idt = bs*i;
159:     s1  = t[idt];  s2 = t[1+idt];  s3 = t[2+idt];  s4 = t[3+idt]; s5 = t[4+idt];
160:     for (j=0; j<nz; j++) {
161:       idx       = bs*vi[j];
162:       t[idx]   -= v[0]*s1  +  v[1]*s2 +  v[2]*s3 +  v[3]*s4 +  v[4]*s5;
163:       t[idx+1] -= v[5]*s1  +  v[6]*s2 +  v[7]*s3 +  v[8]*s4 +  v[9]*s5;
164:       t[idx+2] -= v[10]*s1 + v[11]*s2 + v[12]*s3 + v[13]*s4 + v[14]*s5;
165:       t[idx+3] -= v[15]*s1 + v[16]*s2 + v[17]*s3 + v[18]*s4 + v[19]*s5;
166:       t[idx+4] -= v[20]*s1 + v[21]*s2 + v[22]*s3 + v[23]*s4 + v[24]*s5;
167:       v        += bs2;
168:     }
169:   }

171:   /* copy t into x according to permutation */
172:   for (i=0; i<n; i++) {
173:     ii      = bs*i;  ir = bs*r[i];
174:     x[ir]   = t[ii];  x[ir+1] = t[ii+1]; x[ir+2] = t[ii+2];  x[ir+3] = t[ii+3];
175:     x[ir+4] = t[ii+4];
176:   }

178:   ISRestoreIndices(isrow,&rout);
179:   ISRestoreIndices(iscol,&cout);
180:   VecRestoreArrayRead(bb,&b);
181:   VecRestoreArray(xx,&x);
182:   PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);
183:   return 0;
184: }