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