xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact6.c (revision 503c0ea9b45bcfbcebbb1ea5341243bbc69f0bea)
1 
2 #include <../src/mat/impls/sbaij/seq/sbaij.h>
3 #include <petsc/private/kernels/blockinvert.h>
4 
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;
10   const PetscInt *ai,*aj,*perm_ptr,mbs=a->mbs,*bi=b->i,*bj=b->j;
11   PetscInt       i,j,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
12   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
13   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
14   PetscBool      pivotinblocks = b->pivotinblocks;
15   PetscReal      shift         = info->shiftamount;
16   PetscBool      allowzeropivot,zeropivotdetected=PETSC_FALSE;
17 
18   PetscFunctionBegin;
19   /* initialization */
20   allowzeropivot = PetscNot(A->erroriffailure);
21   PetscCall(PetscCalloc1(16*mbs,&rtmp));
22   PetscCall(PetscMalloc2(mbs,&il,mbs,&jl));
23   il[0] = 0;
24   for (i=0; i<mbs; i++) jl[i] = mbs;
25 
26   PetscCall(PetscMalloc2(16,&dk,16,&uik));
27   PetscCall(ISGetIndices(perm,&perm_ptr));
28 
29   /* check permutation */
30   if (!a->permute) {
31     ai = a->i; aj = a->j; aa = a->a;
32   } else {
33     ai   = a->inew; aj = a->jnew;
34     PetscCall(PetscMalloc1(16*ai[mbs],&aa));
35     PetscCall(PetscArraycpy(aa,a->a,16*ai[mbs]));
36     PetscCall(PetscMalloc1(ai[mbs],&a2anew));
37     PetscCall(PetscArraycpy(a2anew,a->a2anew,ai[mbs]));
38 
39     for (i=0; i<mbs; i++) {
40       jmin = ai[i]; jmax = ai[i+1];
41       for (j=jmin; j<jmax; j++) {
42         while (a2anew[j] != j) {
43           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
44           for (k1=0; k1<16; k1++) {
45             dk[k1]      = aa[k*16+k1];
46             aa[k*16+k1] = aa[j*16+k1];
47             aa[j*16+k1] = dk[k1];
48           }
49         }
50         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
51         if (i > aj[j]) {
52           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
53           ap = aa + j*16;                     /* ptr to the beginning of j-th block of aa */
54           for (k=0; k<16; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
55           for (k=0; k<4; k++) {               /* j-th block of aa <- dk^T */
56             for (k1=0; k1<4; k1++) *ap++ = dk[k + 4*k1];
57           }
58         }
59       }
60     }
61     PetscCall(PetscFree(a2anew));
62   }
63 
64   /* for each row k */
65   for (k = 0; k<mbs; k++) {
66 
67     /*initialize k-th row with elements nonzero in row perm(k) of A */
68     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
69     if (jmin < jmax) {
70       ap = aa + jmin*16;
71       for (j = jmin; j < jmax; j++) {
72         vj       = perm_ptr[aj[j]];   /* block col. index */
73         rtmp_ptr = rtmp + vj*16;
74         for (i=0; i<16; i++) *rtmp_ptr++ = *ap++;
75       }
76     }
77 
78     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
79     PetscCall(PetscArraycpy(dk,rtmp+k*16,16));
80     i    = jl[k]; /* first row to be added to k_th row  */
81 
82     while (i < mbs) {
83       nexti = jl[i]; /* next row to be added to k_th row */
84 
85       /* compute multiplier */
86       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
87 
88       /* uik = -inv(Di)*U_bar(i,k) */
89       diag = ba + i*16;
90       u    = ba + ili*16;
91 
92       uik[0] = -(diag[0]*u[0] + diag[4]*u[1] + diag[8]*u[2] + diag[12]*u[3]);
93       uik[1] = -(diag[1]*u[0] + diag[5]*u[1] + diag[9]*u[2] + diag[13]*u[3]);
94       uik[2] = -(diag[2]*u[0] + diag[6]*u[1] + diag[10]*u[2]+ diag[14]*u[3]);
95       uik[3] = -(diag[3]*u[0] + diag[7]*u[1] + diag[11]*u[2]+ diag[15]*u[3]);
96 
97       uik[4] = -(diag[0]*u[4] + diag[4]*u[5] + diag[8]*u[6] + diag[12]*u[7]);
98       uik[5] = -(diag[1]*u[4] + diag[5]*u[5] + diag[9]*u[6] + diag[13]*u[7]);
99       uik[6] = -(diag[2]*u[4] + diag[6]*u[5] + diag[10]*u[6]+ diag[14]*u[7]);
100       uik[7] = -(diag[3]*u[4] + diag[7]*u[5] + diag[11]*u[6]+ diag[15]*u[7]);
101 
102       uik[8] = -(diag[0]*u[8] + diag[4]*u[9] + diag[8]*u[10] + diag[12]*u[11]);
103       uik[9] = -(diag[1]*u[8] + diag[5]*u[9] + diag[9]*u[10] + diag[13]*u[11]);
104       uik[10]= -(diag[2]*u[8] + diag[6]*u[9] + diag[10]*u[10]+ diag[14]*u[11]);
105       uik[11]= -(diag[3]*u[8] + diag[7]*u[9] + diag[11]*u[10]+ diag[15]*u[11]);
106 
107       uik[12]= -(diag[0]*u[12] + diag[4]*u[13] + diag[8]*u[14] + diag[12]*u[15]);
108       uik[13]= -(diag[1]*u[12] + diag[5]*u[13] + diag[9]*u[14] + diag[13]*u[15]);
109       uik[14]= -(diag[2]*u[12] + diag[6]*u[13] + diag[10]*u[14]+ diag[14]*u[15]);
110       uik[15]= -(diag[3]*u[12] + diag[7]*u[13] + diag[11]*u[14]+ diag[15]*u[15]);
111 
112       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
113       dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3];
114       dk[1] += uik[4]*u[0] + uik[5]*u[1] + uik[6]*u[2] + uik[7]*u[3];
115       dk[2] += uik[8]*u[0] + uik[9]*u[1] + uik[10]*u[2]+ uik[11]*u[3];
116       dk[3] += uik[12]*u[0]+ uik[13]*u[1]+ uik[14]*u[2]+ uik[15]*u[3];
117 
118       dk[4] += uik[0]*u[4] + uik[1]*u[5] + uik[2]*u[6] + uik[3]*u[7];
119       dk[5] += uik[4]*u[4] + uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7];
120       dk[6] += uik[8]*u[4] + uik[9]*u[5] + uik[10]*u[6]+ uik[11]*u[7];
121       dk[7] += uik[12]*u[4]+ uik[13]*u[5]+ uik[14]*u[6]+ uik[15]*u[7];
122 
123       dk[8] += uik[0]*u[8] + uik[1]*u[9] + uik[2]*u[10] + uik[3]*u[11];
124       dk[9] += uik[4]*u[8] + uik[5]*u[9] + uik[6]*u[10] + uik[7]*u[11];
125       dk[10]+= uik[8]*u[8] + uik[9]*u[9] + uik[10]*u[10]+ uik[11]*u[11];
126       dk[11]+= uik[12]*u[8]+ uik[13]*u[9]+ uik[14]*u[10]+ uik[15]*u[11];
127 
128       dk[12]+= uik[0]*u[12] + uik[1]*u[13] + uik[2]*u[14] + uik[3]*u[15];
129       dk[13]+= uik[4]*u[12] + uik[5]*u[13] + uik[6]*u[14] + uik[7]*u[15];
130       dk[14]+= uik[8]*u[12] + uik[9]*u[13] + uik[10]*u[14]+ uik[11]*u[15];
131       dk[15]+= uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14]+ uik[15]*u[15];
132 
133       PetscCall(PetscLogFlops(64.0*4.0));
134 
135       /* update -U(i,k) */
136       PetscCall(PetscArraycpy(ba+ili*16,uik,16));
137 
138       /* add multiple of row i to k-th row ... */
139       jmin = ili + 1; jmax = bi[i+1];
140       if (jmin < jmax) {
141         for (j=jmin; j<jmax; j++) {
142           /* rtmp += -U(i,k)^T * U_bar(i,j) */
143           rtmp_ptr     = rtmp + bj[j]*16;
144           u            = ba + j*16;
145           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3];
146           rtmp_ptr[1] += uik[4]*u[0] + uik[5]*u[1] + uik[6]*u[2] + uik[7]*u[3];
147           rtmp_ptr[2] += uik[8]*u[0] + uik[9]*u[1] + uik[10]*u[2]+ uik[11]*u[3];
148           rtmp_ptr[3] += uik[12]*u[0]+ uik[13]*u[1]+ uik[14]*u[2]+ uik[15]*u[3];
149 
150           rtmp_ptr[4] += uik[0]*u[4] + uik[1]*u[5] + uik[2]*u[6] + uik[3]*u[7];
151           rtmp_ptr[5] += uik[4]*u[4] + uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7];
152           rtmp_ptr[6] += uik[8]*u[4] + uik[9]*u[5] + uik[10]*u[6]+ uik[11]*u[7];
153           rtmp_ptr[7] += uik[12]*u[4]+ uik[13]*u[5]+ uik[14]*u[6]+ uik[15]*u[7];
154 
155           rtmp_ptr[8] += uik[0]*u[8] + uik[1]*u[9] + uik[2]*u[10] + uik[3]*u[11];
156           rtmp_ptr[9] += uik[4]*u[8] + uik[5]*u[9] + uik[6]*u[10] + uik[7]*u[11];
157           rtmp_ptr[10]+= uik[8]*u[8] + uik[9]*u[9] + uik[10]*u[10]+ uik[11]*u[11];
158           rtmp_ptr[11]+= uik[12]*u[8]+ uik[13]*u[9]+ uik[14]*u[10]+ uik[15]*u[11];
159 
160           rtmp_ptr[12]+= uik[0]*u[12] + uik[1]*u[13] + uik[2]*u[14] + uik[3]*u[15];
161           rtmp_ptr[13]+= uik[4]*u[12] + uik[5]*u[13] + uik[6]*u[14] + uik[7]*u[15];
162           rtmp_ptr[14]+= uik[8]*u[12] + uik[9]*u[13] + uik[10]*u[14]+ uik[11]*u[15];
163           rtmp_ptr[15]+= uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14]+ uik[15]*u[15];
164         }
165         PetscCall(PetscLogFlops(2.0*64.0*(jmax-jmin)));
166 
167         /* ... add i to row list for next nonzero entry */
168         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
169         j     = bj[jmin];
170         jl[i] = jl[j]; jl[j] = i; /* update jl */
171       }
172       i = nexti;
173     }
174 
175     /* save nonzero entries in k-th row of U ... */
176 
177     /* invert diagonal block */
178     diag = ba+k*16;
179     PetscCall(PetscArraycpy(diag,dk,16));
180 
181     if (pivotinblocks) {
182       PetscCall(PetscKernel_A_gets_inverse_A_4(diag,shift, allowzeropivot,&zeropivotdetected));
183       if (zeropivotdetected) C->factorerrortype = MAT_FACTOR_NUMERIC_ZEROPIVOT;
184     } else {
185       PetscCall(PetscKernel_A_gets_inverse_A_4_nopivot(diag));
186     }
187 
188     jmin = bi[k]; jmax = bi[k+1];
189     if (jmin < jmax) {
190       for (j=jmin; j<jmax; j++) {
191         vj       = bj[j];      /* block col. index of U */
192         u        = ba + j*16;
193         rtmp_ptr = rtmp + vj*16;
194         for (k1=0; k1<16; k1++) {
195           *u++        = *rtmp_ptr;
196           *rtmp_ptr++ = 0.0;
197         }
198       }
199 
200       /* ... add k to row list for first nonzero entry in k-th row */
201       il[k] = jmin;
202       i     = bj[jmin];
203       jl[k] = jl[i]; jl[i] = k;
204     }
205   }
206 
207   PetscCall(PetscFree(rtmp));
208   PetscCall(PetscFree2(il,jl));
209   PetscCall(PetscFree2(dk,uik));
210   if (a->permute) {
211     PetscCall(PetscFree(aa));
212   }
213 
214   PetscCall(ISRestoreIndices(perm,&perm_ptr));
215 
216   C->ops->solve          = MatSolve_SeqSBAIJ_4_inplace;
217   C->ops->solvetranspose = MatSolve_SeqSBAIJ_4_inplace;
218   C->assembled           = PETSC_TRUE;
219   C->preallocated        = PETSC_TRUE;
220 
221   PetscCall(PetscLogFlops(1.3333*64*b->mbs)); /* from inverting diagonal blocks */
222   PetscFunctionReturn(0);
223 }
224