1 #include <../src/mat/impls/baij/seq/baij.h> 2 3 PetscErrorCode MatSolveTranspose_SeqBAIJ_6_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,s6,x1,x2,x3,x4,x5,x6,*x; 10 11 PetscFunctionBegin; 12 PetscCall(VecCopy(bb,xx)); 13 PetscCall(VecGetArray(xx,&x)); 14 15 /* forward solve the U^T */ 16 idx = 0; 17 for (i=0; i<n; i++) { 18 19 v = aa + 36*diag[i]; 20 /* multiply by the inverse of the block diagonal */ 21 x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx]; x4 = x[3+idx]; x5 = x[4+idx]; 22 x6 = x[5+idx]; 23 s1 = v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5 + v[5]*x6; 24 s2 = v[6]*x1 + v[7]*x2 + v[8]*x3 + v[9]*x4 + v[10]*x5 + v[11]*x6; 25 s3 = v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4 + v[16]*x5 + v[17]*x6; 26 s4 = v[18]*x1 + v[19]*x2 + v[20]*x3 + v[21]*x4 + v[22]*x5 + v[23]*x6; 27 s5 = v[24]*x1 + v[25]*x2 + v[26]*x3 + v[27]*x4 + v[28]*x5 + v[29]*x6; 28 s6 = v[30]*x1 + v[31]*x2 + v[32]*x3 + v[33]*x4 + v[34]*x5 + v[35]*x6; 29 v += 36; 30 31 vi = aj + diag[i] + 1; 32 nz = ai[i+1] - diag[i] - 1; 33 while (nz--) { 34 oidx = 6*(*vi++); 35 x[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5 + v[5]*s6; 36 x[oidx+1] -= v[6]*s1 + v[7]*s2 + v[8]*s3 + v[9]*s4 + v[10]*s5 + v[11]*s6; 37 x[oidx+2] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4 + v[16]*s5 + v[17]*s6; 38 x[oidx+3] -= v[18]*s1 + v[19]*s2 + v[20]*s3 + v[21]*s4 + v[22]*s5 + v[23]*s6; 39 x[oidx+4] -= v[24]*s1 + v[25]*s2 + v[26]*s3 + v[27]*s4 + v[28]*s5 + v[29]*s6; 40 x[oidx+5] -= v[30]*s1 + v[31]*s2 + v[32]*s3 + v[33]*s4 + v[34]*s5 + v[35]*s6; 41 v += 36; 42 } 43 x[idx] = s1;x[1+idx] = s2; x[2+idx] = s3;x[3+idx] = s4; x[4+idx] = s5; 44 x[5+idx] = s6; 45 idx += 6; 46 } 47 /* backward solve the L^T */ 48 for (i=n-1; i>=0; i--) { 49 v = aa + 36*diag[i] - 36; 50 vi = aj + diag[i] - 1; 51 nz = diag[i] - ai[i]; 52 idt = 6*i; 53 s1 = x[idt]; s2 = x[1+idt]; s3 = x[2+idt];s4 = x[3+idt]; s5 = x[4+idt]; 54 s6 = x[5+idt]; 55 while (nz--) { 56 idx = 6*(*vi--); 57 x[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5 + v[5]*s6; 58 x[idx+1] -= v[6]*s1 + v[7]*s2 + v[8]*s3 + v[9]*s4 + v[10]*s5 + v[11]*s6; 59 x[idx+2] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4 + v[16]*s5 + v[17]*s6; 60 x[idx+3] -= v[18]*s1 + v[19]*s2 + v[20]*s3 + v[21]*s4 + v[22]*s5 + v[23]*s6; 61 x[idx+4] -= v[24]*s1 + v[25]*s2 + v[26]*s3 + v[27]*s4 + v[28]*s5 + v[29]*s6; 62 x[idx+5] -= v[30]*s1 + v[31]*s2 + v[32]*s3 + v[33]*s4 + v[34]*s5 + v[35]*s6; 63 v -= 36; 64 } 65 } 66 PetscCall(VecRestoreArray(xx,&x)); 67 PetscCall(PetscLogFlops(2.0*36*(a->nz) - 6.0*A->cmap->n)); 68 PetscFunctionReturn(0); 69 } 70 71 PetscErrorCode MatSolveTranspose_SeqBAIJ_6_NaturalOrdering(Mat A,Vec bb,Vec xx) 72 { 73 Mat_SeqBAIJ *a=(Mat_SeqBAIJ*)A->data; 74 const PetscInt n=a->mbs,*vi,*ai=a->i,*aj=a->j,*diag=a->diag; 75 PetscInt nz,idx,idt,j,i,oidx; 76 const PetscInt bs =A->rmap->bs,bs2=a->bs2; 77 const MatScalar *aa=a->a,*v; 78 PetscScalar s1,s2,s3,s4,s5,s6,x1,x2,x3,x4,x5,x6,*x; 79 80 PetscFunctionBegin; 81 PetscCall(VecCopy(bb,xx)); 82 PetscCall(VecGetArray(xx,&x)); 83 84 /* forward solve the U^T */ 85 idx = 0; 86 for (i=0; i<n; i++) { 87 v = aa + bs2*diag[i]; 88 /* multiply by the inverse of the block diagonal */ 89 x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx]; x4 = x[3+idx]; 90 x5 = x[4+idx]; x6 = x[5+idx]; 91 s1 = v[0]*x1 + v[1]*x2 + v[2]*x3 + v[3]*x4 + v[4]*x5 + v[5]*x6; 92 s2 = v[6]*x1 + v[7]*x2 + v[8]*x3 + v[9]*x4 + v[10]*x5 + v[11]*x6; 93 s3 = v[12]*x1 + v[13]*x2 + v[14]*x3 + v[15]*x4 + v[16]*x5 + v[17]*x6; 94 s4 = v[18]*x1 + v[19]*x2 + v[20]*x3 + v[21]*x4 + v[22]*x5 + v[23]*x6; 95 s5 = v[24]*x1 + v[25]*x2 + v[26]*x3 + v[27]*x4 + v[28]*x5 + v[29]*x6; 96 s6 = v[30]*x1 + v[31]*x2 + v[32]*x3 + v[33]*x4 + v[34]*x5 + v[35]*x6; 97 v -= bs2; 98 99 vi = aj + diag[i] - 1; 100 nz = diag[i] - diag[i+1] - 1; 101 for (j=0; j>-nz; j--) { 102 oidx = bs*vi[j]; 103 x[oidx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5 + v[5]*s6; 104 x[oidx+1] -= v[6]*s1 + v[7]*s2 + v[8]*s3 + v[9]*s4 + v[10]*s5 + v[11]*s6; 105 x[oidx+2] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4 + v[16]*s5 + v[17]*s6; 106 x[oidx+3] -= v[18]*s1 + v[19]*s2 + v[20]*s3 + v[21]*s4 + v[22]*s5 + v[23]*s6; 107 x[oidx+4] -= v[24]*s1 + v[25]*s2 + v[26]*s3 + v[27]*s4 + v[28]*s5 + v[29]*s6; 108 x[oidx+5] -= v[30]*s1 + v[31]*s2 + v[32]*s3 + v[33]*s4 + v[34]*s5 + v[35]*s6; 109 v -= bs2; 110 } 111 x[idx] = s1;x[1+idx] = s2; x[2+idx] = s3; x[3+idx] = s4; x[4+idx] = s5; 112 x[5+idx] = s6; 113 idx += bs; 114 } 115 /* backward solve the L^T */ 116 for (i=n-1; i>=0; i--) { 117 v = aa + bs2*ai[i]; 118 vi = aj + ai[i]; 119 nz = ai[i+1] - ai[i]; 120 idt = bs*i; 121 s1 = x[idt]; s2 = x[1+idt]; s3 = x[2+idt]; s4 = x[3+idt]; s5 = x[4+idt]; 122 s6 = x[5+idt]; 123 for (j=0; j<nz; j++) { 124 idx = bs*vi[j]; 125 x[idx] -= v[0]*s1 + v[1]*s2 + v[2]*s3 + v[3]*s4 + v[4]*s5 + v[5]*s6; 126 x[idx+1] -= v[6]*s1 + v[7]*s2 + v[8]*s3 + v[9]*s4 + v[10]*s5 + v[11]*s6; 127 x[idx+2] -= v[12]*s1 + v[13]*s2 + v[14]*s3 + v[15]*s4 + v[16]*s5 + v[17]*s6; 128 x[idx+3] -= v[18]*s1 + v[19]*s2 + v[20]*s3 + v[21]*s4 + v[22]*s5 + v[23]*s6; 129 x[idx+4] -= v[24]*s1 + v[25]*s2 + v[26]*s3 + v[27]*s4 + v[28]*s5 + v[29]*s6; 130 x[idx+5] -= v[30]*s1 + v[31]*s2 + v[32]*s3 + v[33]*s4 + v[34]*s5 + v[35]*s6; 131 v += bs2; 132 } 133 } 134 PetscCall(VecRestoreArray(xx,&x)); 135 PetscCall(PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n)); 136 PetscFunctionReturn(0); 137 } 138