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 3 by 3 */
MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat C,Mat A,const MatFactorInfo * info)5d71ae5a4SJacob Faibussowitsch PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(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 *a2anew, i, j, 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;
13182b8fbaSHong Zhang PetscReal shift = info->shiftamount;
14a455e926SHong Zhang PetscBool allowzeropivot, zeropivotdetected;
1581278733SSatish Balay
1681278733SSatish Balay PetscFunctionBegin;
1781278733SSatish Balay /* initialization */
180164db54SHong Zhang allowzeropivot = PetscNot(A->erroriffailure);
199566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(9 * mbs, &rtmp));
209566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(mbs, &il, mbs, &jl));
216df5ee2eSHong Zhang il[0] = 0;
226df5ee2eSHong Zhang for (i = 0; i < mbs; i++) jl[i] = mbs;
236df5ee2eSHong Zhang
249566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(9, &dk, 9, &uik));
259566063dSJacob Faibussowitsch PetscCall(ISGetIndices(perm, &perm_ptr));
2681278733SSatish Balay
2781278733SSatish Balay /* check permutation */
2881278733SSatish Balay if (!a->permute) {
299371c9d4SSatish Balay ai = a->i;
309371c9d4SSatish Balay aj = a->j;
319371c9d4SSatish Balay aa = a->a;
3281278733SSatish Balay } else {
339371c9d4SSatish Balay ai = a->inew;
349371c9d4SSatish Balay aj = a->jnew;
359566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(9 * ai[mbs], &aa));
369566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(aa, a->a, 9 * ai[mbs]));
379566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ai[mbs], &a2anew));
389566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(a2anew, a->a2anew, ai[mbs]));
3981278733SSatish Balay
4081278733SSatish Balay for (i = 0; i < mbs; i++) {
419371c9d4SSatish Balay jmin = ai[i];
429371c9d4SSatish Balay jmax = ai[i + 1];
4381278733SSatish Balay for (j = jmin; j < jmax; j++) {
4481278733SSatish Balay while (a2anew[j] != j) {
459371c9d4SSatish Balay k = a2anew[j];
469371c9d4SSatish Balay a2anew[j] = a2anew[k];
479371c9d4SSatish Balay a2anew[k] = k;
4881278733SSatish Balay for (k1 = 0; k1 < 9; k1++) {
4981278733SSatish Balay dk[k1] = aa[k * 9 + k1];
5081278733SSatish Balay aa[k * 9 + k1] = aa[j * 9 + k1];
5181278733SSatish Balay aa[j * 9 + k1] = dk[k1];
5281278733SSatish Balay }
5381278733SSatish Balay }
54*5e116b59SBarry Smith /* transform column-oriented blocks that lie in the lower triangle to row-oriented blocks */
5581278733SSatish Balay if (i > aj[j]) {
5681278733SSatish Balay /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
5781278733SSatish Balay ap = aa + j * 9; /* ptr to the beginning of j-th block of aa */
5881278733SSatish Balay for (k = 0; k < 9; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
5981278733SSatish Balay for (k = 0; k < 3; k++) { /* j-th block of aa <- dk^T */
6081278733SSatish Balay for (k1 = 0; k1 < 3; k1++) *ap++ = dk[k + 3 * k1];
6181278733SSatish Balay }
6281278733SSatish Balay }
6381278733SSatish Balay }
6481278733SSatish Balay }
659566063dSJacob Faibussowitsch PetscCall(PetscFree(a2anew));
6681278733SSatish Balay }
6781278733SSatish Balay
6881278733SSatish Balay /* for each row k */
6981278733SSatish Balay for (k = 0; k < mbs; k++) {
7081278733SSatish Balay /*initialize k-th row with elements nonzero in row perm(k) of A */
719371c9d4SSatish Balay jmin = ai[perm_ptr[k]];
729371c9d4SSatish Balay jmax = ai[perm_ptr[k] + 1];
7381278733SSatish Balay if (jmin < jmax) {
7481278733SSatish Balay ap = aa + jmin * 9;
7581278733SSatish Balay for (j = jmin; j < jmax; j++) {
7681278733SSatish Balay vj = perm_ptr[aj[j]]; /* block col. index */
7781278733SSatish Balay rtmp_ptr = rtmp + vj * 9;
7881278733SSatish Balay for (i = 0; i < 9; i++) *rtmp_ptr++ = *ap++;
7981278733SSatish Balay }
8081278733SSatish Balay }
8181278733SSatish Balay
8281278733SSatish Balay /* modify k-th row by adding in those rows i with U(i,k) != 0 */
839566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(dk, rtmp + k * 9, 9));
8481278733SSatish Balay i = jl[k]; /* first row to be added to k_th row */
8581278733SSatish Balay
8681278733SSatish Balay while (i < mbs) {
8781278733SSatish Balay nexti = jl[i]; /* next row to be added to k_th row */
8881278733SSatish Balay
8981278733SSatish Balay /* compute multiplier */
9081278733SSatish Balay ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
9181278733SSatish Balay
9281278733SSatish Balay /* uik = -inv(Di)*U_bar(i,k) */
9381278733SSatish Balay diag = ba + i * 9;
9481278733SSatish Balay u = ba + ili * 9;
9581278733SSatish Balay
9681278733SSatish Balay uik[0] = -(diag[0] * u[0] + diag[3] * u[1] + diag[6] * u[2]);
9781278733SSatish Balay uik[1] = -(diag[1] * u[0] + diag[4] * u[1] + diag[7] * u[2]);
9881278733SSatish Balay uik[2] = -(diag[2] * u[0] + diag[5] * u[1] + diag[8] * u[2]);
9981278733SSatish Balay
10081278733SSatish Balay uik[3] = -(diag[0] * u[3] + diag[3] * u[4] + diag[6] * u[5]);
10181278733SSatish Balay uik[4] = -(diag[1] * u[3] + diag[4] * u[4] + diag[7] * u[5]);
10281278733SSatish Balay uik[5] = -(diag[2] * u[3] + diag[5] * u[4] + diag[8] * u[5]);
10381278733SSatish Balay
10481278733SSatish Balay uik[6] = -(diag[0] * u[6] + diag[3] * u[7] + diag[6] * u[8]);
10581278733SSatish Balay uik[7] = -(diag[1] * u[6] + diag[4] * u[7] + diag[7] * u[8]);
10681278733SSatish Balay uik[8] = -(diag[2] * u[6] + diag[5] * u[7] + diag[8] * u[8]);
10781278733SSatish Balay
10881278733SSatish Balay /* update D(k) += -U(i,k)^T * U_bar(i,k) */
10981278733SSatish Balay dk[0] += uik[0] * u[0] + uik[1] * u[1] + uik[2] * u[2];
11081278733SSatish Balay dk[1] += uik[3] * u[0] + uik[4] * u[1] + uik[5] * u[2];
11181278733SSatish Balay dk[2] += uik[6] * u[0] + uik[7] * u[1] + uik[8] * u[2];
11281278733SSatish Balay
11381278733SSatish Balay dk[3] += uik[0] * u[3] + uik[1] * u[4] + uik[2] * u[5];
11481278733SSatish Balay dk[4] += uik[3] * u[3] + uik[4] * u[4] + uik[5] * u[5];
11581278733SSatish Balay dk[5] += uik[6] * u[3] + uik[7] * u[4] + uik[8] * u[5];
11681278733SSatish Balay
11781278733SSatish Balay dk[6] += uik[0] * u[6] + uik[1] * u[7] + uik[2] * u[8];
11881278733SSatish Balay dk[7] += uik[3] * u[6] + uik[4] * u[7] + uik[5] * u[8];
11981278733SSatish Balay dk[8] += uik[6] * u[6] + uik[7] * u[7] + uik[8] * u[8];
12081278733SSatish Balay
1219566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(27.0 * 4.0));
122187a9f4bSHong Zhang
12381278733SSatish Balay /* update -U(i,k) */
1249566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(ba + ili * 9, uik, 9));
12581278733SSatish Balay
12681278733SSatish Balay /* add multiple of row i to k-th row ... */
1279371c9d4SSatish Balay jmin = ili + 1;
1289371c9d4SSatish Balay jmax = bi[i + 1];
12981278733SSatish Balay if (jmin < jmax) {
13081278733SSatish Balay for (j = jmin; j < jmax; j++) {
13181278733SSatish Balay /* rtmp += -U(i,k)^T * U_bar(i,j) */
13281278733SSatish Balay rtmp_ptr = rtmp + bj[j] * 9;
13381278733SSatish Balay u = ba + j * 9;
13481278733SSatish Balay rtmp_ptr[0] += uik[0] * u[0] + uik[1] * u[1] + uik[2] * u[2];
13581278733SSatish Balay rtmp_ptr[1] += uik[3] * u[0] + uik[4] * u[1] + uik[5] * u[2];
13681278733SSatish Balay rtmp_ptr[2] += uik[6] * u[0] + uik[7] * u[1] + uik[8] * u[2];
13781278733SSatish Balay
13881278733SSatish Balay rtmp_ptr[3] += uik[0] * u[3] + uik[1] * u[4] + uik[2] * u[5];
13981278733SSatish Balay rtmp_ptr[4] += uik[3] * u[3] + uik[4] * u[4] + uik[5] * u[5];
14081278733SSatish Balay rtmp_ptr[5] += uik[6] * u[3] + uik[7] * u[4] + uik[8] * u[5];
14181278733SSatish Balay
14281278733SSatish Balay rtmp_ptr[6] += uik[0] * u[6] + uik[1] * u[7] + uik[2] * u[8];
14381278733SSatish Balay rtmp_ptr[7] += uik[3] * u[6] + uik[4] * u[7] + uik[5] * u[8];
14481278733SSatish Balay rtmp_ptr[8] += uik[6] * u[6] + uik[7] * u[7] + uik[8] * u[8];
14581278733SSatish Balay }
1469566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(2.0 * 27.0 * (jmax - jmin)));
14781278733SSatish Balay
14881278733SSatish Balay /* ... add i to row list for next nonzero entry */
14981278733SSatish Balay il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
15081278733SSatish Balay j = bj[jmin];
1519371c9d4SSatish Balay jl[i] = jl[j];
1529371c9d4SSatish Balay jl[j] = i; /* update jl */
15381278733SSatish Balay }
15481278733SSatish Balay i = nexti;
15581278733SSatish Balay }
15681278733SSatish Balay
15781278733SSatish Balay /* save nonzero entries in k-th row of U ... */
15881278733SSatish Balay
15981278733SSatish Balay /* invert diagonal block */
16081278733SSatish Balay diag = ba + k * 9;
1619566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(diag, dk, 9));
1629566063dSJacob Faibussowitsch PetscCall(PetscKernel_A_gets_inverse_A_3(diag, shift, allowzeropivot, &zeropivotdetected));
1637b6c816cSBarry Smith if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
16481278733SSatish Balay
1659371c9d4SSatish Balay jmin = bi[k];
1669371c9d4SSatish Balay jmax = bi[k + 1];
16781278733SSatish Balay if (jmin < jmax) {
16881278733SSatish Balay for (j = jmin; j < jmax; j++) {
16981278733SSatish Balay vj = bj[j]; /* block col. index of U */
17081278733SSatish Balay u = ba + j * 9;
17181278733SSatish Balay rtmp_ptr = rtmp + vj * 9;
17281278733SSatish Balay for (k1 = 0; k1 < 9; k1++) {
17381278733SSatish Balay *u++ = *rtmp_ptr;
17481278733SSatish Balay *rtmp_ptr++ = 0.0;
17581278733SSatish Balay }
17681278733SSatish Balay }
17781278733SSatish Balay
17881278733SSatish Balay /* ... add k to row list for first nonzero entry in k-th row */
17981278733SSatish Balay il[k] = jmin;
18081278733SSatish Balay i = bj[jmin];
1819371c9d4SSatish Balay jl[k] = jl[i];
1829371c9d4SSatish Balay jl[i] = k;
18381278733SSatish Balay }
18481278733SSatish Balay }
18581278733SSatish Balay
1869566063dSJacob Faibussowitsch PetscCall(PetscFree(rtmp));
1879566063dSJacob Faibussowitsch PetscCall(PetscFree2(il, jl));
1889566063dSJacob Faibussowitsch PetscCall(PetscFree2(dk, uik));
1891baa6e33SBarry Smith if (a->permute) PetscCall(PetscFree(aa));
19081278733SSatish Balay
1919566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(perm, &perm_ptr));
19226fbe8dcSKarl Rupp
1934f79d315SHong Zhang C->ops->solve = MatSolve_SeqSBAIJ_3_inplace;
1944f79d315SHong Zhang C->ops->solvetranspose = MatSolve_SeqSBAIJ_3_inplace;
19581278733SSatish Balay C->assembled = PETSC_TRUE;
19681278733SSatish Balay C->preallocated = PETSC_TRUE;
19726fbe8dcSKarl Rupp
1989566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(1.3333 * 27 * b->mbs)); /* from inverting diagonal blocks */
1993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
20081278733SSatish Balay }
201