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