Actual source code: sbaijfact6.c
petsc-3.12.5 2020-03-29
2: #include <../src/mat/impls/sbaij/seq/sbaij.h>
3: #include <petsc/private/kernels/blockinvert.h>
5: /* Version for when blocks are 4 by 4 */
6: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat C,Mat A,const MatFactorInfo *info)
7: {
8: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
9: IS perm = b->row;
11: const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
12: PetscInt i,j,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
13: MatScalar *ba = b->a,*aa,*ap,*dk,*uik;
14: MatScalar *u,*diag,*rtmp,*rtmp_ptr;
15: PetscBool pivotinblocks = b->pivotinblocks;
16: PetscReal shift = info->shiftamount;
17: PetscBool allowzeropivot,zeropivotdetected=PETSC_FALSE;
20: /* initialization */
21: allowzeropivot = PetscNot(A->erroriffailure);
22: PetscCalloc1(16*mbs,&rtmp);
23: PetscMalloc2(mbs,&il,mbs,&jl);
24: il[0] = 0;
25: for (i=0; i<mbs; i++) jl[i] = mbs;
27: PetscMalloc2(16,&dk,16,&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: PetscMalloc1(16*ai[mbs],&aa);
36: PetscArraycpy(aa,a->a,16*ai[mbs]);
37: PetscMalloc1(ai[mbs],&a2anew);
38: PetscArraycpy(a2anew,a->a2anew,ai[mbs]);
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<16; k1++) {
46: dk[k1] = aa[k*16+k1];
47: aa[k*16+k1] = aa[j*16+k1];
48: aa[j*16+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*16; /* ptr to the beginning of j-th block of aa */
55: for (k=0; k<16; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
56: for (k=0; k<4; k++) { /* j-th block of aa <- dk^T */
57: for (k1=0; k1<4; k1++) *ap++ = dk[k + 4*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*16;
72: for (j = jmin; j < jmax; j++) {
73: vj = perm_ptr[aj[j]]; /* block col. index */
74: rtmp_ptr = rtmp + vj*16;
75: for (i=0; i<16; i++) *rtmp_ptr++ = *ap++;
76: }
77: }
79: /* modify k-th row by adding in those rows i with U(i,k) != 0 */
80: PetscArraycpy(dk,rtmp+k*16,16);
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*16;
91: u = ba + ili*16;
93: uik[0] = -(diag[0]*u[0] + diag[4]*u[1] + diag[8]*u[2] + diag[12]*u[3]);
94: uik[1] = -(diag[1]*u[0] + diag[5]*u[1] + diag[9]*u[2] + diag[13]*u[3]);
95: uik[2] = -(diag[2]*u[0] + diag[6]*u[1] + diag[10]*u[2]+ diag[14]*u[3]);
96: uik[3] = -(diag[3]*u[0] + diag[7]*u[1] + diag[11]*u[2]+ diag[15]*u[3]);
98: uik[4] = -(diag[0]*u[4] + diag[4]*u[5] + diag[8]*u[6] + diag[12]*u[7]);
99: uik[5] = -(diag[1]*u[4] + diag[5]*u[5] + diag[9]*u[6] + diag[13]*u[7]);
100: uik[6] = -(diag[2]*u[4] + diag[6]*u[5] + diag[10]*u[6]+ diag[14]*u[7]);
101: uik[7] = -(diag[3]*u[4] + diag[7]*u[5] + diag[11]*u[6]+ diag[15]*u[7]);
103: uik[8] = -(diag[0]*u[8] + diag[4]*u[9] + diag[8]*u[10] + diag[12]*u[11]);
104: uik[9] = -(diag[1]*u[8] + diag[5]*u[9] + diag[9]*u[10] + diag[13]*u[11]);
105: uik[10]= -(diag[2]*u[8] + diag[6]*u[9] + diag[10]*u[10]+ diag[14]*u[11]);
106: uik[11]= -(diag[3]*u[8] + diag[7]*u[9] + diag[11]*u[10]+ diag[15]*u[11]);
108: uik[12]= -(diag[0]*u[12] + diag[4]*u[13] + diag[8]*u[14] + diag[12]*u[15]);
109: uik[13]= -(diag[1]*u[12] + diag[5]*u[13] + diag[9]*u[14] + diag[13]*u[15]);
110: uik[14]= -(diag[2]*u[12] + diag[6]*u[13] + diag[10]*u[14]+ diag[14]*u[15]);
111: uik[15]= -(diag[3]*u[12] + diag[7]*u[13] + diag[11]*u[14]+ diag[15]*u[15]);
113: /* update D(k) += -U(i,k)^T * U_bar(i,k) */
114: dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3];
115: dk[1] += uik[4]*u[0] + uik[5]*u[1] + uik[6]*u[2] + uik[7]*u[3];
116: dk[2] += uik[8]*u[0] + uik[9]*u[1] + uik[10]*u[2]+ uik[11]*u[3];
117: dk[3] += uik[12]*u[0]+ uik[13]*u[1]+ uik[14]*u[2]+ uik[15]*u[3];
119: dk[4] += uik[0]*u[4] + uik[1]*u[5] + uik[2]*u[6] + uik[3]*u[7];
120: dk[5] += uik[4]*u[4] + uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7];
121: dk[6] += uik[8]*u[4] + uik[9]*u[5] + uik[10]*u[6]+ uik[11]*u[7];
122: dk[7] += uik[12]*u[4]+ uik[13]*u[5]+ uik[14]*u[6]+ uik[15]*u[7];
124: dk[8] += uik[0]*u[8] + uik[1]*u[9] + uik[2]*u[10] + uik[3]*u[11];
125: dk[9] += uik[4]*u[8] + uik[5]*u[9] + uik[6]*u[10] + uik[7]*u[11];
126: dk[10]+= uik[8]*u[8] + uik[9]*u[9] + uik[10]*u[10]+ uik[11]*u[11];
127: dk[11]+= uik[12]*u[8]+ uik[13]*u[9]+ uik[14]*u[10]+ uik[15]*u[11];
129: dk[12]+= uik[0]*u[12] + uik[1]*u[13] + uik[2]*u[14] + uik[3]*u[15];
130: dk[13]+= uik[4]*u[12] + uik[5]*u[13] + uik[6]*u[14] + uik[7]*u[15];
131: dk[14]+= uik[8]*u[12] + uik[9]*u[13] + uik[10]*u[14]+ uik[11]*u[15];
132: dk[15]+= uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14]+ uik[15]*u[15];
134: PetscLogFlops(64.0*4.0);
136: /* update -U(i,k) */
137: PetscArraycpy(ba+ili*16,uik,16);
139: /* add multiple of row i to k-th row ... */
140: jmin = ili + 1; jmax = bi[i+1];
141: if (jmin < jmax) {
142: for (j=jmin; j<jmax; j++) {
143: /* rtmp += -U(i,k)^T * U_bar(i,j) */
144: rtmp_ptr = rtmp + bj[j]*16;
145: u = ba + j*16;
146: rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3];
147: rtmp_ptr[1] += uik[4]*u[0] + uik[5]*u[1] + uik[6]*u[2] + uik[7]*u[3];
148: rtmp_ptr[2] += uik[8]*u[0] + uik[9]*u[1] + uik[10]*u[2]+ uik[11]*u[3];
149: rtmp_ptr[3] += uik[12]*u[0]+ uik[13]*u[1]+ uik[14]*u[2]+ uik[15]*u[3];
151: rtmp_ptr[4] += uik[0]*u[4] + uik[1]*u[5] + uik[2]*u[6] + uik[3]*u[7];
152: rtmp_ptr[5] += uik[4]*u[4] + uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7];
153: rtmp_ptr[6] += uik[8]*u[4] + uik[9]*u[5] + uik[10]*u[6]+ uik[11]*u[7];
154: rtmp_ptr[7] += uik[12]*u[4]+ uik[13]*u[5]+ uik[14]*u[6]+ uik[15]*u[7];
156: rtmp_ptr[8] += uik[0]*u[8] + uik[1]*u[9] + uik[2]*u[10] + uik[3]*u[11];
157: rtmp_ptr[9] += uik[4]*u[8] + uik[5]*u[9] + uik[6]*u[10] + uik[7]*u[11];
158: rtmp_ptr[10]+= uik[8]*u[8] + uik[9]*u[9] + uik[10]*u[10]+ uik[11]*u[11];
159: rtmp_ptr[11]+= uik[12]*u[8]+ uik[13]*u[9]+ uik[14]*u[10]+ uik[15]*u[11];
161: rtmp_ptr[12]+= uik[0]*u[12] + uik[1]*u[13] + uik[2]*u[14] + uik[3]*u[15];
162: rtmp_ptr[13]+= uik[4]*u[12] + uik[5]*u[13] + uik[6]*u[14] + uik[7]*u[15];
163: rtmp_ptr[14]+= uik[8]*u[12] + uik[9]*u[13] + uik[10]*u[14]+ uik[11]*u[15];
164: rtmp_ptr[15]+= uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14]+ uik[15]*u[15];
165: }
166: PetscLogFlops(2.0*64.0*(jmax-jmin));
168: /* ... add i to row list for next nonzero entry */
169: il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
170: j = bj[jmin];
171: jl[i] = jl[j]; jl[j] = i; /* update jl */
172: }
173: i = nexti;
174: }
176: /* save nonzero entries in k-th row of U ... */
178: /* invert diagonal block */
179: diag = ba+k*16;
180: PetscArraycpy(diag,dk,16);
182: if (pivotinblocks) {
183: PetscKernel_A_gets_inverse_A_4(diag,shift, allowzeropivot,&zeropivotdetected);
184: if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
185: } else {
186: PetscKernel_A_gets_inverse_A_4_nopivot(diag);
187: }
189: jmin = bi[k]; jmax = bi[k+1];
190: if (jmin < jmax) {
191: for (j=jmin; j<jmax; j++) {
192: vj = bj[j]; /* block col. index of U */
193: u = ba + j*16;
194: rtmp_ptr = rtmp + vj*16;
195: for (k1=0; k1<16; k1++) {
196: *u++ = *rtmp_ptr;
197: *rtmp_ptr++ = 0.0;
198: }
199: }
201: /* ... add k to row list for first nonzero entry in k-th row */
202: il[k] = jmin;
203: i = bj[jmin];
204: jl[k] = jl[i]; jl[i] = k;
205: }
206: }
208: PetscFree(rtmp);
209: PetscFree2(il,jl);
210: PetscFree2(dk,uik);
211: if (a->permute) {
212: PetscFree(aa);
213: }
215: ISRestoreIndices(perm,&perm_ptr);
217: C->ops->solve = MatSolve_SeqSBAIJ_4_inplace;
218: C->ops->solvetranspose = MatSolve_SeqSBAIJ_4_inplace;
219: C->assembled = PETSC_TRUE;
220: C->preallocated = PETSC_TRUE;
222: PetscLogFlops(1.3333*64*b->mbs); /* from inverting diagonal blocks */
223: return(0);
224: }