1c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h> 2af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h> 381278733SSatish Balay 481278733SSatish Balay /* Version for when blocks are 4 by 4 */ 5d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat C, Mat A, const MatFactorInfo *info) 6d71ae5a4SJacob Faibussowitsch { 781278733SSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data, *b = (Mat_SeqSBAIJ *)C->data; 881278733SSatish Balay IS perm = b->row; 95d0c19d7SBarry Smith const PetscInt *ai, *aj, *perm_ptr, mbs = a->mbs, *bi = b->i, *bj = b->j; 105d0c19d7SBarry Smith PetscInt i, j, *a2anew, k, k1, jmin, jmax, *jl, *il, vj, nexti, ili; 1181278733SSatish Balay MatScalar *ba = b->a, *aa, *ap, *dk, *uik; 1281278733SSatish Balay MatScalar *u, *diag, *rtmp, *rtmp_ptr; 13ace3abfcSBarry Smith PetscBool pivotinblocks = b->pivotinblocks; 14182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 150164db54SHong Zhang PetscBool allowzeropivot, zeropivotdetected = PETSC_FALSE; 1681278733SSatish Balay 1781278733SSatish Balay PetscFunctionBegin; 1881278733SSatish Balay /* initialization */ 190164db54SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 209566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(16 * mbs, &rtmp)); 219566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(mbs, &il, mbs, &jl)); 226df5ee2eSHong Zhang il[0] = 0; 236df5ee2eSHong Zhang for (i = 0; i < mbs; i++) jl[i] = mbs; 246df5ee2eSHong Zhang 259566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(16, &dk, 16, &uik)); 269566063dSJacob Faibussowitsch PetscCall(ISGetIndices(perm, &perm_ptr)); 2781278733SSatish Balay 2881278733SSatish Balay /* check permutation */ 2981278733SSatish Balay if (!a->permute) { 309371c9d4SSatish Balay ai = a->i; 319371c9d4SSatish Balay aj = a->j; 329371c9d4SSatish Balay aa = a->a; 3381278733SSatish Balay } else { 349371c9d4SSatish Balay ai = a->inew; 359371c9d4SSatish Balay aj = a->jnew; 369566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(16 * ai[mbs], &aa)); 379566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(aa, a->a, 16 * ai[mbs])); 389566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ai[mbs], &a2anew)); 399566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(a2anew, a->a2anew, ai[mbs])); 4081278733SSatish Balay 4181278733SSatish Balay for (i = 0; i < mbs; i++) { 429371c9d4SSatish Balay jmin = ai[i]; 439371c9d4SSatish Balay jmax = ai[i + 1]; 4481278733SSatish Balay for (j = jmin; j < jmax; j++) { 4581278733SSatish Balay while (a2anew[j] != j) { 469371c9d4SSatish Balay k = a2anew[j]; 479371c9d4SSatish Balay a2anew[j] = a2anew[k]; 489371c9d4SSatish Balay a2anew[k] = k; 4981278733SSatish Balay for (k1 = 0; k1 < 16; k1++) { 5081278733SSatish Balay dk[k1] = aa[k * 16 + k1]; 5181278733SSatish Balay aa[k * 16 + k1] = aa[j * 16 + k1]; 5281278733SSatish Balay aa[j * 16 + k1] = dk[k1]; 5381278733SSatish Balay } 5481278733SSatish Balay } 55*5e116b59SBarry Smith /* transform column-oriented blocks that lie in the lower triangle to row-oriented blocks */ 5681278733SSatish Balay if (i > aj[j]) { 5781278733SSatish Balay /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */ 5881278733SSatish Balay ap = aa + j * 16; /* ptr to the beginning of j-th block of aa */ 5981278733SSatish Balay for (k = 0; k < 16; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */ 6081278733SSatish Balay for (k = 0; k < 4; k++) { /* j-th block of aa <- dk^T */ 6181278733SSatish Balay for (k1 = 0; k1 < 4; k1++) *ap++ = dk[k + 4 * k1]; 6281278733SSatish Balay } 6381278733SSatish Balay } 6481278733SSatish Balay } 6581278733SSatish Balay } 669566063dSJacob Faibussowitsch PetscCall(PetscFree(a2anew)); 6781278733SSatish Balay } 6881278733SSatish Balay 6981278733SSatish Balay /* for each row k */ 7081278733SSatish Balay for (k = 0; k < mbs; k++) { 7181278733SSatish Balay /*initialize k-th row with elements nonzero in row perm(k) of A */ 729371c9d4SSatish Balay jmin = ai[perm_ptr[k]]; 739371c9d4SSatish Balay jmax = ai[perm_ptr[k] + 1]; 7481278733SSatish Balay if (jmin < jmax) { 7581278733SSatish Balay ap = aa + jmin * 16; 7681278733SSatish Balay for (j = jmin; j < jmax; j++) { 7781278733SSatish Balay vj = perm_ptr[aj[j]]; /* block col. index */ 7881278733SSatish Balay rtmp_ptr = rtmp + vj * 16; 7981278733SSatish Balay for (i = 0; i < 16; i++) *rtmp_ptr++ = *ap++; 8081278733SSatish Balay } 8181278733SSatish Balay } 8281278733SSatish Balay 8381278733SSatish Balay /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 849566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dk, rtmp + k * 16, 16)); 8581278733SSatish Balay i = jl[k]; /* first row to be added to k_th row */ 8681278733SSatish Balay 8781278733SSatish Balay while (i < mbs) { 8881278733SSatish Balay nexti = jl[i]; /* next row to be added to k_th row */ 8981278733SSatish Balay 9081278733SSatish Balay /* compute multiplier */ 9181278733SSatish Balay ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 9281278733SSatish Balay 9381278733SSatish Balay /* uik = -inv(Di)*U_bar(i,k) */ 9481278733SSatish Balay diag = ba + i * 16; 9581278733SSatish Balay u = ba + ili * 16; 9681278733SSatish Balay 9781278733SSatish Balay uik[0] = -(diag[0] * u[0] + diag[4] * u[1] + diag[8] * u[2] + diag[12] * u[3]); 9881278733SSatish Balay uik[1] = -(diag[1] * u[0] + diag[5] * u[1] + diag[9] * u[2] + diag[13] * u[3]); 9981278733SSatish Balay uik[2] = -(diag[2] * u[0] + diag[6] * u[1] + diag[10] * u[2] + diag[14] * u[3]); 10081278733SSatish Balay uik[3] = -(diag[3] * u[0] + diag[7] * u[1] + diag[11] * u[2] + diag[15] * u[3]); 10181278733SSatish Balay 10281278733SSatish Balay uik[4] = -(diag[0] * u[4] + diag[4] * u[5] + diag[8] * u[6] + diag[12] * u[7]); 10381278733SSatish Balay uik[5] = -(diag[1] * u[4] + diag[5] * u[5] + diag[9] * u[6] + diag[13] * u[7]); 10481278733SSatish Balay uik[6] = -(diag[2] * u[4] + diag[6] * u[5] + diag[10] * u[6] + diag[14] * u[7]); 10581278733SSatish Balay uik[7] = -(diag[3] * u[4] + diag[7] * u[5] + diag[11] * u[6] + diag[15] * u[7]); 10681278733SSatish Balay 10781278733SSatish Balay uik[8] = -(diag[0] * u[8] + diag[4] * u[9] + diag[8] * u[10] + diag[12] * u[11]); 10881278733SSatish Balay uik[9] = -(diag[1] * u[8] + diag[5] * u[9] + diag[9] * u[10] + diag[13] * u[11]); 10981278733SSatish Balay uik[10] = -(diag[2] * u[8] + diag[6] * u[9] + diag[10] * u[10] + diag[14] * u[11]); 11081278733SSatish Balay uik[11] = -(diag[3] * u[8] + diag[7] * u[9] + diag[11] * u[10] + diag[15] * u[11]); 11181278733SSatish Balay 11281278733SSatish Balay uik[12] = -(diag[0] * u[12] + diag[4] * u[13] + diag[8] * u[14] + diag[12] * u[15]); 11381278733SSatish Balay uik[13] = -(diag[1] * u[12] + diag[5] * u[13] + diag[9] * u[14] + diag[13] * u[15]); 11481278733SSatish Balay uik[14] = -(diag[2] * u[12] + diag[6] * u[13] + diag[10] * u[14] + diag[14] * u[15]); 11581278733SSatish Balay uik[15] = -(diag[3] * u[12] + diag[7] * u[13] + diag[11] * u[14] + diag[15] * u[15]); 11681278733SSatish Balay 11781278733SSatish Balay /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 11881278733SSatish Balay dk[0] += uik[0] * u[0] + uik[1] * u[1] + uik[2] * u[2] + uik[3] * u[3]; 11981278733SSatish Balay dk[1] += uik[4] * u[0] + uik[5] * u[1] + uik[6] * u[2] + uik[7] * u[3]; 12081278733SSatish Balay dk[2] += uik[8] * u[0] + uik[9] * u[1] + uik[10] * u[2] + uik[11] * u[3]; 12181278733SSatish Balay dk[3] += uik[12] * u[0] + uik[13] * u[1] + uik[14] * u[2] + uik[15] * u[3]; 12281278733SSatish Balay 12381278733SSatish Balay dk[4] += uik[0] * u[4] + uik[1] * u[5] + uik[2] * u[6] + uik[3] * u[7]; 12481278733SSatish Balay dk[5] += uik[4] * u[4] + uik[5] * u[5] + uik[6] * u[6] + uik[7] * u[7]; 12581278733SSatish Balay dk[6] += uik[8] * u[4] + uik[9] * u[5] + uik[10] * u[6] + uik[11] * u[7]; 12681278733SSatish Balay dk[7] += uik[12] * u[4] + uik[13] * u[5] + uik[14] * u[6] + uik[15] * u[7]; 12781278733SSatish Balay 12881278733SSatish Balay dk[8] += uik[0] * u[8] + uik[1] * u[9] + uik[2] * u[10] + uik[3] * u[11]; 12981278733SSatish Balay dk[9] += uik[4] * u[8] + uik[5] * u[9] + uik[6] * u[10] + uik[7] * u[11]; 13081278733SSatish Balay dk[10] += uik[8] * u[8] + uik[9] * u[9] + uik[10] * u[10] + uik[11] * u[11]; 13181278733SSatish Balay dk[11] += uik[12] * u[8] + uik[13] * u[9] + uik[14] * u[10] + uik[15] * u[11]; 13281278733SSatish Balay 13381278733SSatish Balay dk[12] += uik[0] * u[12] + uik[1] * u[13] + uik[2] * u[14] + uik[3] * u[15]; 13481278733SSatish Balay dk[13] += uik[4] * u[12] + uik[5] * u[13] + uik[6] * u[14] + uik[7] * u[15]; 13581278733SSatish Balay dk[14] += uik[8] * u[12] + uik[9] * u[13] + uik[10] * u[14] + uik[11] * u[15]; 13681278733SSatish Balay dk[15] += uik[12] * u[12] + uik[13] * u[13] + uik[14] * u[14] + uik[15] * u[15]; 13781278733SSatish Balay 1389566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(64.0 * 4.0)); 139187a9f4bSHong Zhang 14081278733SSatish Balay /* update -U(i,k) */ 1419566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ba + ili * 16, uik, 16)); 14281278733SSatish Balay 14381278733SSatish Balay /* add multiple of row i to k-th row ... */ 1449371c9d4SSatish Balay jmin = ili + 1; 1459371c9d4SSatish Balay jmax = bi[i + 1]; 14681278733SSatish Balay if (jmin < jmax) { 14781278733SSatish Balay for (j = jmin; j < jmax; j++) { 14881278733SSatish Balay /* rtmp += -U(i,k)^T * U_bar(i,j) */ 14981278733SSatish Balay rtmp_ptr = rtmp + bj[j] * 16; 15081278733SSatish Balay u = ba + j * 16; 15181278733SSatish Balay rtmp_ptr[0] += uik[0] * u[0] + uik[1] * u[1] + uik[2] * u[2] + uik[3] * u[3]; 15281278733SSatish Balay rtmp_ptr[1] += uik[4] * u[0] + uik[5] * u[1] + uik[6] * u[2] + uik[7] * u[3]; 15381278733SSatish Balay rtmp_ptr[2] += uik[8] * u[0] + uik[9] * u[1] + uik[10] * u[2] + uik[11] * u[3]; 15481278733SSatish Balay rtmp_ptr[3] += uik[12] * u[0] + uik[13] * u[1] + uik[14] * u[2] + uik[15] * u[3]; 15581278733SSatish Balay 15681278733SSatish Balay rtmp_ptr[4] += uik[0] * u[4] + uik[1] * u[5] + uik[2] * u[6] + uik[3] * u[7]; 15781278733SSatish Balay rtmp_ptr[5] += uik[4] * u[4] + uik[5] * u[5] + uik[6] * u[6] + uik[7] * u[7]; 15881278733SSatish Balay rtmp_ptr[6] += uik[8] * u[4] + uik[9] * u[5] + uik[10] * u[6] + uik[11] * u[7]; 15981278733SSatish Balay rtmp_ptr[7] += uik[12] * u[4] + uik[13] * u[5] + uik[14] * u[6] + uik[15] * u[7]; 16081278733SSatish Balay 16181278733SSatish Balay rtmp_ptr[8] += uik[0] * u[8] + uik[1] * u[9] + uik[2] * u[10] + uik[3] * u[11]; 16281278733SSatish Balay rtmp_ptr[9] += uik[4] * u[8] + uik[5] * u[9] + uik[6] * u[10] + uik[7] * u[11]; 16381278733SSatish Balay rtmp_ptr[10] += uik[8] * u[8] + uik[9] * u[9] + uik[10] * u[10] + uik[11] * u[11]; 16481278733SSatish Balay rtmp_ptr[11] += uik[12] * u[8] + uik[13] * u[9] + uik[14] * u[10] + uik[15] * u[11]; 16581278733SSatish Balay 16681278733SSatish Balay rtmp_ptr[12] += uik[0] * u[12] + uik[1] * u[13] + uik[2] * u[14] + uik[3] * u[15]; 16781278733SSatish Balay rtmp_ptr[13] += uik[4] * u[12] + uik[5] * u[13] + uik[6] * u[14] + uik[7] * u[15]; 16881278733SSatish Balay rtmp_ptr[14] += uik[8] * u[12] + uik[9] * u[13] + uik[10] * u[14] + uik[11] * u[15]; 16981278733SSatish Balay rtmp_ptr[15] += uik[12] * u[12] + uik[13] * u[13] + uik[14] * u[14] + uik[15] * u[15]; 17081278733SSatish Balay } 1719566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * 64.0 * (jmax - jmin))); 17281278733SSatish Balay 17381278733SSatish Balay /* ... add i to row list for next nonzero entry */ 17481278733SSatish Balay il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 17581278733SSatish Balay j = bj[jmin]; 1769371c9d4SSatish Balay jl[i] = jl[j]; 1779371c9d4SSatish Balay jl[j] = i; /* update jl */ 17881278733SSatish Balay } 17981278733SSatish Balay i = nexti; 18081278733SSatish Balay } 18181278733SSatish Balay 18281278733SSatish Balay /* save nonzero entries in k-th row of U ... */ 18381278733SSatish Balay 18481278733SSatish Balay /* invert diagonal block */ 18581278733SSatish Balay diag = ba + k * 16; 1869566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(diag, dk, 16)); 187bcd9e38bSBarry Smith 188bcd9e38bSBarry Smith if (pivotinblocks) { 1899566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_4(diag, shift, allowzeropivot, &zeropivotdetected)); 1907b6c816cSBarry Smith if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 191bcd9e38bSBarry Smith } else { 1929566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_4_nopivot(diag)); 193bcd9e38bSBarry Smith } 19481278733SSatish Balay 1959371c9d4SSatish Balay jmin = bi[k]; 1969371c9d4SSatish Balay jmax = bi[k + 1]; 19781278733SSatish Balay if (jmin < jmax) { 19881278733SSatish Balay for (j = jmin; j < jmax; j++) { 19981278733SSatish Balay vj = bj[j]; /* block col. index of U */ 20081278733SSatish Balay u = ba + j * 16; 20181278733SSatish Balay rtmp_ptr = rtmp + vj * 16; 20281278733SSatish Balay for (k1 = 0; k1 < 16; k1++) { 20381278733SSatish Balay *u++ = *rtmp_ptr; 20481278733SSatish Balay *rtmp_ptr++ = 0.0; 20581278733SSatish Balay } 20681278733SSatish Balay } 20781278733SSatish Balay 20881278733SSatish Balay /* ... add k to row list for first nonzero entry in k-th row */ 20981278733SSatish Balay il[k] = jmin; 21081278733SSatish Balay i = bj[jmin]; 2119371c9d4SSatish Balay jl[k] = jl[i]; 2129371c9d4SSatish Balay jl[i] = k; 21381278733SSatish Balay } 21481278733SSatish Balay } 21581278733SSatish Balay 2169566063dSJacob Faibussowitsch PetscCall(PetscFree(rtmp)); 2179566063dSJacob Faibussowitsch PetscCall(PetscFree2(il, jl)); 2189566063dSJacob Faibussowitsch PetscCall(PetscFree2(dk, uik)); 2191baa6e33SBarry Smith if (a->permute) PetscCall(PetscFree(aa)); 22081278733SSatish Balay 2219566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(perm, &perm_ptr)); 22226fbe8dcSKarl Rupp 2234f79d315SHong Zhang C->ops->solve = MatSolve_SeqSBAIJ_4_inplace; 2244f79d315SHong Zhang C->ops->solvetranspose = MatSolve_SeqSBAIJ_4_inplace; 22581278733SSatish Balay C->assembled = PETSC_TRUE; 22681278733SSatish Balay C->preallocated = PETSC_TRUE; 22726fbe8dcSKarl Rupp 2289566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.3333 * 64 * b->mbs)); /* from inverting diagonal blocks */ 2293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 23081278733SSatish Balay } 231