124c02a0fSBarry Smith 2af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 3c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> 424c02a0fSBarry Smith 524c02a0fSBarry Smith /* 624c02a0fSBarry Smith Private context (data structure) for the CP preconditioner. 724c02a0fSBarry Smith */ 824c02a0fSBarry Smith typedef struct { 924c02a0fSBarry Smith PetscInt n,m; 1024c02a0fSBarry Smith Vec work; 1124c02a0fSBarry Smith PetscScalar *d; /* sum of squares of each column */ 1224c02a0fSBarry Smith PetscScalar *a; /* non-zeros by column */ 1324c02a0fSBarry Smith PetscInt *i,*j; /* offsets of nonzeros by column, non-zero indices by column */ 1424c02a0fSBarry Smith } PC_CP; 1524c02a0fSBarry Smith 1624c02a0fSBarry Smith static PetscErrorCode PCSetUp_CP(PC pc) 1724c02a0fSBarry Smith { 1824c02a0fSBarry Smith PC_CP *cp = (PC_CP*)pc->data; 1924c02a0fSBarry Smith PetscInt i,j,*colcnt; 20ace3abfcSBarry Smith PetscBool flg; 2124c02a0fSBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ*)pc->pmat->data; 2224c02a0fSBarry Smith 2324c02a0fSBarry Smith PetscFunctionBegin; 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg)); 252c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Currently only handles SeqAIJ matrices"); 2624c02a0fSBarry Smith 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetLocalSize(pc->pmat,&cp->m,&cp->n)); 282c71b3e2SJacob Faibussowitsch PetscCheckFalse(cp->m != cp->n,PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently only for square matrices"); 2924c02a0fSBarry Smith 30*5f80ce2aSJacob Faibussowitsch if (!cp->work) CHKERRQ(MatCreateVecs(pc->pmat,&cp->work,NULL)); 31*5f80ce2aSJacob Faibussowitsch if (!cp->d) CHKERRQ(PetscMalloc1(cp->n,&cp->d)); 3224c02a0fSBarry Smith if (cp->a && pc->flag != SAME_NONZERO_PATTERN) { 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(cp->a,cp->i,cp->j)); 340a545947SLisandro Dalcin cp->a = NULL; 3524c02a0fSBarry Smith } 3624c02a0fSBarry Smith 3724c02a0fSBarry Smith /* convert to column format */ 3824c02a0fSBarry Smith if (!cp->a) { 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc3(aij->nz,&cp->a,cp->n+1,&cp->i,aij->nz,&cp->j)); 4024c02a0fSBarry Smith } 41*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(cp->n,&colcnt)); 4224c02a0fSBarry Smith 432fa5cd67SKarl Rupp for (i=0; i<aij->nz; i++) colcnt[aij->j[i]]++; 44e60cf9a0SBarry Smith cp->i[0] = 0; 452fa5cd67SKarl Rupp for (i=0; i<cp->n; i++) cp->i[i+1] = cp->i[i] + colcnt[i]; 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArrayzero(colcnt,cp->n)); 4724c02a0fSBarry Smith for (i=0; i<cp->m; i++) { /* over rows */ 4824c02a0fSBarry Smith for (j=aij->i[i]; j<aij->i[i+1]; j++) { /* over columns in row */ 4924c02a0fSBarry Smith cp->j[cp->i[aij->j[j]]+colcnt[aij->j[j]]] = i; 5024c02a0fSBarry Smith cp->a[cp->i[aij->j[j]]+colcnt[aij->j[j]]++] = aij->a[j]; 5124c02a0fSBarry Smith } 5224c02a0fSBarry Smith } 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(colcnt)); 5424c02a0fSBarry Smith 5524c02a0fSBarry Smith /* compute sum of squares of each column d[] */ 5624c02a0fSBarry Smith for (i=0; i<cp->n; i++) { /* over columns */ 5775567043SBarry Smith cp->d[i] = 0.; 582fa5cd67SKarl Rupp for (j=cp->i[i]; j<cp->i[i+1]; j++) cp->d[i] += cp->a[j]*cp->a[j]; /* over rows in column */ 5924c02a0fSBarry Smith cp->d[i] = 1.0/cp->d[i]; 6024c02a0fSBarry Smith } 6124c02a0fSBarry Smith PetscFunctionReturn(0); 6224c02a0fSBarry Smith } 6324c02a0fSBarry Smith /* -------------------------------------------------------------------------- */ 6424c02a0fSBarry Smith static PetscErrorCode PCApply_CP(PC pc,Vec bb,Vec xx) 6524c02a0fSBarry Smith { 6624c02a0fSBarry Smith PC_CP *cp = (PC_CP*)pc->data; 6724c02a0fSBarry Smith PetscScalar *b,*x,xt; 6824c02a0fSBarry Smith PetscInt i,j; 6924c02a0fSBarry Smith 7024c02a0fSBarry Smith PetscFunctionBegin; 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(bb,cp->work)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(cp->work,&b)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(xx,&x)); 7424c02a0fSBarry Smith 7524c02a0fSBarry Smith for (i=0; i<cp->n; i++) { /* over columns */ 7675567043SBarry Smith xt = 0.; 772fa5cd67SKarl Rupp for (j=cp->i[i]; j<cp->i[i+1]; j++) xt += cp->a[j]*b[cp->j[j]]; /* over rows in column */ 7824c02a0fSBarry Smith xt *= cp->d[i]; 7924c02a0fSBarry Smith x[i] = xt; 802fa5cd67SKarl Rupp for (j=cp->i[i]; j<cp->i[i+1]; j++) b[cp->j[j]] -= xt*cp->a[j]; /* over rows in column updating b*/ 8124c02a0fSBarry Smith } 8224c02a0fSBarry Smith for (i=cp->n-1; i>-1; i--) { /* over columns */ 8375567043SBarry Smith xt = 0.; 842fa5cd67SKarl Rupp for (j=cp->i[i]; j<cp->i[i+1]; j++) xt += cp->a[j]*b[cp->j[j]]; /* over rows in column */ 8524c02a0fSBarry Smith xt *= cp->d[i]; 8624c02a0fSBarry Smith x[i] = xt; 872fa5cd67SKarl Rupp for (j=cp->i[i]; j<cp->i[i+1]; j++) b[cp->j[j]] -= xt*cp->a[j]; /* over rows in column updating b*/ 8824c02a0fSBarry Smith } 8924c02a0fSBarry Smith 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(cp->work,&b)); 91*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(xx,&x)); 9224c02a0fSBarry Smith PetscFunctionReturn(0); 9324c02a0fSBarry Smith } 9424c02a0fSBarry Smith /* -------------------------------------------------------------------------- */ 9569d2c0f9SBarry Smith static PetscErrorCode PCReset_CP(PC pc) 9624c02a0fSBarry Smith { 9724c02a0fSBarry Smith PC_CP *cp = (PC_CP*)pc->data; 9824c02a0fSBarry Smith 9924c02a0fSBarry Smith PetscFunctionBegin; 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(cp->d)); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&cp->work)); 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(cp->a,cp->i,cp->j)); 10369d2c0f9SBarry Smith PetscFunctionReturn(0); 10469d2c0f9SBarry Smith } 10569d2c0f9SBarry Smith 10669d2c0f9SBarry Smith static PetscErrorCode PCDestroy_CP(PC pc) 10769d2c0f9SBarry Smith { 10869d2c0f9SBarry Smith PC_CP *cp = (PC_CP*)pc->data; 10969d2c0f9SBarry Smith 11069d2c0f9SBarry Smith PetscFunctionBegin; 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(PCReset_CP(pc)); 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(cp->d)); 113*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree3(cp->a,cp->i,cp->j)); 114*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pc->data)); 11524c02a0fSBarry Smith PetscFunctionReturn(0); 11624c02a0fSBarry Smith } 11724c02a0fSBarry Smith 1184416b707SBarry Smith static PetscErrorCode PCSetFromOptions_CP(PetscOptionItems *PetscOptionsObject,PC pc) 11924c02a0fSBarry Smith { 12024c02a0fSBarry Smith PetscFunctionBegin; 12124c02a0fSBarry Smith PetscFunctionReturn(0); 12224c02a0fSBarry Smith } 12324c02a0fSBarry Smith 12424c02a0fSBarry Smith /*MC 12524c02a0fSBarry Smith PCCP - a "column-projection" preconditioner 12624c02a0fSBarry Smith 12724c02a0fSBarry Smith This is a terrible preconditioner and is not recommended, ever! 12824c02a0fSBarry Smith 12973a58da7SBarry Smith Loops over the entries of x computing dx_i (e_i is the unit vector in the ith direction) to 13024c02a0fSBarry Smith $ 13124c02a0fSBarry Smith $ min || b - A(x + dx_i e_i ||_2 13224c02a0fSBarry Smith $ dx_i 13324c02a0fSBarry Smith $ 13473a58da7SBarry Smith $ That is, it changes a single entry of x to minimize the new residual norm. 13524c02a0fSBarry Smith $ Let A_i represent the ith column of A, then the minimization can be written as 13624c02a0fSBarry Smith $ 13724c02a0fSBarry Smith $ min || r - (dx_i) A e_i ||_2 13824c02a0fSBarry Smith $ dx_i 13924c02a0fSBarry Smith $ or min || r - (dx_i) A_i ||_2 14024c02a0fSBarry Smith $ dx_i 14124c02a0fSBarry Smith $ 14224c02a0fSBarry Smith $ take the derivative with respect to dx_i to obtain 14324c02a0fSBarry Smith $ dx_i = (A_i^T A_i)^(-1) A_i^T r 14424c02a0fSBarry Smith $ 14524c02a0fSBarry Smith $ This algorithm can be thought of as Gauss-Seidel on the normal equations 14624c02a0fSBarry Smith 14795452b02SPatrick Sanan Notes: 14895452b02SPatrick Sanan This proceedure can also be done with block columns or any groups of columns 14924c02a0fSBarry Smith but this is not coded. 15024c02a0fSBarry Smith 15124c02a0fSBarry Smith These "projections" can be done simultaneously for all columns (similar to Jacobi) 15224c02a0fSBarry Smith or sequentially (similar to Gauss-Seidel/SOR). This is only coded for SOR type. 15324c02a0fSBarry Smith 15424c02a0fSBarry Smith This is related to, but not the same as "row projection" methods. 15524c02a0fSBarry Smith 15624c02a0fSBarry Smith This is currently coded only for SeqAIJ matrices in sequential (SOR) form. 15728529972SSatish Balay 15828529972SSatish Balay Level: intermediate 15924c02a0fSBarry Smith 16024c02a0fSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCJACOBI, PCSOR 16124c02a0fSBarry Smith 16224c02a0fSBarry Smith M*/ 16324c02a0fSBarry Smith 1648cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_CP(PC pc) 16524c02a0fSBarry Smith { 16624c02a0fSBarry Smith PC_CP *cp; 16724c02a0fSBarry Smith 16824c02a0fSBarry Smith PetscFunctionBegin; 169*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(pc,&cp)); 17024c02a0fSBarry Smith pc->data = (void*)cp; 17124c02a0fSBarry Smith 17224c02a0fSBarry Smith pc->ops->apply = PCApply_CP; 17324c02a0fSBarry Smith pc->ops->applytranspose = PCApply_CP; 17424c02a0fSBarry Smith pc->ops->setup = PCSetUp_CP; 17569d2c0f9SBarry Smith pc->ops->reset = PCReset_CP; 17624c02a0fSBarry Smith pc->ops->destroy = PCDestroy_CP; 17724c02a0fSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_CP; 1780a545947SLisandro Dalcin pc->ops->view = NULL; 1790a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 18024c02a0fSBarry Smith PetscFunctionReturn(0); 18124c02a0fSBarry Smith } 182