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