Actual source code: sbaijfact5.c
1: #define PETSCMAT_DLL
3: #include ../src/mat/impls/sbaij/seq/sbaij.h
4: #include ../src/mat/blockinvert.h
6: /*
7: Version for when blocks are 4 by 4 Using natural ordering
8: */
11: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
12: {
13: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
15: PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
16: PetscInt *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
17: MatScalar *ba = b->a,*aa,*ap,*dk,*uik;
18: MatScalar *u,*diag,*rtmp,*rtmp_ptr;
19: PetscTruth pivotinblocks = b->pivotinblocks;
20: PetscReal shift = info->shiftamount;
23:
24: /* initialization */
25: PetscMalloc(16*mbs*sizeof(MatScalar),&rtmp);
26: PetscMemzero(rtmp,16*mbs*sizeof(MatScalar));
27: PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);
28: for (i=0; i<mbs; i++) {
29: jl[i] = mbs; il[0] = 0;
30: }
31: PetscMalloc2(16,MatScalar,&dk,16,MatScalar,&uik);
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*16;
41: for (j = jmin; j < jmax; j++){
42: vj = aj[j]; /* block col. index */
43: rtmp_ptr = rtmp + vj*16;
44: for (i=0; i<16; 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*16,16*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*16;
60: u = ba + ili*16;
62: uik[0] = -(diag[0]*u[0] + diag[4]*u[1] + diag[8]*u[2] + diag[12]*u[3]);
63: uik[1] = -(diag[1]*u[0] + diag[5]*u[1] + diag[9]*u[2] + diag[13]*u[3]);
64: uik[2] = -(diag[2]*u[0] + diag[6]*u[1] + diag[10]*u[2]+ diag[14]*u[3]);
65: uik[3] = -(diag[3]*u[0] + diag[7]*u[1] + diag[11]*u[2]+ diag[15]*u[3]);
67: uik[4] = -(diag[0]*u[4] + diag[4]*u[5] + diag[8]*u[6] + diag[12]*u[7]);
68: uik[5] = -(diag[1]*u[4] + diag[5]*u[5] + diag[9]*u[6] + diag[13]*u[7]);
69: uik[6] = -(diag[2]*u[4] + diag[6]*u[5] + diag[10]*u[6]+ diag[14]*u[7]);
70: uik[7] = -(diag[3]*u[4] + diag[7]*u[5] + diag[11]*u[6]+ diag[15]*u[7]);
72: uik[8] = -(diag[0]*u[8] + diag[4]*u[9] + diag[8]*u[10] + diag[12]*u[11]);
73: uik[9] = -(diag[1]*u[8] + diag[5]*u[9] + diag[9]*u[10] + diag[13]*u[11]);
74: uik[10]= -(diag[2]*u[8] + diag[6]*u[9] + diag[10]*u[10]+ diag[14]*u[11]);
75: uik[11]= -(diag[3]*u[8] + diag[7]*u[9] + diag[11]*u[10]+ diag[15]*u[11]);
77: uik[12]= -(diag[0]*u[12] + diag[4]*u[13] + diag[8]*u[14] + diag[12]*u[15]);
78: uik[13]= -(diag[1]*u[12] + diag[5]*u[13] + diag[9]*u[14] + diag[13]*u[15]);
79: uik[14]= -(diag[2]*u[12] + diag[6]*u[13] + diag[10]*u[14]+ diag[14]*u[15]);
80: uik[15]= -(diag[3]*u[12] + diag[7]*u[13] + diag[11]*u[14]+ diag[15]*u[15]);
82: /* update D(k) += -U(i,k)^T * U_bar(i,k) */
83: dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3];
84: dk[1] += uik[4]*u[0] + uik[5]*u[1] + uik[6]*u[2] + uik[7]*u[3];
85: dk[2] += uik[8]*u[0] + uik[9]*u[1] + uik[10]*u[2]+ uik[11]*u[3];
86: dk[3] += uik[12]*u[0]+ uik[13]*u[1]+ uik[14]*u[2]+ uik[15]*u[3];
88: dk[4] += uik[0]*u[4] + uik[1]*u[5] + uik[2]*u[6] + uik[3]*u[7];
89: dk[5] += uik[4]*u[4] + uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7];
90: dk[6] += uik[8]*u[4] + uik[9]*u[5] + uik[10]*u[6]+ uik[11]*u[7];
91: dk[7] += uik[12]*u[4]+ uik[13]*u[5]+ uik[14]*u[6]+ uik[15]*u[7];
93: dk[8] += uik[0]*u[8] + uik[1]*u[9] + uik[2]*u[10] + uik[3]*u[11];
94: dk[9] += uik[4]*u[8] + uik[5]*u[9] + uik[6]*u[10] + uik[7]*u[11];
95: dk[10]+= uik[8]*u[8] + uik[9]*u[9] + uik[10]*u[10]+ uik[11]*u[11];
96: dk[11]+= uik[12]*u[8]+ uik[13]*u[9]+ uik[14]*u[10]+ uik[15]*u[11];
98: dk[12]+= uik[0]*u[12] + uik[1]*u[13] + uik[2]*u[14] + uik[3]*u[15];
99: dk[13]+= uik[4]*u[12] + uik[5]*u[13] + uik[6]*u[14] + uik[7]*u[15];
100: dk[14]+= uik[8]*u[12] + uik[9]*u[13] + uik[10]*u[14]+ uik[11]*u[15];
101: dk[15]+= uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14]+ uik[15]*u[15];
103: PetscLogFlops(64.0*4.0);
104:
105: /* update -U(i,k) */
106: PetscMemcpy(ba+ili*16,uik,16*sizeof(MatScalar));
108: /* add multiple of row i to k-th row ... */
109: jmin = ili + 1; jmax = bi[i+1];
110: if (jmin < jmax){
111: for (j=jmin; j<jmax; j++) {
112: /* rtmp += -U(i,k)^T * U_bar(i,j) */
113: rtmp_ptr = rtmp + bj[j]*16;
114: u = ba + j*16;
115: rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3];
116: rtmp_ptr[1] += uik[4]*u[0] + uik[5]*u[1] + uik[6]*u[2] + uik[7]*u[3];
117: rtmp_ptr[2] += uik[8]*u[0] + uik[9]*u[1] + uik[10]*u[2]+ uik[11]*u[3];
118: rtmp_ptr[3] += uik[12]*u[0]+ uik[13]*u[1]+ uik[14]*u[2]+ uik[15]*u[3];
120: rtmp_ptr[4] += uik[0]*u[4] + uik[1]*u[5] + uik[2]*u[6] + uik[3]*u[7];
121: rtmp_ptr[5] += uik[4]*u[4] + uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7];
122: rtmp_ptr[6] += uik[8]*u[4] + uik[9]*u[5] + uik[10]*u[6]+ uik[11]*u[7];
123: rtmp_ptr[7] += uik[12]*u[4]+ uik[13]*u[5]+ uik[14]*u[6]+ uik[15]*u[7];
125: rtmp_ptr[8] += uik[0]*u[8] + uik[1]*u[9] + uik[2]*u[10] + uik[3]*u[11];
126: rtmp_ptr[9] += uik[4]*u[8] + uik[5]*u[9] + uik[6]*u[10] + uik[7]*u[11];
127: rtmp_ptr[10]+= uik[8]*u[8] + uik[9]*u[9] + uik[10]*u[10]+ uik[11]*u[11];
128: rtmp_ptr[11]+= uik[12]*u[8]+ uik[13]*u[9]+ uik[14]*u[10]+ uik[15]*u[11];
130: rtmp_ptr[12]+= uik[0]*u[12] + uik[1]*u[13] + uik[2]*u[14] + uik[3]*u[15];
131: rtmp_ptr[13]+= uik[4]*u[12] + uik[5]*u[13] + uik[6]*u[14] + uik[7]*u[15];
132: rtmp_ptr[14]+= uik[8]*u[12] + uik[9]*u[13] + uik[10]*u[14]+ uik[11]*u[15];
133: rtmp_ptr[15]+= uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14]+ uik[15]*u[15];
134: }
135: PetscLogFlops(2.0*64.0*(jmax-jmin));
136:
137: /* ... add i to row list for next nonzero entry */
138: il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
139: j = bj[jmin];
140: jl[i] = jl[j]; jl[j] = i; /* update jl */
141: }
142: i = nexti;
143: }
145: /* save nonzero entries in k-th row of U ... */
147: /* invert diagonal block */
148: diag = ba+k*16;
149: PetscMemcpy(diag,dk,16*sizeof(MatScalar));
150: if (pivotinblocks) {
151: Kernel_A_gets_inverse_A_4(diag,shift);
152: } else {
153: Kernel_A_gets_inverse_A_4_nopivot(diag);
154: }
155:
156: jmin = bi[k]; jmax = bi[k+1];
157: if (jmin < jmax) {
158: for (j=jmin; j<jmax; j++){
159: vj = bj[j]; /* block col. index of U */
160: u = ba + j*16;
161: rtmp_ptr = rtmp + vj*16;
162: for (k1=0; k1<16; k1++){
163: *u++ = *rtmp_ptr;
164: *rtmp_ptr++ = 0.0;
165: }
166: }
167:
168: /* ... add k to row list for first nonzero entry in k-th row */
169: il[k] = jmin;
170: i = bj[jmin];
171: jl[k] = jl[i]; jl[i] = k;
172: }
173: }
175: PetscFree(rtmp);
176: PetscFree2(il,jl);
177: PetscFree2(dk,uik);
178:
179: C->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering_inplace;
180: C->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering_inplace;
181: C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_4_NaturalOrdering_inplace;
182: C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_4_NaturalOrdering_inplace;
184: C->assembled = PETSC_TRUE;
185: C->preallocated = PETSC_TRUE;
186: PetscLogFlops(1.3333*64*b->mbs); /* from inverting diagonal blocks */
187: return(0);
188: }