xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact6.c (revision bebe2cf65d55febe21a5af8db2bd2e168caaa2e7)
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 #undef __FUNCT__
7 #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_4"
8 PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(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;
12   PetscErrorCode ierr;
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;
15   MatScalar      *ba = b->a,*aa,*ap,*dk,*uik;
16   MatScalar      *u,*diag,*rtmp,*rtmp_ptr;
17   PetscBool      pivotinblocks = b->pivotinblocks;
18   PetscReal      shift         = info->shiftamount;
19 
20   PetscFunctionBegin;
21   /* initialization */
22   ierr = PetscCalloc1(16*mbs,&rtmp);CHKERRQ(ierr);
23   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
24   for (i=0; i<mbs; i++) {
25     jl[i] = mbs; il[0] = 0;
26   }
27   ierr = PetscMalloc2(16,&dk,16,&uik);CHKERRQ(ierr);
28   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
29 
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     ierr = PetscMalloc1(16*ai[mbs],&aa);CHKERRQ(ierr);
36     ierr = PetscMemcpy(aa,a->a,16*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
37     ierr = PetscMalloc1(ai[mbs],&a2anew);CHKERRQ(ierr);
38     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
39 
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     ierr = PetscFree(a2anew);CHKERRQ(ierr);
63   }
64 
65   /* for each row k */
66   for (k = 0; k<mbs; k++) {
67 
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     }
78 
79     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
80     ierr = PetscMemcpy(dk,rtmp+k*16,16*sizeof(MatScalar));CHKERRQ(ierr);
81     i    = jl[k]; /* first row to be added to k_th row  */
82 
83     while (i < mbs) {
84       nexti = jl[i]; /* next row to be added to k_th row */
85 
86       /* compute multiplier */
87       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
88 
89       /* uik = -inv(Di)*U_bar(i,k) */
90       diag = ba + i*16;
91       u    = ba + ili*16;
92 
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]);
97 
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]);
102 
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]);
107 
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]);
112 
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];
118 
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];
123 
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];
128 
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];
133 
134       ierr = PetscLogFlops(64.0*4.0);CHKERRQ(ierr);
135 
136       /* update -U(i,k) */
137       ierr = PetscMemcpy(ba+ili*16,uik,16*sizeof(MatScalar));CHKERRQ(ierr);
138 
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];
150 
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];
155 
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];
160 
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         ierr = PetscLogFlops(2.0*64.0*(jmax-jmin));CHKERRQ(ierr);
167 
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     }
175 
176     /* save nonzero entries in k-th row of U ... */
177 
178     /* invert diagonal block */
179     diag = ba+k*16;
180     ierr = PetscMemcpy(diag,dk,16*sizeof(MatScalar));CHKERRQ(ierr);
181 
182     if (pivotinblocks) {
183       ierr = PetscKernel_A_gets_inverse_A_4(diag,shift);CHKERRQ(ierr);
184     } else {
185       ierr = PetscKernel_A_gets_inverse_A_4_nopivot(diag);CHKERRQ(ierr);
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   ierr = PetscFree(rtmp);CHKERRQ(ierr);
208   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
209   ierr = PetscFree2(dk,uik);CHKERRQ(ierr);
210   if (a->permute) {
211     ierr = PetscFree(aa);CHKERRQ(ierr);
212   }
213 
214   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
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   ierr = PetscLogFlops(1.3333*64*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
222   PetscFunctionReturn(0);
223 }
224