1 #include <../src/mat/impls/baij/seq/baij.h> 2 #include <petsc/private/kernels/blockinvert.h> 3 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_3_NaturalOrdering_inplace(Mat A,Vec bb,Vec xx) 9 { 10 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 11 const PetscInt n =a->mbs,*ai=a->i,*aj=a->j; 12 const PetscInt *diag = a->diag,*vi; 13 const MatScalar *aa =a->a,*v; 14 PetscScalar *x,s1,s2,s3,x1,x2,x3; 15 const PetscScalar *b; 16 PetscInt jdx,idt,idx,nz,i; 17 18 PetscFunctionBegin; 19 PetscCall(VecGetArrayRead(bb,&b)); 20 PetscCall(VecGetArray(xx,&x)); 21 22 /* forward solve the lower triangular */ 23 idx = 0; 24 x[0] = b[0]; x[1] = b[1]; x[2] = b[2]; 25 for (i=1; i<n; i++) { 26 v = aa + 9*ai[i]; 27 vi = aj + ai[i]; 28 nz = diag[i] - ai[i]; 29 idx += 3; 30 s1 = b[idx];s2 = b[1+idx];s3 = b[2+idx]; 31 while (nz--) { 32 jdx = 3*(*vi++); 33 x1 = x[jdx];x2 = x[1+jdx];x3 = x[2+jdx]; 34 s1 -= v[0]*x1 + v[3]*x2 + v[6]*x3; 35 s2 -= v[1]*x1 + v[4]*x2 + v[7]*x3; 36 s3 -= v[2]*x1 + v[5]*x2 + v[8]*x3; 37 v += 9; 38 } 39 x[idx] = s1; 40 x[1+idx] = s2; 41 x[2+idx] = s3; 42 } 43 /* backward solve the upper triangular */ 44 for (i=n-1; i>=0; i--) { 45 v = aa + 9*diag[i] + 9; 46 vi = aj + diag[i] + 1; 47 nz = ai[i+1] - diag[i] - 1; 48 idt = 3*i; 49 s1 = x[idt]; s2 = x[1+idt]; 50 s3 = x[2+idt]; 51 while (nz--) { 52 idx = 3*(*vi++); 53 x1 = x[idx]; x2 = x[1+idx];x3 = x[2+idx]; 54 s1 -= v[0]*x1 + v[3]*x2 + v[6]*x3; 55 s2 -= v[1]*x1 + v[4]*x2 + v[7]*x3; 56 s3 -= v[2]*x1 + v[5]*x2 + v[8]*x3; 57 v += 9; 58 } 59 v = aa + 9*diag[i]; 60 x[idt] = v[0]*s1 + v[3]*s2 + v[6]*s3; 61 x[1+idt] = v[1]*s1 + v[4]*s2 + v[7]*s3; 62 x[2+idt] = v[2]*s1 + v[5]*s2 + v[8]*s3; 63 } 64 65 PetscCall(VecRestoreArrayRead(bb,&b)); 66 PetscCall(VecRestoreArray(xx,&x)); 67 PetscCall(PetscLogFlops(2.0*9*(a->nz) - 3.0*A->cmap->n)); 68 PetscFunctionReturn(0); 69 } 70 71 PetscErrorCode MatSolve_SeqBAIJ_3_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,*adiag=a->diag; 75 PetscInt i,k,nz,idx,jdx,idt; 76 const PetscInt bs = A->rmap->bs,bs2 = a->bs2; 77 const MatScalar *aa=a->a,*v; 78 PetscScalar *x; 79 const PetscScalar *b; 80 PetscScalar s1,s2,s3,x1,x2,x3; 81 82 PetscFunctionBegin; 83 PetscCall(VecGetArrayRead(bb,&b)); 84 PetscCall(VecGetArray(xx,&x)); 85 /* forward solve the lower triangular */ 86 idx = 0; 87 x[0] = b[idx]; x[1] = b[1+idx];x[2] = b[2+idx]; 88 for (i=1; i<n; i++) { 89 v = aa + bs2*ai[i]; 90 vi = aj + ai[i]; 91 nz = ai[i+1] - ai[i]; 92 idx = bs*i; 93 s1 = b[idx];s2 = b[1+idx];s3 = b[2+idx]; 94 for (k=0; k<nz; k++) { 95 jdx = bs*vi[k]; 96 x1 = x[jdx];x2 = x[1+jdx]; x3 =x[2+jdx]; 97 s1 -= v[0]*x1 + v[3]*x2 + v[6]*x3; 98 s2 -= v[1]*x1 + v[4]*x2 + v[7]*x3; 99 s3 -= v[2]*x1 + v[5]*x2 + v[8]*x3; 100 101 v += bs2; 102 } 103 104 x[idx] = s1; 105 x[1+idx] = s2; 106 x[2+idx] = s3; 107 } 108 109 /* backward solve the upper triangular */ 110 for (i=n-1; i>=0; i--) { 111 v = aa + bs2*(adiag[i+1]+1); 112 vi = aj + adiag[i+1]+1; 113 nz = adiag[i] - adiag[i+1]-1; 114 idt = bs*i; 115 s1 = x[idt]; s2 = x[1+idt];s3 = x[2+idt]; 116 117 for (k=0; k<nz; k++) { 118 idx = bs*vi[k]; 119 x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx]; 120 s1 -= v[0]*x1 + v[3]*x2 + v[6]*x3; 121 s2 -= v[1]*x1 + v[4]*x2 + v[7]*x3; 122 s3 -= v[2]*x1 + v[5]*x2 + v[8]*x3; 123 124 v += bs2; 125 } 126 /* x = inv_diagonal*x */ 127 x[idt] = v[0]*s1 + v[3]*s2 + v[6]*s3; 128 x[1+idt] = v[1]*s1 + v[4]*s2 + v[7]*s3; 129 x[2+idt] = v[2]*s1 + v[5]*s2 + v[8]*s3; 130 131 } 132 133 PetscCall(VecRestoreArrayRead(bb,&b)); 134 PetscCall(VecRestoreArray(xx,&x)); 135 PetscCall(PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n)); 136 PetscFunctionReturn(0); 137 } 138 139 PetscErrorCode MatForwardSolve_SeqBAIJ_3_NaturalOrdering(Mat A,Vec bb,Vec xx) 140 { 141 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 142 const PetscInt n =a->mbs,*vi,*ai=a->i,*aj=a->j; 143 PetscInt i,k,nz,idx,jdx; 144 const PetscInt bs = A->rmap->bs,bs2 = a->bs2; 145 const MatScalar *aa=a->a,*v; 146 PetscScalar *x; 147 const PetscScalar *b; 148 PetscScalar s1,s2,s3,x1,x2,x3; 149 150 PetscFunctionBegin; 151 PetscCall(VecGetArrayRead(bb,&b)); 152 PetscCall(VecGetArray(xx,&x)); 153 /* forward solve the lower triangular */ 154 idx = 0; 155 x[0] = b[idx]; x[1] = b[1+idx];x[2] = b[2+idx]; 156 for (i=1; i<n; i++) { 157 v = aa + bs2*ai[i]; 158 vi = aj + ai[i]; 159 nz = ai[i+1] - ai[i]; 160 idx = bs*i; 161 s1 = b[idx];s2 = b[1+idx];s3 = b[2+idx]; 162 for (k=0; k<nz; k++) { 163 jdx = bs*vi[k]; 164 x1 = x[jdx];x2 = x[1+jdx]; x3 =x[2+jdx]; 165 s1 -= v[0]*x1 + v[3]*x2 + v[6]*x3; 166 s2 -= v[1]*x1 + v[4]*x2 + v[7]*x3; 167 s3 -= v[2]*x1 + v[5]*x2 + v[8]*x3; 168 169 v += bs2; 170 } 171 172 x[idx] = s1; 173 x[1+idx] = s2; 174 x[2+idx] = s3; 175 } 176 177 PetscCall(VecRestoreArrayRead(bb,&b)); 178 PetscCall(VecRestoreArray(xx,&x)); 179 PetscCall(PetscLogFlops(1.0*bs2*(a->nz) - bs*A->cmap->n)); 180 PetscFunctionReturn(0); 181 } 182 183 PetscErrorCode MatBackwardSolve_SeqBAIJ_3_NaturalOrdering(Mat A,Vec bb,Vec xx) 184 { 185 Mat_SeqBAIJ *a = (Mat_SeqBAIJ*)A->data; 186 const PetscInt n =a->mbs,*vi,*aj=a->j,*adiag=a->diag; 187 PetscInt i,k,nz,idx,idt; 188 const PetscInt bs = A->rmap->bs,bs2 = a->bs2; 189 const MatScalar *aa=a->a,*v; 190 PetscScalar *x; 191 const PetscScalar *b; 192 PetscScalar s1,s2,s3,x1,x2,x3; 193 194 PetscFunctionBegin; 195 PetscCall(VecGetArrayRead(bb,&b)); 196 PetscCall(VecGetArray(xx,&x)); 197 198 /* backward solve the upper triangular */ 199 for (i=n-1; i>=0; i--) { 200 v = aa + bs2*(adiag[i+1]+1); 201 vi = aj + adiag[i+1]+1; 202 nz = adiag[i] - adiag[i+1]-1; 203 idt = bs*i; 204 s1 = b[idt]; s2 = b[1+idt];s3 = b[2+idt]; 205 206 for (k=0; k<nz; k++) { 207 idx = bs*vi[k]; 208 x1 = x[idx]; x2 = x[1+idx]; x3 = x[2+idx]; 209 s1 -= v[0]*x1 + v[3]*x2 + v[6]*x3; 210 s2 -= v[1]*x1 + v[4]*x2 + v[7]*x3; 211 s3 -= v[2]*x1 + v[5]*x2 + v[8]*x3; 212 213 v += bs2; 214 } 215 /* x = inv_diagonal*x */ 216 x[idt] = v[0]*s1 + v[3]*s2 + v[6]*s3; 217 x[1+idt] = v[1]*s1 + v[4]*s2 + v[7]*s3; 218 x[2+idt] = v[2]*s1 + v[5]*s2 + v[8]*s3; 219 220 } 221 222 PetscCall(VecRestoreArrayRead(bb,&b)); 223 PetscCall(VecRestoreArray(xx,&x)); 224 PetscCall(PetscLogFlops(2.0*bs2*(a->nz) - bs*A->cmap->n)); 225 PetscFunctionReturn(0); 226 } 227