xref: /petsc/src/mat/impls/baij/seq/baijsolvnat6.c (revision feb237ba1401b79cfa9a99e81307a265c4a68462)
12c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h>
22c733ed4SBarry Smith #include <petsc/private/kernels/blockinvert.h>
32c733ed4SBarry Smith 
42c733ed4SBarry Smith 
52c733ed4SBarry Smith PetscErrorCode MatSolve_SeqBAIJ_6_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx)
62c733ed4SBarry Smith {
72c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
82c733ed4SBarry Smith   PetscInt          i,nz,idx,idt,jdx;
92c733ed4SBarry Smith   PetscErrorCode    ierr;
102c733ed4SBarry Smith   const PetscInt    *diag = a->diag,*vi,n=a->mbs,*ai=a->i,*aj=a->j;
112c733ed4SBarry Smith   const MatScalar   *aa   =a->a,*v;
122c733ed4SBarry Smith   PetscScalar       *x,s1,s2,s3,s4,s5,s6,x1,x2,x3,x4,x5,x6;
132c733ed4SBarry Smith   const PetscScalar *b;
142c733ed4SBarry Smith 
152c733ed4SBarry Smith   PetscFunctionBegin;
162c733ed4SBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
172c733ed4SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
182c733ed4SBarry Smith   /* forward solve the lower triangular */
192c733ed4SBarry Smith   idx  = 0;
202c733ed4SBarry Smith   x[0] = b[idx];   x[1] = b[1+idx]; x[2] = b[2+idx];
212c733ed4SBarry Smith   x[3] = b[3+idx]; x[4] = b[4+idx]; x[5] = b[5+idx];
222c733ed4SBarry Smith   for (i=1; i<n; i++) {
232c733ed4SBarry Smith     v   =  aa + 36*ai[i];
242c733ed4SBarry Smith     vi  =  aj + ai[i];
252c733ed4SBarry Smith     nz  =  diag[i] - ai[i];
262c733ed4SBarry Smith     idx =  6*i;
272c733ed4SBarry Smith     s1  =  b[idx];   s2 = b[1+idx]; s3 = b[2+idx];
282c733ed4SBarry Smith     s4  =  b[3+idx]; s5 = b[4+idx]; s6 = b[5+idx];
292c733ed4SBarry Smith     while (nz--) {
302c733ed4SBarry Smith       jdx = 6*(*vi++);
312c733ed4SBarry Smith       x1  = x[jdx];   x2 = x[1+jdx]; x3 = x[2+jdx];
322c733ed4SBarry Smith       x4  = x[3+jdx]; x5 = x[4+jdx]; x6 = x[5+jdx];
332c733ed4SBarry Smith       s1 -= v[0]*x1 + v[6]*x2  + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6;
342c733ed4SBarry Smith       s2 -= v[1]*x1 + v[7]*x2  + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6;
352c733ed4SBarry Smith       s3 -= v[2]*x1 + v[8]*x2  + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6;
362c733ed4SBarry Smith       s4 -= v[3]*x1 + v[9]*x2  + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6;
372c733ed4SBarry Smith       s5 -= v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6;
382c733ed4SBarry Smith       s6 -= v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6;
392c733ed4SBarry Smith       v  += 36;
402c733ed4SBarry Smith     }
412c733ed4SBarry Smith     x[idx]   = s1;
422c733ed4SBarry Smith     x[1+idx] = s2;
432c733ed4SBarry Smith     x[2+idx] = s3;
442c733ed4SBarry Smith     x[3+idx] = s4;
452c733ed4SBarry Smith     x[4+idx] = s5;
462c733ed4SBarry Smith     x[5+idx] = s6;
472c733ed4SBarry Smith   }
482c733ed4SBarry Smith   /* backward solve the upper triangular */
492c733ed4SBarry Smith   for (i=n-1; i>=0; i--) {
502c733ed4SBarry Smith     v   = aa + 36*diag[i] + 36;
512c733ed4SBarry Smith     vi  = aj + diag[i] + 1;
522c733ed4SBarry Smith     nz  = ai[i+1] - diag[i] - 1;
532c733ed4SBarry Smith     idt = 6*i;
542c733ed4SBarry Smith     s1  = x[idt];   s2 = x[1+idt];
552c733ed4SBarry Smith     s3  = x[2+idt]; s4 = x[3+idt];
562c733ed4SBarry Smith     s5  = x[4+idt]; s6 = x[5+idt];
572c733ed4SBarry Smith     while (nz--) {
582c733ed4SBarry Smith       idx = 6*(*vi++);
592c733ed4SBarry Smith       x1  = x[idx];   x2 = x[1+idx]; x3 = x[2+idx];
602c733ed4SBarry Smith       x4  = x[3+idx]; x5 = x[4+idx]; x6 = x[5+idx];
612c733ed4SBarry Smith       s1 -= v[0]*x1 + v[6]*x2  + v[12]*x3 + v[18]*x4 + v[24]*x5 + v[30]*x6;
622c733ed4SBarry Smith       s2 -= v[1]*x1 + v[7]*x2  + v[13]*x3 + v[19]*x4 + v[25]*x5 + v[31]*x6;
632c733ed4SBarry Smith       s3 -= v[2]*x1 + v[8]*x2  + v[14]*x3 + v[20]*x4 + v[26]*x5 + v[32]*x6;
642c733ed4SBarry Smith       s4 -= v[3]*x1 + v[9]*x2  + v[15]*x3 + v[21]*x4 + v[27]*x5 + v[33]*x6;
652c733ed4SBarry Smith       s5 -= v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4 + v[28]*x5 + v[34]*x6;
662c733ed4SBarry Smith       s6 -= v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4 + v[29]*x5 + v[35]*x6;
672c733ed4SBarry Smith       v  += 36;
682c733ed4SBarry Smith     }
692c733ed4SBarry Smith     v        = aa + 36*diag[i];
702c733ed4SBarry Smith     x[idt]   = v[0]*s1 + v[6]*s2  + v[12]*s3 + v[18]*s4 + v[24]*s5 + v[30]*s6;
712c733ed4SBarry Smith     x[1+idt] = v[1]*s1 + v[7]*s2  + v[13]*s3 + v[19]*s4 + v[25]*s5 + v[31]*s6;
722c733ed4SBarry Smith     x[2+idt] = v[2]*s1 + v[8]*s2  + v[14]*s3 + v[20]*s4 + v[26]*s5 + v[32]*s6;
732c733ed4SBarry Smith     x[3+idt] = v[3]*s1 + v[9]*s2  + v[15]*s3 + v[21]*s4 + v[27]*s5 + v[33]*s6;
742c733ed4SBarry Smith     x[4+idt] = v[4]*s1 + v[10]*s2 + v[16]*s3 + v[22]*s4 + v[28]*s5 + v[34]*s6;
752c733ed4SBarry Smith     x[5+idt] = v[5]*s1 + v[11]*s2 + v[17]*s3 + v[23]*s4 + v[29]*s5 + v[35]*s6;
762c733ed4SBarry Smith   }
772c733ed4SBarry Smith 
782c733ed4SBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
792c733ed4SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
802c733ed4SBarry Smith   ierr = PetscLogFlops(2.0*36*(a->nz) - 6.0*A->cmap->n);CHKERRQ(ierr);
812c733ed4SBarry Smith   PetscFunctionReturn(0);
822c733ed4SBarry Smith }
832c733ed4SBarry Smith 
842c733ed4SBarry Smith PetscErrorCode MatSolve_SeqBAIJ_6_NaturalOrdering(Mat A,Vec bb,Vec xx)
852c733ed4SBarry Smith {
862c733ed4SBarry Smith   Mat_SeqBAIJ       *a = (Mat_SeqBAIJ*)A->data;
872c733ed4SBarry Smith   const PetscInt    n  =a->mbs,*vi,*ai=a->i,*aj=a->j,*adiag=a->diag;
882c733ed4SBarry Smith   PetscErrorCode    ierr;
892c733ed4SBarry Smith   PetscInt          i,k,nz,idx,jdx,idt;
902c733ed4SBarry Smith   const PetscInt    bs = A->rmap->bs,bs2 = a->bs2;
912c733ed4SBarry Smith   const MatScalar   *aa=a->a,*v;
922c733ed4SBarry Smith   PetscScalar       *x;
932c733ed4SBarry Smith   const PetscScalar *b;
942c733ed4SBarry Smith   PetscScalar       s1,s2,s3,s4,s5,s6,x1,x2,x3,x4,x5,x6;
952c733ed4SBarry Smith 
962c733ed4SBarry Smith   PetscFunctionBegin;
972c733ed4SBarry Smith   ierr = VecGetArrayRead(bb,&b);CHKERRQ(ierr);
982c733ed4SBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
992c733ed4SBarry Smith   /* forward solve the lower triangular */
1002c733ed4SBarry Smith   idx  = 0;
1012c733ed4SBarry Smith   x[0] = b[idx]; x[1] = b[1+idx];x[2] = b[2+idx];x[3] = b[3+idx];
1022c733ed4SBarry Smith   x[4] = b[4+idx];x[5] = b[5+idx];
1032c733ed4SBarry Smith   for (i=1; i<n; i++) {
1042c733ed4SBarry Smith     v   = aa + bs2*ai[i];
1052c733ed4SBarry Smith     vi  = aj + ai[i];
1062c733ed4SBarry Smith     nz  = ai[i+1] - ai[i];
1072c733ed4SBarry Smith     idx = bs*i;
1082c733ed4SBarry Smith     s1  = b[idx];s2 = b[1+idx];s3 = b[2+idx];s4 = b[3+idx];
1092c733ed4SBarry Smith     s5  = b[4+idx];s6 = b[5+idx];
1102c733ed4SBarry Smith     for (k=0; k<nz; k++) {
1112c733ed4SBarry Smith       jdx = bs*vi[k];
1122c733ed4SBarry Smith       x1  = x[jdx];x2 = x[1+jdx]; x3 =x[2+jdx];x4 =x[3+jdx];
1132c733ed4SBarry Smith       x5  = x[4+jdx]; x6 = x[5+jdx];
1142c733ed4SBarry Smith       s1 -= v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4  + v[24]*x5 + v[30]*x6;
115*feb237baSPierre Jolivet       s2 -= v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4  + v[25]*x5 + v[31]*x6;
1162c733ed4SBarry Smith       s3 -= v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4  + v[26]*x5 + v[32]*x6;
1172c733ed4SBarry Smith       s4 -= v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4  + v[27]*x5 + v[33]*x6;
1182c733ed4SBarry Smith       s5 -= v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4  + v[28]*x5 + v[34]*x6;
1192c733ed4SBarry Smith       s6 -= v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4  + v[29]*x5 + v[35]*x6;
1202c733ed4SBarry Smith       v  +=  bs2;
1212c733ed4SBarry Smith     }
1222c733ed4SBarry Smith 
1232c733ed4SBarry Smith     x[idx]   = s1;
1242c733ed4SBarry Smith     x[1+idx] = s2;
1252c733ed4SBarry Smith     x[2+idx] = s3;
1262c733ed4SBarry Smith     x[3+idx] = s4;
1272c733ed4SBarry Smith     x[4+idx] = s5;
1282c733ed4SBarry Smith     x[5+idx] = s6;
1292c733ed4SBarry Smith   }
1302c733ed4SBarry Smith 
1312c733ed4SBarry Smith   /* backward solve the upper triangular */
1322c733ed4SBarry Smith   for (i=n-1; i>=0; i--) {
1332c733ed4SBarry Smith     v   = aa + bs2*(adiag[i+1]+1);
1342c733ed4SBarry Smith     vi  = aj + adiag[i+1]+1;
1352c733ed4SBarry Smith     nz  = adiag[i] - adiag[i+1]-1;
1362c733ed4SBarry Smith     idt = bs*i;
1372c733ed4SBarry Smith     s1  = x[idt];  s2 = x[1+idt];s3 = x[2+idt];s4 = x[3+idt];
1382c733ed4SBarry Smith     s5  = x[4+idt];s6 = x[5+idt];
1392c733ed4SBarry Smith     for (k=0; k<nz; k++) {
1402c733ed4SBarry Smith       idx = bs*vi[k];
1412c733ed4SBarry Smith       x1  = x[idx];   x2 = x[1+idx]; x3 = x[2+idx];x4 = x[3+idx];
1422c733ed4SBarry Smith       x5  = x[4+idx];x6 = x[5+idx];
1432c733ed4SBarry Smith       s1 -= v[0]*x1 + v[6]*x2 + v[12]*x3 + v[18]*x4  + v[24]*x5 + v[30]*x6;
144*feb237baSPierre Jolivet       s2 -= v[1]*x1 + v[7]*x2 + v[13]*x3 + v[19]*x4  + v[25]*x5 + v[31]*x6;
1452c733ed4SBarry Smith       s3 -= v[2]*x1 + v[8]*x2 + v[14]*x3 + v[20]*x4  + v[26]*x5 + v[32]*x6;
1462c733ed4SBarry Smith       s4 -= v[3]*x1 + v[9]*x2 + v[15]*x3 + v[21]*x4  + v[27]*x5 + v[33]*x6;
1472c733ed4SBarry Smith       s5 -= v[4]*x1 + v[10]*x2 + v[16]*x3 + v[22]*x4  + v[28]*x5 + v[34]*x6;
1482c733ed4SBarry Smith       s6 -= v[5]*x1 + v[11]*x2 + v[17]*x3 + v[23]*x4  + v[29]*x5 + v[35]*x6;
1492c733ed4SBarry Smith       v  +=  bs2;
1502c733ed4SBarry Smith     }
1512c733ed4SBarry Smith     /* x = inv_diagonal*x */
1522c733ed4SBarry Smith     x[idt]   = v[0]*s1 + v[6]*s2 + v[12]*s3 + v[18]*s4 + v[24]*s5 + v[30]*s6;
1532c733ed4SBarry Smith     x[1+idt] = v[1]*s1 + v[7]*s2 + v[13]*s3 + v[19]*s4 + v[25]*s5 + v[31]*s6;
1542c733ed4SBarry Smith     x[2+idt] = v[2]*s1 + v[8]*s2 + v[14]*s3 + v[20]*s4 + v[26]*s5 + v[32]*s6;
1552c733ed4SBarry Smith     x[3+idt] = v[3]*s1 + v[9]*s2 + v[15]*s3 + v[21]*s4 + v[27]*s5 + v[33]*s6;
1562c733ed4SBarry Smith     x[4+idt] = v[4]*s1 + v[10]*s2 + v[16]*s3 + v[22]*s4 + v[28]*s5 + v[34]*s6;
1572c733ed4SBarry Smith     x[5+idt] = v[5]*s1 + v[11]*s2 + v[17]*s3 + v[23]*s4 + v[29]*s5 + v[35]*s6;
1582c733ed4SBarry Smith   }
1592c733ed4SBarry Smith 
1602c733ed4SBarry Smith   ierr = VecRestoreArrayRead(bb,&b);CHKERRQ(ierr);
1612c733ed4SBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
1622c733ed4SBarry Smith   ierr = PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n);CHKERRQ(ierr);
1632c733ed4SBarry Smith   PetscFunctionReturn(0);
1642c733ed4SBarry Smith }
165