1be1d678aSKris Buschelman 283287d42SBarry Smith /* 383287d42SBarry Smith Factorization code for BAIJ format. 483287d42SBarry Smith */ 5c6db04a5SJed Brown #include <../src/mat/impls/baij/seq/baij.h> 6af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h> 783287d42SBarry Smith 883287d42SBarry Smith /* ------------------------------------------------------------*/ 983287d42SBarry Smith /* 1083287d42SBarry Smith Version for when blocks are 5 by 5 1183287d42SBarry Smith */ 129371c9d4SSatish Balay PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_inplace(Mat C, Mat A, const MatFactorInfo *info) { 1383287d42SBarry Smith Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data; 1483287d42SBarry Smith IS isrow = b->row, isicol = b->icol; 155d0c19d7SBarry Smith const PetscInt *r, *ic, *bi = b->i, *bj = b->j, *ajtmpold, *ajtmp; 16766f9fbaSBarry Smith PetscInt i, j, n = a->mbs, nz, row, idx, ipvt[5]; 175d0c19d7SBarry Smith const PetscInt *diag_offset = b->diag, *ai = a->i, *aj = a->j, *pj; 188397fe1aSBarry Smith MatScalar *w, *pv, *rtmp, *x, *pc; 198397fe1aSBarry Smith const MatScalar *v, *aa = a->a; 2083287d42SBarry Smith MatScalar p1, p2, p3, p4, m1, m2, m3, m4, m5, m6, m7, m8, m9, x1, x2, x3, x4; 2183287d42SBarry Smith MatScalar p5, p6, p7, p8, p9, x5, x6, x7, x8, x9, x10, x11, x12, x13, x14, x15, x16; 2283287d42SBarry Smith MatScalar x17, x18, x19, x20, x21, x22, x23, x24, x25, p10, p11, p12, p13, p14; 2383287d42SBarry Smith MatScalar p15, p16, p17, p18, p19, p20, p21, p22, p23, p24, p25, m10, m11, m12; 2483287d42SBarry Smith MatScalar m13, m14, m15, m16, m17, m18, m19, m20, m21, m22, m23, m24, m25; 25766f9fbaSBarry Smith MatScalar *ba = b->a, work[25]; 26182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 27a455e926SHong Zhang PetscBool allowzeropivot, zeropivotdetected; 2883287d42SBarry Smith 2983287d42SBarry Smith PetscFunctionBegin; 300164db54SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 319566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrow, &r)); 329566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isicol, &ic)); 339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(25 * (n + 1), &rtmp)); 3483287d42SBarry Smith 358397fe1aSBarry Smith #define PETSC_USE_MEMZERO 1 368397fe1aSBarry Smith #define PETSC_USE_MEMCPY 1 378397fe1aSBarry Smith 3883287d42SBarry Smith for (i = 0; i < n; i++) { 3983287d42SBarry Smith nz = bi[i + 1] - bi[i]; 4083287d42SBarry Smith ajtmp = bj + bi[i]; 4183287d42SBarry Smith for (j = 0; j < nz; j++) { 428397fe1aSBarry Smith #if defined(PETSC_USE_MEMZERO) 439566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(rtmp + 25 * ajtmp[j], 25)); 448397fe1aSBarry Smith #else 4583287d42SBarry Smith x = rtmp + 25 * ajtmp[j]; 4683287d42SBarry Smith x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0; 4783287d42SBarry Smith x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = x[16] = x[17] = 0.0; 4883287d42SBarry Smith x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = 0.0; 498397fe1aSBarry Smith #endif 5083287d42SBarry Smith } 5183287d42SBarry Smith /* load in initial (unfactored row) */ 5283287d42SBarry Smith idx = r[i]; 5383287d42SBarry Smith nz = ai[idx + 1] - ai[idx]; 5483287d42SBarry Smith ajtmpold = aj + ai[idx]; 5583287d42SBarry Smith v = aa + 25 * ai[idx]; 5683287d42SBarry Smith for (j = 0; j < nz; j++) { 578397fe1aSBarry Smith #if defined(PETSC_USE_MEMCPY) 589566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(rtmp + 25 * ic[ajtmpold[j]], v, 25)); 598397fe1aSBarry Smith #else 6083287d42SBarry Smith x = rtmp + 25 * ic[ajtmpold[j]]; 619371c9d4SSatish Balay x[0] = v[0]; 629371c9d4SSatish Balay x[1] = v[1]; 639371c9d4SSatish Balay x[2] = v[2]; 649371c9d4SSatish Balay x[3] = v[3]; 659371c9d4SSatish Balay x[4] = v[4]; 669371c9d4SSatish Balay x[5] = v[5]; 679371c9d4SSatish Balay x[6] = v[6]; 689371c9d4SSatish Balay x[7] = v[7]; 699371c9d4SSatish Balay x[8] = v[8]; 709371c9d4SSatish Balay x[9] = v[9]; 719371c9d4SSatish Balay x[10] = v[10]; 729371c9d4SSatish Balay x[11] = v[11]; 739371c9d4SSatish Balay x[12] = v[12]; 749371c9d4SSatish Balay x[13] = v[13]; 759371c9d4SSatish Balay x[14] = v[14]; 769371c9d4SSatish Balay x[15] = v[15]; 779371c9d4SSatish Balay x[16] = v[16]; 789371c9d4SSatish Balay x[17] = v[17]; 799371c9d4SSatish Balay x[18] = v[18]; 809371c9d4SSatish Balay x[19] = v[19]; 819371c9d4SSatish Balay x[20] = v[20]; 829371c9d4SSatish Balay x[21] = v[21]; 839371c9d4SSatish Balay x[22] = v[22]; 849371c9d4SSatish Balay x[23] = v[23]; 859371c9d4SSatish Balay x[24] = v[24]; 868397fe1aSBarry Smith #endif 8783287d42SBarry Smith v += 25; 8883287d42SBarry Smith } 8983287d42SBarry Smith row = *ajtmp++; 9083287d42SBarry Smith while (row < i) { 9183287d42SBarry Smith pc = rtmp + 25 * row; 929371c9d4SSatish Balay p1 = pc[0]; 939371c9d4SSatish Balay p2 = pc[1]; 949371c9d4SSatish Balay p3 = pc[2]; 959371c9d4SSatish Balay p4 = pc[3]; 969371c9d4SSatish Balay p5 = pc[4]; 979371c9d4SSatish Balay p6 = pc[5]; 989371c9d4SSatish Balay p7 = pc[6]; 999371c9d4SSatish Balay p8 = pc[7]; 1009371c9d4SSatish Balay p9 = pc[8]; 1019371c9d4SSatish Balay p10 = pc[9]; 1029371c9d4SSatish Balay p11 = pc[10]; 1039371c9d4SSatish Balay p12 = pc[11]; 1049371c9d4SSatish Balay p13 = pc[12]; 1059371c9d4SSatish Balay p14 = pc[13]; 1069371c9d4SSatish Balay p15 = pc[14]; 1079371c9d4SSatish Balay p16 = pc[15]; 1089371c9d4SSatish Balay p17 = pc[16]; 1099371c9d4SSatish Balay p18 = pc[17]; 1109371c9d4SSatish Balay p19 = pc[18]; 1119371c9d4SSatish Balay p20 = pc[19]; 1129371c9d4SSatish Balay p21 = pc[20]; 1139371c9d4SSatish Balay p22 = pc[21]; 1149371c9d4SSatish Balay p23 = pc[22]; 1159371c9d4SSatish Balay p24 = pc[23]; 11683287d42SBarry Smith p25 = pc[24]; 1179371c9d4SSatish Balay if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 || p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 || p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 || p25 != 0.0) { 11883287d42SBarry Smith pv = ba + 25 * diag_offset[row]; 11983287d42SBarry Smith pj = bj + diag_offset[row] + 1; 1209371c9d4SSatish Balay x1 = pv[0]; 1219371c9d4SSatish Balay x2 = pv[1]; 1229371c9d4SSatish Balay x3 = pv[2]; 1239371c9d4SSatish Balay x4 = pv[3]; 1249371c9d4SSatish Balay x5 = pv[4]; 1259371c9d4SSatish Balay x6 = pv[5]; 1269371c9d4SSatish Balay x7 = pv[6]; 1279371c9d4SSatish Balay x8 = pv[7]; 1289371c9d4SSatish Balay x9 = pv[8]; 1299371c9d4SSatish Balay x10 = pv[9]; 1309371c9d4SSatish Balay x11 = pv[10]; 1319371c9d4SSatish Balay x12 = pv[11]; 1329371c9d4SSatish Balay x13 = pv[12]; 1339371c9d4SSatish Balay x14 = pv[13]; 1349371c9d4SSatish Balay x15 = pv[14]; 1359371c9d4SSatish Balay x16 = pv[15]; 1369371c9d4SSatish Balay x17 = pv[16]; 1379371c9d4SSatish Balay x18 = pv[17]; 1389371c9d4SSatish Balay x19 = pv[18]; 1399371c9d4SSatish Balay x20 = pv[19]; 1409371c9d4SSatish Balay x21 = pv[20]; 1419371c9d4SSatish Balay x22 = pv[21]; 1429371c9d4SSatish Balay x23 = pv[22]; 1439371c9d4SSatish Balay x24 = pv[23]; 1449371c9d4SSatish Balay x25 = pv[24]; 14583287d42SBarry Smith pc[0] = m1 = p1 * x1 + p6 * x2 + p11 * x3 + p16 * x4 + p21 * x5; 14683287d42SBarry Smith pc[1] = m2 = p2 * x1 + p7 * x2 + p12 * x3 + p17 * x4 + p22 * x5; 14783287d42SBarry Smith pc[2] = m3 = p3 * x1 + p8 * x2 + p13 * x3 + p18 * x4 + p23 * x5; 14883287d42SBarry Smith pc[3] = m4 = p4 * x1 + p9 * x2 + p14 * x3 + p19 * x4 + p24 * x5; 14983287d42SBarry Smith pc[4] = m5 = p5 * x1 + p10 * x2 + p15 * x3 + p20 * x4 + p25 * x5; 15083287d42SBarry Smith 15183287d42SBarry Smith pc[5] = m6 = p1 * x6 + p6 * x7 + p11 * x8 + p16 * x9 + p21 * x10; 15283287d42SBarry Smith pc[6] = m7 = p2 * x6 + p7 * x7 + p12 * x8 + p17 * x9 + p22 * x10; 15383287d42SBarry Smith pc[7] = m8 = p3 * x6 + p8 * x7 + p13 * x8 + p18 * x9 + p23 * x10; 15483287d42SBarry Smith pc[8] = m9 = p4 * x6 + p9 * x7 + p14 * x8 + p19 * x9 + p24 * x10; 15583287d42SBarry Smith pc[9] = m10 = p5 * x6 + p10 * x7 + p15 * x8 + p20 * x9 + p25 * x10; 15683287d42SBarry Smith 15783287d42SBarry Smith pc[10] = m11 = p1 * x11 + p6 * x12 + p11 * x13 + p16 * x14 + p21 * x15; 15883287d42SBarry Smith pc[11] = m12 = p2 * x11 + p7 * x12 + p12 * x13 + p17 * x14 + p22 * x15; 15983287d42SBarry Smith pc[12] = m13 = p3 * x11 + p8 * x12 + p13 * x13 + p18 * x14 + p23 * x15; 16083287d42SBarry Smith pc[13] = m14 = p4 * x11 + p9 * x12 + p14 * x13 + p19 * x14 + p24 * x15; 16183287d42SBarry Smith pc[14] = m15 = p5 * x11 + p10 * x12 + p15 * x13 + p20 * x14 + p25 * x15; 16283287d42SBarry Smith 16383287d42SBarry Smith pc[15] = m16 = p1 * x16 + p6 * x17 + p11 * x18 + p16 * x19 + p21 * x20; 16483287d42SBarry Smith pc[16] = m17 = p2 * x16 + p7 * x17 + p12 * x18 + p17 * x19 + p22 * x20; 16583287d42SBarry Smith pc[17] = m18 = p3 * x16 + p8 * x17 + p13 * x18 + p18 * x19 + p23 * x20; 16683287d42SBarry Smith pc[18] = m19 = p4 * x16 + p9 * x17 + p14 * x18 + p19 * x19 + p24 * x20; 16783287d42SBarry Smith pc[19] = m20 = p5 * x16 + p10 * x17 + p15 * x18 + p20 * x19 + p25 * x20; 16883287d42SBarry Smith 16983287d42SBarry Smith pc[20] = m21 = p1 * x21 + p6 * x22 + p11 * x23 + p16 * x24 + p21 * x25; 17083287d42SBarry Smith pc[21] = m22 = p2 * x21 + p7 * x22 + p12 * x23 + p17 * x24 + p22 * x25; 17183287d42SBarry Smith pc[22] = m23 = p3 * x21 + p8 * x22 + p13 * x23 + p18 * x24 + p23 * x25; 17283287d42SBarry Smith pc[23] = m24 = p4 * x21 + p9 * x22 + p14 * x23 + p19 * x24 + p24 * x25; 17383287d42SBarry Smith pc[24] = m25 = p5 * x21 + p10 * x22 + p15 * x23 + p20 * x24 + p25 * x25; 17483287d42SBarry Smith 17583287d42SBarry Smith nz = bi[row + 1] - diag_offset[row] - 1; 17683287d42SBarry Smith pv += 25; 17783287d42SBarry Smith for (j = 0; j < nz; j++) { 1789371c9d4SSatish Balay x1 = pv[0]; 1799371c9d4SSatish Balay x2 = pv[1]; 1809371c9d4SSatish Balay x3 = pv[2]; 1819371c9d4SSatish Balay x4 = pv[3]; 1829371c9d4SSatish Balay x5 = pv[4]; 1839371c9d4SSatish Balay x6 = pv[5]; 1849371c9d4SSatish Balay x7 = pv[6]; 1859371c9d4SSatish Balay x8 = pv[7]; 1869371c9d4SSatish Balay x9 = pv[8]; 1879371c9d4SSatish Balay x10 = pv[9]; 1889371c9d4SSatish Balay x11 = pv[10]; 1899371c9d4SSatish Balay x12 = pv[11]; 1909371c9d4SSatish Balay x13 = pv[12]; 1919371c9d4SSatish Balay x14 = pv[13]; 1929371c9d4SSatish Balay x15 = pv[14]; 1939371c9d4SSatish Balay x16 = pv[15]; 1949371c9d4SSatish Balay x17 = pv[16]; 1959371c9d4SSatish Balay x18 = pv[17]; 1969371c9d4SSatish Balay x19 = pv[18]; 1979371c9d4SSatish Balay x20 = pv[19]; 1989371c9d4SSatish Balay x21 = pv[20]; 1999371c9d4SSatish Balay x22 = pv[21]; 2009371c9d4SSatish Balay x23 = pv[22]; 2019371c9d4SSatish Balay x24 = pv[23]; 2029371c9d4SSatish Balay x25 = pv[24]; 20383287d42SBarry Smith x = rtmp + 25 * pj[j]; 20483287d42SBarry Smith x[0] -= m1 * x1 + m6 * x2 + m11 * x3 + m16 * x4 + m21 * x5; 20583287d42SBarry Smith x[1] -= m2 * x1 + m7 * x2 + m12 * x3 + m17 * x4 + m22 * x5; 20683287d42SBarry Smith x[2] -= m3 * x1 + m8 * x2 + m13 * x3 + m18 * x4 + m23 * x5; 20783287d42SBarry Smith x[3] -= m4 * x1 + m9 * x2 + m14 * x3 + m19 * x4 + m24 * x5; 20883287d42SBarry Smith x[4] -= m5 * x1 + m10 * x2 + m15 * x3 + m20 * x4 + m25 * x5; 20983287d42SBarry Smith 21083287d42SBarry Smith x[5] -= m1 * x6 + m6 * x7 + m11 * x8 + m16 * x9 + m21 * x10; 21183287d42SBarry Smith x[6] -= m2 * x6 + m7 * x7 + m12 * x8 + m17 * x9 + m22 * x10; 21283287d42SBarry Smith x[7] -= m3 * x6 + m8 * x7 + m13 * x8 + m18 * x9 + m23 * x10; 21383287d42SBarry Smith x[8] -= m4 * x6 + m9 * x7 + m14 * x8 + m19 * x9 + m24 * x10; 21483287d42SBarry Smith x[9] -= m5 * x6 + m10 * x7 + m15 * x8 + m20 * x9 + m25 * x10; 21583287d42SBarry Smith 21683287d42SBarry Smith x[10] -= m1 * x11 + m6 * x12 + m11 * x13 + m16 * x14 + m21 * x15; 21783287d42SBarry Smith x[11] -= m2 * x11 + m7 * x12 + m12 * x13 + m17 * x14 + m22 * x15; 21883287d42SBarry Smith x[12] -= m3 * x11 + m8 * x12 + m13 * x13 + m18 * x14 + m23 * x15; 21983287d42SBarry Smith x[13] -= m4 * x11 + m9 * x12 + m14 * x13 + m19 * x14 + m24 * x15; 22083287d42SBarry Smith x[14] -= m5 * x11 + m10 * x12 + m15 * x13 + m20 * x14 + m25 * x15; 22183287d42SBarry Smith 22283287d42SBarry Smith x[15] -= m1 * x16 + m6 * x17 + m11 * x18 + m16 * x19 + m21 * x20; 22383287d42SBarry Smith x[16] -= m2 * x16 + m7 * x17 + m12 * x18 + m17 * x19 + m22 * x20; 22483287d42SBarry Smith x[17] -= m3 * x16 + m8 * x17 + m13 * x18 + m18 * x19 + m23 * x20; 22583287d42SBarry Smith x[18] -= m4 * x16 + m9 * x17 + m14 * x18 + m19 * x19 + m24 * x20; 22683287d42SBarry Smith x[19] -= m5 * x16 + m10 * x17 + m15 * x18 + m20 * x19 + m25 * x20; 22783287d42SBarry Smith 22883287d42SBarry Smith x[20] -= m1 * x21 + m6 * x22 + m11 * x23 + m16 * x24 + m21 * x25; 22983287d42SBarry Smith x[21] -= m2 * x21 + m7 * x22 + m12 * x23 + m17 * x24 + m22 * x25; 23083287d42SBarry Smith x[22] -= m3 * x21 + m8 * x22 + m13 * x23 + m18 * x24 + m23 * x25; 23183287d42SBarry Smith x[23] -= m4 * x21 + m9 * x22 + m14 * x23 + m19 * x24 + m24 * x25; 23283287d42SBarry Smith x[24] -= m5 * x21 + m10 * x22 + m15 * x23 + m20 * x24 + m25 * x25; 23383287d42SBarry Smith 23483287d42SBarry Smith pv += 25; 23583287d42SBarry Smith } 2369566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(250.0 * nz + 225.0)); 23783287d42SBarry Smith } 23883287d42SBarry Smith row = *ajtmp++; 23983287d42SBarry Smith } 24083287d42SBarry Smith /* finished row so stick it into b->a */ 24183287d42SBarry Smith pv = ba + 25 * bi[i]; 24283287d42SBarry Smith pj = bj + bi[i]; 24383287d42SBarry Smith nz = bi[i + 1] - bi[i]; 24483287d42SBarry Smith for (j = 0; j < nz; j++) { 2458397fe1aSBarry Smith #if defined(PETSC_USE_MEMCPY) 2469566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(pv, rtmp + 25 * pj[j], 25)); 2478397fe1aSBarry Smith #else 24883287d42SBarry Smith x = rtmp + 25 * pj[j]; 2499371c9d4SSatish Balay pv[0] = x[0]; 2509371c9d4SSatish Balay pv[1] = x[1]; 2519371c9d4SSatish Balay pv[2] = x[2]; 2529371c9d4SSatish Balay pv[3] = x[3]; 2539371c9d4SSatish Balay pv[4] = x[4]; 2549371c9d4SSatish Balay pv[5] = x[5]; 2559371c9d4SSatish Balay pv[6] = x[6]; 2569371c9d4SSatish Balay pv[7] = x[7]; 2579371c9d4SSatish Balay pv[8] = x[8]; 2589371c9d4SSatish Balay pv[9] = x[9]; 2599371c9d4SSatish Balay pv[10] = x[10]; 2609371c9d4SSatish Balay pv[11] = x[11]; 2619371c9d4SSatish Balay pv[12] = x[12]; 2629371c9d4SSatish Balay pv[13] = x[13]; 2639371c9d4SSatish Balay pv[14] = x[14]; 2649371c9d4SSatish Balay pv[15] = x[15]; 2659371c9d4SSatish Balay pv[16] = x[16]; 2669371c9d4SSatish Balay pv[17] = x[17]; 2679371c9d4SSatish Balay pv[18] = x[18]; 2689371c9d4SSatish Balay pv[19] = x[19]; 2699371c9d4SSatish Balay pv[20] = x[20]; 2709371c9d4SSatish Balay pv[21] = x[21]; 2719371c9d4SSatish Balay pv[22] = x[22]; 2729371c9d4SSatish Balay pv[23] = x[23]; 2739371c9d4SSatish Balay pv[24] = x[24]; 2748397fe1aSBarry Smith #endif 27583287d42SBarry Smith pv += 25; 27683287d42SBarry Smith } 27783287d42SBarry Smith /* invert diagonal block */ 27883287d42SBarry Smith w = ba + 25 * diag_offset[i]; 2799566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_5(w, ipvt, work, shift, allowzeropivot, &zeropivotdetected)); 2807b6c816cSBarry Smith if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 28183287d42SBarry Smith } 28283287d42SBarry Smith 2839566063dSJacob Faibussowitsch PetscCall(PetscFree(rtmp)); 2849566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isicol, &ic)); 2859566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &r)); 28626fbe8dcSKarl Rupp 28706e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_5_inplace; 28806e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_inplace; 28983287d42SBarry Smith C->assembled = PETSC_TRUE; 29026fbe8dcSKarl Rupp 2919566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.333333333333 * 5 * 5 * 5 * b->mbs)); /* from inverting diagonal blocks */ 29283287d42SBarry Smith PetscFunctionReturn(0); 29383287d42SBarry Smith } 2940deeaf61SShri Abhyankar 2954dd39f65SShri Abhyankar /* MatLUFactorNumeric_SeqBAIJ_5 - 2964dd39f65SShri Abhyankar copied from MatLUFactorNumeric_SeqBAIJ_N_inplace() and manually re-implemented 29796b95a6bSBarry Smith PetscKernel_A_gets_A_times_B() 29896b95a6bSBarry Smith PetscKernel_A_gets_A_minus_B_times_C() 29996b95a6bSBarry Smith PetscKernel_A_gets_inverse_A() 3000deeaf61SShri Abhyankar */ 301c0c7eb62SShri Abhyankar 3029371c9d4SSatish Balay PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5(Mat B, Mat A, const MatFactorInfo *info) { 3030deeaf61SShri Abhyankar Mat C = B; 3040deeaf61SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data; 3050deeaf61SShri Abhyankar IS isrow = b->row, isicol = b->icol; 3065a586d82SBarry Smith const PetscInt *r, *ic; 307bbd65245SShri Abhyankar PetscInt i, j, k, nz, nzL, row; 308bbd65245SShri Abhyankar const PetscInt n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j; 309bbd65245SShri Abhyankar const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag, *pj, bs2 = a->bs2; 310766f9fbaSBarry Smith MatScalar *rtmp, *pc, *mwork, *v, *pv, *aa = a->a, work[25]; 311bbd65245SShri Abhyankar PetscInt flg, ipvt[5]; 312182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 313a455e926SHong Zhang PetscBool allowzeropivot, zeropivotdetected; 3140deeaf61SShri Abhyankar 3150deeaf61SShri Abhyankar PetscFunctionBegin; 3160164db54SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 3179566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrow, &r)); 3189566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isicol, &ic)); 3190deeaf61SShri Abhyankar 3200deeaf61SShri Abhyankar /* generate work space needed by the factorization */ 3219566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork)); 3229566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(rtmp, bs2 * n)); 3230deeaf61SShri Abhyankar 3240deeaf61SShri Abhyankar for (i = 0; i < n; i++) { 3250deeaf61SShri Abhyankar /* zero rtmp */ 3260deeaf61SShri Abhyankar /* L part */ 3270deeaf61SShri Abhyankar nz = bi[i + 1] - bi[i]; 3280deeaf61SShri Abhyankar bjtmp = bj + bi[i]; 329*48a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2)); 3300deeaf61SShri Abhyankar 3310deeaf61SShri Abhyankar /* U part */ 33278bb4007SShri Abhyankar nz = bdiag[i] - bdiag[i + 1]; 33378bb4007SShri Abhyankar bjtmp = bj + bdiag[i + 1] + 1; 334*48a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2)); 33578bb4007SShri Abhyankar 33678bb4007SShri Abhyankar /* load in initial (unfactored row) */ 33778bb4007SShri Abhyankar nz = ai[r[i] + 1] - ai[r[i]]; 33878bb4007SShri Abhyankar ajtmp = aj + ai[r[i]]; 33978bb4007SShri Abhyankar v = aa + bs2 * ai[r[i]]; 340*48a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ic[ajtmp[j]], v + bs2 * j, bs2)); 34178bb4007SShri Abhyankar 34278bb4007SShri Abhyankar /* elimination */ 34378bb4007SShri Abhyankar bjtmp = bj + bi[i]; 34478bb4007SShri Abhyankar nzL = bi[i + 1] - bi[i]; 34578bb4007SShri Abhyankar for (k = 0; k < nzL; k++) { 34678bb4007SShri Abhyankar row = bjtmp[k]; 34778bb4007SShri Abhyankar pc = rtmp + bs2 * row; 348c35f09e5SBarry Smith for (flg = 0, j = 0; j < bs2; j++) { 349c35f09e5SBarry Smith if (pc[j] != 0.0) { 350c35f09e5SBarry Smith flg = 1; 351c35f09e5SBarry Smith break; 352c35f09e5SBarry Smith } 353c35f09e5SBarry Smith } 35478bb4007SShri Abhyankar if (flg) { 35578bb4007SShri Abhyankar pv = b->a + bs2 * bdiag[row]; 35696b95a6bSBarry Smith /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 3579566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_A_times_B_5(pc, pv, mwork)); 35878bb4007SShri Abhyankar 359a5b23f4aSJose E. Roman pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */ 36078bb4007SShri Abhyankar pv = b->a + bs2 * (bdiag[row + 1] + 1); 36178bb4007SShri Abhyankar nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */ 36278bb4007SShri Abhyankar for (j = 0; j < nz; j++) { 36396b95a6bSBarry Smith /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 36478bb4007SShri Abhyankar /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 36578bb4007SShri Abhyankar v = rtmp + bs2 * pj[j]; 3669566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_A_minus_B_times_C_5(v, pc, pv)); 36778bb4007SShri Abhyankar pv += bs2; 36878bb4007SShri Abhyankar } 3699566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(250.0 * nz + 225)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 37078bb4007SShri Abhyankar } 37178bb4007SShri Abhyankar } 37278bb4007SShri Abhyankar 37378bb4007SShri Abhyankar /* finished row so stick it into b->a */ 37478bb4007SShri Abhyankar /* L part */ 37578bb4007SShri Abhyankar pv = b->a + bs2 * bi[i]; 37678bb4007SShri Abhyankar pj = b->j + bi[i]; 37778bb4007SShri Abhyankar nz = bi[i + 1] - bi[i]; 378*48a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2)); 37978bb4007SShri Abhyankar 380a5b23f4aSJose E. Roman /* Mark diagonal and invert diagonal for simpler triangular solves */ 38178bb4007SShri Abhyankar pv = b->a + bs2 * bdiag[i]; 38278bb4007SShri Abhyankar pj = b->j + bdiag[i]; 3839566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2)); 3849566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_5(pv, ipvt, work, shift, allowzeropivot, &zeropivotdetected)); 3857b6c816cSBarry Smith if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 38678bb4007SShri Abhyankar 38778bb4007SShri Abhyankar /* U part */ 38878bb4007SShri Abhyankar pv = b->a + bs2 * (bdiag[i + 1] + 1); 38978bb4007SShri Abhyankar pj = b->j + bdiag[i + 1] + 1; 39078bb4007SShri Abhyankar nz = bdiag[i] - bdiag[i + 1] - 1; 391*48a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2)); 39278bb4007SShri Abhyankar } 39378bb4007SShri Abhyankar 3949566063dSJacob Faibussowitsch PetscCall(PetscFree2(rtmp, mwork)); 3959566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isicol, &ic)); 3969566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrow, &r)); 39726fbe8dcSKarl Rupp 3984dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_5; 3994dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5; 40078bb4007SShri Abhyankar C->assembled = PETSC_TRUE; 40126fbe8dcSKarl Rupp 4029566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.333333333333 * 5 * 5 * 5 * n)); /* from inverting diagonal blocks */ 40378bb4007SShri Abhyankar PetscFunctionReturn(0); 40478bb4007SShri Abhyankar } 40578bb4007SShri Abhyankar 4060deeaf61SShri Abhyankar /* 4070deeaf61SShri Abhyankar Version for when blocks are 5 by 5 Using natural ordering 4080deeaf61SShri Abhyankar */ 4099371c9d4SSatish Balay PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering_inplace(Mat C, Mat A, const MatFactorInfo *info) { 4100deeaf61SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data; 411766f9fbaSBarry Smith PetscInt i, j, n = a->mbs, *bi = b->i, *bj = b->j, ipvt[5]; 4120deeaf61SShri Abhyankar PetscInt *ajtmpold, *ajtmp, nz, row; 4130deeaf61SShri Abhyankar PetscInt *diag_offset = b->diag, *ai = a->i, *aj = a->j, *pj; 4140deeaf61SShri Abhyankar MatScalar *pv, *v, *rtmp, *pc, *w, *x; 4150deeaf61SShri Abhyankar MatScalar x1, x2, x3, x4, x5, x6, x7, x8, x9, x10, x11, x12, x13, x14, x15; 4160deeaf61SShri Abhyankar MatScalar x16, x17, x18, x19, x20, x21, x22, x23, x24, x25; 4170deeaf61SShri Abhyankar MatScalar p1, p2, p3, p4, p5, p6, p7, p8, p9, p10, p11, p12, p13, p14, p15; 4180deeaf61SShri Abhyankar MatScalar p16, p17, p18, p19, p20, p21, p22, p23, p24, p25; 4190deeaf61SShri Abhyankar MatScalar m1, m2, m3, m4, m5, m6, m7, m8, m9, m10, m11, m12, m13, m14, m15; 4200deeaf61SShri Abhyankar MatScalar m16, m17, m18, m19, m20, m21, m22, m23, m24, m25; 421766f9fbaSBarry Smith MatScalar *ba = b->a, *aa = a->a, work[25]; 422182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 423a455e926SHong Zhang PetscBool allowzeropivot, zeropivotdetected; 4240deeaf61SShri Abhyankar 4250deeaf61SShri Abhyankar PetscFunctionBegin; 4260164db54SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 4279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(25 * (n + 1), &rtmp)); 4280deeaf61SShri Abhyankar for (i = 0; i < n; i++) { 4290deeaf61SShri Abhyankar nz = bi[i + 1] - bi[i]; 4300deeaf61SShri Abhyankar ajtmp = bj + bi[i]; 4310deeaf61SShri Abhyankar for (j = 0; j < nz; j++) { 4320deeaf61SShri Abhyankar x = rtmp + 25 * ajtmp[j]; 4330deeaf61SShri Abhyankar x[0] = x[1] = x[2] = x[3] = x[4] = x[5] = x[6] = x[7] = x[8] = x[9] = 0.0; 4340deeaf61SShri Abhyankar x[10] = x[11] = x[12] = x[13] = x[14] = x[15] = 0.0; 4350deeaf61SShri Abhyankar x[16] = x[17] = x[18] = x[19] = x[20] = x[21] = x[22] = x[23] = x[24] = 0.0; 4360deeaf61SShri Abhyankar } 4370deeaf61SShri Abhyankar /* load in initial (unfactored row) */ 4380deeaf61SShri Abhyankar nz = ai[i + 1] - ai[i]; 4390deeaf61SShri Abhyankar ajtmpold = aj + ai[i]; 4400deeaf61SShri Abhyankar v = aa + 25 * ai[i]; 4410deeaf61SShri Abhyankar for (j = 0; j < nz; j++) { 4420deeaf61SShri Abhyankar x = rtmp + 25 * ajtmpold[j]; 4439371c9d4SSatish Balay x[0] = v[0]; 4449371c9d4SSatish Balay x[1] = v[1]; 4459371c9d4SSatish Balay x[2] = v[2]; 4469371c9d4SSatish Balay x[3] = v[3]; 4479371c9d4SSatish Balay x[4] = v[4]; 4489371c9d4SSatish Balay x[5] = v[5]; 4499371c9d4SSatish Balay x[6] = v[6]; 4509371c9d4SSatish Balay x[7] = v[7]; 4519371c9d4SSatish Balay x[8] = v[8]; 4529371c9d4SSatish Balay x[9] = v[9]; 4539371c9d4SSatish Balay x[10] = v[10]; 4549371c9d4SSatish Balay x[11] = v[11]; 4559371c9d4SSatish Balay x[12] = v[12]; 4569371c9d4SSatish Balay x[13] = v[13]; 4579371c9d4SSatish Balay x[14] = v[14]; 4589371c9d4SSatish Balay x[15] = v[15]; 4599371c9d4SSatish Balay x[16] = v[16]; 4609371c9d4SSatish Balay x[17] = v[17]; 4619371c9d4SSatish Balay x[18] = v[18]; 4629371c9d4SSatish Balay x[19] = v[19]; 4639371c9d4SSatish Balay x[20] = v[20]; 4649371c9d4SSatish Balay x[21] = v[21]; 4659371c9d4SSatish Balay x[22] = v[22]; 4669371c9d4SSatish Balay x[23] = v[23]; 4670deeaf61SShri Abhyankar x[24] = v[24]; 4680deeaf61SShri Abhyankar v += 25; 4690deeaf61SShri Abhyankar } 4700deeaf61SShri Abhyankar row = *ajtmp++; 4710deeaf61SShri Abhyankar while (row < i) { 4720deeaf61SShri Abhyankar pc = rtmp + 25 * row; 4739371c9d4SSatish Balay p1 = pc[0]; 4749371c9d4SSatish Balay p2 = pc[1]; 4759371c9d4SSatish Balay p3 = pc[2]; 4769371c9d4SSatish Balay p4 = pc[3]; 4779371c9d4SSatish Balay p5 = pc[4]; 4789371c9d4SSatish Balay p6 = pc[5]; 4799371c9d4SSatish Balay p7 = pc[6]; 4809371c9d4SSatish Balay p8 = pc[7]; 4819371c9d4SSatish Balay p9 = pc[8]; 4829371c9d4SSatish Balay p10 = pc[9]; 4839371c9d4SSatish Balay p11 = pc[10]; 4849371c9d4SSatish Balay p12 = pc[11]; 4859371c9d4SSatish Balay p13 = pc[12]; 4869371c9d4SSatish Balay p14 = pc[13]; 4879371c9d4SSatish Balay p15 = pc[14]; 4889371c9d4SSatish Balay p16 = pc[15]; 4899371c9d4SSatish Balay p17 = pc[16]; 4909371c9d4SSatish Balay p18 = pc[17]; 4919371c9d4SSatish Balay p19 = pc[18]; 4929371c9d4SSatish Balay p20 = pc[19]; 4939371c9d4SSatish Balay p21 = pc[20]; 4949371c9d4SSatish Balay p22 = pc[21]; 4959371c9d4SSatish Balay p23 = pc[22]; 4969371c9d4SSatish Balay p24 = pc[23]; 4979371c9d4SSatish Balay p25 = pc[24]; 4989371c9d4SSatish Balay if (p1 != 0.0 || p2 != 0.0 || p3 != 0.0 || p4 != 0.0 || p5 != 0.0 || p6 != 0.0 || p7 != 0.0 || p8 != 0.0 || p9 != 0.0 || p10 != 0.0 || p11 != 0.0 || p12 != 0.0 || p13 != 0.0 || p14 != 0.0 || p15 != 0.0 || p16 != 0.0 || p17 != 0.0 || p18 != 0.0 || p19 != 0.0 || p20 != 0.0 || p21 != 0.0 || p22 != 0.0 || p23 != 0.0 || p24 != 0.0 || p25 != 0.0) { 4990deeaf61SShri Abhyankar pv = ba + 25 * diag_offset[row]; 5000deeaf61SShri Abhyankar pj = bj + diag_offset[row] + 1; 5019371c9d4SSatish Balay x1 = pv[0]; 5029371c9d4SSatish Balay x2 = pv[1]; 5039371c9d4SSatish Balay x3 = pv[2]; 5049371c9d4SSatish Balay x4 = pv[3]; 5059371c9d4SSatish Balay x5 = pv[4]; 5069371c9d4SSatish Balay x6 = pv[5]; 5079371c9d4SSatish Balay x7 = pv[6]; 5089371c9d4SSatish Balay x8 = pv[7]; 5099371c9d4SSatish Balay x9 = pv[8]; 5109371c9d4SSatish Balay x10 = pv[9]; 5119371c9d4SSatish Balay x11 = pv[10]; 5129371c9d4SSatish Balay x12 = pv[11]; 5139371c9d4SSatish Balay x13 = pv[12]; 5149371c9d4SSatish Balay x14 = pv[13]; 5159371c9d4SSatish Balay x15 = pv[14]; 5169371c9d4SSatish Balay x16 = pv[15]; 5179371c9d4SSatish Balay x17 = pv[16]; 5189371c9d4SSatish Balay x18 = pv[17]; 5199371c9d4SSatish Balay x19 = pv[18]; 5209371c9d4SSatish Balay x20 = pv[19]; 5219371c9d4SSatish Balay x21 = pv[20]; 5229371c9d4SSatish Balay x22 = pv[21]; 5239371c9d4SSatish Balay x23 = pv[22]; 5249371c9d4SSatish Balay x24 = pv[23]; 5250deeaf61SShri Abhyankar x25 = pv[24]; 5260deeaf61SShri Abhyankar pc[0] = m1 = p1 * x1 + p6 * x2 + p11 * x3 + p16 * x4 + p21 * x5; 5270deeaf61SShri Abhyankar pc[1] = m2 = p2 * x1 + p7 * x2 + p12 * x3 + p17 * x4 + p22 * x5; 5280deeaf61SShri Abhyankar pc[2] = m3 = p3 * x1 + p8 * x2 + p13 * x3 + p18 * x4 + p23 * x5; 5290deeaf61SShri Abhyankar pc[3] = m4 = p4 * x1 + p9 * x2 + p14 * x3 + p19 * x4 + p24 * x5; 5300deeaf61SShri Abhyankar pc[4] = m5 = p5 * x1 + p10 * x2 + p15 * x3 + p20 * x4 + p25 * x5; 5310deeaf61SShri Abhyankar 5320deeaf61SShri Abhyankar pc[5] = m6 = p1 * x6 + p6 * x7 + p11 * x8 + p16 * x9 + p21 * x10; 5330deeaf61SShri Abhyankar pc[6] = m7 = p2 * x6 + p7 * x7 + p12 * x8 + p17 * x9 + p22 * x10; 5340deeaf61SShri Abhyankar pc[7] = m8 = p3 * x6 + p8 * x7 + p13 * x8 + p18 * x9 + p23 * x10; 5350deeaf61SShri Abhyankar pc[8] = m9 = p4 * x6 + p9 * x7 + p14 * x8 + p19 * x9 + p24 * x10; 5360deeaf61SShri Abhyankar pc[9] = m10 = p5 * x6 + p10 * x7 + p15 * x8 + p20 * x9 + p25 * x10; 5370deeaf61SShri Abhyankar 5380deeaf61SShri Abhyankar pc[10] = m11 = p1 * x11 + p6 * x12 + p11 * x13 + p16 * x14 + p21 * x15; 5390deeaf61SShri Abhyankar pc[11] = m12 = p2 * x11 + p7 * x12 + p12 * x13 + p17 * x14 + p22 * x15; 5400deeaf61SShri Abhyankar pc[12] = m13 = p3 * x11 + p8 * x12 + p13 * x13 + p18 * x14 + p23 * x15; 5410deeaf61SShri Abhyankar pc[13] = m14 = p4 * x11 + p9 * x12 + p14 * x13 + p19 * x14 + p24 * x15; 5420deeaf61SShri Abhyankar pc[14] = m15 = p5 * x11 + p10 * x12 + p15 * x13 + p20 * x14 + p25 * x15; 5430deeaf61SShri Abhyankar 5440deeaf61SShri Abhyankar pc[15] = m16 = p1 * x16 + p6 * x17 + p11 * x18 + p16 * x19 + p21 * x20; 5450deeaf61SShri Abhyankar pc[16] = m17 = p2 * x16 + p7 * x17 + p12 * x18 + p17 * x19 + p22 * x20; 5460deeaf61SShri Abhyankar pc[17] = m18 = p3 * x16 + p8 * x17 + p13 * x18 + p18 * x19 + p23 * x20; 5470deeaf61SShri Abhyankar pc[18] = m19 = p4 * x16 + p9 * x17 + p14 * x18 + p19 * x19 + p24 * x20; 5480deeaf61SShri Abhyankar pc[19] = m20 = p5 * x16 + p10 * x17 + p15 * x18 + p20 * x19 + p25 * x20; 5490deeaf61SShri Abhyankar 5500deeaf61SShri Abhyankar pc[20] = m21 = p1 * x21 + p6 * x22 + p11 * x23 + p16 * x24 + p21 * x25; 5510deeaf61SShri Abhyankar pc[21] = m22 = p2 * x21 + p7 * x22 + p12 * x23 + p17 * x24 + p22 * x25; 5520deeaf61SShri Abhyankar pc[22] = m23 = p3 * x21 + p8 * x22 + p13 * x23 + p18 * x24 + p23 * x25; 5530deeaf61SShri Abhyankar pc[23] = m24 = p4 * x21 + p9 * x22 + p14 * x23 + p19 * x24 + p24 * x25; 5540deeaf61SShri Abhyankar pc[24] = m25 = p5 * x21 + p10 * x22 + p15 * x23 + p20 * x24 + p25 * x25; 5550deeaf61SShri Abhyankar 5560deeaf61SShri Abhyankar nz = bi[row + 1] - diag_offset[row] - 1; 5570deeaf61SShri Abhyankar pv += 25; 5580deeaf61SShri Abhyankar for (j = 0; j < nz; j++) { 5599371c9d4SSatish Balay x1 = pv[0]; 5609371c9d4SSatish Balay x2 = pv[1]; 5619371c9d4SSatish Balay x3 = pv[2]; 5629371c9d4SSatish Balay x4 = pv[3]; 5639371c9d4SSatish Balay x5 = pv[4]; 5649371c9d4SSatish Balay x6 = pv[5]; 5659371c9d4SSatish Balay x7 = pv[6]; 5669371c9d4SSatish Balay x8 = pv[7]; 5679371c9d4SSatish Balay x9 = pv[8]; 5689371c9d4SSatish Balay x10 = pv[9]; 5699371c9d4SSatish Balay x11 = pv[10]; 5709371c9d4SSatish Balay x12 = pv[11]; 5719371c9d4SSatish Balay x13 = pv[12]; 5729371c9d4SSatish Balay x14 = pv[13]; 5739371c9d4SSatish Balay x15 = pv[14]; 5749371c9d4SSatish Balay x16 = pv[15]; 5759371c9d4SSatish Balay x17 = pv[16]; 5769371c9d4SSatish Balay x18 = pv[17]; 5779371c9d4SSatish Balay x19 = pv[18]; 5789371c9d4SSatish Balay x20 = pv[19]; 5799371c9d4SSatish Balay x21 = pv[20]; 5809371c9d4SSatish Balay x22 = pv[21]; 5819371c9d4SSatish Balay x23 = pv[22]; 5829371c9d4SSatish Balay x24 = pv[23]; 5839371c9d4SSatish Balay x25 = pv[24]; 5840deeaf61SShri Abhyankar x = rtmp + 25 * pj[j]; 5850deeaf61SShri Abhyankar x[0] -= m1 * x1 + m6 * x2 + m11 * x3 + m16 * x4 + m21 * x5; 5860deeaf61SShri Abhyankar x[1] -= m2 * x1 + m7 * x2 + m12 * x3 + m17 * x4 + m22 * x5; 5870deeaf61SShri Abhyankar x[2] -= m3 * x1 + m8 * x2 + m13 * x3 + m18 * x4 + m23 * x5; 5880deeaf61SShri Abhyankar x[3] -= m4 * x1 + m9 * x2 + m14 * x3 + m19 * x4 + m24 * x5; 5890deeaf61SShri Abhyankar x[4] -= m5 * x1 + m10 * x2 + m15 * x3 + m20 * x4 + m25 * x5; 5900deeaf61SShri Abhyankar 5910deeaf61SShri Abhyankar x[5] -= m1 * x6 + m6 * x7 + m11 * x8 + m16 * x9 + m21 * x10; 5920deeaf61SShri Abhyankar x[6] -= m2 * x6 + m7 * x7 + m12 * x8 + m17 * x9 + m22 * x10; 5930deeaf61SShri Abhyankar x[7] -= m3 * x6 + m8 * x7 + m13 * x8 + m18 * x9 + m23 * x10; 5940deeaf61SShri Abhyankar x[8] -= m4 * x6 + m9 * x7 + m14 * x8 + m19 * x9 + m24 * x10; 5950deeaf61SShri Abhyankar x[9] -= m5 * x6 + m10 * x7 + m15 * x8 + m20 * x9 + m25 * x10; 5960deeaf61SShri Abhyankar 5970deeaf61SShri Abhyankar x[10] -= m1 * x11 + m6 * x12 + m11 * x13 + m16 * x14 + m21 * x15; 5980deeaf61SShri Abhyankar x[11] -= m2 * x11 + m7 * x12 + m12 * x13 + m17 * x14 + m22 * x15; 5990deeaf61SShri Abhyankar x[12] -= m3 * x11 + m8 * x12 + m13 * x13 + m18 * x14 + m23 * x15; 6000deeaf61SShri Abhyankar x[13] -= m4 * x11 + m9 * x12 + m14 * x13 + m19 * x14 + m24 * x15; 6010deeaf61SShri Abhyankar x[14] -= m5 * x11 + m10 * x12 + m15 * x13 + m20 * x14 + m25 * x15; 6020deeaf61SShri Abhyankar 6030deeaf61SShri Abhyankar x[15] -= m1 * x16 + m6 * x17 + m11 * x18 + m16 * x19 + m21 * x20; 6040deeaf61SShri Abhyankar x[16] -= m2 * x16 + m7 * x17 + m12 * x18 + m17 * x19 + m22 * x20; 6050deeaf61SShri Abhyankar x[17] -= m3 * x16 + m8 * x17 + m13 * x18 + m18 * x19 + m23 * x20; 6060deeaf61SShri Abhyankar x[18] -= m4 * x16 + m9 * x17 + m14 * x18 + m19 * x19 + m24 * x20; 6070deeaf61SShri Abhyankar x[19] -= m5 * x16 + m10 * x17 + m15 * x18 + m20 * x19 + m25 * x20; 6080deeaf61SShri Abhyankar 6090deeaf61SShri Abhyankar x[20] -= m1 * x21 + m6 * x22 + m11 * x23 + m16 * x24 + m21 * x25; 6100deeaf61SShri Abhyankar x[21] -= m2 * x21 + m7 * x22 + m12 * x23 + m17 * x24 + m22 * x25; 6110deeaf61SShri Abhyankar x[22] -= m3 * x21 + m8 * x22 + m13 * x23 + m18 * x24 + m23 * x25; 6120deeaf61SShri Abhyankar x[23] -= m4 * x21 + m9 * x22 + m14 * x23 + m19 * x24 + m24 * x25; 6130deeaf61SShri Abhyankar x[24] -= m5 * x21 + m10 * x22 + m15 * x23 + m20 * x24 + m25 * x25; 6140deeaf61SShri Abhyankar pv += 25; 6150deeaf61SShri Abhyankar } 6169566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(250.0 * nz + 225.0)); 6170deeaf61SShri Abhyankar } 6180deeaf61SShri Abhyankar row = *ajtmp++; 6190deeaf61SShri Abhyankar } 6200deeaf61SShri Abhyankar /* finished row so stick it into b->a */ 6210deeaf61SShri Abhyankar pv = ba + 25 * bi[i]; 6220deeaf61SShri Abhyankar pj = bj + bi[i]; 6230deeaf61SShri Abhyankar nz = bi[i + 1] - bi[i]; 6240deeaf61SShri Abhyankar for (j = 0; j < nz; j++) { 6250deeaf61SShri Abhyankar x = rtmp + 25 * pj[j]; 6269371c9d4SSatish Balay pv[0] = x[0]; 6279371c9d4SSatish Balay pv[1] = x[1]; 6289371c9d4SSatish Balay pv[2] = x[2]; 6299371c9d4SSatish Balay pv[3] = x[3]; 6309371c9d4SSatish Balay pv[4] = x[4]; 6319371c9d4SSatish Balay pv[5] = x[5]; 6329371c9d4SSatish Balay pv[6] = x[6]; 6339371c9d4SSatish Balay pv[7] = x[7]; 6349371c9d4SSatish Balay pv[8] = x[8]; 6359371c9d4SSatish Balay pv[9] = x[9]; 6369371c9d4SSatish Balay pv[10] = x[10]; 6379371c9d4SSatish Balay pv[11] = x[11]; 6389371c9d4SSatish Balay pv[12] = x[12]; 6399371c9d4SSatish Balay pv[13] = x[13]; 6409371c9d4SSatish Balay pv[14] = x[14]; 6419371c9d4SSatish Balay pv[15] = x[15]; 6429371c9d4SSatish Balay pv[16] = x[16]; 6439371c9d4SSatish Balay pv[17] = x[17]; 6449371c9d4SSatish Balay pv[18] = x[18]; 6459371c9d4SSatish Balay pv[19] = x[19]; 6469371c9d4SSatish Balay pv[20] = x[20]; 6479371c9d4SSatish Balay pv[21] = x[21]; 6489371c9d4SSatish Balay pv[22] = x[22]; 6499371c9d4SSatish Balay pv[23] = x[23]; 6509371c9d4SSatish Balay pv[24] = x[24]; 6510deeaf61SShri Abhyankar pv += 25; 6520deeaf61SShri Abhyankar } 6530deeaf61SShri Abhyankar /* invert diagonal block */ 6540deeaf61SShri Abhyankar w = ba + 25 * diag_offset[i]; 6559566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_5(w, ipvt, work, shift, allowzeropivot, &zeropivotdetected)); 6567b6c816cSBarry Smith if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 6570deeaf61SShri Abhyankar } 6580deeaf61SShri Abhyankar 6599566063dSJacob Faibussowitsch PetscCall(PetscFree(rtmp)); 66026fbe8dcSKarl Rupp 66106e38f1dSHong Zhang C->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering_inplace; 66206e38f1dSHong Zhang C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering_inplace; 6630deeaf61SShri Abhyankar C->assembled = PETSC_TRUE; 66426fbe8dcSKarl Rupp 6659566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.333333333333 * 5 * 5 * 5 * b->mbs)); /* from inverting diagonal blocks */ 6660deeaf61SShri Abhyankar PetscFunctionReturn(0); 6670deeaf61SShri Abhyankar } 6680deeaf61SShri Abhyankar 6699371c9d4SSatish Balay PetscErrorCode MatLUFactorNumeric_SeqBAIJ_5_NaturalOrdering(Mat B, Mat A, const MatFactorInfo *info) { 6700deeaf61SShri Abhyankar Mat C = B; 6710deeaf61SShri Abhyankar Mat_SeqBAIJ *a = (Mat_SeqBAIJ *)A->data, *b = (Mat_SeqBAIJ *)C->data; 672bbd65245SShri Abhyankar PetscInt i, j, k, nz, nzL, row; 673766f9fbaSBarry Smith const PetscInt n = a->mbs, *ai = a->i, *aj = a->j, *bi = b->i, *bj = b->j; 674bbd65245SShri Abhyankar const PetscInt *ajtmp, *bjtmp, *bdiag = b->diag, *pj, bs2 = a->bs2; 675bbd65245SShri Abhyankar MatScalar *rtmp, *pc, *mwork, *v, *vv, *pv, *aa = a->a, work[25]; 676bbd65245SShri Abhyankar PetscInt flg, ipvt[5]; 677182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 678a455e926SHong Zhang PetscBool allowzeropivot, zeropivotdetected; 6790deeaf61SShri Abhyankar 6800deeaf61SShri Abhyankar PetscFunctionBegin; 6810164db54SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 6820164db54SHong Zhang 6830deeaf61SShri Abhyankar /* generate work space needed by the factorization */ 6849566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(bs2 * n, &rtmp, bs2, &mwork)); 6859566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(rtmp, bs2 * n)); 6860deeaf61SShri Abhyankar 6870deeaf61SShri Abhyankar for (i = 0; i < n; i++) { 6880deeaf61SShri Abhyankar /* zero rtmp */ 6890deeaf61SShri Abhyankar /* L part */ 6900deeaf61SShri Abhyankar nz = bi[i + 1] - bi[i]; 6910deeaf61SShri Abhyankar bjtmp = bj + bi[i]; 692*48a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2)); 6930deeaf61SShri Abhyankar 6940deeaf61SShri Abhyankar /* U part */ 69553cca76cSShri Abhyankar nz = bdiag[i] - bdiag[i + 1]; 69653cca76cSShri Abhyankar bjtmp = bj + bdiag[i + 1] + 1; 697*48a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArrayzero(rtmp + bs2 * bjtmp[j], bs2)); 69853cca76cSShri Abhyankar 69953cca76cSShri Abhyankar /* load in initial (unfactored row) */ 70053cca76cSShri Abhyankar nz = ai[i + 1] - ai[i]; 70153cca76cSShri Abhyankar ajtmp = aj + ai[i]; 70253cca76cSShri Abhyankar v = aa + bs2 * ai[i]; 703*48a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(rtmp + bs2 * ajtmp[j], v + bs2 * j, bs2)); 70453cca76cSShri Abhyankar 70553cca76cSShri Abhyankar /* elimination */ 70653cca76cSShri Abhyankar bjtmp = bj + bi[i]; 70753cca76cSShri Abhyankar nzL = bi[i + 1] - bi[i]; 70853cca76cSShri Abhyankar for (k = 0; k < nzL; k++) { 70953cca76cSShri Abhyankar row = bjtmp[k]; 71053cca76cSShri Abhyankar pc = rtmp + bs2 * row; 711c35f09e5SBarry Smith for (flg = 0, j = 0; j < bs2; j++) { 712c35f09e5SBarry Smith if (pc[j] != 0.0) { 713c35f09e5SBarry Smith flg = 1; 714c35f09e5SBarry Smith break; 715c35f09e5SBarry Smith } 716c35f09e5SBarry Smith } 71753cca76cSShri Abhyankar if (flg) { 71853cca76cSShri Abhyankar pv = b->a + bs2 * bdiag[row]; 71996b95a6bSBarry Smith /* PetscKernel_A_gets_A_times_B(bs,pc,pv,mwork); *pc = *pc * (*pv); */ 7209566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_A_times_B_5(pc, pv, mwork)); 72153cca76cSShri Abhyankar 722a5b23f4aSJose E. Roman pj = b->j + bdiag[row + 1] + 1; /* beginning of U(row,:) */ 72353cca76cSShri Abhyankar pv = b->a + bs2 * (bdiag[row + 1] + 1); 72453cca76cSShri Abhyankar nz = bdiag[row] - bdiag[row + 1] - 1; /* num of entries inU(row,:), excluding diag */ 72553cca76cSShri Abhyankar for (j = 0; j < nz; j++) { 72696b95a6bSBarry Smith /* PetscKernel_A_gets_A_minus_B_times_C(bs,rtmp+bs2*pj[j],pc,pv+bs2*j); */ 72753cca76cSShri Abhyankar /* rtmp+bs2*pj[j] = rtmp+bs2*pj[j] - (*pc)*(pv+bs2*j) */ 728766f9fbaSBarry Smith vv = rtmp + bs2 * pj[j]; 7299566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_A_minus_B_times_C_5(vv, pc, pv)); 73053cca76cSShri Abhyankar pv += bs2; 73153cca76cSShri Abhyankar } 7329566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(250.0 * nz + 225)); /* flops = 2*bs^3*nz + 2*bs^3 - bs2) */ 73353cca76cSShri Abhyankar } 73453cca76cSShri Abhyankar } 73553cca76cSShri Abhyankar 73653cca76cSShri Abhyankar /* finished row so stick it into b->a */ 73753cca76cSShri Abhyankar /* L part */ 73853cca76cSShri Abhyankar pv = b->a + bs2 * bi[i]; 73953cca76cSShri Abhyankar pj = b->j + bi[i]; 74053cca76cSShri Abhyankar nz = bi[i + 1] - bi[i]; 741*48a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2)); 74253cca76cSShri Abhyankar 743a5b23f4aSJose E. Roman /* Mark diagonal and invert diagonal for simpler triangular solves */ 74453cca76cSShri Abhyankar pv = b->a + bs2 * bdiag[i]; 74553cca76cSShri Abhyankar pj = b->j + bdiag[i]; 7469566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(pv, rtmp + bs2 * pj[0], bs2)); 7479566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_5(pv, ipvt, work, shift, allowzeropivot, &zeropivotdetected)); 7487b6c816cSBarry Smith if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 74953cca76cSShri Abhyankar 75053cca76cSShri Abhyankar /* U part */ 75153cca76cSShri Abhyankar pv = b->a + bs2 * (bdiag[i + 1] + 1); 75253cca76cSShri Abhyankar pj = b->j + bdiag[i + 1] + 1; 75353cca76cSShri Abhyankar nz = bdiag[i] - bdiag[i + 1] - 1; 754*48a46eb9SPierre Jolivet for (j = 0; j < nz; j++) PetscCall(PetscArraycpy(pv + bs2 * j, rtmp + bs2 * pj[j], bs2)); 75553cca76cSShri Abhyankar } 7569566063dSJacob Faibussowitsch PetscCall(PetscFree2(rtmp, mwork)); 75726fbe8dcSKarl Rupp 7584dd39f65SShri Abhyankar C->ops->solve = MatSolve_SeqBAIJ_5_NaturalOrdering; 7594dd39f65SShri Abhyankar C->ops->solvetranspose = MatSolveTranspose_SeqBAIJ_5_NaturalOrdering; 76053cca76cSShri Abhyankar C->assembled = PETSC_TRUE; 76126fbe8dcSKarl Rupp 7629566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.333333333333 * 5 * 5 * 5 * n)); /* from inverting diagonal blocks */ 76353cca76cSShri Abhyankar PetscFunctionReturn(0); 76453cca76cSShri Abhyankar } 765