xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact6.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
11c2a3de1SBarry Smith 
23a7fca6bSBarry Smith #include "src/mat/impls/sbaij/seq/sbaij.h"
381278733SSatish Balay #include "src/inline/ilu.h"
481278733SSatish Balay 
581278733SSatish Balay /* Version for when blocks are 4 by 4  */
64a2ae208SSatish Balay #undef __FUNCT__
74a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_4"
8dfbe8321SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_4(Mat A,Mat *B)
981278733SSatish Balay {
1081278733SSatish Balay   Mat                C = *B;
1181278733SSatish Balay   Mat_SeqSBAIJ       *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
1281278733SSatish Balay   IS                 perm = b->row;
13*6849ba73SBarry Smith   PetscErrorCode ierr;
14*6849ba73SBarry Smith   int                *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
1581278733SSatish Balay   int                *ai,*aj,*a2anew,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili;
1681278733SSatish Balay   MatScalar          *ba = b->a,*aa,*ap,*dk,*uik;
1781278733SSatish Balay   MatScalar          *u,*diag,*rtmp,*rtmp_ptr;
18bcd9e38bSBarry Smith   PetscTruth         pivotinblocks = b->pivotinblocks;
1981278733SSatish Balay 
2081278733SSatish Balay   PetscFunctionBegin;
2181278733SSatish Balay   /* initialization */
2281278733SSatish Balay   ierr = PetscMalloc(16*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
2381278733SSatish Balay   ierr = PetscMemzero(rtmp,16*mbs*sizeof(MatScalar));CHKERRQ(ierr);
2481278733SSatish Balay   ierr = PetscMalloc(2*mbs*sizeof(int),&il);CHKERRQ(ierr);
2581278733SSatish Balay   jl   = il + mbs;
2681278733SSatish Balay   for (i=0; i<mbs; i++) {
2781278733SSatish Balay     jl[i] = mbs; il[0] = 0;
2881278733SSatish Balay   }
2981278733SSatish Balay   ierr = PetscMalloc(32*sizeof(MatScalar),&dk);CHKERRQ(ierr);
3081278733SSatish Balay   uik  = dk + 16;
3181278733SSatish Balay   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
3281278733SSatish Balay 
3381278733SSatish Balay   /* check permutation */
3481278733SSatish Balay   if (!a->permute){
3581278733SSatish Balay     ai = a->i; aj = a->j; aa = a->a;
3681278733SSatish Balay   } else {
3781278733SSatish Balay     ai   = a->inew; aj = a->jnew;
3881278733SSatish Balay     ierr = PetscMalloc(16*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
3981278733SSatish Balay     ierr = PetscMemcpy(aa,a->a,16*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
4081278733SSatish Balay     ierr = PetscMalloc(ai[mbs]*sizeof(int),&a2anew);CHKERRQ(ierr);
4181278733SSatish Balay     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(int));CHKERRQ(ierr);
4281278733SSatish Balay 
4381278733SSatish Balay     for (i=0; i<mbs; i++){
4481278733SSatish Balay       jmin = ai[i]; jmax = ai[i+1];
4581278733SSatish Balay       for (j=jmin; j<jmax; j++){
4681278733SSatish Balay         while (a2anew[j] != j){
4781278733SSatish Balay           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
4881278733SSatish Balay           for (k1=0; k1<16; k1++){
4981278733SSatish Balay             dk[k1]       = aa[k*16+k1];
5081278733SSatish Balay             aa[k*16+k1] = aa[j*16+k1];
5181278733SSatish Balay             aa[j*16+k1] = dk[k1];
5281278733SSatish Balay           }
5381278733SSatish Balay         }
5481278733SSatish Balay         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
5581278733SSatish Balay         if (i > aj[j]){
5681278733SSatish Balay           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
5781278733SSatish Balay           ap = aa + j*16;                     /* ptr to the beginning of j-th block of aa */
5881278733SSatish Balay           for (k=0; k<16; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
5981278733SSatish Balay           for (k=0; k<4; k++){               /* j-th block of aa <- dk^T */
6081278733SSatish Balay             for (k1=0; k1<4; k1++) *ap++ = dk[k + 4*k1];
6181278733SSatish Balay           }
6281278733SSatish Balay         }
6381278733SSatish Balay       }
6481278733SSatish Balay     }
6581278733SSatish Balay     ierr = PetscFree(a2anew);CHKERRQ(ierr);
6681278733SSatish Balay   }
6781278733SSatish Balay 
6881278733SSatish Balay   /* for each row k */
6981278733SSatish Balay   for (k = 0; k<mbs; k++){
7081278733SSatish Balay 
7181278733SSatish Balay     /*initialize k-th row with elements nonzero in row perm(k) of A */
7281278733SSatish Balay     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
7381278733SSatish Balay     if (jmin < jmax) {
7481278733SSatish Balay       ap = aa + jmin*16;
7581278733SSatish Balay       for (j = jmin; j < jmax; j++){
7681278733SSatish Balay         vj = perm_ptr[aj[j]];         /* block col. index */
7781278733SSatish Balay         rtmp_ptr = rtmp + vj*16;
7881278733SSatish Balay         for (i=0; i<16; i++) *rtmp_ptr++ = *ap++;
7981278733SSatish Balay       }
8081278733SSatish Balay     }
8181278733SSatish Balay 
8281278733SSatish Balay     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
8381278733SSatish Balay     ierr = PetscMemcpy(dk,rtmp+k*16,16*sizeof(MatScalar));CHKERRQ(ierr);
8481278733SSatish Balay     i = jl[k]; /* first row to be added to k_th row  */
8581278733SSatish Balay 
8681278733SSatish Balay     while (i < mbs){
8781278733SSatish Balay       nexti = jl[i]; /* next row to be added to k_th row */
8881278733SSatish Balay 
8981278733SSatish Balay       /* compute multiplier */
9081278733SSatish Balay       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
9181278733SSatish Balay 
9281278733SSatish Balay       /* uik = -inv(Di)*U_bar(i,k) */
9381278733SSatish Balay       diag = ba + i*16;
9481278733SSatish Balay       u    = ba + ili*16;
9581278733SSatish Balay 
9681278733SSatish Balay       uik[0] = -(diag[0]*u[0] + diag[4]*u[1] + diag[8]*u[2] + diag[12]*u[3]);
9781278733SSatish Balay       uik[1] = -(diag[1]*u[0] + diag[5]*u[1] + diag[9]*u[2] + diag[13]*u[3]);
9881278733SSatish Balay       uik[2] = -(diag[2]*u[0] + diag[6]*u[1] + diag[10]*u[2]+ diag[14]*u[3]);
9981278733SSatish Balay       uik[3] = -(diag[3]*u[0] + diag[7]*u[1] + diag[11]*u[2]+ diag[15]*u[3]);
10081278733SSatish Balay 
10181278733SSatish Balay       uik[4] = -(diag[0]*u[4] + diag[4]*u[5] + diag[8]*u[6] + diag[12]*u[7]);
10281278733SSatish Balay       uik[5] = -(diag[1]*u[4] + diag[5]*u[5] + diag[9]*u[6] + diag[13]*u[7]);
10381278733SSatish Balay       uik[6] = -(diag[2]*u[4] + diag[6]*u[5] + diag[10]*u[6]+ diag[14]*u[7]);
10481278733SSatish Balay       uik[7] = -(diag[3]*u[4] + diag[7]*u[5] + diag[11]*u[6]+ diag[15]*u[7]);
10581278733SSatish Balay 
10681278733SSatish Balay       uik[8] = -(diag[0]*u[8] + diag[4]*u[9] + diag[8]*u[10] + diag[12]*u[11]);
10781278733SSatish Balay       uik[9] = -(diag[1]*u[8] + diag[5]*u[9] + diag[9]*u[10] + diag[13]*u[11]);
10881278733SSatish Balay       uik[10]= -(diag[2]*u[8] + diag[6]*u[9] + diag[10]*u[10]+ diag[14]*u[11]);
10981278733SSatish Balay       uik[11]= -(diag[3]*u[8] + diag[7]*u[9] + diag[11]*u[10]+ diag[15]*u[11]);
11081278733SSatish Balay 
11181278733SSatish Balay       uik[12]= -(diag[0]*u[12] + diag[4]*u[13] + diag[8]*u[14] + diag[12]*u[15]);
11281278733SSatish Balay       uik[13]= -(diag[1]*u[12] + diag[5]*u[13] + diag[9]*u[14] + diag[13]*u[15]);
11381278733SSatish Balay       uik[14]= -(diag[2]*u[12] + diag[6]*u[13] + diag[10]*u[14]+ diag[14]*u[15]);
11481278733SSatish Balay       uik[15]= -(diag[3]*u[12] + diag[7]*u[13] + diag[11]*u[14]+ diag[15]*u[15]);
11581278733SSatish Balay 
11681278733SSatish Balay       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
11781278733SSatish Balay       dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3];
11881278733SSatish Balay       dk[1] += uik[4]*u[0] + uik[5]*u[1] + uik[6]*u[2] + uik[7]*u[3];
11981278733SSatish Balay       dk[2] += uik[8]*u[0] + uik[9]*u[1] + uik[10]*u[2]+ uik[11]*u[3];
12081278733SSatish Balay       dk[3] += uik[12]*u[0]+ uik[13]*u[1]+ uik[14]*u[2]+ uik[15]*u[3];
12181278733SSatish Balay 
12281278733SSatish Balay       dk[4] += uik[0]*u[4] + uik[1]*u[5] + uik[2]*u[6] + uik[3]*u[7];
12381278733SSatish Balay       dk[5] += uik[4]*u[4] + uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7];
12481278733SSatish Balay       dk[6] += uik[8]*u[4] + uik[9]*u[5] + uik[10]*u[6]+ uik[11]*u[7];
12581278733SSatish Balay       dk[7] += uik[12]*u[4]+ uik[13]*u[5]+ uik[14]*u[6]+ uik[15]*u[7];
12681278733SSatish Balay 
12781278733SSatish Balay       dk[8] += uik[0]*u[8] + uik[1]*u[9] + uik[2]*u[10] + uik[3]*u[11];
12881278733SSatish Balay       dk[9] += uik[4]*u[8] + uik[5]*u[9] + uik[6]*u[10] + uik[7]*u[11];
12981278733SSatish Balay       dk[10]+= uik[8]*u[8] + uik[9]*u[9] + uik[10]*u[10]+ uik[11]*u[11];
13081278733SSatish Balay       dk[11]+= uik[12]*u[8]+ uik[13]*u[9]+ uik[14]*u[10]+ uik[15]*u[11];
13181278733SSatish Balay 
13281278733SSatish Balay       dk[12]+= uik[0]*u[12] + uik[1]*u[13] + uik[2]*u[14] + uik[3]*u[15];
13381278733SSatish Balay       dk[13]+= uik[4]*u[12] + uik[5]*u[13] + uik[6]*u[14] + uik[7]*u[15];
13481278733SSatish Balay       dk[14]+= uik[8]*u[12] + uik[9]*u[13] + uik[10]*u[14]+ uik[11]*u[15];
13581278733SSatish Balay       dk[15]+= uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14]+ uik[15]*u[15];
13681278733SSatish Balay 
13781278733SSatish Balay       /* update -U(i,k) */
13881278733SSatish Balay       ierr = PetscMemcpy(ba+ili*16,uik,16*sizeof(MatScalar));CHKERRQ(ierr);
13981278733SSatish Balay 
14081278733SSatish Balay       /* add multiple of row i to k-th row ... */
14181278733SSatish Balay       jmin = ili + 1; jmax = bi[i+1];
14281278733SSatish Balay       if (jmin < jmax){
14381278733SSatish Balay         for (j=jmin; j<jmax; j++) {
14481278733SSatish Balay           /* rtmp += -U(i,k)^T * U_bar(i,j) */
14581278733SSatish Balay           rtmp_ptr = rtmp + bj[j]*16;
14681278733SSatish Balay           u = ba + j*16;
14781278733SSatish Balay           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2] + uik[3]*u[3];
14881278733SSatish Balay           rtmp_ptr[1] += uik[4]*u[0] + uik[5]*u[1] + uik[6]*u[2] + uik[7]*u[3];
14981278733SSatish Balay           rtmp_ptr[2] += uik[8]*u[0] + uik[9]*u[1] + uik[10]*u[2]+ uik[11]*u[3];
15081278733SSatish Balay           rtmp_ptr[3] += uik[12]*u[0]+ uik[13]*u[1]+ uik[14]*u[2]+ uik[15]*u[3];
15181278733SSatish Balay 
15281278733SSatish Balay           rtmp_ptr[4] += uik[0]*u[4] + uik[1]*u[5] + uik[2]*u[6] + uik[3]*u[7];
15381278733SSatish Balay           rtmp_ptr[5] += uik[4]*u[4] + uik[5]*u[5] + uik[6]*u[6] + uik[7]*u[7];
15481278733SSatish Balay           rtmp_ptr[6] += uik[8]*u[4] + uik[9]*u[5] + uik[10]*u[6]+ uik[11]*u[7];
15581278733SSatish Balay           rtmp_ptr[7] += uik[12]*u[4]+ uik[13]*u[5]+ uik[14]*u[6]+ uik[15]*u[7];
15681278733SSatish Balay 
15781278733SSatish Balay           rtmp_ptr[8] += uik[0]*u[8] + uik[1]*u[9] + uik[2]*u[10] + uik[3]*u[11];
15881278733SSatish Balay           rtmp_ptr[9] += uik[4]*u[8] + uik[5]*u[9] + uik[6]*u[10] + uik[7]*u[11];
15981278733SSatish Balay           rtmp_ptr[10]+= uik[8]*u[8] + uik[9]*u[9] + uik[10]*u[10]+ uik[11]*u[11];
16081278733SSatish Balay           rtmp_ptr[11]+= uik[12]*u[8]+ uik[13]*u[9]+ uik[14]*u[10]+ uik[15]*u[11];
16181278733SSatish Balay 
16281278733SSatish Balay           rtmp_ptr[12]+= uik[0]*u[12] + uik[1]*u[13] + uik[2]*u[14] + uik[3]*u[15];
16381278733SSatish Balay           rtmp_ptr[13]+= uik[4]*u[12] + uik[5]*u[13] + uik[6]*u[14] + uik[7]*u[15];
16481278733SSatish Balay           rtmp_ptr[14]+= uik[8]*u[12] + uik[9]*u[13] + uik[10]*u[14]+ uik[11]*u[15];
16581278733SSatish Balay           rtmp_ptr[15]+= uik[12]*u[12]+ uik[13]*u[13]+ uik[14]*u[14]+ uik[15]*u[15];
16681278733SSatish Balay         }
16781278733SSatish Balay 
16881278733SSatish Balay         /* ... add i to row list for next nonzero entry */
16981278733SSatish Balay         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
17081278733SSatish Balay         j     = bj[jmin];
17181278733SSatish Balay         jl[i] = jl[j]; jl[j] = i; /* update jl */
17281278733SSatish Balay       }
17381278733SSatish Balay       i = nexti;
17481278733SSatish Balay     }
17581278733SSatish Balay 
17681278733SSatish Balay     /* save nonzero entries in k-th row of U ... */
17781278733SSatish Balay 
17881278733SSatish Balay     /* invert diagonal block */
17981278733SSatish Balay     diag = ba+k*16;
18081278733SSatish Balay     ierr = PetscMemcpy(diag,dk,16*sizeof(MatScalar));CHKERRQ(ierr);
181bcd9e38bSBarry Smith 
182bcd9e38bSBarry Smith     if (pivotinblocks) {
18381278733SSatish Balay       ierr = Kernel_A_gets_inverse_A_4(diag);CHKERRQ(ierr);
184bcd9e38bSBarry Smith     } else {
185bcd9e38bSBarry Smith       ierr = Kernel_A_gets_inverse_A_4_nopivot(diag);CHKERRQ(ierr);
186bcd9e38bSBarry Smith     }
18781278733SSatish Balay 
18881278733SSatish Balay     jmin = bi[k]; jmax = bi[k+1];
18981278733SSatish Balay     if (jmin < jmax) {
19081278733SSatish Balay       for (j=jmin; j<jmax; j++){
19181278733SSatish Balay          vj = bj[j];           /* block col. index of U */
19281278733SSatish Balay          u   = ba + j*16;
19381278733SSatish Balay          rtmp_ptr = rtmp + vj*16;
19481278733SSatish Balay          for (k1=0; k1<16; k1++){
19581278733SSatish Balay            *u++        = *rtmp_ptr;
19681278733SSatish Balay            *rtmp_ptr++ = 0.0;
19781278733SSatish Balay          }
19881278733SSatish Balay       }
19981278733SSatish Balay 
20081278733SSatish Balay       /* ... add k to row list for first nonzero entry in k-th row */
20181278733SSatish Balay       il[k] = jmin;
20281278733SSatish Balay       i     = bj[jmin];
20381278733SSatish Balay       jl[k] = jl[i]; jl[i] = k;
20481278733SSatish Balay     }
20581278733SSatish Balay   }
20681278733SSatish Balay 
20781278733SSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
20881278733SSatish Balay   ierr = PetscFree(il);CHKERRQ(ierr);
20981278733SSatish Balay   ierr = PetscFree(dk);CHKERRQ(ierr);
21081278733SSatish Balay   if (a->permute){
21181278733SSatish Balay     ierr = PetscFree(aa);CHKERRQ(ierr);
21281278733SSatish Balay   }
21381278733SSatish Balay 
21481278733SSatish Balay   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
21581278733SSatish Balay   C->factor    = FACTOR_CHOLESKY;
21681278733SSatish Balay   C->assembled = PETSC_TRUE;
21781278733SSatish Balay   C->preallocated = PETSC_TRUE;
21881278733SSatish Balay   PetscLogFlops(1.3333*64*b->mbs); /* from inverting diagonal blocks */
21981278733SSatish Balay   PetscFunctionReturn(0);
22081278733SSatish Balay }
221