124c02a0fSBarry Smith #define PETSCKSP_DLL 224c02a0fSBarry Smith 324c02a0fSBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 47c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h" 524c02a0fSBarry Smith 624c02a0fSBarry Smith /* 724c02a0fSBarry Smith Private context (data structure) for the CP preconditioner. 824c02a0fSBarry Smith */ 924c02a0fSBarry Smith typedef struct { 1024c02a0fSBarry Smith PetscInt n,m; 1124c02a0fSBarry Smith Vec work; 1224c02a0fSBarry Smith PetscScalar *d; /* sum of squares of each column */ 1324c02a0fSBarry Smith PetscScalar *a; /* non-zeros by column */ 1424c02a0fSBarry Smith PetscInt *i,*j; /* offsets of nonzeros by column, non-zero indices by column */ 1524c02a0fSBarry Smith } PC_CP; 1624c02a0fSBarry Smith 1724c02a0fSBarry Smith 1824c02a0fSBarry Smith #undef __FUNCT__ 1924c02a0fSBarry Smith #define __FUNCT__ "PCSetUp_CP" 2024c02a0fSBarry Smith static PetscErrorCode PCSetUp_CP(PC pc) 2124c02a0fSBarry Smith { 2224c02a0fSBarry Smith PC_CP *cp = (PC_CP*)pc->data; 2324c02a0fSBarry Smith PetscInt i,j,*colcnt; 2424c02a0fSBarry Smith PetscErrorCode ierr; 2524c02a0fSBarry Smith PetscTruth flg; 2624c02a0fSBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ*)pc->pmat->data; 2724c02a0fSBarry Smith 2824c02a0fSBarry Smith PetscFunctionBegin; 2924c02a0fSBarry Smith ierr = PetscTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg);CHKERRQ(ierr); 3024c02a0fSBarry Smith if (!flg) SETERRQ(PETSC_ERR_SUP,"Currently only handles SeqAIJ matrices"); 3124c02a0fSBarry Smith 3224c02a0fSBarry Smith ierr = MatGetLocalSize(pc->pmat,&cp->m,&cp->n);CHKERRQ(ierr); 3324c02a0fSBarry Smith if (cp->m != cp->n) SETERRQ(PETSC_ERR_SUP,"Currently only for square matrices"); 3424c02a0fSBarry Smith 3524c02a0fSBarry Smith if (!cp->work) {ierr = MatGetVecs(pc->pmat,&cp->work,PETSC_NULL);CHKERRQ(ierr);} 3624c02a0fSBarry Smith if (!cp->d) {ierr = PetscMalloc(cp->n*sizeof(PetscScalar),&cp->d);CHKERRQ(ierr);} 3724c02a0fSBarry Smith if (cp->a && pc->flag != SAME_NONZERO_PATTERN) { 3824c02a0fSBarry Smith ierr = PetscFree3(cp->a,cp->i,cp->j);CHKERRQ(ierr); 3924c02a0fSBarry Smith cp->a = 0; 4024c02a0fSBarry Smith } 4124c02a0fSBarry Smith 4224c02a0fSBarry Smith /* convert to column format */ 4324c02a0fSBarry Smith if (!cp->a) { 4424c02a0fSBarry Smith ierr = PetscMalloc3(aij->nz,PetscScalar,&cp->a,cp->n+1,PetscInt,&cp->i,aij->nz,PetscInt,&cp->j);CHKERRQ(ierr); 4524c02a0fSBarry Smith } 4624c02a0fSBarry Smith ierr = PetscMalloc(cp->n*sizeof(PetscInt),&colcnt);CHKERRQ(ierr); 4724c02a0fSBarry Smith ierr = PetscMemzero(colcnt,cp->n*sizeof(PetscInt));CHKERRQ(ierr); 4824c02a0fSBarry Smith 4924c02a0fSBarry Smith for (i=0; i<aij->nz; i++) { 5024c02a0fSBarry Smith colcnt[aij->j[i]]++; 5124c02a0fSBarry Smith } 52*e60cf9a0SBarry Smith cp->i[0] = 0; 5324c02a0fSBarry Smith for (i=0; i<cp->n; i++) { 5424c02a0fSBarry Smith cp->i[i+1] = cp->i[i] + colcnt[i]; 5524c02a0fSBarry Smith } 5624c02a0fSBarry Smith ierr = PetscMemzero(colcnt,cp->n*sizeof(PetscInt));CHKERRQ(ierr); 5724c02a0fSBarry Smith for (i=0; i<cp->m; i++) { /* over rows */ 5824c02a0fSBarry Smith for (j=aij->i[i]; j<aij->i[i+1]; j++) { /* over columns in row */ 5924c02a0fSBarry Smith cp->j[cp->i[aij->j[j]]+colcnt[aij->j[j]]] = i; 6024c02a0fSBarry Smith cp->a[cp->i[aij->j[j]]+colcnt[aij->j[j]]++] = aij->a[j]; 6124c02a0fSBarry Smith } 6224c02a0fSBarry Smith } 6324c02a0fSBarry Smith ierr = PetscFree(colcnt);CHKERRQ(ierr); 6424c02a0fSBarry Smith 6524c02a0fSBarry Smith /* compute sum of squares of each column d[] */ 6624c02a0fSBarry Smith for (i=0; i<cp->n; i++) { /* over columns */ 6775567043SBarry Smith cp->d[i] = 0.; 6824c02a0fSBarry Smith for (j=cp->i[i]; j<cp->i[i+1]; j++) { /* over rows in column */ 6924c02a0fSBarry Smith cp->d[i] += cp->a[j]*cp->a[j]; 7024c02a0fSBarry Smith } 7124c02a0fSBarry Smith cp->d[i] = 1.0/cp->d[i]; 7224c02a0fSBarry Smith } 7324c02a0fSBarry Smith PetscFunctionReturn(0); 7424c02a0fSBarry Smith } 7524c02a0fSBarry Smith /* -------------------------------------------------------------------------- */ 7624c02a0fSBarry Smith #undef __FUNCT__ 7724c02a0fSBarry Smith #define __FUNCT__ "PCApply_CP" 7824c02a0fSBarry Smith static PetscErrorCode PCApply_CP(PC pc,Vec bb,Vec xx) 7924c02a0fSBarry Smith { 8024c02a0fSBarry Smith PC_CP *cp = (PC_CP*)pc->data; 8124c02a0fSBarry Smith PetscErrorCode ierr; 8224c02a0fSBarry Smith PetscScalar *b,*x,xt; 8324c02a0fSBarry Smith PetscInt i,j; 8424c02a0fSBarry Smith 8524c02a0fSBarry Smith PetscFunctionBegin; 8624c02a0fSBarry Smith ierr = VecCopy(bb,cp->work);CHKERRQ(ierr); 8724c02a0fSBarry Smith ierr = VecGetArray(cp->work,&b);CHKERRQ(ierr); 8824c02a0fSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 8924c02a0fSBarry Smith 9024c02a0fSBarry Smith for (i=0; i<cp->n; i++) { /* over columns */ 9175567043SBarry Smith xt = 0.; 9224c02a0fSBarry Smith for (j=cp->i[i]; j<cp->i[i+1]; j++) { /* over rows in column */ 9324c02a0fSBarry Smith xt += cp->a[j]*b[cp->j[j]]; 9424c02a0fSBarry Smith } 9524c02a0fSBarry Smith xt *= cp->d[i]; 9624c02a0fSBarry Smith x[i] = xt; 9724c02a0fSBarry Smith for (j=cp->i[i]; j<cp->i[i+1]; j++) { /* over rows in column updating b*/ 9824c02a0fSBarry Smith b[cp->j[j]] -= xt*cp->a[j]; 9924c02a0fSBarry Smith } 10024c02a0fSBarry Smith } 10124c02a0fSBarry Smith for (i=cp->n-1; i>-1; i--) { /* over columns */ 10275567043SBarry Smith xt = 0.; 10324c02a0fSBarry Smith for (j=cp->i[i]; j<cp->i[i+1]; j++) { /* over rows in column */ 10424c02a0fSBarry Smith xt += cp->a[j]*b[cp->j[j]]; 10524c02a0fSBarry Smith } 10624c02a0fSBarry Smith xt *= cp->d[i]; 10724c02a0fSBarry Smith x[i] = xt; 10824c02a0fSBarry Smith for (j=cp->i[i]; j<cp->i[i+1]; j++) { /* over rows in column updating b*/ 10924c02a0fSBarry Smith b[cp->j[j]] -= xt*cp->a[j]; 11024c02a0fSBarry Smith } 11124c02a0fSBarry Smith } 11224c02a0fSBarry Smith 11324c02a0fSBarry Smith ierr = VecRestoreArray(cp->work,&b);CHKERRQ(ierr); 11424c02a0fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 11524c02a0fSBarry Smith PetscFunctionReturn(0); 11624c02a0fSBarry Smith } 11724c02a0fSBarry Smith /* -------------------------------------------------------------------------- */ 11824c02a0fSBarry Smith #undef __FUNCT__ 11924c02a0fSBarry Smith #define __FUNCT__ "PCDestroy_CP" 12024c02a0fSBarry Smith static PetscErrorCode PCDestroy_CP(PC pc) 12124c02a0fSBarry Smith { 12224c02a0fSBarry Smith PC_CP *cp = (PC_CP*)pc->data; 12324c02a0fSBarry Smith PetscErrorCode ierr; 12424c02a0fSBarry Smith 12524c02a0fSBarry Smith PetscFunctionBegin; 12624c02a0fSBarry Smith ierr = PetscFree(cp->d);CHKERRQ(ierr); 12724c02a0fSBarry Smith if (cp->work) {ierr = VecDestroy(cp->work);CHKERRQ(ierr);} 12824c02a0fSBarry Smith ierr = PetscFree3(cp->a,cp->i,cp->j);CHKERRQ(ierr); 12924c02a0fSBarry Smith ierr = PetscFree(cp);CHKERRQ(ierr); 13024c02a0fSBarry Smith PetscFunctionReturn(0); 13124c02a0fSBarry Smith } 13224c02a0fSBarry Smith 13324c02a0fSBarry Smith #undef __FUNCT__ 13424c02a0fSBarry Smith #define __FUNCT__ "PCSetFromOptions_CP" 13524c02a0fSBarry Smith static PetscErrorCode PCSetFromOptions_CP(PC pc) 13624c02a0fSBarry Smith { 13724c02a0fSBarry Smith PetscFunctionBegin; 13824c02a0fSBarry Smith PetscFunctionReturn(0); 13924c02a0fSBarry Smith } 14024c02a0fSBarry Smith 14124c02a0fSBarry Smith 14224c02a0fSBarry Smith /*MC 14324c02a0fSBarry Smith PCCP - a "column-projection" preconditioner 14424c02a0fSBarry Smith 14524c02a0fSBarry Smith This is a terrible preconditioner and is not recommended, ever! 14624c02a0fSBarry Smith 14724c02a0fSBarry Smith Loops over the entries of x computing dx_i to 14824c02a0fSBarry Smith $ 14924c02a0fSBarry Smith $ min || b - A(x + dx_i e_i ||_2 15024c02a0fSBarry Smith $ dx_i 15124c02a0fSBarry Smith $ 15224c02a0fSBarry Smith $ That is, it changes a single entry of x to minimize the new residual. 15324c02a0fSBarry Smith $ Let A_i represent the ith column of A, then the minimization can be written as 15424c02a0fSBarry Smith $ 15524c02a0fSBarry Smith $ min || r - (dx_i) A e_i ||_2 15624c02a0fSBarry Smith $ dx_i 15724c02a0fSBarry Smith $ or min || r - (dx_i) A_i ||_2 15824c02a0fSBarry Smith $ dx_i 15924c02a0fSBarry Smith $ 16024c02a0fSBarry Smith $ take the derivative with respect to dx_i to obtain 16124c02a0fSBarry Smith $ dx_i = (A_i^T A_i)^(-1) A_i^T r 16224c02a0fSBarry Smith $ 16324c02a0fSBarry Smith $ This algorithm can be thought of as Gauss-Seidel on the normal equations 16424c02a0fSBarry Smith 16524c02a0fSBarry Smith Notes: This proceedure can also be done with block columns or any groups of columns 16624c02a0fSBarry Smith but this is not coded. 16724c02a0fSBarry Smith 16824c02a0fSBarry Smith These "projections" can be done simultaneously for all columns (similar to Jacobi) 16924c02a0fSBarry Smith or sequentially (similar to Gauss-Seidel/SOR). This is only coded for SOR type. 17024c02a0fSBarry Smith 17124c02a0fSBarry Smith This is related to, but not the same as "row projection" methods. 17224c02a0fSBarry Smith 17324c02a0fSBarry Smith This is currently coded only for SeqAIJ matrices in sequential (SOR) form. 17428529972SSatish Balay 17528529972SSatish Balay Level: intermediate 17624c02a0fSBarry Smith 17724c02a0fSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCJACOBI, PCSOR 17824c02a0fSBarry Smith 17924c02a0fSBarry Smith M*/ 18024c02a0fSBarry Smith 18124c02a0fSBarry Smith EXTERN_C_BEGIN 18224c02a0fSBarry Smith #undef __FUNCT__ 18324c02a0fSBarry Smith #define __FUNCT__ "PCCreate_CP" 18424c02a0fSBarry Smith PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_CP(PC pc) 18524c02a0fSBarry Smith { 18624c02a0fSBarry Smith PC_CP *cp; 18724c02a0fSBarry Smith PetscErrorCode ierr; 18824c02a0fSBarry Smith 18924c02a0fSBarry Smith PetscFunctionBegin; 19024c02a0fSBarry Smith ierr = PetscNewLog(pc,PC_CP,&cp);CHKERRQ(ierr); 19124c02a0fSBarry Smith pc->data = (void*)cp; 19224c02a0fSBarry Smith 19324c02a0fSBarry Smith pc->ops->apply = PCApply_CP; 19424c02a0fSBarry Smith pc->ops->applytranspose = PCApply_CP; 19524c02a0fSBarry Smith pc->ops->setup = PCSetUp_CP; 19624c02a0fSBarry Smith pc->ops->destroy = PCDestroy_CP; 19724c02a0fSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_CP; 19824c02a0fSBarry Smith pc->ops->view = 0; 19924c02a0fSBarry Smith pc->ops->applyrichardson = 0; 20024c02a0fSBarry Smith PetscFunctionReturn(0); 20124c02a0fSBarry Smith } 20224c02a0fSBarry Smith EXTERN_C_END 20324c02a0fSBarry Smith 20424c02a0fSBarry Smith 205