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