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: */
8: int MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat A,Mat *B)
9: {
10: Mat C = *B;
11: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
12: int ierr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
13: int *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
14: MatScalar *ba = b->a,*aa,*ap,*dk,*uik;
15: MatScalar *u,*diag,*rtmp,*rtmp_ptr;
18:
19: /* initialization */
20: PetscMalloc(9*mbs*sizeof(MatScalar),&rtmp);
21: PetscMemzero(rtmp,9*mbs*sizeof(MatScalar));
22: PetscMalloc(2*mbs*sizeof(int),&il);
23: jl = il + mbs;
24: for (i=0; i<mbs; i++) {
25: jl[i] = mbs; il[0] = 0;
26: }
27: PetscMalloc(18*sizeof(MatScalar),&dk);
28: uik = dk + 9;
30: ai = a->i; aj = a->j; aa = a->a;
32: /* for each row k */
33: for (k = 0; k<mbs; k++){
35: /*initialize k-th row with elements nonzero in row k of A */
36: jmin = ai[k]; jmax = ai[k+1];
37: if (jmin < jmax) {
38: ap = aa + jmin*9;
39: for (j = jmin; j < jmax; j++){
40: vj = aj[j]; /* block col. index */
41: rtmp_ptr = rtmp + vj*9;
42: for (i=0; i<9; i++) *rtmp_ptr++ = *ap++;
43: }
44: }
46: /* modify k-th row by adding in those rows i with U(i,k) != 0 */
47: PetscMemcpy(dk,rtmp+k*9,9*sizeof(MatScalar));
48: i = jl[k]; /* first row to be added to k_th row */
50: while (i < mbs){
51: nexti = jl[i]; /* next row to be added to k_th row */
53: /* compute multiplier */
54: ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
56: /* uik = -inv(Di)*U_bar(i,k) */
57: diag = ba + i*9;
58: u = ba + ili*9;
60: uik[0] = -(diag[0]*u[0] + diag[3]*u[1] + diag[6]*u[2]);
61: uik[1] = -(diag[1]*u[0] + diag[4]*u[1] + diag[7]*u[2]);
62: uik[2] = -(diag[2]*u[0] + diag[5]*u[1] + diag[8]*u[2]);
64: uik[3] = -(diag[0]*u[3] + diag[3]*u[4] + diag[6]*u[5]);
65: uik[4] = -(diag[1]*u[3] + diag[4]*u[4] + diag[7]*u[5]);
66: uik[5] = -(diag[2]*u[3] + diag[5]*u[4] + diag[8]*u[5]);
68: uik[6] = -(diag[0]*u[6] + diag[3]*u[7] + diag[6]*u[8]);
69: uik[7] = -(diag[1]*u[6] + diag[4]*u[7] + diag[7]*u[8]);
70: uik[8] = -(diag[2]*u[6] + diag[5]*u[7] + diag[8]*u[8]);
72: /* update D(k) += -U(i,k)^T * U_bar(i,k) */
73: dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
74: dk[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
75: dk[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];
77: dk[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
78: dk[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
79: dk[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];
81: dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
82: dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
83: dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
85: /* update -U(i,k) */
86: PetscMemcpy(ba+ili*9,uik,9*sizeof(MatScalar));
88: /* add multiple of row i to k-th row ... */
89: jmin = ili + 1; jmax = bi[i+1];
90: if (jmin < jmax){
91: for (j=jmin; j<jmax; j++) {
92: /* rtmp += -U(i,k)^T * U_bar(i,j) */
93: rtmp_ptr = rtmp + bj[j]*9;
94: u = ba + j*9;
95: rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
96: rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
97: rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];
99: rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
100: rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
101: rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];
103: rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
104: rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
105: rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
106: }
107:
108: /* ... add i to row list for next nonzero entry */
109: il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
110: j = bj[jmin];
111: jl[i] = jl[j]; jl[j] = i; /* update jl */
112: }
113: i = nexti;
114: }
116: /* save nonzero entries in k-th row of U ... */
118: /* invert diagonal block */
119: diag = ba+k*9;
120: PetscMemcpy(diag,dk,9*sizeof(MatScalar));
121: Kernel_A_gets_inverse_A_3(diag);
122:
123: jmin = bi[k]; jmax = bi[k+1];
124: if (jmin < jmax) {
125: for (j=jmin; j<jmax; j++){
126: vj = bj[j]; /* block col. index of U */
127: u = ba + j*9;
128: rtmp_ptr = rtmp + vj*9;
129: for (k1=0; k1<9; k1++){
130: *u++ = *rtmp_ptr;
131: *rtmp_ptr++ = 0.0;
132: }
133: }
134:
135: /* ... add k to row list for first nonzero entry in k-th row */
136: il[k] = jmin;
137: i = bj[jmin];
138: jl[k] = jl[i]; jl[i] = k;
139: }
140: }
142: PetscFree(rtmp);
143: PetscFree(il);
144: PetscFree(dk);
146: C->factor = FACTOR_CHOLESKY;
147: C->assembled = PETSC_TRUE;
148: C->preallocated = PETSC_TRUE;
149: PetscLogFlops(1.3333*27*b->mbs); /* from inverting diagonal blocks */
150: return(0);
151: }