12c733ed4SBarry Smith #include <../src/mat/impls/baij/seq/baij.h> 22c733ed4SBarry Smith #include <petsc/private/kernels/blockinvert.h> 32c733ed4SBarry Smith 4d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SeqBAIJ_5_NaturalOrdering_inplace(Mat A, Vec bb, Vec xx) 5d71ae5a4SJacob Faibussowitsch { 62c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 72c733ed4SBarry Smith const PetscInt *diag = a->diag, n = a->mbs, *vi, *ai = a->i, *aj = a->j; 82c733ed4SBarry Smith PetscInt i, nz, idx, idt, jdx; 92c733ed4SBarry Smith const MatScalar *aa = a->a, *v; 102c733ed4SBarry Smith PetscScalar *x, s1, s2, s3, s4, s5, x1, x2, x3, x4, x5; 112c733ed4SBarry Smith const PetscScalar *b; 122c733ed4SBarry Smith 132c733ed4SBarry Smith PetscFunctionBegin; 149566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 159566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x)); 162c733ed4SBarry Smith /* forward solve the lower triangular */ 172c733ed4SBarry Smith idx = 0; 189371c9d4SSatish Balay x[0] = b[idx]; 199371c9d4SSatish Balay x[1] = b[1 + idx]; 209371c9d4SSatish Balay x[2] = b[2 + idx]; 219371c9d4SSatish Balay x[3] = b[3 + idx]; 229371c9d4SSatish Balay x[4] = b[4 + idx]; 232c733ed4SBarry Smith for (i = 1; i < n; i++) { 242c733ed4SBarry Smith v = aa + 25 * ai[i]; 252c733ed4SBarry Smith vi = aj + ai[i]; 262c733ed4SBarry Smith nz = diag[i] - ai[i]; 272c733ed4SBarry Smith idx = 5 * i; 289371c9d4SSatish Balay s1 = b[idx]; 299371c9d4SSatish Balay s2 = b[1 + idx]; 309371c9d4SSatish Balay s3 = b[2 + idx]; 319371c9d4SSatish Balay s4 = b[3 + idx]; 329371c9d4SSatish Balay s5 = b[4 + idx]; 332c733ed4SBarry Smith while (nz--) { 342c733ed4SBarry Smith jdx = 5 * (*vi++); 359371c9d4SSatish Balay x1 = x[jdx]; 369371c9d4SSatish Balay x2 = x[1 + jdx]; 379371c9d4SSatish Balay x3 = x[2 + jdx]; 389371c9d4SSatish Balay x4 = x[3 + jdx]; 399371c9d4SSatish Balay x5 = x[4 + jdx]; 402c733ed4SBarry Smith s1 -= v[0] * x1 + v[5] * x2 + v[10] * x3 + v[15] * x4 + v[20] * x5; 412c733ed4SBarry Smith s2 -= v[1] * x1 + v[6] * x2 + v[11] * x3 + v[16] * x4 + v[21] * x5; 422c733ed4SBarry Smith s3 -= v[2] * x1 + v[7] * x2 + v[12] * x3 + v[17] * x4 + v[22] * x5; 432c733ed4SBarry Smith s4 -= v[3] * x1 + v[8] * x2 + v[13] * x3 + v[18] * x4 + v[23] * x5; 442c733ed4SBarry Smith s5 -= v[4] * x1 + v[9] * x2 + v[14] * x3 + v[19] * x4 + v[24] * x5; 452c733ed4SBarry Smith v += 25; 462c733ed4SBarry Smith } 472c733ed4SBarry Smith x[idx] = s1; 482c733ed4SBarry Smith x[1 + idx] = s2; 492c733ed4SBarry Smith x[2 + idx] = s3; 502c733ed4SBarry Smith x[3 + idx] = s4; 512c733ed4SBarry Smith x[4 + idx] = s5; 522c733ed4SBarry Smith } 532c733ed4SBarry Smith /* backward solve the upper triangular */ 542c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) { 552c733ed4SBarry Smith v = aa + 25 * diag[i] + 25; 562c733ed4SBarry Smith vi = aj + diag[i] + 1; 572c733ed4SBarry Smith nz = ai[i + 1] - diag[i] - 1; 582c733ed4SBarry Smith idt = 5 * i; 599371c9d4SSatish Balay s1 = x[idt]; 609371c9d4SSatish Balay s2 = x[1 + idt]; 619371c9d4SSatish Balay s3 = x[2 + idt]; 629371c9d4SSatish Balay s4 = x[3 + idt]; 639371c9d4SSatish Balay s5 = x[4 + idt]; 642c733ed4SBarry Smith while (nz--) { 652c733ed4SBarry Smith idx = 5 * (*vi++); 669371c9d4SSatish Balay x1 = x[idx]; 679371c9d4SSatish Balay x2 = x[1 + idx]; 689371c9d4SSatish Balay x3 = x[2 + idx]; 699371c9d4SSatish Balay x4 = x[3 + idx]; 709371c9d4SSatish Balay x5 = x[4 + idx]; 712c733ed4SBarry Smith s1 -= v[0] * x1 + v[5] * x2 + v[10] * x3 + v[15] * x4 + v[20] * x5; 722c733ed4SBarry Smith s2 -= v[1] * x1 + v[6] * x2 + v[11] * x3 + v[16] * x4 + v[21] * x5; 732c733ed4SBarry Smith s3 -= v[2] * x1 + v[7] * x2 + v[12] * x3 + v[17] * x4 + v[22] * x5; 742c733ed4SBarry Smith s4 -= v[3] * x1 + v[8] * x2 + v[13] * x3 + v[18] * x4 + v[23] * x5; 752c733ed4SBarry Smith s5 -= v[4] * x1 + v[9] * x2 + v[14] * x3 + v[19] * x4 + v[24] * x5; 762c733ed4SBarry Smith v += 25; 772c733ed4SBarry Smith } 782c733ed4SBarry Smith v = aa + 25 * diag[i]; 792c733ed4SBarry Smith x[idt] = v[0] * s1 + v[5] * s2 + v[10] * s3 + v[15] * s4 + v[20] * s5; 802c733ed4SBarry Smith x[1 + idt] = v[1] * s1 + v[6] * s2 + v[11] * s3 + v[16] * s4 + v[21] * s5; 812c733ed4SBarry Smith x[2 + idt] = v[2] * s1 + v[7] * s2 + v[12] * s3 + v[17] * s4 + v[22] * s5; 822c733ed4SBarry Smith x[3 + idt] = v[3] * s1 + v[8] * s2 + v[13] * s3 + v[18] * s4 + v[23] * s5; 832c733ed4SBarry Smith x[4 + idt] = v[4] * s1 + v[9] * s2 + v[14] * s3 + v[19] * s4 + v[24] * s5; 842c733ed4SBarry Smith } 852c733ed4SBarry Smith 869566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 879566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 889566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * 25 * (a->nz) - 5.0 * A->cmap->n)); 89*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 902c733ed4SBarry Smith } 912c733ed4SBarry Smith 92d71ae5a4SJacob Faibussowitsch PetscErrorCode MatSolve_SeqBAIJ_5_NaturalOrdering(Mat A, Vec bb, Vec xx) 93d71ae5a4SJacob Faibussowitsch { 942c733ed4SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data; 952c733ed4SBarry Smith const PetscInt n = a->mbs, *vi, *ai = a->i, *aj = a->j, *adiag = a->diag; 962c733ed4SBarry Smith PetscInt i, k, nz, idx, idt, jdx; 972c733ed4SBarry Smith const MatScalar *aa = a->a, *v; 982c733ed4SBarry Smith PetscScalar *x, s1, s2, s3, s4, s5, x1, x2, x3, x4, x5; 992c733ed4SBarry Smith const PetscScalar *b; 1002c733ed4SBarry Smith 1012c733ed4SBarry Smith PetscFunctionBegin; 1029566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(bb, &b)); 1039566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x)); 1042c733ed4SBarry Smith /* forward solve the lower triangular */ 1052c733ed4SBarry Smith idx = 0; 1069371c9d4SSatish Balay x[0] = b[idx]; 1079371c9d4SSatish Balay x[1] = b[1 + idx]; 1089371c9d4SSatish Balay x[2] = b[2 + idx]; 1099371c9d4SSatish Balay x[3] = b[3 + idx]; 1109371c9d4SSatish Balay x[4] = b[4 + idx]; 1112c733ed4SBarry Smith for (i = 1; i < n; i++) { 1122c733ed4SBarry Smith v = aa + 25 * ai[i]; 1132c733ed4SBarry Smith vi = aj + ai[i]; 1142c733ed4SBarry Smith nz = ai[i + 1] - ai[i]; 1152c733ed4SBarry Smith idx = 5 * i; 1169371c9d4SSatish Balay s1 = b[idx]; 1179371c9d4SSatish Balay s2 = b[1 + idx]; 1189371c9d4SSatish Balay s3 = b[2 + idx]; 1199371c9d4SSatish Balay s4 = b[3 + idx]; 1209371c9d4SSatish Balay s5 = b[4 + idx]; 1212c733ed4SBarry Smith for (k = 0; k < nz; k++) { 1222c733ed4SBarry Smith jdx = 5 * vi[k]; 1239371c9d4SSatish Balay x1 = x[jdx]; 1249371c9d4SSatish Balay x2 = x[1 + jdx]; 1259371c9d4SSatish Balay x3 = x[2 + jdx]; 1269371c9d4SSatish Balay x4 = x[3 + jdx]; 1279371c9d4SSatish Balay x5 = x[4 + jdx]; 1282c733ed4SBarry Smith s1 -= v[0] * x1 + v[5] * x2 + v[10] * x3 + v[15] * x4 + v[20] * x5; 1292c733ed4SBarry Smith s2 -= v[1] * x1 + v[6] * x2 + v[11] * x3 + v[16] * x4 + v[21] * x5; 1302c733ed4SBarry Smith s3 -= v[2] * x1 + v[7] * x2 + v[12] * x3 + v[17] * x4 + v[22] * x5; 1312c733ed4SBarry Smith s4 -= v[3] * x1 + v[8] * x2 + v[13] * x3 + v[18] * x4 + v[23] * x5; 1322c733ed4SBarry Smith s5 -= v[4] * x1 + v[9] * x2 + v[14] * x3 + v[19] * x4 + v[24] * x5; 1332c733ed4SBarry Smith v += 25; 1342c733ed4SBarry Smith } 1352c733ed4SBarry Smith x[idx] = s1; 1362c733ed4SBarry Smith x[1 + idx] = s2; 1372c733ed4SBarry Smith x[2 + idx] = s3; 1382c733ed4SBarry Smith x[3 + idx] = s4; 1392c733ed4SBarry Smith x[4 + idx] = s5; 1402c733ed4SBarry Smith } 1412c733ed4SBarry Smith 1422c733ed4SBarry Smith /* backward solve the upper triangular */ 1432c733ed4SBarry Smith for (i = n - 1; i >= 0; i--) { 1442c733ed4SBarry Smith v = aa + 25 * (adiag[i + 1] + 1); 1452c733ed4SBarry Smith vi = aj + adiag[i + 1] + 1; 1462c733ed4SBarry Smith nz = adiag[i] - adiag[i + 1] - 1; 1472c733ed4SBarry Smith idt = 5 * i; 1489371c9d4SSatish Balay s1 = x[idt]; 1499371c9d4SSatish Balay s2 = x[1 + idt]; 1509371c9d4SSatish Balay s3 = x[2 + idt]; 1519371c9d4SSatish Balay s4 = x[3 + idt]; 1529371c9d4SSatish Balay s5 = x[4 + idt]; 1532c733ed4SBarry Smith for (k = 0; k < nz; k++) { 1542c733ed4SBarry Smith idx = 5 * vi[k]; 1559371c9d4SSatish Balay x1 = x[idx]; 1569371c9d4SSatish Balay x2 = x[1 + idx]; 1579371c9d4SSatish Balay x3 = x[2 + idx]; 1589371c9d4SSatish Balay x4 = x[3 + idx]; 1599371c9d4SSatish Balay x5 = x[4 + idx]; 1602c733ed4SBarry Smith s1 -= v[0] * x1 + v[5] * x2 + v[10] * x3 + v[15] * x4 + v[20] * x5; 1612c733ed4SBarry Smith s2 -= v[1] * x1 + v[6] * x2 + v[11] * x3 + v[16] * x4 + v[21] * x5; 1622c733ed4SBarry Smith s3 -= v[2] * x1 + v[7] * x2 + v[12] * x3 + v[17] * x4 + v[22] * x5; 1632c733ed4SBarry Smith s4 -= v[3] * x1 + v[8] * x2 + v[13] * x3 + v[18] * x4 + v[23] * x5; 1642c733ed4SBarry Smith s5 -= v[4] * x1 + v[9] * x2 + v[14] * x3 + v[19] * x4 + v[24] * x5; 1652c733ed4SBarry Smith v += 25; 1662c733ed4SBarry Smith } 1672c733ed4SBarry Smith /* x = inv_diagonal*x */ 1682c733ed4SBarry Smith x[idt] = v[0] * s1 + v[5] * s2 + v[10] * s3 + v[15] * s4 + v[20] * s5; 1692c733ed4SBarry Smith x[1 + idt] = v[1] * s1 + v[6] * s2 + v[11] * s3 + v[16] * s4 + v[21] * s5; 1702c733ed4SBarry Smith x[2 + idt] = v[2] * s1 + v[7] * s2 + v[12] * s3 + v[17] * s4 + v[22] * s5; 1712c733ed4SBarry Smith x[3 + idt] = v[3] * s1 + v[8] * s2 + v[13] * s3 + v[18] * s4 + v[23] * s5; 1722c733ed4SBarry Smith x[4 + idt] = v[4] * s1 + v[9] * s2 + v[14] * s3 + v[19] * s4 + v[24] * s5; 1732c733ed4SBarry Smith } 1742c733ed4SBarry Smith 1759566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(bb, &b)); 1769566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 1779566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * 25 * (a->nz) - 5.0 * A->cmap->n)); 178*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1792c733ed4SBarry Smith } 180