Actual source code: baijsolvnat2.c

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

  4: /*
  5:       Special case where the matrix was ILU(0) factored in the natural
  6:    ordering. This eliminates the need for the column and row permutation.
  7: */
  8: PetscErrorCode MatSolve_SeqBAIJ_2_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)
  9: {
 10:   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
 11:   const PetscInt    n  =a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag;
 12:   const MatScalar   *aa=a->a,*v;
 13:   PetscScalar       *x,s1,s2,x1,x2;
 14:   const PetscScalar *b;
 15:   PetscInt          jdx,idt,idx,nz,i;

 17:   VecGetArrayRead(bb,&b);
 18:   VecGetArray(xx,&x);

 20:   /* forward solve the lower triangular */
 21:   idx  = 0;
 22:   x[0] = b[0]; x[1] = b[1];
 23:   for (i=1; i<n; i++) {
 24:     v    =  aa      + 4*ai[i];
 25:     vi   =  aj      + ai[i];
 26:     nz   =  diag[i] - ai[i];
 27:     idx +=  2;
 28:     s1   =  b[idx];s2 = b[1+idx];
 29:     while (nz--) {
 30:       jdx = 2*(*vi++);
 31:       x1  = x[jdx];x2 = x[1+jdx];
 32:       s1 -= v[0]*x1 + v[2]*x2;
 33:       s2 -= v[1]*x1 + v[3]*x2;
 34:       v  += 4;
 35:     }
 36:     x[idx]   = s1;
 37:     x[1+idx] = s2;
 38:   }
 39:   /* backward solve the upper triangular */
 40:   for (i=n-1; i>=0; i--) {
 41:     v   = aa + 4*diag[i] + 4;
 42:     vi  = aj + diag[i] + 1;
 43:     nz  = ai[i+1] - diag[i] - 1;
 44:     idt = 2*i;
 45:     s1  = x[idt];  s2 = x[1+idt];
 46:     while (nz--) {
 47:       idx = 2*(*vi++);
 48:       x1  = x[idx];   x2 = x[1+idx];
 49:       s1 -= v[0]*x1 + v[2]*x2;
 50:       s2 -= v[1]*x1 + v[3]*x2;
 51:       v  += 4;
 52:     }
 53:     v        = aa +  4*diag[i];
 54:     x[idt]   = v[0]*s1 + v[2]*s2;
 55:     x[1+idt] = v[1]*s1 + v[3]*s2;
 56:   }

 58:   VecRestoreArrayRead(bb,&b);
 59:   VecRestoreArray(xx,&x);
 60:   PetscLogFlops(2.0*4*(a->nz) - 2.0*A->cmap->n);
 61:   return 0;
 62: }

 64: PetscErrorCode MatSolve_SeqBAIJ_2_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,*adiag=a->diag;
 68:   PetscInt          i,k,nz,idx,idt,jdx;
 69:   const MatScalar   *aa=a->a,*v;
 70:   PetscScalar       *x,s1,s2,x1,x2;
 71:   const PetscScalar *b;

 73:   VecGetArrayRead(bb,&b);
 74:   VecGetArray(xx,&x);
 75:   /* forward solve the lower triangular */
 76:   idx  = 0;
 77:   x[0] = b[idx]; x[1] = b[1+idx];
 78:   for (i=1; i<n; i++) {
 79:     v   = aa + 4*ai[i];
 80:     vi  = aj + ai[i];
 81:     nz  = ai[i+1] - ai[i];
 82:     idx = 2*i;
 83:     s1  = b[idx];s2 = b[1+idx];
 84:     PetscPrefetchBlock(vi+nz,nz,0,PETSC_PREFETCH_HINT_NTA);
 85:     PetscPrefetchBlock(v+4*nz,4*nz,0,PETSC_PREFETCH_HINT_NTA);
 86:     for (k=0; k<nz; k++) {
 87:       jdx = 2*vi[k];
 88:       x1  = x[jdx];x2 = x[1+jdx];
 89:       s1 -= v[0]*x1 + v[2]*x2;
 90:       s2 -= v[1]*x1 + v[3]*x2;
 91:       v  +=  4;
 92:     }
 93:     x[idx]   = s1;
 94:     x[1+idx] = s2;
 95:   }

 97:   /* backward solve the upper triangular */
 98:   for (i=n-1; i>=0; i--) {
 99:     v   = aa + 4*(adiag[i+1]+1);
100:     vi  = aj + adiag[i+1]+1;
101:     nz  = adiag[i] - adiag[i+1]-1;
102:     idt = 2*i;
103:     s1  = x[idt];  s2 = x[1+idt];
104:     PetscPrefetchBlock(vi+nz,nz,0,PETSC_PREFETCH_HINT_NTA);
105:     PetscPrefetchBlock(v+4*nz,4*nz,0,PETSC_PREFETCH_HINT_NTA);
106:     for (k=0; k<nz; k++) {
107:       idx = 2*vi[k];
108:       x1  = x[idx];   x2 = x[1+idx];
109:       s1 -= v[0]*x1 + v[2]*x2;
110:       s2 -= v[1]*x1 + v[3]*x2;
111:       v  += 4;
112:     }
113:     /* x = inv_diagonal*x */
114:     x[idt]   = v[0]*s1 + v[2]*s2;
115:     x[1+idt] = v[1]*s1 + v[3]*s2;
116:   }

118:   VecRestoreArrayRead(bb,&b);
119:   VecRestoreArray(xx,&x);
120:   PetscLogFlops(2.0*4*(a->nz) - 2.0*A->cmap->n);
121:   return 0;
122: }

124: PetscErrorCode MatForwardSolve_SeqBAIJ_2_NaturalOrdering(Mat A,Vec bb,Vec xx)
125: {
126:   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
127:   const PetscInt    n  = a->mbs,*vi,*ai=a->i,*aj=a->j;
128:   PetscInt          i,k,nz,idx,jdx;
129:   const MatScalar   *aa=a->a,*v;
130:   PetscScalar       *x,s1,s2,x1,x2;
131:   const PetscScalar *b;

133:   VecGetArrayRead(bb,&b);
134:   VecGetArray(xx,&x);
135:   /* forward solve the lower triangular */
136:   idx  = 0;
137:   x[0] = b[idx]; x[1] = b[1+idx];
138:   for (i=1; i<n; i++) {
139:     v   = aa + 4*ai[i];
140:     vi  = aj + ai[i];
141:     nz  = ai[i+1] - ai[i];
142:     idx = 2*i;
143:     s1  = b[idx];s2 = b[1+idx];
144:     PetscPrefetchBlock(vi+nz,nz,0,PETSC_PREFETCH_HINT_NTA);
145:     PetscPrefetchBlock(v+4*nz,4*nz,0,PETSC_PREFETCH_HINT_NTA);
146:     for (k=0; k<nz; k++) {
147:       jdx = 2*vi[k];
148:       x1  = x[jdx];x2 = x[1+jdx];
149:       s1 -= v[0]*x1 + v[2]*x2;
150:       s2 -= v[1]*x1 + v[3]*x2;
151:       v  +=  4;
152:     }
153:     x[idx]   = s1;
154:     x[1+idx] = s2;
155:   }

157:   VecRestoreArrayRead(bb,&b);
158:   VecRestoreArray(xx,&x);
159:   PetscLogFlops(4.0*(a->nz) - A->cmap->n);
160:   return 0;
161: }

163: PetscErrorCode MatBackwardSolve_SeqBAIJ_2_NaturalOrdering(Mat A,Vec bb,Vec xx)
164: {
165:   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
166:   const PetscInt    n  = a->mbs,*vi,*aj=a->j,*adiag=a->diag;
167:   PetscInt          i,k,nz,idx,idt;
168:   const MatScalar   *aa=a->a,*v;
169:   PetscScalar       *x,s1,s2,x1,x2;
170:   const PetscScalar *b;

172:   VecGetArrayRead(bb,&b);
173:   VecGetArray(xx,&x);

175:   /* backward solve the upper triangular */
176:   for (i=n-1; i>=0; i--) {
177:     v   = aa + 4*(adiag[i+1]+1);
178:     vi  = aj + adiag[i+1]+1;
179:     nz  = adiag[i] - adiag[i+1]-1;
180:     idt = 2*i;
181:     s1  = b[idt];  s2 = b[1+idt];
182:     PetscPrefetchBlock(vi+nz,nz,0,PETSC_PREFETCH_HINT_NTA);
183:     PetscPrefetchBlock(v+4*nz,4*nz,0,PETSC_PREFETCH_HINT_NTA);
184:     for (k=0; k<nz; k++) {
185:       idx = 2*vi[k];
186:       x1  = x[idx];   x2 = x[1+idx];
187:       s1 -= v[0]*x1 + v[2]*x2;
188:       s2 -= v[1]*x1 + v[3]*x2;
189:       v  += 4;
190:     }
191:     /* x = inv_diagonal*x */
192:     x[idt]   = v[0]*s1 + v[2]*s2;
193:     x[1+idt] = v[1]*s1 + v[3]*s2;
194:   }

196:   VecRestoreArrayRead(bb,&b);
197:   VecRestoreArray(xx,&x);
198:   PetscLogFlops(4.0*a->nz - A->cmap->n);
199:   return 0;
200: }