Actual source code: sbaijfact5.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: /*
6: Version for when blocks are 4 by 4 Using natural ordering
7: */
10: PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
11: {
12: Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
14: PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
15: PetscInt *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
16: MatScalar *ba = b->a,*aa,*ap,*dk,*uik;
17: MatScalar *u,*diag,*rtmp,*rtmp_ptr;
18: PetscBool pivotinblocks = b->pivotinblocks;
19: PetscReal shift = info->shiftamount;
20: PetscBool allowzeropivot,zeropivotdetected=PETSC_FALSE;
23: /* initialization */
24: allowzeropivot = PetscNot(A->erroriffailure);
25: PetscCalloc1(16*mbs,&rtmp);
26: PetscMalloc2(mbs,&il,mbs,&jl);
27: il[0] = 0;
28: for (i=0; i<mbs; i++) jl[i] = mbs;
29:
30: PetscMalloc2(16,&dk,16,&uik);
31: ai = a->i; aj = a->j; aa = a->a;
33: /* for each row k */
34: for (k = 0; k<mbs; k++) {
36: /*initialize k-th row with elements nonzero in row k of A */
37: jmin = ai[k]; jmax = ai[k+1];
38: if (jmin < jmax) {
39: ap = aa + jmin*16;
40: for (j = jmin; j < jmax; j++) {
41: vj = aj[j]; /* block col. index */
42: rtmp_ptr = rtmp + vj*16;
43: for (i=0; i<16; i++) *rtmp_ptr++ = *ap++;
44: }
45: }
47: /* modify k-th row by adding in those rows i with U(i,k) != 0 */
48: PetscMemcpy(dk,rtmp+k*16,16*sizeof(MatScalar));
49: i = jl[k]; /* first row to be added to k_th row */
51: while (i < mbs) {
52: nexti = jl[i]; /* next row to be added to k_th row */
54: /* compute multiplier */
55: ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */
57: /* uik = -inv(Di)*U_bar(i,k) */
58: diag = ba + i*16;
59: u = ba + ili*16;
61: uik[0] = -(diag[0]*u[0] + diag[4]*u[1] + diag[8]*u[2] + diag[12]*u[3]);
62: uik[1] = -(diag[1]*u[0] + diag[5]*u[1] + diag[9]*u[2] + diag[13]*u[3]);
63: uik[2] = -(diag[2]*u[0] + diag[6]*u[1] + diag[10]*u[2]+ diag[14]*u[3]);
64: uik[3] = -(diag[3]*u[0] + diag[7]*u[1] + diag[11]*u[2]+ diag[15]*u[3]);
66: uik[4] = -(diag[0]*u[4] + diag[4]*u[5] + diag[8]*u[6] + diag[12]*u[7]);
67: uik[5] = -(diag[1]*u[4] + diag[5]*u[5] + diag[9]*u[6] + diag[13]*u[7]);
68: uik[6] = -(diag[2]*u[4] + diag[6]*u[5] + diag[10]*u[6]+ diag[14]*u[7]);
69: uik[7] = -(diag[3]*u[4] + diag[7]*u[5] + diag[11]*u[6]+ diag[15]*u[7]);
71: uik[8] = -(diag[0]*u[8] + diag[4]*u[9] + diag[8]*u[10] + diag[12]*u[11]);
72: uik[9] = -(diag[1]*u[8] + diag[5]*u[9] + diag[9]*u[10] + diag[13]*u[11]);
73: uik[10]= -(diag[2]*u[8] + diag[6]*u[9] + diag[10]*u[10]+ diag[14]*u[11]);
74: uik[11]= -(diag[3]*u[8] + diag[7]*u[9] + diag[11]*u[10]+ diag[15]*u[11]);
76: uik[12]= -(diag[0]*u[12] + diag[4]*u[13] + diag[8]*u[14] + diag[12]*u[15]);
77: uik[13]= -(diag[1]*u[12] + diag[5]*u[13] + diag[9]*u[14] + diag[13]*u[15]);
78: uik[14]= -(diag[2]*u[12] + diag[6]*u[13] + diag[10]*u[14]+ diag[14]*u[15]);
79: uik[15]= -(diag[3]*u[12] + diag[7]*u[13] + diag[11]*u[14]+ diag[15]*u[15]);
81: /* update D(k) += -U(i,k)^T * U_bar(i,k) */
82: dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3];
83: dk[1] += uik[4]*u[0] + uik[5]*u[1] + uik[6]*u[2] + uik[7]*u[3];
84: dk[2] += uik[8]*u[0] + uik[9]*u[1] + uik[10]*u[2]+ uik[11]*u[3];
85: dk[3] += uik[12]*u[0]+ uik[13]*u[1]+ uik[14]*u[2]+ uik[15]*u[3];
87: dk[4] += uik[0]*u[4] + uik[1]*u[5] + uik[2]*u[6] + uik[3]*u[7];
88: dk[5] += uik[4]*u[4] + uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7];
89: dk[6] += uik[8]*u[4] + uik[9]*u[5] + uik[10]*u[6]+ uik[11]*u[7];
90: dk[7] += uik[12]*u[4]+ uik[13]*u[5]+ uik[14]*u[6]+ uik[15]*u[7];
92: dk[8] += uik[0]*u[8] + uik[1]*u[9] + uik[2]*u[10] + uik[3]*u[11];
93: dk[9] += uik[4]*u[8] + uik[5]*u[9] + uik[6]*u[10] + uik[7]*u[11];
94: dk[10]+= uik[8]*u[8] + uik[9]*u[9] + uik[10]*u[10]+ uik[11]*u[11];
95: dk[11]+= uik[12]*u[8]+ uik[13]*u[9]+ uik[14]*u[10]+ uik[15]*u[11];
97: dk[12]+= uik[0]*u[12] + uik[1]*u[13] + uik[2]*u[14] + uik[3]*u[15];
98: dk[13]+= uik[4]*u[12] + uik[5]*u[13] + uik[6]*u[14] + uik[7]*u[15];
99: dk[14]+= uik[8]*u[12] + uik[9]*u[13] + uik[10]*u[14]+ uik[11]*u[15];
100: dk[15]+= uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14]+ uik[15]*u[15];
102: PetscLogFlops(64.0*4.0);
104: /* update -U(i,k) */
105: PetscMemcpy(ba+ili*16,uik,16*sizeof(MatScalar));
107: /* add multiple of row i to k-th row ... */
108: jmin = ili + 1; jmax = bi[i+1];
109: if (jmin < jmax) {
110: for (j=jmin; j<jmax; j++) {
111: /* rtmp += -U(i,k)^T * U_bar(i,j) */
112: rtmp_ptr = rtmp + bj[j]*16;
113: u = ba + j*16;
114: rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3];
115: rtmp_ptr[1] += uik[4]*u[0] + uik[5]*u[1] + uik[6]*u[2] + uik[7]*u[3];
116: rtmp_ptr[2] += uik[8]*u[0] + uik[9]*u[1] + uik[10]*u[2]+ uik[11]*u[3];
117: rtmp_ptr[3] += uik[12]*u[0]+ uik[13]*u[1]+ uik[14]*u[2]+ uik[15]*u[3];
119: rtmp_ptr[4] += uik[0]*u[4] + uik[1]*u[5] + uik[2]*u[6] + uik[3]*u[7];
120: rtmp_ptr[5] += uik[4]*u[4] + uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7];
121: rtmp_ptr[6] += uik[8]*u[4] + uik[9]*u[5] + uik[10]*u[6]+ uik[11]*u[7];
122: rtmp_ptr[7] += uik[12]*u[4]+ uik[13]*u[5]+ uik[14]*u[6]+ uik[15]*u[7];
124: rtmp_ptr[8] += uik[0]*u[8] + uik[1]*u[9] + uik[2]*u[10] + uik[3]*u[11];
125: rtmp_ptr[9] += uik[4]*u[8] + uik[5]*u[9] + uik[6]*u[10] + uik[7]*u[11];
126: rtmp_ptr[10]+= uik[8]*u[8] + uik[9]*u[9] + uik[10]*u[10]+ uik[11]*u[11];
127: rtmp_ptr[11]+= uik[12]*u[8]+ uik[13]*u[9]+ uik[14]*u[10]+ uik[15]*u[11];
129: rtmp_ptr[12]+= uik[0]*u[12] + uik[1]*u[13] + uik[2]*u[14] + uik[3]*u[15];
130: rtmp_ptr[13]+= uik[4]*u[12] + uik[5]*u[13] + uik[6]*u[14] + uik[7]*u[15];
131: rtmp_ptr[14]+= uik[8]*u[12] + uik[9]*u[13] + uik[10]*u[14]+ uik[11]*u[15];
132: rtmp_ptr[15]+= uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14]+ uik[15]*u[15];
133: }
134: PetscLogFlops(2.0*64.0*(jmax-jmin));
136: /* ... add i to row list for next nonzero entry */
137: il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */
138: j = bj[jmin];
139: jl[i] = jl[j]; jl[j] = i; /* update jl */
140: }
141: i = nexti;
142: }
144: /* save nonzero entries in k-th row of U ... */
146: /* invert diagonal block */
147: diag = ba+k*16;
148: PetscMemcpy(diag,dk,16*sizeof(MatScalar));
149: if (pivotinblocks) {
150: PetscKernel_A_gets_inverse_A_4(diag,shift,allowzeropivot,&zeropivotdetected);
151: if (zeropivotdetected) C->errortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
152: } else {
153: PetscKernel_A_gets_inverse_A_4_nopivot(diag);
154: }
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: }
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);
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;
187: PetscLogFlops(1.3333*64*b->mbs); /* from inverting diagonal blocks */
188: return(0);
189: }