Actual source code: sbaijfact3.c
petsc-3.4.5 2014-06-29
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 */
8: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat C,Mat A,const MatFactorInfo *info)
9: {
10: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
11: IS perm = b->row;
13: const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
14: PetscInt *a2anew,i,j,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
15: MatScalar *ba = b->a,*aa,*ap,*dk,*uik;
16: MatScalar *u,*diag,*rtmp,*rtmp_ptr;
17: PetscReal shift = info->shiftamount;
20: /* initialization */
21: PetscMalloc(9*mbs*sizeof(MatScalar),&rtmp);
22: PetscMemzero(rtmp,9*mbs*sizeof(MatScalar));
23: PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);
24: for (i=0; i<mbs; i++) {
25: jl[i] = mbs; il[0] = 0;
26: }
27: PetscMalloc2(9,MatScalar,&dk,9,MatScalar,&uik);
28: ISGetIndices(perm,&perm_ptr);
30: /* check permutation */
31: if (!a->permute) {
32: ai = a->i; aj = a->j; aa = a->a;
33: } else {
34: ai = a->inew; aj = a->jnew;
35: PetscMalloc(9*ai[mbs]*sizeof(MatScalar),&aa);
36: PetscMemcpy(aa,a->a,9*ai[mbs]*sizeof(MatScalar));
37: PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);
38: PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));
40: for (i=0; i<mbs; i++) {
41: jmin = ai[i]; jmax = ai[i+1];
42: for (j=jmin; j<jmax; j++) {
43: while (a2anew[j] != j) {
44: k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
45: for (k1=0; k1<9; k1++) {
46: dk[k1] = aa[k*9+k1];
47: aa[k*9+k1] = aa[j*9+k1];
48: aa[j*9+k1] = dk[k1];
49: }
50: }
51: /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
52: if (i > aj[j]) {
53: /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
54: ap = aa + j*9; /* ptr to the beginning of j-th block of aa */
55: for (k=0; k<9; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
56: for (k=0; k<3; k++) { /* j-th block of aa <- dk^T */
57: for (k1=0; k1<3; k1++) *ap++ = dk[k + 3*k1];
58: }
59: }
60: }
61: }
62: PetscFree(a2anew);
63: }
65: /* for each row k */
66: for (k = 0; k<mbs; k++) {
68: /*initialize k-th row with elements nonzero in row perm(k) of A */
69: jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
70: if (jmin < jmax) {
71: ap = aa + jmin*9;
72: for (j = jmin; j < jmax; j++) {
73: vj = perm_ptr[aj[j]]; /* block col. index */
74: rtmp_ptr = rtmp + vj*9;
75: for (i=0; i<9; i++) *rtmp_ptr++ = *ap++;
76: }
77: }
79: /* modify k-th row by adding in those rows i with U(i,k) != 0 */
80: PetscMemcpy(dk,rtmp+k*9,9*sizeof(MatScalar));
81: i = jl[k]; /* first row to be added to k_th row */
83: while (i < mbs) {
84: nexti = jl[i]; /* next row to be added to k_th row */
86: /* compute multiplier */
87: ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
89: /* uik = -inv(Di)*U_bar(i,k) */
90: diag = ba + i*9;
91: u = ba + ili*9;
93: uik[0] = -(diag[0]*u[0] + diag[3]*u[1] + diag[6]*u[2]);
94: uik[1] = -(diag[1]*u[0] + diag[4]*u[1] + diag[7]*u[2]);
95: uik[2] = -(diag[2]*u[0] + diag[5]*u[1] + diag[8]*u[2]);
97: uik[3] = -(diag[0]*u[3] + diag[3]*u[4] + diag[6]*u[5]);
98: uik[4] = -(diag[1]*u[3] + diag[4]*u[4] + diag[7]*u[5]);
99: uik[5] = -(diag[2]*u[3] + diag[5]*u[4] + diag[8]*u[5]);
101: uik[6] = -(diag[0]*u[6] + diag[3]*u[7] + diag[6]*u[8]);
102: uik[7] = -(diag[1]*u[6] + diag[4]*u[7] + diag[7]*u[8]);
103: uik[8] = -(diag[2]*u[6] + diag[5]*u[7] + diag[8]*u[8]);
105: /* update D(k) += -U(i,k)^T * U_bar(i,k) */
106: dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
107: dk[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
108: dk[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];
110: dk[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
111: dk[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
112: dk[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];
114: dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
115: dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
116: dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
118: PetscLogFlops(27.0*4.0);
120: /* update -U(i,k) */
121: PetscMemcpy(ba+ili*9,uik,9*sizeof(MatScalar));
123: /* add multiple of row i to k-th row ... */
124: jmin = ili + 1; jmax = bi[i+1];
125: if (jmin < jmax) {
126: for (j=jmin; j<jmax; j++) {
127: /* rtmp += -U(i,k)^T * U_bar(i,j) */
128: rtmp_ptr = rtmp + bj[j]*9;
129: u = ba + j*9;
130: rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
131: rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
132: rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];
134: rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
135: rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
136: rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];
138: rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
139: rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
140: rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
141: }
142: PetscLogFlops(2.0*27.0*(jmax-jmin));
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: PetscKernel_A_gets_inverse_A_3(diag,shift);
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: }
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: PetscFree2(il,jl);
180: PetscFree2(dk,uik);
181: if (a->permute) {
182: PetscFree(aa);
183: }
185: ISRestoreIndices(perm,&perm_ptr);
187: C->ops->solve = MatSolve_SeqSBAIJ_3_inplace;
188: C->ops->solvetranspose = MatSolve_SeqSBAIJ_3_inplace;
189: C->assembled = PETSC_TRUE;
190: C->preallocated = PETSC_TRUE;
192: PetscLogFlops(1.3333*27*b->mbs); /* from inverting diagonal blocks */
193: return(0);
194: }