xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact3.c (revision 719d5645761d844e4357b7ee00a3296dffe0b787)
1be1d678aSKris Buschelman #define PETSCMAT_DLL
21c2a3de1SBarry Smith 
33a7fca6bSBarry Smith #include "src/mat/impls/sbaij/seq/sbaij.h"
481278733SSatish Balay #include "src/inline/ilu.h"
581278733SSatish Balay 
681278733SSatish Balay /* Version for when blocks are 3 by 3  */
74a2ae208SSatish Balay #undef __FUNCT__
84a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_3"
9*719d5645SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3(Mat C,Mat A,MatFactorInfo *info)
1081278733SSatish Balay {
1181278733SSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data;
1281278733SSatish Balay   IS             perm = b->row;
136849ba73SBarry Smith   PetscErrorCode ierr;
1413f74950SBarry Smith   PetscInt       *perm_ptr,i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
1513f74950SBarry Smith   PetscInt       *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;
1862bba022SBarry Smith   PetscReal      shift = info->shiftinblocks;
1981278733SSatish Balay 
2081278733SSatish Balay   PetscFunctionBegin;
2181278733SSatish Balay 
2281278733SSatish Balay   /* initialization */
2381278733SSatish Balay   ierr = PetscMalloc(9*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr);
2481278733SSatish Balay   ierr = PetscMemzero(rtmp,9*mbs*sizeof(MatScalar));CHKERRQ(ierr);
2513f74950SBarry Smith   ierr = PetscMalloc(2*mbs*sizeof(PetscInt),&il);CHKERRQ(ierr);
2681278733SSatish Balay   jl   = il + mbs;
2781278733SSatish Balay   for (i=0; i<mbs; i++) {
2881278733SSatish Balay     jl[i] = mbs; il[0] = 0;
2981278733SSatish Balay   }
3081278733SSatish Balay   ierr = PetscMalloc(18*sizeof(MatScalar),&dk);CHKERRQ(ierr);
3181278733SSatish Balay   uik  = dk + 9;
3281278733SSatish Balay   ierr = ISGetIndices(perm,&perm_ptr);CHKERRQ(ierr);
3381278733SSatish Balay 
3481278733SSatish Balay   /* check permutation */
3581278733SSatish Balay   if (!a->permute){
3681278733SSatish Balay     ai = a->i; aj = a->j; aa = a->a;
3781278733SSatish Balay   } else {
3881278733SSatish Balay     ai   = a->inew; aj = a->jnew;
3981278733SSatish Balay     ierr = PetscMalloc(9*ai[mbs]*sizeof(MatScalar),&aa);CHKERRQ(ierr);
4081278733SSatish Balay     ierr = PetscMemcpy(aa,a->a,9*ai[mbs]*sizeof(MatScalar));CHKERRQ(ierr);
4113f74950SBarry Smith     ierr = PetscMalloc(ai[mbs]*sizeof(PetscInt),&a2anew);CHKERRQ(ierr);
4213f74950SBarry Smith     ierr = PetscMemcpy(a2anew,a->a2anew,(ai[mbs])*sizeof(PetscInt));CHKERRQ(ierr);
4381278733SSatish Balay 
4481278733SSatish Balay     for (i=0; i<mbs; i++){
4581278733SSatish Balay       jmin = ai[i]; jmax = ai[i+1];
4681278733SSatish Balay       for (j=jmin; j<jmax; j++){
4781278733SSatish Balay         while (a2anew[j] != j){
4881278733SSatish Balay           k = a2anew[j]; a2anew[j] = a2anew[k]; a2anew[k] = k;
4981278733SSatish Balay           for (k1=0; k1<9; k1++){
5081278733SSatish Balay             dk[k1]       = aa[k*9+k1];
5181278733SSatish Balay             aa[k*9+k1] = aa[j*9+k1];
5281278733SSatish Balay             aa[j*9+k1] = dk[k1];
5381278733SSatish Balay           }
5481278733SSatish Balay         }
5581278733SSatish Balay         /* transform columnoriented blocks that lie in the lower triangle to roworiented blocks */
5681278733SSatish Balay         if (i > aj[j]){
5781278733SSatish Balay           /* printf("change orientation, row: %d, col: %d\n",i,aj[j]); */
5881278733SSatish Balay           ap = aa + j*9;                     /* ptr to the beginning of j-th block of aa */
5981278733SSatish Balay           for (k=0; k<9; k++) dk[k] = ap[k]; /* dk <- j-th block of aa */
6081278733SSatish Balay           for (k=0; k<3; k++){               /* j-th block of aa <- dk^T */
6181278733SSatish Balay             for (k1=0; k1<3; k1++) *ap++ = dk[k + 3*k1];
6281278733SSatish Balay           }
6381278733SSatish Balay         }
6481278733SSatish Balay       }
6581278733SSatish Balay     }
6681278733SSatish Balay     ierr = PetscFree(a2anew);CHKERRQ(ierr);
6781278733SSatish Balay   }
6881278733SSatish Balay 
6981278733SSatish Balay   /* for each row k */
7081278733SSatish Balay   for (k = 0; k<mbs; k++){
7181278733SSatish Balay 
7281278733SSatish Balay     /*initialize k-th row with elements nonzero in row perm(k) of A */
7381278733SSatish Balay     jmin = ai[perm_ptr[k]]; jmax = ai[perm_ptr[k]+1];
7481278733SSatish Balay     if (jmin < jmax) {
7581278733SSatish Balay       ap = aa + jmin*9;
7681278733SSatish Balay       for (j = jmin; j < jmax; j++){
7781278733SSatish Balay         vj = perm_ptr[aj[j]];         /* block col. index */
7881278733SSatish Balay         rtmp_ptr = rtmp + vj*9;
7981278733SSatish Balay         for (i=0; i<9; i++) *rtmp_ptr++ = *ap++;
8081278733SSatish Balay       }
8181278733SSatish Balay     }
8281278733SSatish Balay 
8381278733SSatish Balay     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
8481278733SSatish Balay     ierr = PetscMemcpy(dk,rtmp+k*9,9*sizeof(MatScalar));CHKERRQ(ierr);
8581278733SSatish Balay     i = jl[k]; /* first row to be added to k_th row  */
8681278733SSatish Balay 
8781278733SSatish Balay     while (i < mbs){
8881278733SSatish Balay       nexti = jl[i]; /* next row to be added to k_th row */
8981278733SSatish Balay 
9081278733SSatish Balay       /* compute multiplier */
9181278733SSatish Balay       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
9281278733SSatish Balay 
9381278733SSatish Balay       /* uik = -inv(Di)*U_bar(i,k) */
9481278733SSatish Balay       diag = ba + i*9;
9581278733SSatish Balay       u    = ba + ili*9;
9681278733SSatish Balay 
9781278733SSatish Balay       uik[0] = -(diag[0]*u[0] + diag[3]*u[1] + diag[6]*u[2]);
9881278733SSatish Balay       uik[1] = -(diag[1]*u[0] + diag[4]*u[1] + diag[7]*u[2]);
9981278733SSatish Balay       uik[2] = -(diag[2]*u[0] + diag[5]*u[1] + diag[8]*u[2]);
10081278733SSatish Balay 
10181278733SSatish Balay       uik[3] = -(diag[0]*u[3] + diag[3]*u[4] + diag[6]*u[5]);
10281278733SSatish Balay       uik[4] = -(diag[1]*u[3] + diag[4]*u[4] + diag[7]*u[5]);
10381278733SSatish Balay       uik[5] = -(diag[2]*u[3] + diag[5]*u[4] + diag[8]*u[5]);
10481278733SSatish Balay 
10581278733SSatish Balay       uik[6] = -(diag[0]*u[6] + diag[3]*u[7] + diag[6]*u[8]);
10681278733SSatish Balay       uik[7] = -(diag[1]*u[6] + diag[4]*u[7] + diag[7]*u[8]);
10781278733SSatish Balay       uik[8] = -(diag[2]*u[6] + diag[5]*u[7] + diag[8]*u[8]);
10881278733SSatish Balay 
10981278733SSatish Balay       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
11081278733SSatish Balay       dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
11181278733SSatish Balay       dk[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
11281278733SSatish Balay       dk[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];
11381278733SSatish Balay 
11481278733SSatish Balay       dk[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
11581278733SSatish Balay       dk[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
11681278733SSatish Balay       dk[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];
11781278733SSatish Balay 
11881278733SSatish Balay       dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
11981278733SSatish Balay       dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
12081278733SSatish Balay       dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
12181278733SSatish Balay 
122187a9f4bSHong Zhang       ierr = PetscLogFlops(27*4);CHKERRQ(ierr);
123187a9f4bSHong Zhang 
12481278733SSatish Balay       /* update -U(i,k) */
12581278733SSatish Balay       ierr = PetscMemcpy(ba+ili*9,uik,9*sizeof(MatScalar));CHKERRQ(ierr);
12681278733SSatish Balay 
12781278733SSatish Balay       /* add multiple of row i to k-th row ... */
12881278733SSatish Balay       jmin = ili + 1; jmax = bi[i+1];
12981278733SSatish Balay       if (jmin < jmax){
13081278733SSatish Balay         for (j=jmin; j<jmax; j++) {
13181278733SSatish Balay           /* rtmp += -U(i,k)^T * U_bar(i,j) */
13281278733SSatish Balay           rtmp_ptr = rtmp + bj[j]*9;
13381278733SSatish Balay           u = ba + j*9;
13481278733SSatish Balay           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
13581278733SSatish Balay           rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
13681278733SSatish Balay           rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];
13781278733SSatish Balay 
13881278733SSatish Balay           rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
13981278733SSatish Balay           rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
14081278733SSatish Balay           rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];
14181278733SSatish Balay 
14281278733SSatish Balay           rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
14381278733SSatish Balay           rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
14481278733SSatish Balay           rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
14581278733SSatish Balay         }
146187a9f4bSHong Zhang         ierr = PetscLogFlops(2*27*(jmax-jmin));CHKERRQ(ierr);
14781278733SSatish Balay 
14881278733SSatish Balay         /* ... add i to row list for next nonzero entry */
14981278733SSatish Balay         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
15081278733SSatish Balay         j     = bj[jmin];
15181278733SSatish Balay         jl[i] = jl[j]; jl[j] = i; /* update jl */
15281278733SSatish Balay       }
15381278733SSatish Balay       i = nexti;
15481278733SSatish Balay     }
15581278733SSatish Balay 
15681278733SSatish Balay     /* save nonzero entries in k-th row of U ... */
15781278733SSatish Balay 
15881278733SSatish Balay     /* invert diagonal block */
15981278733SSatish Balay     diag = ba+k*9;
16081278733SSatish Balay     ierr = PetscMemcpy(diag,dk,9*sizeof(MatScalar));CHKERRQ(ierr);
16162bba022SBarry Smith     ierr = Kernel_A_gets_inverse_A_3(diag,shift);CHKERRQ(ierr);
16281278733SSatish Balay 
16381278733SSatish Balay     jmin = bi[k]; jmax = bi[k+1];
16481278733SSatish Balay     if (jmin < jmax) {
16581278733SSatish Balay       for (j=jmin; j<jmax; j++){
16681278733SSatish Balay          vj = bj[j];           /* block col. index of U */
16781278733SSatish Balay          u   = ba + j*9;
16881278733SSatish Balay          rtmp_ptr = rtmp + vj*9;
16981278733SSatish Balay          for (k1=0; k1<9; k1++){
17081278733SSatish Balay            *u++        = *rtmp_ptr;
17181278733SSatish Balay            *rtmp_ptr++ = 0.0;
17281278733SSatish Balay          }
17381278733SSatish Balay       }
17481278733SSatish Balay 
17581278733SSatish Balay       /* ... add k to row list for first nonzero entry in k-th row */
17681278733SSatish Balay       il[k] = jmin;
17781278733SSatish Balay       i     = bj[jmin];
17881278733SSatish Balay       jl[k] = jl[i]; jl[i] = k;
17981278733SSatish Balay     }
18081278733SSatish Balay   }
18181278733SSatish Balay 
18281278733SSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
18381278733SSatish Balay   ierr = PetscFree(il);CHKERRQ(ierr);
18481278733SSatish Balay   ierr = PetscFree(dk);CHKERRQ(ierr);
18581278733SSatish Balay   if (a->permute){
18681278733SSatish Balay     ierr = PetscFree(aa);CHKERRQ(ierr);
18781278733SSatish Balay   }
18881278733SSatish Balay 
18981278733SSatish Balay   ierr = ISRestoreIndices(perm,&perm_ptr);CHKERRQ(ierr);
190db4efbfdSBarry Smith   C->ops->solve          = MatSolve_SeqSBAIJ_3;
191db4efbfdSBarry Smith   C->ops->solvetranspose = MatSolve_SeqSBAIJ_3;
19281278733SSatish Balay   C->assembled = PETSC_TRUE;
19381278733SSatish Balay   C->preallocated = PETSC_TRUE;
194efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*27*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
19581278733SSatish Balay   PetscFunctionReturn(0);
19681278733SSatish Balay }
197