Actual source code: sbaijfact4.c
1: /*$Id: sbaijfact4.c,v 1.4 2001/06/21 21:17:00 bsmith Exp $*/
2: #include src/mat/impls/sbaij/seq/sbaij.h
3: #include src/inline/ilu.h
5: /*
6: Version for when blocks are 3 by 3 Using natural ordering
7: */
10: int MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat A,Mat *B)
11: {
12: Mat C = *B;
13: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
14: int ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
15: int *ai,*aj,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(int),&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;
32: ai = a->i; aj = a->j; aa = a->a;
34: /* for each row k */
35: for (k = 0; k<mbs; k++){
37: /*initialize k-th row with elements nonzero in row k of A */
38: jmin = ai[k]; jmax = ai[k+1];
39: if (jmin < jmax) {
40: ap = aa + jmin*9;
41: for (j = jmin; j < jmax; j++){
42: vj = aj[j]; /* block col. index */
43: rtmp_ptr = rtmp + vj*9;
44: for (i=0; i<9; i++) *rtmp_ptr++ = *ap++;
45: }
46: }
48: /* modify k-th row by adding in those rows i with U(i,k) != 0 */
49: PetscMemcpy(dk,rtmp+k*9,9*sizeof(MatScalar));
50: i = jl[k]; /* first row to be added to k_th row */
52: while (i < mbs){
53: nexti = jl[i]; /* next row to be added to k_th row */
55: /* compute multiplier */
56: ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
58: /* uik = -inv(Di)*U_bar(i,k) */
59: diag = ba + i*9;
60: u = ba + ili*9;
62: uik[0] = -(diag[0]*u[0] + diag[3]*u[1] + diag[6]*u[2]);
63: uik[1] = -(diag[1]*u[0] + diag[4]*u[1] + diag[7]*u[2]);
64: uik[2] = -(diag[2]*u[0] + diag[5]*u[1] + diag[8]*u[2]);
66: uik[3] = -(diag[0]*u[3] + diag[3]*u[4] + diag[6]*u[5]);
67: uik[4] = -(diag[1]*u[3] + diag[4]*u[4] + diag[7]*u[5]);
68: uik[5] = -(diag[2]*u[3] + diag[5]*u[4] + diag[8]*u[5]);
70: uik[6] = -(diag[0]*u[6] + diag[3]*u[7] + diag[6]*u[8]);
71: uik[7] = -(diag[1]*u[6] + diag[4]*u[7] + diag[7]*u[8]);
72: uik[8] = -(diag[2]*u[6] + diag[5]*u[7] + diag[8]*u[8]);
74: /* update D(k) += -U(i,k)^T * U_bar(i,k) */
75: dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
76: dk[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
77: dk[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];
79: dk[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
80: dk[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
81: dk[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];
83: dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
84: dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
85: dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
87: /* update -U(i,k) */
88: PetscMemcpy(ba+ili*9,uik,9*sizeof(MatScalar));
90: /* add multiple of row i to k-th row ... */
91: jmin = ili + 1; jmax = bi[i+1];
92: if (jmin < jmax){
93: for (j=jmin; j<jmax; j++) {
94: /* rtmp += -U(i,k)^T * U_bar(i,j) */
95: rtmp_ptr = rtmp + bj[j]*9;
96: u = ba + j*9;
97: rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
98: rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
99: rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];
101: rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
102: rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
103: rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];
105: rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
106: rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
107: rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
108: }
109:
110: /* ... add i to row list for next nonzero entry */
111: il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
112: j = bj[jmin];
113: jl[i] = jl[j]; jl[j] = i; /* update jl */
114: }
115: i = nexti;
116: }
118: /* save nonzero entries in k-th row of U ... */
120: /* invert diagonal block */
121: diag = ba+k*9;
122: PetscMemcpy(diag,dk,9*sizeof(MatScalar));
123: Kernel_A_gets_inverse_A_3(diag);
124:
125: jmin = bi[k]; jmax = bi[k+1];
126: if (jmin < jmax) {
127: for (j=jmin; j<jmax; j++){
128: vj = bj[j]; /* block col. index of U */
129: u = ba + j*9;
130: rtmp_ptr = rtmp + vj*9;
131: for (k1=0; k1<9; k1++){
132: *u++ = *rtmp_ptr;
133: *rtmp_ptr++ = 0.0;
134: }
135: }
136:
137: /* ... add k to row list for first nonzero entry in k-th row */
138: il[k] = jmin;
139: i = bj[jmin];
140: jl[k] = jl[i]; jl[i] = k;
141: }
142: }
144: PetscFree(rtmp);
145: PetscFree(il);
146: PetscFree(dk);
148: C->factor = FACTOR_CHOLESKY;
149: C->assembled = PETSC_TRUE;
150: C->preallocated = PETSC_TRUE;
151: PetscLogFlops(1.3333*27*b->mbs); /* from inverting diagonal blocks */
152: return(0);
153: }