Actual source code: baijsolvtrannat5.c

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

  3: PetscErrorCode MatSolveTranspose_SeqBAIJ_5_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,x1,x2,x3,x4,x5,*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++) {

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

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

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

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

 76:   /* forward solve the U^T */
 77:   idx = 0;
 78:   for (i=0; i<n; i++) {
 79:     v = aa + bs2*diag[i];
 80:     /* multiply by the inverse of the block diagonal */
 81:     x1 = x[idx];   x2 = x[1+idx];  x3 = x[2+idx];  x4 = x[3+idx];
 82:     x5 = x[4+idx];
 83:     s1 =  v[0]*x1   +  v[1]*x2   + v[2]*x3   + v[3]*x4   + v[4]*x5;
 84:     s2 =  v[5]*x1   +  v[6]*x2   + v[7]*x3   + v[8]*x4   + v[9]*x5;
 85:     s3 =  v[10]*x1  +  v[11]*x2  + v[12]*x3  + v[13]*x4  + v[14]*x5;
 86:     s4 =  v[15]*x1  +  v[16]*x2  + v[17]*x3  + v[18]*x4  + v[19]*x5;
 87:     s5 =  v[20]*x1  +  v[21]*x2  + v[22]*x3  + v[23]*x4   + v[24]*x5;
 88:     v -= bs2;

 90:     vi = aj + diag[i] - 1;
 91:     nz = diag[i] - diag[i+1] - 1;
 92:     for (j=0; j>-nz; j--) {
 93:       oidx       = bs*vi[j];
 94:       x[oidx]   -=  v[0]*s1   +  v[1]*s2   + v[2]*s3   + v[3]*s4   + v[4]*s5;
 95:       x[oidx+1] -=  v[5]*s1   +  v[6]*s2   + v[7]*s3   + v[8]*s4   + v[9]*s5;
 96:       x[oidx+2] -=  v[10]*s1  +  v[11]*s2  + v[12]*s3  + v[13]*s4  + v[14]*s5;
 97:       x[oidx+3] -=  v[15]*s1  +  v[16]*s2  + v[17]*s3  + v[18]*s4  + v[19]*s5;
 98:       x[oidx+4] -=  v[20]*s1  +  v[21]*s2  + v[22]*s3  + v[23]*s4   + v[24]*s5;
 99:       v         -= bs2;
100:     }
101:     x[idx] = s1;x[1+idx] = s2;  x[2+idx] = s3;  x[3+idx] = s4; x[4+idx] = s5;
102:     idx   += bs;
103:   }
104:   /* backward solve the L^T */
105:   for (i=n-1; i>=0; i--) {
106:     v   = aa + bs2*ai[i];
107:     vi  = aj + ai[i];
108:     nz  = ai[i+1] - ai[i];
109:     idt = bs*i;
110:     s1  = x[idt];  s2 = x[1+idt];  s3 = x[2+idt];  s4 = x[3+idt];  s5 = x[4+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   + v[3]*s4   + v[4]*s5;
114:       x[idx+1] -=  v[5]*s1   +  v[6]*s2   + v[7]*s3   + v[8]*s4   + v[9]*s5;
115:       x[idx+2] -=  v[10]*s1  +  v[11]*s2  + v[12]*s3  + v[13]*s4  + v[14]*s5;
116:       x[idx+3] -=  v[15]*s1  +  v[16]*s2  + v[17]*s3  + v[18]*s4  + v[19]*s5;
117:       x[idx+4] -=  v[20]*s1  +  v[21]*s2  + v[22]*s3  + v[23]*s4   + v[24]*s5;
118:       v        += bs2;
119:     }
120:   }
121:   VecRestoreArray(xx,&x);
122:   PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);
123:   return 0;
124: }