Actual source code: sbaijfact5.c
petsc-3.11.4 2019-09-28
2: #include <../src/mat/impls/sbaij/seq/sbaij.h>
3: #include <petsc/private/kernels/blockinvert.h>
5: /*
6: Version for when blocks are 4 by 4 Using natural ordering
7: */
8: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
9: {
10: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
12: PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
13: PetscInt *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
14: MatScalar *ba = b->a,*aa,*ap,*dk,*uik;
15: MatScalar *u,*diag,*rtmp,*rtmp_ptr;
16: PetscBool pivotinblocks = b->pivotinblocks;
17: PetscReal shift = info->shiftamount;
18: PetscBool allowzeropivot,zeropivotdetected=PETSC_FALSE;
21: /* initialization */
22: allowzeropivot = PetscNot(A->erroriffailure);
23: PetscCalloc1(16*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(16,&dk,16,&uik);
29: ai = a->i; aj = a->j; aa = a->a;
31: /* for each row k */
32: for (k = 0; k<mbs; k++) {
34: /*initialize k-th row with elements nonzero in row k of A */
35: jmin = ai[k]; jmax = ai[k+1];
36: if (jmin < jmax) {
37: ap = aa + jmin*16;
38: for (j = jmin; j < jmax; j++) {
39: vj = aj[j]; /* block col. index */
40: rtmp_ptr = rtmp + vj*16;
41: for (i=0; i<16; i++) *rtmp_ptr++ = *ap++;
42: }
43: }
45: /* modify k-th row by adding in those rows i with U(i,k) != 0 */
46: PetscMemcpy(dk,rtmp+k*16,16*sizeof(MatScalar));
47: i = jl[k]; /* first row to be added to k_th row */
49: while (i < mbs) {
50: nexti = jl[i]; /* next row to be added to k_th row */
52: /* compute multiplier */
53: ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
55: /* uik = -inv(Di)*U_bar(i,k) */
56: diag = ba + i*16;
57: u = ba + ili*16;
59: uik[0] = -(diag[0]*u[0] + diag[4]*u[1] + diag[8]*u[2] + diag[12]*u[3]);
60: uik[1] = -(diag[1]*u[0] + diag[5]*u[1] + diag[9]*u[2] + diag[13]*u[3]);
61: uik[2] = -(diag[2]*u[0] + diag[6]*u[1] + diag[10]*u[2]+ diag[14]*u[3]);
62: uik[3] = -(diag[3]*u[0] + diag[7]*u[1] + diag[11]*u[2]+ diag[15]*u[3]);
64: uik[4] = -(diag[0]*u[4] + diag[4]*u[5] + diag[8]*u[6] + diag[12]*u[7]);
65: uik[5] = -(diag[1]*u[4] + diag[5]*u[5] + diag[9]*u[6] + diag[13]*u[7]);
66: uik[6] = -(diag[2]*u[4] + diag[6]*u[5] + diag[10]*u[6]+ diag[14]*u[7]);
67: uik[7] = -(diag[3]*u[4] + diag[7]*u[5] + diag[11]*u[6]+ diag[15]*u[7]);
69: uik[8] = -(diag[0]*u[8] + diag[4]*u[9] + diag[8]*u[10] + diag[12]*u[11]);
70: uik[9] = -(diag[1]*u[8] + diag[5]*u[9] + diag[9]*u[10] + diag[13]*u[11]);
71: uik[10]= -(diag[2]*u[8] + diag[6]*u[9] + diag[10]*u[10]+ diag[14]*u[11]);
72: uik[11]= -(diag[3]*u[8] + diag[7]*u[9] + diag[11]*u[10]+ diag[15]*u[11]);
74: uik[12]= -(diag[0]*u[12] + diag[4]*u[13] + diag[8]*u[14] + diag[12]*u[15]);
75: uik[13]= -(diag[1]*u[12] + diag[5]*u[13] + diag[9]*u[14] + diag[13]*u[15]);
76: uik[14]= -(diag[2]*u[12] + diag[6]*u[13] + diag[10]*u[14]+ diag[14]*u[15]);
77: uik[15]= -(diag[3]*u[12] + diag[7]*u[13] + diag[11]*u[14]+ diag[15]*u[15]);
79: /* update D(k) += -U(i,k)^T * U_bar(i,k) */
80: dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3];
81: dk[1] += uik[4]*u[0] + uik[5]*u[1] + uik[6]*u[2] + uik[7]*u[3];
82: dk[2] += uik[8]*u[0] + uik[9]*u[1] + uik[10]*u[2]+ uik[11]*u[3];
83: dk[3] += uik[12]*u[0]+ uik[13]*u[1]+ uik[14]*u[2]+ uik[15]*u[3];
85: dk[4] += uik[0]*u[4] + uik[1]*u[5] + uik[2]*u[6] + uik[3]*u[7];
86: dk[5] += uik[4]*u[4] + uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7];
87: dk[6] += uik[8]*u[4] + uik[9]*u[5] + uik[10]*u[6]+ uik[11]*u[7];
88: dk[7] += uik[12]*u[4]+ uik[13]*u[5]+ uik[14]*u[6]+ uik[15]*u[7];
90: dk[8] += uik[0]*u[8] + uik[1]*u[9] + uik[2]*u[10] + uik[3]*u[11];
91: dk[9] += uik[4]*u[8] + uik[5]*u[9] + uik[6]*u[10] + uik[7]*u[11];
92: dk[10]+= uik[8]*u[8] + uik[9]*u[9] + uik[10]*u[10]+ uik[11]*u[11];
93: dk[11]+= uik[12]*u[8]+ uik[13]*u[9]+ uik[14]*u[10]+ uik[15]*u[11];
95: dk[12]+= uik[0]*u[12] + uik[1]*u[13] + uik[2]*u[14] + uik[3]*u[15];
96: dk[13]+= uik[4]*u[12] + uik[5]*u[13] + uik[6]*u[14] + uik[7]*u[15];
97: dk[14]+= uik[8]*u[12] + uik[9]*u[13] + uik[10]*u[14]+ uik[11]*u[15];
98: dk[15]+= uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14]+ uik[15]*u[15];
100: PetscLogFlops(64.0*4.0);
102: /* update -U(i,k) */
103: PetscMemcpy(ba+ili*16,uik,16*sizeof(MatScalar));
105: /* add multiple of row i to k-th row ... */
106: jmin = ili + 1; jmax = bi[i+1];
107: if (jmin < jmax) {
108: for (j=jmin; j<jmax; j++) {
109: /* rtmp += -U(i,k)^T * U_bar(i,j) */
110: rtmp_ptr = rtmp + bj[j]*16;
111: u = ba + j*16;
112: rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3];
113: rtmp_ptr[1] += uik[4]*u[0] + uik[5]*u[1] + uik[6]*u[2] + uik[7]*u[3];
114: rtmp_ptr[2] += uik[8]*u[0] + uik[9]*u[1] + uik[10]*u[2]+ uik[11]*u[3];
115: rtmp_ptr[3] += uik[12]*u[0]+ uik[13]*u[1]+ uik[14]*u[2]+ uik[15]*u[3];
117: rtmp_ptr[4] += uik[0]*u[4] + uik[1]*u[5] + uik[2]*u[6] + uik[3]*u[7];
118: rtmp_ptr[5] += uik[4]*u[4] + uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7];
119: rtmp_ptr[6] += uik[8]*u[4] + uik[9]*u[5] + uik[10]*u[6]+ uik[11]*u[7];
120: rtmp_ptr[7] += uik[12]*u[4]+ uik[13]*u[5]+ uik[14]*u[6]+ uik[15]*u[7];
122: rtmp_ptr[8] += uik[0]*u[8] + uik[1]*u[9] + uik[2]*u[10] + uik[3]*u[11];
123: rtmp_ptr[9] += uik[4]*u[8] + uik[5]*u[9] + uik[6]*u[10] + uik[7]*u[11];
124: rtmp_ptr[10]+= uik[8]*u[8] + uik[9]*u[9] + uik[10]*u[10]+ uik[11]*u[11];
125: rtmp_ptr[11]+= uik[12]*u[8]+ uik[13]*u[9]+ uik[14]*u[10]+ uik[15]*u[11];
127: rtmp_ptr[12]+= uik[0]*u[12] + uik[1]*u[13] + uik[2]*u[14] + uik[3]*u[15];
128: rtmp_ptr[13]+= uik[4]*u[12] + uik[5]*u[13] + uik[6]*u[14] + uik[7]*u[15];
129: rtmp_ptr[14]+= uik[8]*u[12] + uik[9]*u[13] + uik[10]*u[14]+ uik[11]*u[15];
130: rtmp_ptr[15]+= uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14]+ uik[15]*u[15];
131: }
132: PetscLogFlops(2.0*64.0*(jmax-jmin));
134: /* ... add i to row list for next nonzero entry */
135: il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
136: j = bj[jmin];
137: jl[i] = jl[j]; jl[j] = i; /* update jl */
138: }
139: i = nexti;
140: }
142: /* save nonzero entries in k-th row of U ... */
144: /* invert diagonal block */
145: diag = ba+k*16;
146: PetscMemcpy(diag,dk,16*sizeof(MatScalar));
147: if (pivotinblocks) {
148: PetscKernel_A_gets_inverse_A_4(diag,shift,allowzeropivot,&zeropivotdetected);
149: if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
150: } else {
151: PetscKernel_A_gets_inverse_A_4_nopivot(diag);
152: }
154: jmin = bi[k]; jmax = bi[k+1];
155: if (jmin < jmax) {
156: for (j=jmin; j<jmax; j++) {
157: vj = bj[j]; /* block col. index of U */
158: u = ba + j*16;
159: rtmp_ptr = rtmp + vj*16;
160: for (k1=0; k1<16; k1++) {
161: *u++ = *rtmp_ptr;
162: *rtmp_ptr++ = 0.0;
163: }
164: }
166: /* ... add k to row list for first nonzero entry in k-th row */
167: il[k] = jmin;
168: i = bj[jmin];
169: jl[k] = jl[i]; jl[i] = k;
170: }
171: }
173: PetscFree(rtmp);
174: PetscFree2(il,jl);
175: PetscFree2(dk,uik);
177: C->ops->solve = MatSolve_SeqSBAIJ_4_NaturalOrdering_inplace;
178: C->ops->solvetranspose = MatSolve_SeqSBAIJ_4_NaturalOrdering_inplace;
179: C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_4_NaturalOrdering_inplace;
180: C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_4_NaturalOrdering_inplace;
182: C->assembled = PETSC_TRUE;
183: C->preallocated = PETSC_TRUE;
185: PetscLogFlops(1.3333*64*b->mbs); /* from inverting diagonal blocks */
186: return(0);
187: }