Actual source code: sbaijfact3.c
2: #include src/mat/impls/sbaij/seq/sbaij.h
3: #include src/inline/ilu.h
5: /* Version for when blocks are 3 by 3 */
8: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat A,Mat *B)
9: {
10: Mat C = *B;
11: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
12: IS perm = b->row;
14: PetscInt *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
15: PetscInt *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
16: MatScalar *ba = b->a,*aa,*ap,*dk,*uik;
17: MatScalar *u,*diag,*rtmp,*rtmp_ptr;
20:
21: /* initialization */
22: PetscMalloc(9*mbs*sizeof(MatScalar),&rtmp);
23: PetscMemzero(rtmp,9*mbs*sizeof(MatScalar));
24: PetscMalloc(2*mbs*sizeof(PetscInt),&il);
25: jl = il + mbs;
26: for (i=0; i<mbs; i++) {
27: jl[i] = mbs; il[0] = 0;
28: }
29: PetscMalloc(18*sizeof(MatScalar),&dk);
30: uik = dk + 9;
31: ISGetIndices(perm,&perm_ptr);
33: /* check permutation */
34: if (!a->permute){
35: ai = a->i; aj = a->j; aa = a->a;
36: } else {
37: ai = a->inew; aj = a->jnew;
38: PetscMalloc(9*ai[mbs]*sizeof(MatScalar),&aa);
39: PetscMemcpy(aa,a->a,9*ai[mbs]*sizeof(MatScalar));
40: PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);
41: PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));
43: for (i=0; i<mbs; i++){
44: jmin = ai[i]; jmax = ai[i+1];
45: for (j=jmin; j<jmax; j++){
46: while (a2anew[j] != j){
47: k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
48: for (k1=0; k1<9; k1++){
49: dk[k1] = aa[k*9+k1];
50: aa[k*9+k1] = aa[j*9+k1];
51: aa[j*9+k1] = dk[k1];
52: }
53: }
54: /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
55: if (i > aj[j]){
56: /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
57: ap = aa + j*9; /* ptr to the beginning of j-th block of aa */
58: for (k=0; k<9; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
59: for (k=0; k<3; k++){ /* j-th block of aa <- dk^T */
60: for (k1=0; k1<3; k1++) *ap++ = dk[k + 3*k1];
61: }
62: }
63: }
64: }
65: PetscFree(a2anew);
66: }
67:
68: /* for each row k */
69: for (k = 0; k<mbs; k++){
71: /*initialize k-th row with elements nonzero in row perm(k) of A */
72: jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
73: if (jmin < jmax) {
74: ap = aa + jmin*9;
75: for (j = jmin; j < jmax; j++){
76: vj = perm_ptr[aj[j]]; /* block col. index */
77: rtmp_ptr = rtmp + vj*9;
78: for (i=0; i<9; i++) *rtmp_ptr++ = *ap++;
79: }
80: }
82: /* modify k-th row by adding in those rows i with U(i,k) != 0 */
83: PetscMemcpy(dk,rtmp+k*9,9*sizeof(MatScalar));
84: i = jl[k]; /* first row to be added to k_th row */
86: while (i < mbs){
87: nexti = jl[i]; /* next row to be added to k_th row */
89: /* compute multiplier */
90: ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
92: /* uik = -inv(Di)*U_bar(i,k) */
93: diag = ba + i*9;
94: u = ba + ili*9;
96: uik[0] = -(diag[0]*u[0] + diag[3]*u[1] + diag[6]*u[2]);
97: uik[1] = -(diag[1]*u[0] + diag[4]*u[1] + diag[7]*u[2]);
98: uik[2] = -(diag[2]*u[0] + diag[5]*u[1] + diag[8]*u[2]);
100: uik[3] = -(diag[0]*u[3] + diag[3]*u[4] + diag[6]*u[5]);
101: uik[4] = -(diag[1]*u[3] + diag[4]*u[4] + diag[7]*u[5]);
102: uik[5] = -(diag[2]*u[3] + diag[5]*u[4] + diag[8]*u[5]);
104: uik[6] = -(diag[0]*u[6] + diag[3]*u[7] + diag[6]*u[8]);
105: uik[7] = -(diag[1]*u[6] + diag[4]*u[7] + diag[7]*u[8]);
106: uik[8] = -(diag[2]*u[6] + diag[5]*u[7] + diag[8]*u[8]);
108: /* update D(k) += -U(i,k)^T * U_bar(i,k) */
109: dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
110: dk[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
111: dk[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];
113: dk[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
114: dk[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
115: dk[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];
117: dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
118: dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
119: dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
121: /* update -U(i,k) */
122: PetscMemcpy(ba+ili*9,uik,9*sizeof(MatScalar));
124: /* add multiple of row i to k-th row ... */
125: jmin = ili + 1; jmax = bi[i+1];
126: if (jmin < jmax){
127: for (j=jmin; j<jmax; j++) {
128: /* rtmp += -U(i,k)^T * U_bar(i,j) */
129: rtmp_ptr = rtmp + bj[j]*9;
130: u = ba + j*9;
131: rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
132: rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
133: rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];
135: rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
136: rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
137: rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];
139: rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
140: rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
141: rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
142: }
143:
144: /* ... add i to row list for next nonzero entry */
145: il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
146: j = bj[jmin];
147: jl[i] = jl[j]; jl[j] = i; /* update jl */
148: }
149: i = nexti;
150: }
152: /* save nonzero entries in k-th row of U ... */
154: /* invert diagonal block */
155: diag = ba+k*9;
156: PetscMemcpy(diag,dk,9*sizeof(MatScalar));
157: Kernel_A_gets_inverse_A_3(diag);
158:
159: jmin = bi[k]; jmax = bi[k+1];
160: if (jmin < jmax) {
161: for (j=jmin; j<jmax; j++){
162: vj = bj[j]; /* block col. index of U */
163: u = ba + j*9;
164: rtmp_ptr = rtmp + vj*9;
165: for (k1=0; k1<9; k1++){
166: *u++ = *rtmp_ptr;
167: *rtmp_ptr++ = 0.0;
168: }
169: }
170:
171: /* ... add k to row list for first nonzero entry in k-th row */
172: il[k] = jmin;
173: i = bj[jmin];
174: jl[k] = jl[i]; jl[i] = k;
175: }
176: }
178: PetscFree(rtmp);
179: PetscFree(il);
180: PetscFree(dk);
181: if (a->permute){
182: PetscFree(aa);
183: }
185: ISRestoreIndices(perm,&perm_ptr);
186: C->factor = FACTOR_CHOLESKY;
187: C->assembled = PETSC_TRUE;
188: C->preallocated = PETSC_TRUE;
189: PetscLogFlops(1.3333*27*b->mbs); /* from inverting diagonal blocks */
190: return(0);
191: }