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; 2024c02a0fSBarry Smith PetscErrorCode ierr; 21ace3abfcSBarry Smith PetscBool flg; 2224c02a0fSBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ*)pc->pmat->data; 2324c02a0fSBarry Smith 2424c02a0fSBarry Smith PetscFunctionBegin; 25251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg);CHKERRQ(ierr); 26*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Currently only handles SeqAIJ matrices"); 2724c02a0fSBarry Smith 2824c02a0fSBarry Smith ierr = MatGetLocalSize(pc->pmat,&cp->m,&cp->n);CHKERRQ(ierr); 29*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(cp->m != cp->n,PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently only for square matrices"); 3024c02a0fSBarry Smith 312a7a6963SBarry Smith if (!cp->work) {ierr = MatCreateVecs(pc->pmat,&cp->work,NULL);CHKERRQ(ierr);} 32785e854fSJed Brown if (!cp->d) {ierr = PetscMalloc1(cp->n,&cp->d);CHKERRQ(ierr);} 3324c02a0fSBarry Smith if (cp->a && pc->flag != SAME_NONZERO_PATTERN) { 3424c02a0fSBarry Smith ierr = PetscFree3(cp->a,cp->i,cp->j);CHKERRQ(ierr); 350a545947SLisandro Dalcin cp->a = NULL; 3624c02a0fSBarry Smith } 3724c02a0fSBarry Smith 3824c02a0fSBarry Smith /* convert to column format */ 3924c02a0fSBarry Smith if (!cp->a) { 40dcca6d9dSJed Brown ierr = PetscMalloc3(aij->nz,&cp->a,cp->n+1,&cp->i,aij->nz,&cp->j);CHKERRQ(ierr); 4124c02a0fSBarry Smith } 421795a4d1SJed Brown ierr = PetscCalloc1(cp->n,&colcnt);CHKERRQ(ierr); 4324c02a0fSBarry Smith 442fa5cd67SKarl Rupp for (i=0; i<aij->nz; i++) colcnt[aij->j[i]]++; 45e60cf9a0SBarry Smith cp->i[0] = 0; 462fa5cd67SKarl Rupp for (i=0; i<cp->n; i++) cp->i[i+1] = cp->i[i] + colcnt[i]; 47580bdb30SBarry Smith ierr = PetscArrayzero(colcnt,cp->n);CHKERRQ(ierr); 4824c02a0fSBarry Smith for (i=0; i<cp->m; i++) { /* over rows */ 4924c02a0fSBarry Smith for (j=aij->i[i]; j<aij->i[i+1]; j++) { /* over columns in row */ 5024c02a0fSBarry Smith cp->j[cp->i[aij->j[j]]+colcnt[aij->j[j]]] = i; 5124c02a0fSBarry Smith cp->a[cp->i[aij->j[j]]+colcnt[aij->j[j]]++] = aij->a[j]; 5224c02a0fSBarry Smith } 5324c02a0fSBarry Smith } 5424c02a0fSBarry Smith ierr = PetscFree(colcnt);CHKERRQ(ierr); 5524c02a0fSBarry Smith 5624c02a0fSBarry Smith /* compute sum of squares of each column d[] */ 5724c02a0fSBarry Smith for (i=0; i<cp->n; i++) { /* over columns */ 5875567043SBarry Smith cp->d[i] = 0.; 592fa5cd67SKarl 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 */ 6024c02a0fSBarry Smith cp->d[i] = 1.0/cp->d[i]; 6124c02a0fSBarry Smith } 6224c02a0fSBarry Smith PetscFunctionReturn(0); 6324c02a0fSBarry Smith } 6424c02a0fSBarry Smith /* -------------------------------------------------------------------------- */ 6524c02a0fSBarry Smith static PetscErrorCode PCApply_CP(PC pc,Vec bb,Vec xx) 6624c02a0fSBarry Smith { 6724c02a0fSBarry Smith PC_CP *cp = (PC_CP*)pc->data; 6824c02a0fSBarry Smith PetscErrorCode ierr; 6924c02a0fSBarry Smith PetscScalar *b,*x,xt; 7024c02a0fSBarry Smith PetscInt i,j; 7124c02a0fSBarry Smith 7224c02a0fSBarry Smith PetscFunctionBegin; 7324c02a0fSBarry Smith ierr = VecCopy(bb,cp->work);CHKERRQ(ierr); 7424c02a0fSBarry Smith ierr = VecGetArray(cp->work,&b);CHKERRQ(ierr); 7524c02a0fSBarry Smith ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 7624c02a0fSBarry Smith 7724c02a0fSBarry Smith for (i=0; i<cp->n; i++) { /* over columns */ 7875567043SBarry Smith xt = 0.; 792fa5cd67SKarl Rupp for (j=cp->i[i]; j<cp->i[i+1]; j++) xt += cp->a[j]*b[cp->j[j]]; /* over rows in column */ 8024c02a0fSBarry Smith xt *= cp->d[i]; 8124c02a0fSBarry Smith x[i] = xt; 822fa5cd67SKarl 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*/ 8324c02a0fSBarry Smith } 8424c02a0fSBarry Smith for (i=cp->n-1; i>-1; i--) { /* over columns */ 8575567043SBarry Smith xt = 0.; 862fa5cd67SKarl Rupp for (j=cp->i[i]; j<cp->i[i+1]; j++) xt += cp->a[j]*b[cp->j[j]]; /* over rows in column */ 8724c02a0fSBarry Smith xt *= cp->d[i]; 8824c02a0fSBarry Smith x[i] = xt; 892fa5cd67SKarl 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*/ 9024c02a0fSBarry Smith } 9124c02a0fSBarry Smith 9224c02a0fSBarry Smith ierr = VecRestoreArray(cp->work,&b);CHKERRQ(ierr); 9324c02a0fSBarry Smith ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 9424c02a0fSBarry Smith PetscFunctionReturn(0); 9524c02a0fSBarry Smith } 9624c02a0fSBarry Smith /* -------------------------------------------------------------------------- */ 9769d2c0f9SBarry Smith static PetscErrorCode PCReset_CP(PC pc) 9824c02a0fSBarry Smith { 9924c02a0fSBarry Smith PC_CP *cp = (PC_CP*)pc->data; 10024c02a0fSBarry Smith PetscErrorCode ierr; 10124c02a0fSBarry Smith 10224c02a0fSBarry Smith PetscFunctionBegin; 10324c02a0fSBarry Smith ierr = PetscFree(cp->d);CHKERRQ(ierr); 1046bf464f9SBarry Smith ierr = VecDestroy(&cp->work);CHKERRQ(ierr); 10524c02a0fSBarry Smith ierr = PetscFree3(cp->a,cp->i,cp->j);CHKERRQ(ierr); 10669d2c0f9SBarry Smith PetscFunctionReturn(0); 10769d2c0f9SBarry Smith } 10869d2c0f9SBarry Smith 10969d2c0f9SBarry Smith static PetscErrorCode PCDestroy_CP(PC pc) 11069d2c0f9SBarry Smith { 11169d2c0f9SBarry Smith PC_CP *cp = (PC_CP*)pc->data; 11269d2c0f9SBarry Smith PetscErrorCode ierr; 11369d2c0f9SBarry Smith 11469d2c0f9SBarry Smith PetscFunctionBegin; 11569d2c0f9SBarry Smith ierr = PCReset_CP(pc);CHKERRQ(ierr); 11669d2c0f9SBarry Smith ierr = PetscFree(cp->d);CHKERRQ(ierr); 11769d2c0f9SBarry Smith ierr = PetscFree3(cp->a,cp->i,cp->j);CHKERRQ(ierr); 118c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 11924c02a0fSBarry Smith PetscFunctionReturn(0); 12024c02a0fSBarry Smith } 12124c02a0fSBarry Smith 1224416b707SBarry Smith static PetscErrorCode PCSetFromOptions_CP(PetscOptionItems *PetscOptionsObject,PC pc) 12324c02a0fSBarry Smith { 12424c02a0fSBarry Smith PetscFunctionBegin; 12524c02a0fSBarry Smith PetscFunctionReturn(0); 12624c02a0fSBarry Smith } 12724c02a0fSBarry Smith 12824c02a0fSBarry Smith /*MC 12924c02a0fSBarry Smith PCCP - a "column-projection" preconditioner 13024c02a0fSBarry Smith 13124c02a0fSBarry Smith This is a terrible preconditioner and is not recommended, ever! 13224c02a0fSBarry Smith 13373a58da7SBarry Smith Loops over the entries of x computing dx_i (e_i is the unit vector in the ith direction) to 13424c02a0fSBarry Smith $ 13524c02a0fSBarry Smith $ min || b - A(x + dx_i e_i ||_2 13624c02a0fSBarry Smith $ dx_i 13724c02a0fSBarry Smith $ 13873a58da7SBarry Smith $ That is, it changes a single entry of x to minimize the new residual norm. 13924c02a0fSBarry Smith $ Let A_i represent the ith column of A, then the minimization can be written as 14024c02a0fSBarry Smith $ 14124c02a0fSBarry Smith $ min || r - (dx_i) A e_i ||_2 14224c02a0fSBarry Smith $ dx_i 14324c02a0fSBarry Smith $ or min || r - (dx_i) A_i ||_2 14424c02a0fSBarry Smith $ dx_i 14524c02a0fSBarry Smith $ 14624c02a0fSBarry Smith $ take the derivative with respect to dx_i to obtain 14724c02a0fSBarry Smith $ dx_i = (A_i^T A_i)^(-1) A_i^T r 14824c02a0fSBarry Smith $ 14924c02a0fSBarry Smith $ This algorithm can be thought of as Gauss-Seidel on the normal equations 15024c02a0fSBarry Smith 15195452b02SPatrick Sanan Notes: 15295452b02SPatrick Sanan This proceedure can also be done with block columns or any groups of columns 15324c02a0fSBarry Smith but this is not coded. 15424c02a0fSBarry Smith 15524c02a0fSBarry Smith These "projections" can be done simultaneously for all columns (similar to Jacobi) 15624c02a0fSBarry Smith or sequentially (similar to Gauss-Seidel/SOR). This is only coded for SOR type. 15724c02a0fSBarry Smith 15824c02a0fSBarry Smith This is related to, but not the same as "row projection" methods. 15924c02a0fSBarry Smith 16024c02a0fSBarry Smith This is currently coded only for SeqAIJ matrices in sequential (SOR) form. 16128529972SSatish Balay 16228529972SSatish Balay Level: intermediate 16324c02a0fSBarry Smith 16424c02a0fSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCJACOBI, PCSOR 16524c02a0fSBarry Smith 16624c02a0fSBarry Smith M*/ 16724c02a0fSBarry Smith 1688cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_CP(PC pc) 16924c02a0fSBarry Smith { 17024c02a0fSBarry Smith PC_CP *cp; 17124c02a0fSBarry Smith PetscErrorCode ierr; 17224c02a0fSBarry Smith 17324c02a0fSBarry Smith PetscFunctionBegin; 174b00a9115SJed Brown ierr = PetscNewLog(pc,&cp);CHKERRQ(ierr); 17524c02a0fSBarry Smith pc->data = (void*)cp; 17624c02a0fSBarry Smith 17724c02a0fSBarry Smith pc->ops->apply = PCApply_CP; 17824c02a0fSBarry Smith pc->ops->applytranspose = PCApply_CP; 17924c02a0fSBarry Smith pc->ops->setup = PCSetUp_CP; 18069d2c0f9SBarry Smith pc->ops->reset = PCReset_CP; 18124c02a0fSBarry Smith pc->ops->destroy = PCDestroy_CP; 18224c02a0fSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_CP; 1830a545947SLisandro Dalcin pc->ops->view = NULL; 1840a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 18524c02a0fSBarry Smith PetscFunctionReturn(0); 18624c02a0fSBarry Smith } 18724c02a0fSBarry Smith 188