xref: /petsc/src/mat/impls/sbaij/seq/sbaijfact4.c (revision 1795a4d16c893ec2fc06bbbc6c5ce592a2de75d4)
11c2a3de1SBarry Smith 
2c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h>
306873bf2SBarry Smith #include <petsc-private/kernels/blockinvert.h>
481278733SSatish Balay 
581278733SSatish Balay /*
681278733SSatish Balay       Version for when blocks are 3 by 3 Using natural ordering
781278733SSatish Balay */
84a2ae208SSatish Balay #undef __FUNCT__
94a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering"
100481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info)
1181278733SSatish Balay {
1281278733SSatish Balay   Mat_SeqSBAIJ   *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ*)C->data;
13dfbe8321SBarry Smith   PetscErrorCode ierr;
1413f74950SBarry Smith   PetscInt       i,j,mbs=a->mbs,*bi=b->i,*bj=b->j;
1513f74950SBarry Smith   PetscInt       *ai,*aj,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;
18182b8fbaSHong Zhang   PetscReal      shift = info->shiftamount;
1981278733SSatish Balay 
2081278733SSatish Balay   PetscFunctionBegin;
2181278733SSatish Balay   /* initialization */
22*1795a4d1SJed Brown   ierr = PetscCalloc1(9*mbs,&rtmp);CHKERRQ(ierr);
23dcca6d9dSJed Brown   ierr = PetscMalloc2(mbs,&il,mbs,&jl);CHKERRQ(ierr);
2481278733SSatish Balay   for (i=0; i<mbs; i++) {
2581278733SSatish Balay     jl[i] = mbs; il[0] = 0;
2681278733SSatish Balay   }
27dcca6d9dSJed Brown   ierr = PetscMalloc2(9,&dk,9,&uik);CHKERRQ(ierr);
2881278733SSatish Balay   ai   = a->i; aj = a->j; aa = a->a;
2981278733SSatish Balay 
3081278733SSatish Balay   /* for each row k */
3181278733SSatish Balay   for (k = 0; k<mbs; k++) {
3281278733SSatish Balay 
3381278733SSatish Balay     /*initialize k-th row with elements nonzero in row k of A */
3481278733SSatish Balay     jmin = ai[k]; jmax = ai[k+1];
3581278733SSatish Balay     if (jmin < jmax) {
3681278733SSatish Balay       ap = aa + jmin*9;
3781278733SSatish Balay       for (j = jmin; j < jmax; j++) {
3881278733SSatish Balay         vj       = aj[j];   /* block col. index */
3981278733SSatish Balay         rtmp_ptr = rtmp + vj*9;
4081278733SSatish Balay         for (i=0; i<9; i++) *rtmp_ptr++ = *ap++;
4181278733SSatish Balay       }
4281278733SSatish Balay     }
4381278733SSatish Balay 
4481278733SSatish Balay     /* modify k-th row by adding in those rows i with U(i,k) != 0 */
4581278733SSatish Balay     ierr = PetscMemcpy(dk,rtmp+k*9,9*sizeof(MatScalar));CHKERRQ(ierr);
4681278733SSatish Balay     i    = jl[k]; /* first row to be added to k_th row  */
4781278733SSatish Balay 
4881278733SSatish Balay     while (i < mbs) {
4981278733SSatish Balay       nexti = jl[i]; /* next row to be added to k_th row */
5081278733SSatish Balay 
5181278733SSatish Balay       /* compute multiplier */
5281278733SSatish Balay       ili = il[i];  /* index of first nonzero element in U(i,k:bms-1) */
5381278733SSatish Balay 
5481278733SSatish Balay       /* uik = -inv(Di)*U_bar(i,k) */
5581278733SSatish Balay       diag = ba + i*9;
5681278733SSatish Balay       u    = ba + ili*9;
5781278733SSatish Balay 
5881278733SSatish Balay       uik[0] = -(diag[0]*u[0] + diag[3]*u[1] + diag[6]*u[2]);
5981278733SSatish Balay       uik[1] = -(diag[1]*u[0] + diag[4]*u[1] + diag[7]*u[2]);
6081278733SSatish Balay       uik[2] = -(diag[2]*u[0] + diag[5]*u[1] + diag[8]*u[2]);
6181278733SSatish Balay 
6281278733SSatish Balay       uik[3] = -(diag[0]*u[3] + diag[3]*u[4] + diag[6]*u[5]);
6381278733SSatish Balay       uik[4] = -(diag[1]*u[3] + diag[4]*u[4] + diag[7]*u[5]);
6481278733SSatish Balay       uik[5] = -(diag[2]*u[3] + diag[5]*u[4] + diag[8]*u[5]);
6581278733SSatish Balay 
6681278733SSatish Balay       uik[6] = -(diag[0]*u[6] + diag[3]*u[7] + diag[6]*u[8]);
6781278733SSatish Balay       uik[7] = -(diag[1]*u[6] + diag[4]*u[7] + diag[7]*u[8]);
6881278733SSatish Balay       uik[8] = -(diag[2]*u[6] + diag[5]*u[7] + diag[8]*u[8]);
6981278733SSatish Balay 
7081278733SSatish Balay       /* update D(k) += -U(i,k)^T * U_bar(i,k) */
7181278733SSatish Balay       dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
7281278733SSatish Balay       dk[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
7381278733SSatish Balay       dk[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];
7481278733SSatish Balay 
7581278733SSatish Balay       dk[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
7681278733SSatish Balay       dk[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
7781278733SSatish Balay       dk[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];
7881278733SSatish Balay 
7981278733SSatish Balay       dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
8081278733SSatish Balay       dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
8181278733SSatish Balay       dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
8281278733SSatish Balay 
83dc0b31edSSatish Balay       ierr = PetscLogFlops(27.0*4.0);CHKERRQ(ierr);
84187a9f4bSHong Zhang 
8581278733SSatish Balay       /* update -U(i,k) */
8681278733SSatish Balay       ierr = PetscMemcpy(ba+ili*9,uik,9*sizeof(MatScalar));CHKERRQ(ierr);
8781278733SSatish Balay 
8881278733SSatish Balay       /* add multiple of row i to k-th row ... */
8981278733SSatish Balay       jmin = ili + 1; jmax = bi[i+1];
9081278733SSatish Balay       if (jmin < jmax) {
9181278733SSatish Balay         for (j=jmin; j<jmax; j++) {
9281278733SSatish Balay           /* rtmp += -U(i,k)^T * U_bar(i,j) */
9381278733SSatish Balay           rtmp_ptr     = rtmp + bj[j]*9;
9481278733SSatish Balay           u            = ba + j*9;
9581278733SSatish Balay           rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2];
9681278733SSatish Balay           rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2];
9781278733SSatish Balay           rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2];
9881278733SSatish Balay 
9981278733SSatish Balay           rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5];
10081278733SSatish Balay           rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5];
10181278733SSatish Balay           rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5];
10281278733SSatish Balay 
10381278733SSatish Balay           rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8];
10481278733SSatish Balay           rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8];
10581278733SSatish Balay           rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8];
10681278733SSatish Balay         }
107dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*27.0*(jmax-jmin));CHKERRQ(ierr);
10881278733SSatish Balay 
10981278733SSatish Balay         /* ... add i to row list for next nonzero entry */
11081278733SSatish Balay         il[i] = jmin;             /* update il(i) in column k+1, ... mbs-1 */
11181278733SSatish Balay         j     = bj[jmin];
11281278733SSatish Balay         jl[i] = jl[j]; jl[j] = i; /* update jl */
11381278733SSatish Balay       }
11481278733SSatish Balay       i = nexti;
11581278733SSatish Balay     }
11681278733SSatish Balay 
11781278733SSatish Balay     /* save nonzero entries in k-th row of U ... */
11881278733SSatish Balay 
11981278733SSatish Balay     /* invert diagonal block */
12081278733SSatish Balay     diag = ba+k*9;
12181278733SSatish Balay     ierr = PetscMemcpy(diag,dk,9*sizeof(MatScalar));CHKERRQ(ierr);
12296b95a6bSBarry Smith     ierr = PetscKernel_A_gets_inverse_A_3(diag,shift);CHKERRQ(ierr);
12381278733SSatish Balay 
12481278733SSatish Balay     jmin = bi[k]; jmax = bi[k+1];
12581278733SSatish Balay     if (jmin < jmax) {
12681278733SSatish Balay       for (j=jmin; j<jmax; j++) {
12781278733SSatish Balay         vj       = bj[j];      /* block col. index of U */
12881278733SSatish Balay         u        = ba + j*9;
12981278733SSatish Balay         rtmp_ptr = rtmp + vj*9;
13081278733SSatish Balay         for (k1=0; k1<9; k1++) {
13181278733SSatish Balay           *u++        = *rtmp_ptr;
13281278733SSatish Balay           *rtmp_ptr++ = 0.0;
13381278733SSatish Balay         }
13481278733SSatish Balay       }
13581278733SSatish Balay 
13681278733SSatish Balay       /* ... add k to row list for first nonzero entry in k-th row */
13781278733SSatish Balay       il[k] = jmin;
13881278733SSatish Balay       i     = bj[jmin];
13981278733SSatish Balay       jl[k] = jl[i]; jl[i] = k;
14081278733SSatish Balay     }
14181278733SSatish Balay   }
14281278733SSatish Balay 
14381278733SSatish Balay   ierr = PetscFree(rtmp);CHKERRQ(ierr);
144d8c74875SBarry Smith   ierr = PetscFree2(il,jl);CHKERRQ(ierr);
145d8c74875SBarry Smith   ierr = PetscFree2(dk,uik);CHKERRQ(ierr);
14681278733SSatish Balay 
1474f79d315SHong Zhang   C->ops->solve          = MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace;
1484f79d315SHong Zhang   C->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace;
1494f79d315SHong Zhang   C->ops->forwardsolve   = MatForwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace;
1504f79d315SHong Zhang   C->ops->backwardsolve  = MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace;
151db4efbfdSBarry Smith 
15281278733SSatish Balay   C->assembled    = PETSC_TRUE;
15381278733SSatish Balay   C->preallocated = PETSC_TRUE;
15426fbe8dcSKarl Rupp 
155efee365bSSatish Balay   ierr = PetscLogFlops(1.3333*27*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */
15681278733SSatish Balay   PetscFunctionReturn(0);
15781278733SSatish Balay }
158