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