Actual source code: sbaijfact3.c
2: #include <../src/mat/impls/sbaij/seq/sbaij.h>
3: #include <petsc/private/kernels/blockinvert.h>
5: /* Version for when blocks are 3 by 3 */
6: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat C,Mat A,const MatFactorInfo *info)
7: {
8: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
9: IS perm = b->row;
10: const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
11: PetscInt *a2anew,i,j,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
12: MatScalar *ba = b->a,*aa,*ap,*dk,*uik;
13: MatScalar *u,*diag,*rtmp,*rtmp_ptr;
14: PetscReal shift = info->shiftamount;
15: PetscBool allowzeropivot,zeropivotdetected;
17: /* initialization */
18: allowzeropivot = PetscNot(A->erroriffailure);
19: PetscCalloc1(9*mbs,&rtmp);
20: PetscMalloc2(mbs,&il,mbs,&jl);
21: il[0] = 0;
22: for (i=0; i<mbs; i++) jl[i] = mbs;
24: PetscMalloc2(9,&dk,9,&uik);
25: ISGetIndices(perm,&perm_ptr);
27: /* check permutation */
28: if (!a->permute) {
29: ai = a->i; aj = a->j; aa = a->a;
30: } else {
31: ai = a->inew; aj = a->jnew;
32: PetscMalloc1(9*ai[mbs],&aa);
33: PetscArraycpy(aa,a->a,9*ai[mbs]);
34: PetscMalloc1(ai[mbs],&a2anew);
35: PetscArraycpy(a2anew,a->a2anew,ai[mbs]);
37: for (i=0; i<mbs; i++) {
38: jmin = ai[i]; jmax = ai[i+1];
39: for (j=jmin; j<jmax; j++) {
40: while (a2anew[j] != j) {
41: k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
42: for (k1=0; k1<9; k1++) {
43: dk[k1] = aa[k*9+k1];
44: aa[k*9+k1] = aa[j*9+k1];
45: aa[j*9+k1] = dk[k1];
46: }
47: }
48: /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
49: if (i > aj[j]) {
50: /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
51: ap = aa + j*9; /* ptr to the beginning of j-th block of aa */
52: for (k=0; k<9; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
53: for (k=0; k<3; k++) { /* j-th block of aa <- dk^T */
54: for (k1=0; k1<3; k1++) *ap++ = dk[k + 3*k1];
55: }
56: }
57: }
58: }
59: PetscFree(a2anew);
60: }
62: /* for each row k */
63: for (k = 0; k<mbs; k++) {
65: /*initialize k-th row with elements nonzero in row perm(k) of A */
66: jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
67: if (jmin < jmax) {
68: ap = aa + jmin*9;
69: for (j = jmin; j < jmax; j++) {
70: vj = perm_ptr[aj[j]]; /* block col. index */
71: rtmp_ptr = rtmp + vj*9;
72: for (i=0; i<9; i++) *rtmp_ptr++ = *ap++;
73: }
74: }
76: /* modify k-th row by adding in those rows i with U(i,k) != 0 */
77: PetscArraycpy(dk,rtmp+k*9,9);
78: i = jl[k]; /* first row to be added to k_th row */
80: while (i < mbs) {
81: nexti = jl[i]; /* next row to be added to k_th row */
83: /* compute multiplier */
84: ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
86: /* uik = -inv(Di)*U_bar(i,k) */
87: diag = ba + i*9;
88: u = ba + ili*9;
90: uik[0] = -(diag[0]*u[0] + diag[3]*u[1] + diag[6]*u[2]);
91: uik[1] = -(diag[1]*u[0] + diag[4]*u[1] + diag[7]*u[2]);
92: uik[2] = -(diag[2]*u[0] + diag[5]*u[1] + diag[8]*u[2]);
94: uik[3] = -(diag[0]*u[3] + diag[3]*u[4] + diag[6]*u[5]);
95: uik[4] = -(diag[1]*u[3] + diag[4]*u[4] + diag[7]*u[5]);
96: uik[5] = -(diag[2]*u[3] + diag[5]*u[4] + diag[8]*u[5]);
98: uik[6] = -(diag[0]*u[6] + diag[3]*u[7] + diag[6]*u[8]);
99: uik[7] = -(diag[1]*u[6] + diag[4]*u[7] + diag[7]*u[8]);
100: uik[8] = -(diag[2]*u[6] + diag[5]*u[7] + diag[8]*u[8]);
102: /* update D(k) += -U(i,k)^T * U_bar(i,k) */
103: dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
104: dk[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
105: dk[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];
107: dk[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
108: dk[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
109: dk[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];
111: dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
112: dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
113: dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
115: PetscLogFlops(27.0*4.0);
117: /* update -U(i,k) */
118: PetscArraycpy(ba+ili*9,uik,9);
120: /* add multiple of row i to k-th row ... */
121: jmin = ili + 1; jmax = bi[i+1];
122: if (jmin < jmax) {
123: for (j=jmin; j<jmax; j++) {
124: /* rtmp += -U(i,k)^T * U_bar(i,j) */
125: rtmp_ptr = rtmp + bj[j]*9;
126: u = ba + j*9;
127: rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
128: rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
129: rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];
131: rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
132: rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
133: rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];
135: rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
136: rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
137: rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
138: }
139: PetscLogFlops(2.0*27.0*(jmax-jmin));
141: /* ... add i to row list for next nonzero entry */
142: il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
143: j = bj[jmin];
144: jl[i] = jl[j]; jl[j] = i; /* update jl */
145: }
146: i = nexti;
147: }
149: /* save nonzero entries in k-th row of U ... */
151: /* invert diagonal block */
152: diag = ba+k*9;
153: PetscArraycpy(diag,dk,9);
154: PetscKernel_A_gets_inverse_A_3(diag,shift,allowzeropivot,&zeropivotdetected);
155: if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
157: jmin = bi[k]; jmax = bi[k+1];
158: if (jmin < jmax) {
159: for (j=jmin; j<jmax; j++) {
160: vj = bj[j]; /* block col. index of U */
161: u = ba + j*9;
162: rtmp_ptr = rtmp + vj*9;
163: for (k1=0; k1<9; k1++) {
164: *u++ = *rtmp_ptr;
165: *rtmp_ptr++ = 0.0;
166: }
167: }
169: /* ... add k to row list for first nonzero entry in k-th row */
170: il[k] = jmin;
171: i = bj[jmin];
172: jl[k] = jl[i]; jl[i] = k;
173: }
174: }
176: PetscFree(rtmp);
177: PetscFree2(il,jl);
178: PetscFree2(dk,uik);
179: if (a->permute) {
180: PetscFree(aa);
181: }
183: ISRestoreIndices(perm,&perm_ptr);
185: C->ops->solve = MatSolve_SeqSBAIJ_3_inplace;
186: C->ops->solvetranspose = MatSolve_SeqSBAIJ_3_inplace;
187: C->assembled = PETSC_TRUE;
188: C->preallocated = PETSC_TRUE;
190: PetscLogFlops(1.3333*27*b->mbs); /* from inverting diagonal blocks */
191: return 0;
192: }