1c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h> 2af0996ceSBarry Smith #include <petsc/private/kernels/blockinvert.h> 381278733SSatish Balay 481278733SSatish Balay /* 581278733SSatish Balay Version for when blocks are 4 by 4 Using natural ordering 681278733SSatish Balay */ 7d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering(Mat C, Mat A, const MatFactorInfo *info) 8d71ae5a4SJacob Faibussowitsch { 981278733SSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ *)A->data, *b = (Mat_SeqSBAIJ *)C->data; 1013f74950SBarry Smith PetscInt i, j, mbs = a->mbs, *bi = b->i, *bj = b->j; 1113f74950SBarry Smith PetscInt *ai, *aj, k, k1, jmin, jmax, *jl, *il, vj, nexti, ili; 1281278733SSatish Balay MatScalar *ba = b->a, *aa, *ap, *dk, *uik; 1381278733SSatish Balay MatScalar *u, *diag, *rtmp, *rtmp_ptr; 14a455e926SHong Zhang PetscBool pivotinblocks = b->pivotinblocks; 15182b8fbaSHong Zhang PetscReal shift = info->shiftamount; 160164db54SHong Zhang PetscBool allowzeropivot, zeropivotdetected = PETSC_FALSE; 1781278733SSatish Balay 1881278733SSatish Balay PetscFunctionBegin; 1981278733SSatish Balay /* initialization */ 200164db54SHong Zhang allowzeropivot = PetscNot(A->erroriffailure); 219566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(16 * mbs, &rtmp)); 229566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(mbs, &il, mbs, &jl)); 236df5ee2eSHong Zhang il[0] = 0; 246df5ee2eSHong Zhang for (i = 0; i < mbs; i++) jl[i] = mbs; 256df5ee2eSHong Zhang 269566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(16, &dk, 16, &uik)); 279371c9d4SSatish Balay ai = a->i; 289371c9d4SSatish Balay aj = a->j; 299371c9d4SSatish Balay aa = a->a; 3081278733SSatish Balay 3181278733SSatish Balay /* for each row k */ 3281278733SSatish Balay for (k = 0; k < mbs; k++) { 3381278733SSatish Balay /*initialize k-th row with elements nonzero in row k of A */ 349371c9d4SSatish Balay jmin = ai[k]; 359371c9d4SSatish Balay jmax = ai[k + 1]; 3681278733SSatish Balay if (jmin < jmax) { 3781278733SSatish Balay ap = aa + jmin * 16; 3881278733SSatish Balay for (j = jmin; j < jmax; j++) { 3981278733SSatish Balay vj = aj[j]; /* block col. index */ 4081278733SSatish Balay rtmp_ptr = rtmp + vj * 16; 4181278733SSatish Balay for (i = 0; i < 16; i++) *rtmp_ptr++ = *ap++; 4281278733SSatish Balay } 4381278733SSatish Balay } 4481278733SSatish Balay 4581278733SSatish Balay /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 469566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dk, rtmp + k * 16, 16)); 4781278733SSatish Balay i = jl[k]; /* first row to be added to k_th row */ 4881278733SSatish Balay 4981278733SSatish Balay while (i < mbs) { 5081278733SSatish Balay nexti = jl[i]; /* next row to be added to k_th row */ 5181278733SSatish Balay 5281278733SSatish Balay /* compute multiplier */ 5381278733SSatish Balay ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 5481278733SSatish Balay 5581278733SSatish Balay /* uik = -inv(Di)*U_bar(i,k) */ 5681278733SSatish Balay diag = ba + i * 16; 5781278733SSatish Balay u = ba + ili * 16; 5881278733SSatish Balay 5981278733SSatish Balay uik[0] = -(diag[0] * u[0] + diag[4] * u[1] + diag[8] * u[2] + diag[12] * u[3]); 6081278733SSatish Balay uik[1] = -(diag[1] * u[0] + diag[5] * u[1] + diag[9] * u[2] + diag[13] * u[3]); 6181278733SSatish Balay uik[2] = -(diag[2] * u[0] + diag[6] * u[1] + diag[10] * u[2] + diag[14] * u[3]); 6281278733SSatish Balay uik[3] = -(diag[3] * u[0] + diag[7] * u[1] + diag[11] * u[2] + diag[15] * u[3]); 6381278733SSatish Balay 6481278733SSatish Balay uik[4] = -(diag[0] * u[4] + diag[4] * u[5] + diag[8] * u[6] + diag[12] * u[7]); 6581278733SSatish Balay uik[5] = -(diag[1] * u[4] + diag[5] * u[5] + diag[9] * u[6] + diag[13] * u[7]); 6681278733SSatish Balay uik[6] = -(diag[2] * u[4] + diag[6] * u[5] + diag[10] * u[6] + diag[14] * u[7]); 6781278733SSatish Balay uik[7] = -(diag[3] * u[4] + diag[7] * u[5] + diag[11] * u[6] + diag[15] * u[7]); 6881278733SSatish Balay 6981278733SSatish Balay uik[8] = -(diag[0] * u[8] + diag[4] * u[9] + diag[8] * u[10] + diag[12] * u[11]); 7081278733SSatish Balay uik[9] = -(diag[1] * u[8] + diag[5] * u[9] + diag[9] * u[10] + diag[13] * u[11]); 7181278733SSatish Balay uik[10] = -(diag[2] * u[8] + diag[6] * u[9] + diag[10] * u[10] + diag[14] * u[11]); 7281278733SSatish Balay uik[11] = -(diag[3] * u[8] + diag[7] * u[9] + diag[11] * u[10] + diag[15] * u[11]); 7381278733SSatish Balay 7481278733SSatish Balay uik[12] = -(diag[0] * u[12] + diag[4] * u[13] + diag[8] * u[14] + diag[12] * u[15]); 7581278733SSatish Balay uik[13] = -(diag[1] * u[12] + diag[5] * u[13] + diag[9] * u[14] + diag[13] * u[15]); 7681278733SSatish Balay uik[14] = -(diag[2] * u[12] + diag[6] * u[13] + diag[10] * u[14] + diag[14] * u[15]); 7781278733SSatish Balay uik[15] = -(diag[3] * u[12] + diag[7] * u[13] + diag[11] * u[14] + diag[15] * u[15]); 7881278733SSatish Balay 7981278733SSatish Balay /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 8081278733SSatish Balay dk[0] += uik[0] * u[0] + uik[1] * u[1] + uik[2] * u[2] + uik[3] * u[3]; 8181278733SSatish Balay dk[1] += uik[4] * u[0] + uik[5] * u[1] + uik[6] * u[2] + uik[7] * u[3]; 8281278733SSatish Balay dk[2] += uik[8] * u[0] + uik[9] * u[1] + uik[10] * u[2] + uik[11] * u[3]; 8381278733SSatish Balay dk[3] += uik[12] * u[0] + uik[13] * u[1] + uik[14] * u[2] + uik[15] * u[3]; 8481278733SSatish Balay 8581278733SSatish Balay dk[4] += uik[0] * u[4] + uik[1] * u[5] + uik[2] * u[6] + uik[3] * u[7]; 8681278733SSatish Balay dk[5] += uik[4] * u[4] + uik[5] * u[5] + uik[6] * u[6] + uik[7] * u[7]; 8781278733SSatish Balay dk[6] += uik[8] * u[4] + uik[9] * u[5] + uik[10] * u[6] + uik[11] * u[7]; 8881278733SSatish Balay dk[7] += uik[12] * u[4] + uik[13] * u[5] + uik[14] * u[6] + uik[15] * u[7]; 8981278733SSatish Balay 9081278733SSatish Balay dk[8] += uik[0] * u[8] + uik[1] * u[9] + uik[2] * u[10] + uik[3] * u[11]; 9181278733SSatish Balay dk[9] += uik[4] * u[8] + uik[5] * u[9] + uik[6] * u[10] + uik[7] * u[11]; 9281278733SSatish Balay dk[10] += uik[8] * u[8] + uik[9] * u[9] + uik[10] * u[10] + uik[11] * u[11]; 9381278733SSatish Balay dk[11] += uik[12] * u[8] + uik[13] * u[9] + uik[14] * u[10] + uik[15] * u[11]; 9481278733SSatish Balay 9581278733SSatish Balay dk[12] += uik[0] * u[12] + uik[1] * u[13] + uik[2] * u[14] + uik[3] * u[15]; 9681278733SSatish Balay dk[13] += uik[4] * u[12] + uik[5] * u[13] + uik[6] * u[14] + uik[7] * u[15]; 9781278733SSatish Balay dk[14] += uik[8] * u[12] + uik[9] * u[13] + uik[10] * u[14] + uik[11] * u[15]; 9881278733SSatish Balay dk[15] += uik[12] * u[12] + uik[13] * u[13] + uik[14] * u[14] + uik[15] * u[15]; 9981278733SSatish Balay 1009566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(64.0 * 4.0)); 101187a9f4bSHong Zhang 10281278733SSatish Balay /* update -U(i,k) */ 1039566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ba + ili * 16, uik, 16)); 10481278733SSatish Balay 10581278733SSatish Balay /* add multiple of row i to k-th row ... */ 1069371c9d4SSatish Balay jmin = ili + 1; 1079371c9d4SSatish Balay jmax = bi[i + 1]; 10881278733SSatish Balay if (jmin < jmax) { 10981278733SSatish Balay for (j = jmin; j < jmax; j++) { 11081278733SSatish Balay /* rtmp += -U(i,k)^T * U_bar(i,j) */ 11181278733SSatish Balay rtmp_ptr = rtmp + bj[j] * 16; 11281278733SSatish Balay u = ba + j * 16; 11381278733SSatish Balay rtmp_ptr[0] += uik[0] * u[0] + uik[1] * u[1] + uik[2] * u[2] + uik[3] * u[3]; 11481278733SSatish Balay rtmp_ptr[1] += uik[4] * u[0] + uik[5] * u[1] + uik[6] * u[2] + uik[7] * u[3]; 11581278733SSatish Balay rtmp_ptr[2] += uik[8] * u[0] + uik[9] * u[1] + uik[10] * u[2] + uik[11] * u[3]; 11681278733SSatish Balay rtmp_ptr[3] += uik[12] * u[0] + uik[13] * u[1] + uik[14] * u[2] + uik[15] * u[3]; 11781278733SSatish Balay 11881278733SSatish Balay rtmp_ptr[4] += uik[0] * u[4] + uik[1] * u[5] + uik[2] * u[6] + uik[3] * u[7]; 11981278733SSatish Balay rtmp_ptr[5] += uik[4] * u[4] + uik[5] * u[5] + uik[6] * u[6] + uik[7] * u[7]; 12081278733SSatish Balay rtmp_ptr[6] += uik[8] * u[4] + uik[9] * u[5] + uik[10] * u[6] + uik[11] * u[7]; 12181278733SSatish Balay rtmp_ptr[7] += uik[12] * u[4] + uik[13] * u[5] + uik[14] * u[6] + uik[15] * u[7]; 12281278733SSatish Balay 12381278733SSatish Balay rtmp_ptr[8] += uik[0] * u[8] + uik[1] * u[9] + uik[2] * u[10] + uik[3] * u[11]; 12481278733SSatish Balay rtmp_ptr[9] += uik[4] * u[8] + uik[5] * u[9] + uik[6] * u[10] + uik[7] * u[11]; 12581278733SSatish Balay rtmp_ptr[10] += uik[8] * u[8] + uik[9] * u[9] + uik[10] * u[10] + uik[11] * u[11]; 12681278733SSatish Balay rtmp_ptr[11] += uik[12] * u[8] + uik[13] * u[9] + uik[14] * u[10] + uik[15] * u[11]; 12781278733SSatish Balay 12881278733SSatish Balay rtmp_ptr[12] += uik[0] * u[12] + uik[1] * u[13] + uik[2] * u[14] + uik[3] * u[15]; 12981278733SSatish Balay rtmp_ptr[13] += uik[4] * u[12] + uik[5] * u[13] + uik[6] * u[14] + uik[7] * u[15]; 13081278733SSatish Balay rtmp_ptr[14] += uik[8] * u[12] + uik[9] * u[13] + uik[10] * u[14] + uik[11] * u[15]; 13181278733SSatish Balay rtmp_ptr[15] += uik[12] * u[12] + uik[13] * u[13] + uik[14] * u[14] + uik[15] * u[15]; 13281278733SSatish Balay } 1339566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * 64.0 * (jmax - jmin))); 13481278733SSatish Balay 13581278733SSatish Balay /* ... add i to row list for next nonzero entry */ 13681278733SSatish Balay il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 13781278733SSatish Balay j = bj[jmin]; 1389371c9d4SSatish Balay jl[i] = jl[j]; 1399371c9d4SSatish Balay jl[j] = i; /* update jl */ 14081278733SSatish Balay } 14181278733SSatish Balay i = nexti; 14281278733SSatish Balay } 14381278733SSatish Balay 14481278733SSatish Balay /* save nonzero entries in k-th row of U ... */ 14581278733SSatish Balay 14681278733SSatish Balay /* invert diagonal block */ 14781278733SSatish Balay diag = ba + k * 16; 1489566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(diag, dk, 16)); 149bcd9e38bSBarry Smith if (pivotinblocks) { 1509566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_4(diag, shift, allowzeropivot, &zeropivotdetected)); 1517b6c816cSBarry Smith if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT; 152bcd9e38bSBarry Smith } else { 1539566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_4_nopivot(diag)); 154bcd9e38bSBarry Smith } 15581278733SSatish Balay 1569371c9d4SSatish Balay jmin = bi[k]; 1579371c9d4SSatish Balay jmax = bi[k + 1]; 15881278733SSatish Balay if (jmin < jmax) { 15981278733SSatish Balay for (j = jmin; j < jmax; j++) { 16081278733SSatish Balay vj = bj[j]; /* block col. index of U */ 16181278733SSatish Balay u = ba + j * 16; 16281278733SSatish Balay rtmp_ptr = rtmp + vj * 16; 16381278733SSatish Balay for (k1 = 0; k1 < 16; k1++) { 16481278733SSatish Balay *u++ = *rtmp_ptr; 16581278733SSatish Balay *rtmp_ptr++ = 0.0; 16681278733SSatish Balay } 16781278733SSatish Balay } 16881278733SSatish Balay 16981278733SSatish Balay /* ... add k to row list for first nonzero entry in k-th row */ 17081278733SSatish Balay il[k] = jmin; 17181278733SSatish Balay i = bj[jmin]; 1729371c9d4SSatish Balay jl[k] = jl[i]; 1739371c9d4SSatish Balay jl[i] = k; 17481278733SSatish Balay } 17581278733SSatish Balay } 17681278733SSatish Balay 1779566063dSJacob Faibussowitsch PetscCall(PetscFree(rtmp)); 1789566063dSJacob Faibussowitsch PetscCall(PetscFree2(il, jl)); 1799566063dSJacob Faibussowitsch PetscCall(PetscFree2(dk, uik)); 18081278733SSatish Balay 1814f79d315SHong Zhang C->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering_inplace; 1824f79d315SHong Zhang C->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering_inplace; 1834f79d315SHong Zhang C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_4_NaturalOrdering_inplace; 1844f79d315SHong Zhang C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_4_NaturalOrdering_inplace; 185db4efbfdSBarry Smith 18681278733SSatish Balay C->assembled = PETSC_TRUE; 18781278733SSatish Balay C->preallocated = PETSC_TRUE; 18826fbe8dcSKarl Rupp 1899566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.3333 * 64 * b->mbs)); /* from inverting diagonal blocks */ 190*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19181278733SSatish Balay } 192