12c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h> 22c733ed4SBarry Smith #include <petsc/private/kernels/blockinvert.h> 32c733ed4SBarry Smith 4*fb56d528SJed Brown #if defined(PETSC_HAVE_XMMINTRIN_H) 5*fb56d528SJed Brown #include <xmmintrin.h> 6*fb56d528SJed Brown #endif 7*fb56d528SJed Brown 82c733ed4SBarry Smith /* 92c733ed4SBarry Smith Special case where the matrix was ILU(0) factored in the natural 102c733ed4SBarry Smith ordering. This eliminates the need for the column and row permutation. 112c733ed4SBarry Smith */ 12d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SeqBAIJ_2_NaturalOrdering_inplace(Mat A, Vec bb, Vec xx) 13d71ae5a4SJacob Faibussowitsch { 142c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 152c733ed4SBarry Smith const PetscInt n = a->mbs, *vi, *ai = a->i, *aj = a->j, *diag = a->diag; 162c733ed4SBarry Smith const MatScalar *aa = a->a, *v; 172c733ed4SBarry Smith PetscScalar *x, s1, s2, x1, x2; 182c733ed4SBarry Smith const PetscScalar *b; 192c733ed4SBarry Smith PetscInt jdx, idt, idx, nz, i; 202c733ed4SBarry Smith 212c733ed4SBarry Smith PetscFunctionBegin; 229566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 239566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x)); 242c733ed4SBarry Smith 252c733ed4SBarry Smith /* forward solve the lower triangular */ 262c733ed4SBarry Smith idx = 0; 279371c9d4SSatish Balay x[0] = b[0]; 289371c9d4SSatish Balay x[1] = b[1]; 292c733ed4SBarry Smith for (i = 1; i < n; i++) { 302c733ed4SBarry Smith v = aa + 4 * ai[i]; 312c733ed4SBarry Smith vi = aj + ai[i]; 322c733ed4SBarry Smith nz = diag[i] - ai[i]; 332c733ed4SBarry Smith idx += 2; 349371c9d4SSatish Balay s1 = b[idx]; 359371c9d4SSatish Balay s2 = b[1 + idx]; 362c733ed4SBarry Smith while (nz--) { 372c733ed4SBarry Smith jdx = 2 * (*vi++); 389371c9d4SSatish Balay x1 = x[jdx]; 399371c9d4SSatish Balay x2 = x[1 + jdx]; 402c733ed4SBarry Smith s1 -= v[0] * x1 + v[2] * x2; 412c733ed4SBarry Smith s2 -= v[1] * x1 + v[3] * x2; 422c733ed4SBarry Smith v += 4; 432c733ed4SBarry Smith } 442c733ed4SBarry Smith x[idx] = s1; 452c733ed4SBarry Smith x[1 + idx] = s2; 462c733ed4SBarry Smith } 472c733ed4SBarry Smith /* backward solve the upper triangular */ 482c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) { 492c733ed4SBarry Smith v = aa + 4 * diag[i] + 4; 502c733ed4SBarry Smith vi = aj + diag[i] + 1; 512c733ed4SBarry Smith nz = ai[i + 1] - diag[i] - 1; 522c733ed4SBarry Smith idt = 2 * i; 539371c9d4SSatish Balay s1 = x[idt]; 549371c9d4SSatish Balay s2 = x[1 + idt]; 552c733ed4SBarry Smith while (nz--) { 562c733ed4SBarry Smith idx = 2 * (*vi++); 579371c9d4SSatish Balay x1 = x[idx]; 589371c9d4SSatish Balay x2 = x[1 + idx]; 592c733ed4SBarry Smith s1 -= v[0] * x1 + v[2] * x2; 602c733ed4SBarry Smith s2 -= v[1] * x1 + v[3] * x2; 612c733ed4SBarry Smith v += 4; 622c733ed4SBarry Smith } 632c733ed4SBarry Smith v = aa + 4 * diag[i]; 642c733ed4SBarry Smith x[idt] = v[0] * s1 + v[2] * s2; 652c733ed4SBarry Smith x[1 + idt] = v[1] * s1 + v[3] * s2; 662c733ed4SBarry Smith } 672c733ed4SBarry Smith 689566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 699566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 709566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * 4 * (a->nz) - 2.0 * A->cmap->n)); 713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 722c733ed4SBarry Smith } 732c733ed4SBarry Smith 74d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SeqBAIJ_2_NaturalOrdering(Mat A, Vec bb, Vec xx) 75d71ae5a4SJacob Faibussowitsch { 762c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 772c733ed4SBarry Smith const PetscInt n = a->mbs, *vi, *ai = a->i, *aj = a->j, *adiag = a->diag; 782c733ed4SBarry Smith PetscInt i, k, nz, idx, idt, jdx; 792c733ed4SBarry Smith const MatScalar *aa = a->a, *v; 802c733ed4SBarry Smith PetscScalar *x, s1, s2, x1, x2; 812c733ed4SBarry Smith const PetscScalar *b; 822c733ed4SBarry Smith 832c733ed4SBarry Smith PetscFunctionBegin; 849566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 859566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x)); 862c733ed4SBarry Smith /* forward solve the lower triangular */ 872c733ed4SBarry Smith idx = 0; 889371c9d4SSatish Balay x[0] = b[idx]; 899371c9d4SSatish Balay x[1] = b[1 + idx]; 902c733ed4SBarry Smith for (i = 1; i < n; i++) { 912c733ed4SBarry Smith v = aa + 4 * ai[i]; 922c733ed4SBarry Smith vi = aj + ai[i]; 932c733ed4SBarry Smith nz = ai[i + 1] - ai[i]; 942c733ed4SBarry Smith idx = 2 * i; 959371c9d4SSatish Balay s1 = b[idx]; 969371c9d4SSatish Balay s2 = b[1 + idx]; 972c733ed4SBarry Smith PetscPrefetchBlock(vi + nz, nz, 0, PETSC_PREFETCH_HINT_NTA); 982c733ed4SBarry Smith PetscPrefetchBlock(v + 4 * nz, 4 * nz, 0, PETSC_PREFETCH_HINT_NTA); 992c733ed4SBarry Smith for (k = 0; k < nz; k++) { 1002c733ed4SBarry Smith jdx = 2 * vi[k]; 1019371c9d4SSatish Balay x1 = x[jdx]; 1029371c9d4SSatish Balay x2 = x[1 + jdx]; 1032c733ed4SBarry Smith s1 -= v[0] * x1 + v[2] * x2; 1042c733ed4SBarry Smith s2 -= v[1] * x1 + v[3] * x2; 1052c733ed4SBarry Smith v += 4; 1062c733ed4SBarry Smith } 1072c733ed4SBarry Smith x[idx] = s1; 1082c733ed4SBarry Smith x[1 + idx] = s2; 1092c733ed4SBarry Smith } 1102c733ed4SBarry Smith 1112c733ed4SBarry Smith /* backward solve the upper triangular */ 1122c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) { 1132c733ed4SBarry Smith v = aa + 4 * (adiag[i + 1] + 1); 1142c733ed4SBarry Smith vi = aj + adiag[i + 1] + 1; 1152c733ed4SBarry Smith nz = adiag[i] - adiag[i + 1] - 1; 1162c733ed4SBarry Smith idt = 2 * i; 1179371c9d4SSatish Balay s1 = x[idt]; 1189371c9d4SSatish Balay s2 = x[1 + idt]; 1192c733ed4SBarry Smith PetscPrefetchBlock(vi + nz, nz, 0, PETSC_PREFETCH_HINT_NTA); 1202c733ed4SBarry Smith PetscPrefetchBlock(v + 4 * nz, 4 * nz, 0, PETSC_PREFETCH_HINT_NTA); 1212c733ed4SBarry Smith for (k = 0; k < nz; k++) { 1222c733ed4SBarry Smith idx = 2 * vi[k]; 1239371c9d4SSatish Balay x1 = x[idx]; 1249371c9d4SSatish Balay x2 = x[1 + idx]; 1252c733ed4SBarry Smith s1 -= v[0] * x1 + v[2] * x2; 1262c733ed4SBarry Smith s2 -= v[1] * x1 + v[3] * x2; 1272c733ed4SBarry Smith v += 4; 1282c733ed4SBarry Smith } 1292c733ed4SBarry Smith /* x = inv_diagonal*x */ 1302c733ed4SBarry Smith x[idt] = v[0] * s1 + v[2] * s2; 1312c733ed4SBarry Smith x[1 + idt] = v[1] * s1 + v[3] * s2; 1322c733ed4SBarry Smith } 1332c733ed4SBarry Smith 1349566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 1359566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 1369566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * 4 * (a->nz) - 2.0 * A->cmap->n)); 1373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1382c733ed4SBarry Smith } 1392c733ed4SBarry Smith 140d71ae5a4SJacob Faibussowitsch PetscErrorCode MatForwardSolve_SeqBAIJ_2_NaturalOrdering(Mat A, Vec bb, Vec xx) 141d71ae5a4SJacob Faibussowitsch { 1422c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 1432c733ed4SBarry Smith const PetscInt n = a->mbs, *vi, *ai = a->i, *aj = a->j; 1442c733ed4SBarry Smith PetscInt i, k, nz, idx, jdx; 1452c733ed4SBarry Smith const MatScalar *aa = a->a, *v; 1462c733ed4SBarry Smith PetscScalar *x, s1, s2, x1, x2; 1472c733ed4SBarry Smith const PetscScalar *b; 1482c733ed4SBarry Smith 1492c733ed4SBarry Smith PetscFunctionBegin; 1509566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 1519566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x)); 1522c733ed4SBarry Smith /* forward solve the lower triangular */ 1532c733ed4SBarry Smith idx = 0; 1549371c9d4SSatish Balay x[0] = b[idx]; 1559371c9d4SSatish Balay x[1] = b[1 + idx]; 1562c733ed4SBarry Smith for (i = 1; i < n; i++) { 1572c733ed4SBarry Smith v = aa + 4 * ai[i]; 1582c733ed4SBarry Smith vi = aj + ai[i]; 1592c733ed4SBarry Smith nz = ai[i + 1] - ai[i]; 1602c733ed4SBarry Smith idx = 2 * i; 1619371c9d4SSatish Balay s1 = b[idx]; 1629371c9d4SSatish Balay s2 = b[1 + idx]; 1632c733ed4SBarry Smith PetscPrefetchBlock(vi + nz, nz, 0, PETSC_PREFETCH_HINT_NTA); 1642c733ed4SBarry Smith PetscPrefetchBlock(v + 4 * nz, 4 * nz, 0, PETSC_PREFETCH_HINT_NTA); 1652c733ed4SBarry Smith for (k = 0; k < nz; k++) { 1662c733ed4SBarry Smith jdx = 2 * vi[k]; 1679371c9d4SSatish Balay x1 = x[jdx]; 1689371c9d4SSatish Balay x2 = x[1 + jdx]; 1692c733ed4SBarry Smith s1 -= v[0] * x1 + v[2] * x2; 1702c733ed4SBarry Smith s2 -= v[1] * x1 + v[3] * x2; 1712c733ed4SBarry Smith v += 4; 1722c733ed4SBarry Smith } 1732c733ed4SBarry Smith x[idx] = s1; 1742c733ed4SBarry Smith x[1 + idx] = s2; 1752c733ed4SBarry Smith } 1762c733ed4SBarry Smith 1779566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 1789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 1799566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(4.0 * (a->nz) - A->cmap->n)); 1803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1812c733ed4SBarry Smith } 1822c733ed4SBarry Smith 183d71ae5a4SJacob Faibussowitsch PetscErrorCode MatBackwardSolve_SeqBAIJ_2_NaturalOrdering(Mat A, Vec bb, Vec xx) 184d71ae5a4SJacob Faibussowitsch { 1852c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 1862c733ed4SBarry Smith const PetscInt n = a->mbs, *vi, *aj = a->j, *adiag = a->diag; 1872c733ed4SBarry Smith PetscInt i, k, nz, idx, idt; 1882c733ed4SBarry Smith const MatScalar *aa = a->a, *v; 1892c733ed4SBarry Smith PetscScalar *x, s1, s2, x1, x2; 1902c733ed4SBarry Smith const PetscScalar *b; 1912c733ed4SBarry Smith 1922c733ed4SBarry Smith PetscFunctionBegin; 1939566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 1949566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x)); 1952c733ed4SBarry Smith 1962c733ed4SBarry Smith /* backward solve the upper triangular */ 1972c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) { 1982c733ed4SBarry Smith v = aa + 4 * (adiag[i + 1] + 1); 1992c733ed4SBarry Smith vi = aj + adiag[i + 1] + 1; 2002c733ed4SBarry Smith nz = adiag[i] - adiag[i + 1] - 1; 2012c733ed4SBarry Smith idt = 2 * i; 2029371c9d4SSatish Balay s1 = b[idt]; 2039371c9d4SSatish Balay s2 = b[1 + idt]; 2042c733ed4SBarry Smith PetscPrefetchBlock(vi + nz, nz, 0, PETSC_PREFETCH_HINT_NTA); 2052c733ed4SBarry Smith PetscPrefetchBlock(v + 4 * nz, 4 * nz, 0, PETSC_PREFETCH_HINT_NTA); 2062c733ed4SBarry Smith for (k = 0; k < nz; k++) { 2072c733ed4SBarry Smith idx = 2 * vi[k]; 2089371c9d4SSatish Balay x1 = x[idx]; 2099371c9d4SSatish Balay x2 = x[1 + idx]; 2102c733ed4SBarry Smith s1 -= v[0] * x1 + v[2] * x2; 2112c733ed4SBarry Smith s2 -= v[1] * x1 + v[3] * x2; 2122c733ed4SBarry Smith v += 4; 2132c733ed4SBarry Smith } 2142c733ed4SBarry Smith /* x = inv_diagonal*x */ 2152c733ed4SBarry Smith x[idt] = v[0] * s1 + v[2] * s2; 2162c733ed4SBarry Smith x[1 + idt] = v[1] * s1 + v[3] * s2; 2172c733ed4SBarry Smith } 2182c733ed4SBarry Smith 2199566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 2209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 2219566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(4.0 * a->nz - A->cmap->n)); 2223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2232c733ed4SBarry Smith } 224