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