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