Actual source code: baijsolvnat5.c

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

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

 13:   VecGetArrayRead(bb,&b);
 14:   VecGetArray(xx,&x);
 15:   /* forward solve the lower triangular */
 16:   idx  = 0;
 17:   x[0] = b[idx]; x[1] = b[1+idx]; x[2] = b[2+idx]; x[3] = b[3+idx];x[4] = b[4+idx];
 18:   for (i=1; i<n; i++) {
 19:     v   =  aa + 25*ai[i];
 20:     vi  =  aj + ai[i];
 21:     nz  =  diag[i] - ai[i];
 22:     idx =  5*i;
 23:     s1  =  b[idx];s2 = b[1+idx];s3 = b[2+idx];s4 = b[3+idx];s5 = b[4+idx];
 24:     while (nz--) {
 25:       jdx = 5*(*vi++);
 26:       x1  = x[jdx];x2 = x[1+jdx];x3 = x[2+jdx];x4 = x[3+jdx];x5 = x[4+jdx];
 27:       s1 -= v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
 28:       s2 -= v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
 29:       s3 -= v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
 30:       s4 -= v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
 31:       s5 -= v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
 32:       v  += 25;
 33:     }
 34:     x[idx]   = s1;
 35:     x[1+idx] = s2;
 36:     x[2+idx] = s3;
 37:     x[3+idx] = s4;
 38:     x[4+idx] = s5;
 39:   }
 40:   /* backward solve the upper triangular */
 41:   for (i=n-1; i>=0; i--) {
 42:     v   = aa + 25*diag[i] + 25;
 43:     vi  = aj + diag[i] + 1;
 44:     nz  = ai[i+1] - diag[i] - 1;
 45:     idt = 5*i;
 46:     s1  = x[idt];  s2 = x[1+idt];
 47:     s3  = x[2+idt];s4 = x[3+idt]; s5 = x[4+idt];
 48:     while (nz--) {
 49:       idx = 5*(*vi++);
 50:       x1  = x[idx];   x2 = x[1+idx];x3    = x[2+idx]; x4 = x[3+idx]; x5 = x[4+idx];
 51:       s1 -= v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
 52:       s2 -= v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
 53:       s3 -= v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
 54:       s4 -= v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
 55:       s5 -= v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
 56:       v  += 25;
 57:     }
 58:     v        = aa + 25*diag[i];
 59:     x[idt]   = v[0]*s1 + v[5]*s2 + v[10]*s3  + v[15]*s4 + v[20]*s5;
 60:     x[1+idt] = v[1]*s1 + v[6]*s2 + v[11]*s3  + v[16]*s4 + v[21]*s5;
 61:     x[2+idt] = v[2]*s1 + v[7]*s2 + v[12]*s3  + v[17]*s4 + v[22]*s5;
 62:     x[3+idt] = v[3]*s1 + v[8]*s2 + v[13]*s3  + v[18]*s4 + v[23]*s5;
 63:     x[4+idt] = v[4]*s1 + v[9]*s2 + v[14]*s3  + v[19]*s4 + v[24]*s5;
 64:   }

 66:   VecRestoreArrayRead(bb,&b);
 67:   VecRestoreArray(xx,&x);
 68:   PetscLogFlops(2.0*25*(a->nz) - 5.0*A->cmap->n);
 69:   return 0;
 70: }

 72: PetscErrorCode MatSolve_SeqBAIJ_5_NaturalOrdering(Mat A,Vec bb,Vec xx)
 73: {
 74:   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
 75:   const PetscInt    n  = a->mbs,*vi,*ai=a->i,*aj=a->j,*adiag=a->diag;
 76:   PetscInt          i,k,nz,idx,idt,jdx;
 77:   const MatScalar   *aa=a->a,*v;
 78:   PetscScalar       *x,s1,s2,s3,s4,s5,x1,x2,x3,x4,x5;
 79:   const PetscScalar *b;

 81:   VecGetArrayRead(bb,&b);
 82:   VecGetArray(xx,&x);
 83:   /* forward solve the lower triangular */
 84:   idx  = 0;
 85:   x[0] = b[idx]; x[1] = b[1+idx]; x[2] = b[2+idx]; x[3] = b[3+idx];x[4] = b[4+idx];
 86:   for (i=1; i<n; i++) {
 87:     v   = aa + 25*ai[i];
 88:     vi  = aj + ai[i];
 89:     nz  = ai[i+1] - ai[i];
 90:     idx = 5*i;
 91:     s1  = b[idx];s2 = b[1+idx];s3 = b[2+idx];s4 = b[3+idx];s5 = b[4+idx];
 92:     for (k=0; k<nz; k++) {
 93:       jdx = 5*vi[k];
 94:       x1  = x[jdx];x2 = x[1+jdx];x3 = x[2+jdx];x4 = x[3+jdx];x5 = x[4+jdx];
 95:       s1 -= v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
 96:       s2 -= v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
 97:       s3 -= v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
 98:       s4 -= v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
 99:       s5 -= v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
100:       v  += 25;
101:     }
102:     x[idx]   = s1;
103:     x[1+idx] = s2;
104:     x[2+idx] = s3;
105:     x[3+idx] = s4;
106:     x[4+idx] = s5;
107:   }

109:   /* backward solve the upper triangular */
110:   for (i=n-1; i>=0; i--) {
111:     v   = aa + 25*(adiag[i+1]+1);
112:     vi  = aj + adiag[i+1]+1;
113:     nz  = adiag[i] - adiag[i+1]-1;
114:     idt = 5*i;
115:     s1  = x[idt];  s2 = x[1+idt];
116:     s3  = x[2+idt];s4 = x[3+idt]; s5 = x[4+idt];
117:     for (k=0; k<nz; k++) {
118:       idx = 5*vi[k];
119:       x1  = x[idx];   x2 = x[1+idx];x3    = x[2+idx]; x4 = x[3+idx]; x5 = x[4+idx];
120:       s1 -= v[0]*x1 + v[5]*x2 + v[10]*x3  + v[15]*x4 + v[20]*x5;
121:       s2 -= v[1]*x1 + v[6]*x2 + v[11]*x3  + v[16]*x4 + v[21]*x5;
122:       s3 -= v[2]*x1 + v[7]*x2 + v[12]*x3  + v[17]*x4 + v[22]*x5;
123:       s4 -= v[3]*x1 + v[8]*x2 + v[13]*x3  + v[18]*x4 + v[23]*x5;
124:       s5 -= v[4]*x1 + v[9]*x2 + v[14]*x3  + v[19]*x4 + v[24]*x5;
125:       v  += 25;
126:     }
127:     /* x = inv_diagonal*x */
128:     x[idt]   = v[0]*s1 + v[5]*s2 + v[10]*s3  + v[15]*s4 + v[20]*s5;
129:     x[1+idt] = v[1]*s1 + v[6]*s2 + v[11]*s3  + v[16]*s4 + v[21]*s5;
130:     x[2+idt] = v[2]*s1 + v[7]*s2 + v[12]*s3  + v[17]*s4 + v[22]*s5;
131:     x[3+idt] = v[3]*s1 + v[8]*s2 + v[13]*s3  + v[18]*s4 + v[23]*s5;
132:     x[4+idt] = v[4]*s1 + v[9]*s2 + v[14]*s3  + v[19]*s4 + v[24]*s5;
133:   }

135:   VecRestoreArrayRead(bb,&b);
136:   VecRestoreArray(xx,&x);
137:   PetscLogFlops(2.0*25*(a->nz) - 5.0*A->cmap->n);
138:   return 0;
139: }