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