1be1d678aSKris Buschelman #define PETSCMAT_DLL 21c2a3de1SBarry Smith 37c4f633dSBarry Smith #include "../src/mat/impls/sbaij/seq/sbaij.h" 4c60f0209SBarry Smith #include "../src/mat/blockinvert.h" 581278733SSatish Balay 681278733SSatish Balay /* 781278733SSatish Balay Version for when blocks are 3 by 3 Using natural ordering 881278733SSatish Balay */ 94a2ae208SSatish Balay #undef __FUNCT__ 104a2ae208SSatish Balay #define __FUNCT__ "MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering" 110481f469SBarry Smith PetscErrorCode MatCholeskyFactorNumeric_SeqSBAIJ_3_NaturalOrdering(Mat C,Mat A,const MatFactorInfo *info) 1281278733SSatish Balay { 1381278733SSatish Balay Mat_SeqSBAIJ *a = (Mat_SeqSBAIJ*)A->data,*b = (Mat_SeqSBAIJ *)C->data; 14dfbe8321SBarry Smith PetscErrorCode ierr; 1513f74950SBarry Smith PetscInt i,j,mbs=a->mbs,*bi=b->i,*bj=b->j; 1613f74950SBarry Smith PetscInt *ai,*aj,k,k1,jmin,jmax,*jl,*il,vj,nexti,ili; 1781278733SSatish Balay MatScalar *ba = b->a,*aa,*ap,*dk,*uik; 1881278733SSatish Balay MatScalar *u,*diag,*rtmp,*rtmp_ptr; 1962bba022SBarry Smith PetscReal shift = info->shiftinblocks; 2081278733SSatish Balay 2181278733SSatish Balay PetscFunctionBegin; 2281278733SSatish Balay 2381278733SSatish Balay /* initialization */ 2481278733SSatish Balay ierr = PetscMalloc(9*mbs*sizeof(MatScalar),&rtmp);CHKERRQ(ierr); 2581278733SSatish Balay ierr = PetscMemzero(rtmp,9*mbs*sizeof(MatScalar));CHKERRQ(ierr); 26d8c74875SBarry Smith ierr = PetscMalloc2(mbs,PetscInt,&il,mbs,PetscInt,&jl);CHKERRQ(ierr); 2781278733SSatish Balay for (i=0; i<mbs; i++) { 2881278733SSatish Balay jl[i] = mbs; il[0] = 0; 2981278733SSatish Balay } 30d8c74875SBarry Smith ierr = PetscMalloc2(9,MatScalar,&dk,9,MatScalar,&uik);CHKERRQ(ierr); 3181278733SSatish Balay ai = a->i; aj = a->j; aa = a->a; 3281278733SSatish Balay 3381278733SSatish Balay /* for each row k */ 3481278733SSatish Balay for (k = 0; k<mbs; k++){ 3581278733SSatish Balay 3681278733SSatish Balay /*initialize k-th row with elements nonzero in row k of A */ 3781278733SSatish Balay jmin = ai[k]; jmax = ai[k+1]; 3881278733SSatish Balay if (jmin < jmax) { 3981278733SSatish Balay ap = aa + jmin*9; 4081278733SSatish Balay for (j = jmin; j < jmax; j++){ 4181278733SSatish Balay vj = aj[j]; /* block col. index */ 4281278733SSatish Balay rtmp_ptr = rtmp + vj*9; 4381278733SSatish Balay for (i=0; i<9; i++) *rtmp_ptr++ = *ap++; 4481278733SSatish Balay } 4581278733SSatish Balay } 4681278733SSatish Balay 4781278733SSatish Balay /* modify k-th row by adding in those rows i with U(i,k) != 0 */ 4881278733SSatish Balay ierr = PetscMemcpy(dk,rtmp+k*9,9*sizeof(MatScalar));CHKERRQ(ierr); 4981278733SSatish Balay i = jl[k]; /* first row to be added to k_th row */ 5081278733SSatish Balay 5181278733SSatish Balay while (i < mbs){ 5281278733SSatish Balay nexti = jl[i]; /* next row to be added to k_th row */ 5381278733SSatish Balay 5481278733SSatish Balay /* compute multiplier */ 5581278733SSatish Balay ili = il[i]; /* index of first nonzero element in U(i,k:bms-1) */ 5681278733SSatish Balay 5781278733SSatish Balay /* uik = -inv(Di)*U_bar(i,k) */ 5881278733SSatish Balay diag = ba + i*9; 5981278733SSatish Balay u = ba + ili*9; 6081278733SSatish Balay 6181278733SSatish Balay uik[0] = -(diag[0]*u[0] + diag[3]*u[1] + diag[6]*u[2]); 6281278733SSatish Balay uik[1] = -(diag[1]*u[0] + diag[4]*u[1] + diag[7]*u[2]); 6381278733SSatish Balay uik[2] = -(diag[2]*u[0] + diag[5]*u[1] + diag[8]*u[2]); 6481278733SSatish Balay 6581278733SSatish Balay uik[3] = -(diag[0]*u[3] + diag[3]*u[4] + diag[6]*u[5]); 6681278733SSatish Balay uik[4] = -(diag[1]*u[3] + diag[4]*u[4] + diag[7]*u[5]); 6781278733SSatish Balay uik[5] = -(diag[2]*u[3] + diag[5]*u[4] + diag[8]*u[5]); 6881278733SSatish Balay 6981278733SSatish Balay uik[6] = -(diag[0]*u[6] + diag[3]*u[7] + diag[6]*u[8]); 7081278733SSatish Balay uik[7] = -(diag[1]*u[6] + diag[4]*u[7] + diag[7]*u[8]); 7181278733SSatish Balay uik[8] = -(diag[2]*u[6] + diag[5]*u[7] + diag[8]*u[8]); 7281278733SSatish Balay 7381278733SSatish Balay /* update D(k) += -U(i,k)^T * U_bar(i,k) */ 7481278733SSatish Balay dk[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2]; 7581278733SSatish Balay dk[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2]; 7681278733SSatish Balay dk[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2]; 7781278733SSatish Balay 7881278733SSatish Balay dk[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5]; 7981278733SSatish Balay dk[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5]; 8081278733SSatish Balay dk[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5]; 8181278733SSatish Balay 8281278733SSatish Balay dk[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8]; 8381278733SSatish Balay dk[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8]; 8481278733SSatish Balay dk[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8]; 8581278733SSatish Balay 86dc0b31edSSatish Balay ierr = PetscLogFlops(27.0*4.0);CHKERRQ(ierr); 87187a9f4bSHong Zhang 8881278733SSatish Balay /* update -U(i,k) */ 8981278733SSatish Balay ierr = PetscMemcpy(ba+ili*9,uik,9*sizeof(MatScalar));CHKERRQ(ierr); 9081278733SSatish Balay 9181278733SSatish Balay /* add multiple of row i to k-th row ... */ 9281278733SSatish Balay jmin = ili + 1; jmax = bi[i+1]; 9381278733SSatish Balay if (jmin < jmax){ 9481278733SSatish Balay for (j=jmin; j<jmax; j++) { 9581278733SSatish Balay /* rtmp += -U(i,k)^T * U_bar(i,j) */ 9681278733SSatish Balay rtmp_ptr = rtmp + bj[j]*9; 9781278733SSatish Balay u = ba + j*9; 9881278733SSatish Balay rtmp_ptr[0] += uik[0]*u[0] + uik[1]*u[1] + uik[2]*u[2]; 9981278733SSatish Balay rtmp_ptr[1] += uik[3]*u[0] + uik[4]*u[1] + uik[5]*u[2]; 10081278733SSatish Balay rtmp_ptr[2] += uik[6]*u[0] + uik[7]*u[1] + uik[8]*u[2]; 10181278733SSatish Balay 10281278733SSatish Balay rtmp_ptr[3] += uik[0]*u[3] + uik[1]*u[4] + uik[2]*u[5]; 10381278733SSatish Balay rtmp_ptr[4] += uik[3]*u[3] + uik[4]*u[4] + uik[5]*u[5]; 10481278733SSatish Balay rtmp_ptr[5] += uik[6]*u[3] + uik[7]*u[4] + uik[8]*u[5]; 10581278733SSatish Balay 10681278733SSatish Balay rtmp_ptr[6] += uik[0]*u[6] + uik[1]*u[7] + uik[2]*u[8]; 10781278733SSatish Balay rtmp_ptr[7] += uik[3]*u[6] + uik[4]*u[7] + uik[5]*u[8]; 10881278733SSatish Balay rtmp_ptr[8] += uik[6]*u[6] + uik[7]*u[7] + uik[8]*u[8]; 10981278733SSatish Balay } 110dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*27.0*(jmax-jmin));CHKERRQ(ierr); 11181278733SSatish Balay 11281278733SSatish Balay /* ... add i to row list for next nonzero entry */ 11381278733SSatish Balay il[i] = jmin; /* update il(i) in column k+1, ... mbs-1 */ 11481278733SSatish Balay j = bj[jmin]; 11581278733SSatish Balay jl[i] = jl[j]; jl[j] = i; /* update jl */ 11681278733SSatish Balay } 11781278733SSatish Balay i = nexti; 11881278733SSatish Balay } 11981278733SSatish Balay 12081278733SSatish Balay /* save nonzero entries in k-th row of U ... */ 12181278733SSatish Balay 12281278733SSatish Balay /* invert diagonal block */ 12381278733SSatish Balay diag = ba+k*9; 12481278733SSatish Balay ierr = PetscMemcpy(diag,dk,9*sizeof(MatScalar));CHKERRQ(ierr); 12562bba022SBarry Smith ierr = Kernel_A_gets_inverse_A_3(diag,shift);CHKERRQ(ierr); 12681278733SSatish Balay 12781278733SSatish Balay jmin = bi[k]; jmax = bi[k+1]; 12881278733SSatish Balay if (jmin < jmax) { 12981278733SSatish Balay for (j=jmin; j<jmax; j++){ 13081278733SSatish Balay vj = bj[j]; /* block col. index of U */ 13181278733SSatish Balay u = ba + j*9; 13281278733SSatish Balay rtmp_ptr = rtmp + vj*9; 13381278733SSatish Balay for (k1=0; k1<9; k1++){ 13481278733SSatish Balay *u++ = *rtmp_ptr; 13581278733SSatish Balay *rtmp_ptr++ = 0.0; 13681278733SSatish Balay } 13781278733SSatish Balay } 13881278733SSatish Balay 13981278733SSatish Balay /* ... add k to row list for first nonzero entry in k-th row */ 14081278733SSatish Balay il[k] = jmin; 14181278733SSatish Balay i = bj[jmin]; 14281278733SSatish Balay jl[k] = jl[i]; jl[i] = k; 14381278733SSatish Balay } 14481278733SSatish Balay } 14581278733SSatish Balay 14681278733SSatish Balay ierr = PetscFree(rtmp);CHKERRQ(ierr); 147d8c74875SBarry Smith ierr = PetscFree2(il,jl);CHKERRQ(ierr); 148d8c74875SBarry Smith ierr = PetscFree2(dk,uik);CHKERRQ(ierr); 14981278733SSatish Balay 150*4f79d315SHong Zhang C->ops->solve = MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace; 151*4f79d315SHong Zhang C->ops->solvetranspose = MatSolve_SeqSBAIJ_3_NaturalOrdering_inplace; 152*4f79d315SHong Zhang C->ops->forwardsolve = MatForwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace; 153*4f79d315SHong Zhang C->ops->backwardsolve = MatBackwardSolve_SeqSBAIJ_3_NaturalOrdering_inplace; 154db4efbfdSBarry Smith 15581278733SSatish Balay C->assembled = PETSC_TRUE; 15681278733SSatish Balay C->preallocated = PETSC_TRUE; 157efee365bSSatish Balay ierr = PetscLogFlops(1.3333*27*b->mbs);CHKERRQ(ierr); /* from inverting diagonal blocks */ 15881278733SSatish Balay PetscFunctionReturn(0); 15981278733SSatish Balay } 160