Actual source code: sbaijfact7.c
petsc-3.7.3 2016-08-01
2: #include <../src/mat/impls/sbaij/seq/sbaij.h>
3: #include <petsc/private/kernels/blockinvert.h>
5: /* Version for when blocks are 5 by 5 */
8: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_5(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 i,j,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili,ipvt[5];
15: MatScalar *ba = b->a,*aa,*ap,*dk,*uik;
16: MatScalar *u,*d,*rtmp,*rtmp_ptr,work[25];
17: PetscReal shift = info->shiftamount;
18: PetscBool allowzeropivot,zeropivotdetected;
21: /* initialization */
22: allowzeropivot = PetscNot(A->erroriffailure);
23: PetscCalloc1(25*mbs,&rtmp);
24: PetscMalloc2(mbs,&il,mbs,&jl);
25: il[0] = 0;
26: for (i=0; i<mbs; i++) jl[i] = mbs;
27:
28: PetscMalloc2(25,&dk,25,&uik);
29: ISGetIndices(perm,&perm_ptr);
31: /* check permutation */
32: if (!a->permute) {
33: ai = a->i; aj = a->j; aa = a->a;
34: } else {
35: ai = a->inew; aj = a->jnew;
36: PetscMalloc1(25*ai[mbs],&aa);
37: PetscMemcpy(aa,a->a,25*ai[mbs]*sizeof(MatScalar));
38: PetscMalloc1(ai[mbs],&a2anew);
39: PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));
41: for (i=0; i<mbs; i++) {
42: jmin = ai[i]; jmax = ai[i+1];
43: for (j=jmin; j<jmax; j++) {
44: while (a2anew[j] != j) {
45: k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
46: for (k1=0; k1<25; k1++) {
47: dk[k1] = aa[k*25+k1];
48: aa[k*25+k1] = aa[j*25+k1];
49: aa[j*25+k1] = dk[k1];
50: }
51: }
52: /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
53: if (i > aj[j]) {
54: /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
55: ap = aa + j*25; /* ptr to the beginning of j-th block of aa */
56: for (k=0; k<25; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
57: for (k=0; k<5; k++) { /* j-th block of aa <- dk^T */
58: for (k1=0; k1<5; k1++) *ap++ = dk[k + 5*k1];
59: }
60: }
61: }
62: }
63: PetscFree(a2anew);
64: }
66: /* for each row k */
67: for (k = 0; k<mbs; k++) {
69: /*initialize k-th row with elements nonzero in row perm(k) of A */
70: jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
71: if (jmin < jmax) {
72: ap = aa + jmin*25;
73: for (j = jmin; j < jmax; j++) {
74: vj = perm_ptr[aj[j]]; /* block col. index */
75: rtmp_ptr = rtmp + vj*25;
76: for (i=0; i<25; i++) *rtmp_ptr++ = *ap++;
77: }
78: }
80: /* modify k-th row by adding in those rows i with U(i,k) != 0 */
81: PetscMemcpy(dk,rtmp+k*25,25*sizeof(MatScalar));
82: i = jl[k]; /* first row to be added to k_th row */
84: while (i < mbs) {
85: nexti = jl[i]; /* next row to be added to k_th row */
87: /* compute multiplier */
88: ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
90: /* uik = -inv(Di)*U_bar(i,k) */
91: d = ba + i*25;
92: u = ba + ili*25;
94: uik[0] = -(d[0]*u[0] + d[5]*u[1] + d[10]*u[2] + d[15]*u[3] + d[20]*u[4]);
95: uik[1] = -(d[1]*u[0] + d[6]*u[1] + d[11]*u[2] + d[16]*u[3] + d[21]*u[4]);
96: uik[2] = -(d[2]*u[0] + d[7]*u[1] + d[12]*u[2] + d[17]*u[3] + d[22]*u[4]);
97: uik[3] = -(d[3]*u[0] + d[8]*u[1] + d[13]*u[2] + d[18]*u[3] + d[23]*u[4]);
98: uik[4] = -(d[4]*u[0] + d[9]*u[1] + d[14]*u[2] + d[19]*u[3] + d[24]*u[4]);
100: uik[5] = -(d[0]*u[5] + d[5]*u[6] + d[10]*u[7] + d[15]*u[8] + d[20]*u[9]);
101: uik[6] = -(d[1]*u[5] + d[6]*u[6] + d[11]*u[7] + d[16]*u[8] + d[21]*u[9]);
102: uik[7] = -(d[2]*u[5] + d[7]*u[6] + d[12]*u[7] + d[17]*u[8] + d[22]*u[9]);
103: uik[8] = -(d[3]*u[5] + d[8]*u[6] + d[13]*u[7] + d[18]*u[8] + d[23]*u[9]);
104: uik[9] = -(d[4]*u[5] + d[9]*u[6] + d[14]*u[7] + d[19]*u[8] + d[24]*u[9]);
106: uik[10]= -(d[0]*u[10] + d[5]*u[11] + d[10]*u[12] + d[15]*u[13] + d[20]*u[14]);
107: uik[11]= -(d[1]*u[10] + d[6]*u[11] + d[11]*u[12] + d[16]*u[13] + d[21]*u[14]);
108: uik[12]= -(d[2]*u[10] + d[7]*u[11] + d[12]*u[12] + d[17]*u[13] + d[22]*u[14]);
109: uik[13]= -(d[3]*u[10] + d[8]*u[11] + d[13]*u[12] + d[18]*u[13] + d[23]*u[14]);
110: uik[14]= -(d[4]*u[10] + d[9]*u[11] + d[14]*u[12] + d[19]*u[13] + d[24]*u[14]);
112: uik[15]= -(d[0]*u[15] + d[5]*u[16] + d[10]*u[17] + d[15]*u[18] + d[20]*u[19]);
113: uik[16]= -(d[1]*u[15] + d[6]*u[16] + d[11]*u[17] + d[16]*u[18] + d[21]*u[19]);
114: uik[17]= -(d[2]*u[15] + d[7]*u[16] + d[12]*u[17] + d[17]*u[18] + d[22]*u[19]);
115: uik[18]= -(d[3]*u[15] + d[8]*u[16] + d[13]*u[17] + d[18]*u[18] + d[23]*u[19]);
116: uik[19]= -(d[4]*u[15] + d[9]*u[16] + d[14]*u[17] + d[19]*u[18] + d[24]*u[19]);
118: uik[20]= -(d[0]*u[20] + d[5]*u[21] + d[10]*u[22] + d[15]*u[23] + d[20]*u[24]);
119: uik[21]= -(d[1]*u[20] + d[6]*u[21] + d[11]*u[22] + d[16]*u[23] + d[21]*u[24]);
120: uik[22]= -(d[2]*u[20] + d[7]*u[21] + d[12]*u[22] + d[17]*u[23] + d[22]*u[24]);
121: uik[23]= -(d[3]*u[20] + d[8]*u[21] + d[13]*u[22] + d[18]*u[23] + d[23]*u[24]);
122: uik[24]= -(d[4]*u[20] + d[9]*u[21] + d[14]*u[22] + d[19]*u[23] + d[24]*u[24]);
125: /* update D(k) += -U(i,k)^T * U_bar(i,k) */
126: dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3] + uik[4]*u[4];
127: dk[1] += uik[5]*u[0] + uik[6]*u[1] + uik[7]*u[2] + uik[8]*u[3] + uik[9]*u[4];
128: dk[2] += uik[10]*u[0]+ uik[11]*u[1]+ uik[12]*u[2]+ uik[13]*u[3]+ uik[14]*u[4];
129: dk[3] += uik[15]*u[0]+ uik[16]*u[1]+ uik[17]*u[2]+ uik[18]*u[3]+ uik[19]*u[4];
130: dk[4] += uik[20]*u[0]+ uik[21]*u[1]+ uik[22]*u[2]+ uik[23]*u[3]+ uik[24]*u[4];
132: dk[5] += uik[0]*u[5] + uik[1]*u[6] + uik[2]*u[7] + uik[3]*u[8] + uik[4]*u[9];
133: dk[6] += uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8] + uik[9]*u[9];
134: dk[7] += uik[10]*u[5]+ uik[11]*u[6]+ uik[12]*u[7]+ uik[13]*u[8]+ uik[14]*u[9];
135: dk[8] += uik[15]*u[5]+ uik[16]*u[6]+ uik[17]*u[7]+ uik[18]*u[8]+ uik[19]*u[9];
136: dk[9] += uik[20]*u[5]+ uik[21]*u[6]+ uik[22]*u[7]+ uik[23]*u[8]+ uik[24]*u[9];
138: dk[10] += uik[0]*u[10] + uik[1]*u[11] + uik[2]*u[12] + uik[3]*u[13] + uik[4]*u[14];
139: dk[11] += uik[5]*u[10] + uik[6]*u[11] + uik[7]*u[12] + uik[8]*u[13] + uik[9]*u[14];
140: dk[12] += uik[10]*u[10]+ uik[11]*u[11]+ uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14];
141: dk[13] += uik[15]*u[10]+ uik[16]*u[11]+ uik[17]*u[12]+ uik[18]*u[13]+ uik[19]*u[14];
142: dk[14] += uik[20]*u[10]+ uik[21]*u[11]+ uik[22]*u[12]+ uik[23]*u[13]+ uik[24]*u[14];
144: dk[15] += uik[0]*u[15] + uik[1]*u[16] + uik[2]*u[17] + uik[3]*u[18] + uik[4]*u[19];
145: dk[16] += uik[5]*u[15] + uik[6]*u[16] + uik[7]*u[17] + uik[8]*u[18] + uik[9]*u[19];
146: dk[17] += uik[10]*u[15]+ uik[11]*u[16]+ uik[12]*u[17]+ uik[13]*u[18]+ uik[14]*u[19];
147: dk[18] += uik[15]*u[15]+ uik[16]*u[16]+ uik[17]*u[17]+ uik[18]*u[18]+ uik[19]*u[19];
148: dk[19] += uik[20]*u[15]+ uik[21]*u[16]+ uik[22]*u[17]+ uik[23]*u[18]+ uik[24]*u[19];
150: dk[20] += uik[0]*u[20] + uik[1]*u[21] + uik[2]*u[22] + uik[3]*u[23] + uik[4]*u[24];
151: dk[21] += uik[5]*u[20] + uik[6]*u[21] + uik[7]*u[22] + uik[8]*u[23] + uik[9]*u[24];
152: dk[22] += uik[10]*u[20]+ uik[11]*u[21]+ uik[12]*u[22]+ uik[13]*u[23]+ uik[14]*u[24];
153: dk[23] += uik[15]*u[20]+ uik[16]*u[21]+ uik[17]*u[22]+ uik[18]*u[23]+ uik[19]*u[24];
154: dk[24] += uik[20]*u[20]+ uik[21]*u[21]+ uik[22]*u[22]+ uik[23]*u[23]+ uik[24]*u[24];
156: PetscLogFlops(125.0*4.0);
158: /* update -U(i,k) */
159: PetscMemcpy(ba+ili*25,uik,25*sizeof(MatScalar));
161: /* add multiple of row i to k-th row ... */
162: jmin = ili + 1; jmax = bi[i+1];
163: if (jmin < jmax) {
164: for (j=jmin; j<jmax; j++) {
165: /* rtmp += -U(i,k)^T * U_bar(i,j) */
166: rtmp_ptr = rtmp + bj[j]*25;
167: u = ba + j*25;
168: rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3] + uik[4]*u[4];
169: rtmp_ptr[1] += uik[5]*u[0] + uik[6]*u[1] + uik[7]*u[2] + uik[8]*u[3] + uik[9]*u[4];
170: rtmp_ptr[2] += uik[10]*u[0]+ uik[11]*u[1]+ uik[12]*u[2]+ uik[13]*u[3]+ uik[14]*u[4];
171: rtmp_ptr[3] += uik[15]*u[0]+ uik[16]*u[1]+ uik[17]*u[2]+ uik[18]*u[3]+ uik[19]*u[4];
172: rtmp_ptr[4] += uik[20]*u[0]+ uik[21]*u[1]+ uik[22]*u[2]+ uik[23]*u[3]+ uik[24]*u[4];
174: rtmp_ptr[5] += uik[0]*u[5] + uik[1]*u[6] + uik[2]*u[7] + uik[3]*u[8] + uik[4]*u[9];
175: rtmp_ptr[6] += uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8] + uik[9]*u[9];
176: rtmp_ptr[7] += uik[10]*u[5]+ uik[11]*u[6]+ uik[12]*u[7]+ uik[13]*u[8]+ uik[14]*u[9];
177: rtmp_ptr[8] += uik[15]*u[5]+ uik[16]*u[6]+ uik[17]*u[7]+ uik[18]*u[8]+ uik[19]*u[9];
178: rtmp_ptr[9] += uik[20]*u[5]+ uik[21]*u[6]+ uik[22]*u[7]+ uik[23]*u[8]+ uik[24]*u[9];
180: rtmp_ptr[10] += uik[0]*u[10] + uik[1]*u[11] + uik[2]*u[12] + uik[3]*u[13] + uik[4]*u[14];
181: rtmp_ptr[11] += uik[5]*u[10] + uik[6]*u[11] + uik[7]*u[12] + uik[8]*u[13] + uik[9]*u[14];
182: rtmp_ptr[12] += uik[10]*u[10]+ uik[11]*u[11]+ uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14];
183: rtmp_ptr[13] += uik[15]*u[10]+ uik[16]*u[11]+ uik[17]*u[12]+ uik[18]*u[13]+ uik[19]*u[14];
184: rtmp_ptr[14] += uik[20]*u[10]+ uik[21]*u[11]+ uik[22]*u[12]+ uik[23]*u[13]+ uik[24]*u[14];
186: rtmp_ptr[15] += uik[0]*u[15] + uik[1]*u[16] + uik[2]*u[17] + uik[3]*u[18] + uik[4]*u[19];
187: rtmp_ptr[16] += uik[5]*u[15] + uik[6]*u[16] + uik[7]*u[17] + uik[8]*u[18] + uik[9]*u[19];
188: rtmp_ptr[17] += uik[10]*u[15]+ uik[11]*u[16]+ uik[12]*u[17]+ uik[13]*u[18]+ uik[14]*u[19];
189: rtmp_ptr[18] += uik[15]*u[15]+ uik[16]*u[16]+ uik[17]*u[17]+ uik[18]*u[18]+ uik[19]*u[19];
190: rtmp_ptr[19] += uik[20]*u[15]+ uik[21]*u[16]+ uik[22]*u[17]+ uik[23]*u[18]+ uik[24]*u[19];
192: rtmp_ptr[20] += uik[0]*u[20] + uik[1]*u[21] + uik[2]*u[22] + uik[3]*u[23] + uik[4]*u[24];
193: rtmp_ptr[21] += uik[5]*u[20] + uik[6]*u[21] + uik[7]*u[22] + uik[8]*u[23] + uik[9]*u[24];
194: rtmp_ptr[22] += uik[10]*u[20]+ uik[11]*u[21]+ uik[12]*u[22]+ uik[13]*u[23]+ uik[14]*u[24];
195: rtmp_ptr[23] += uik[15]*u[20]+ uik[16]*u[21]+ uik[17]*u[22]+ uik[18]*u[23]+ uik[19]*u[24];
196: rtmp_ptr[24] += uik[20]*u[20]+ uik[21]*u[21]+ uik[22]*u[22]+ uik[23]*u[23]+ uik[24]*u[24];
197: }
198: PetscLogFlops(2.0*125.0*(jmax-jmin));
200: /* ... add i to row list for next nonzero entry */
201: il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
202: j = bj[jmin];
203: jl[i] = jl[j]; jl[j] = i; /* update jl */
204: }
205: i = nexti;
206: }
208: /* save nonzero entries in k-th row of U ... */
210: /* invert diagonal block */
211: d = ba+k*25;
212: PetscMemcpy(d,dk,25*sizeof(MatScalar));
213: PetscKernel_A_gets_inverse_A_5(d,ipvt,work,shift,allowzeropivot,&zeropivotdetected);
214: if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
216: jmin = bi[k]; jmax = bi[k+1];
217: if (jmin < jmax) {
218: for (j=jmin; j<jmax; j++) {
219: vj = bj[j]; /* block col. index of U */
220: u = ba + j*25;
221: rtmp_ptr = rtmp + vj*25;
222: for (k1=0; k1<25; k1++) {
223: *u++ = *rtmp_ptr;
224: *rtmp_ptr++ = 0.0;
225: }
226: }
228: /* ... add k to row list for first nonzero entry in k-th row */
229: il[k] = jmin;
230: i = bj[jmin];
231: jl[k] = jl[i]; jl[i] = k;
232: }
233: }
235: PetscFree(rtmp);
236: PetscFree2(il,jl);
237: PetscFree2(dk,uik);
238: if (a->permute) {
239: PetscFree(aa);
240: }
242: ISRestoreIndices(perm,&perm_ptr);
244: C->ops->solve = MatSolve_SeqSBAIJ_5_inplace;
245: C->ops->solvetranspose = MatSolve_SeqSBAIJ_5_inplace;
246: C->assembled = PETSC_TRUE;
247: C->preallocated = PETSC_TRUE;
249: PetscLogFlops(1.3333*125*b->mbs); /* from inverting diagonal blocks */
250: return(0);
251: }