xref: /petsc/src/ksp/pc/impls/cp/cp.c (revision 2fa5cd679192b9b390e47ae2d0650965e6b1d9fa)
124c02a0fSBarry Smith 
2b45d2f2cSJed Brown #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 
1724c02a0fSBarry Smith #undef __FUNCT__
1824c02a0fSBarry Smith #define __FUNCT__ "PCSetUp_CP"
1924c02a0fSBarry Smith static PetscErrorCode PCSetUp_CP(PC pc)
2024c02a0fSBarry Smith {
2124c02a0fSBarry Smith   PC_CP          *cp = (PC_CP*)pc->data;
2224c02a0fSBarry Smith   PetscInt       i,j,*colcnt;
2324c02a0fSBarry Smith   PetscErrorCode ierr;
24ace3abfcSBarry Smith   PetscBool      flg;
2524c02a0fSBarry Smith   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)pc->pmat->data;
2624c02a0fSBarry Smith 
2724c02a0fSBarry Smith   PetscFunctionBegin;
28251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg);CHKERRQ(ierr);
29e7e72b3dSBarry Smith   if (!flg) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Currently only handles SeqAIJ matrices");
3024c02a0fSBarry Smith 
3124c02a0fSBarry Smith   ierr = MatGetLocalSize(pc->pmat,&cp->m,&cp->n);CHKERRQ(ierr);
32e32f2f54SBarry Smith   if (cp->m != cp->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently only for square matrices");
3324c02a0fSBarry Smith 
3424c02a0fSBarry Smith   if (!cp->work) {ierr = MatGetVecs(pc->pmat,&cp->work,PETSC_NULL);CHKERRQ(ierr);}
3524c02a0fSBarry Smith   if (!cp->d) {ierr = PetscMalloc(cp->n*sizeof(PetscScalar),&cp->d);CHKERRQ(ierr);}
3624c02a0fSBarry Smith   if (cp->a && pc->flag != SAME_NONZERO_PATTERN) {
3724c02a0fSBarry Smith     ierr  = PetscFree3(cp->a,cp->i,cp->j);CHKERRQ(ierr);
3824c02a0fSBarry Smith     cp->a = 0;
3924c02a0fSBarry Smith   }
4024c02a0fSBarry Smith 
4124c02a0fSBarry Smith   /* convert to column format */
4224c02a0fSBarry Smith   if (!cp->a) {
4324c02a0fSBarry Smith     ierr = PetscMalloc3(aij->nz,PetscScalar,&cp->a,cp->n+1,PetscInt,&cp->i,aij->nz,PetscInt,&cp->j);CHKERRQ(ierr);
4424c02a0fSBarry Smith   }
4524c02a0fSBarry Smith   ierr = PetscMalloc(cp->n*sizeof(PetscInt),&colcnt);CHKERRQ(ierr);
4624c02a0fSBarry Smith   ierr = PetscMemzero(colcnt,cp->n*sizeof(PetscInt));CHKERRQ(ierr);
4724c02a0fSBarry Smith 
48*2fa5cd67SKarl Rupp   for (i=0; i<aij->nz; i++) colcnt[aij->j[i]]++;
49e60cf9a0SBarry Smith   cp->i[0] = 0;
50*2fa5cd67SKarl Rupp   for (i=0; i<cp->n; i++) cp->i[i+1] = cp->i[i] + colcnt[i];
5124c02a0fSBarry Smith   ierr = PetscMemzero(colcnt,cp->n*sizeof(PetscInt));CHKERRQ(ierr);
5224c02a0fSBarry Smith   for (i=0; i<cp->m; i++) {  /* over rows */
5324c02a0fSBarry Smith     for (j=aij->i[i]; j<aij->i[i+1]; j++) {  /* over columns in row */
5424c02a0fSBarry Smith       cp->j[cp->i[aij->j[j]]+colcnt[aij->j[j]]]   = i;
5524c02a0fSBarry Smith       cp->a[cp->i[aij->j[j]]+colcnt[aij->j[j]]++] = aij->a[j];
5624c02a0fSBarry Smith     }
5724c02a0fSBarry Smith   }
5824c02a0fSBarry Smith   ierr = PetscFree(colcnt);CHKERRQ(ierr);
5924c02a0fSBarry Smith 
6024c02a0fSBarry Smith   /* compute sum of squares of each column d[] */
6124c02a0fSBarry Smith   for (i=0; i<cp->n; i++) {  /* over columns */
6275567043SBarry Smith     cp->d[i] = 0.;
63*2fa5cd67SKarl 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 */
6424c02a0fSBarry Smith     cp->d[i] = 1.0/cp->d[i];
6524c02a0fSBarry Smith   }
6624c02a0fSBarry Smith   PetscFunctionReturn(0);
6724c02a0fSBarry Smith }
6824c02a0fSBarry Smith /* -------------------------------------------------------------------------- */
6924c02a0fSBarry Smith #undef __FUNCT__
7024c02a0fSBarry Smith #define __FUNCT__ "PCApply_CP"
7124c02a0fSBarry Smith static PetscErrorCode PCApply_CP(PC pc,Vec bb,Vec xx)
7224c02a0fSBarry Smith {
7324c02a0fSBarry Smith   PC_CP          *cp = (PC_CP*)pc->data;
7424c02a0fSBarry Smith   PetscErrorCode ierr;
7524c02a0fSBarry Smith   PetscScalar    *b,*x,xt;
7624c02a0fSBarry Smith   PetscInt       i,j;
7724c02a0fSBarry Smith 
7824c02a0fSBarry Smith   PetscFunctionBegin;
7924c02a0fSBarry Smith   ierr = VecCopy(bb,cp->work);CHKERRQ(ierr);
8024c02a0fSBarry Smith   ierr = VecGetArray(cp->work,&b);CHKERRQ(ierr);
8124c02a0fSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
8224c02a0fSBarry Smith 
8324c02a0fSBarry Smith   for (i=0; i<cp->n; i++) {  /* over columns */
8475567043SBarry Smith     xt = 0.;
85*2fa5cd67SKarl Rupp     for (j=cp->i[i]; j<cp->i[i+1]; j++) xt += cp->a[j]*b[cp->j[j]]; /* over rows in column */
8624c02a0fSBarry Smith     xt  *= cp->d[i];
8724c02a0fSBarry Smith     x[i] = xt;
88*2fa5cd67SKarl 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*/
8924c02a0fSBarry Smith   }
9024c02a0fSBarry Smith   for (i=cp->n-1; i>-1; i--) {  /* over columns */
9175567043SBarry Smith     xt = 0.;
92*2fa5cd67SKarl Rupp     for (j=cp->i[i]; j<cp->i[i+1]; j++) xt += cp->a[j]*b[cp->j[j]]; /* over rows in column */
9324c02a0fSBarry Smith     xt  *= cp->d[i];
9424c02a0fSBarry Smith     x[i] = xt;
95*2fa5cd67SKarl 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*/
9624c02a0fSBarry Smith   }
9724c02a0fSBarry Smith 
9824c02a0fSBarry Smith   ierr = VecRestoreArray(cp->work,&b);CHKERRQ(ierr);
9924c02a0fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
10024c02a0fSBarry Smith   PetscFunctionReturn(0);
10124c02a0fSBarry Smith }
10224c02a0fSBarry Smith /* -------------------------------------------------------------------------- */
10324c02a0fSBarry Smith #undef __FUNCT__
10469d2c0f9SBarry Smith #define __FUNCT__ "PCReset_CP"
10569d2c0f9SBarry Smith static PetscErrorCode PCReset_CP(PC pc)
10624c02a0fSBarry Smith {
10724c02a0fSBarry Smith   PC_CP          *cp = (PC_CP*)pc->data;
10824c02a0fSBarry Smith   PetscErrorCode ierr;
10924c02a0fSBarry Smith 
11024c02a0fSBarry Smith   PetscFunctionBegin;
11124c02a0fSBarry Smith   ierr = PetscFree(cp->d);CHKERRQ(ierr);
1126bf464f9SBarry Smith   ierr = VecDestroy(&cp->work);CHKERRQ(ierr);
11324c02a0fSBarry Smith   ierr = PetscFree3(cp->a,cp->i,cp->j);CHKERRQ(ierr);
11469d2c0f9SBarry Smith   PetscFunctionReturn(0);
11569d2c0f9SBarry Smith }
11669d2c0f9SBarry Smith 
11769d2c0f9SBarry Smith #undef __FUNCT__
11869d2c0f9SBarry Smith #define __FUNCT__ "PCDestroy_CP"
11969d2c0f9SBarry Smith static PetscErrorCode PCDestroy_CP(PC pc)
12069d2c0f9SBarry Smith {
12169d2c0f9SBarry Smith   PC_CP          *cp = (PC_CP*)pc->data;
12269d2c0f9SBarry Smith   PetscErrorCode ierr;
12369d2c0f9SBarry Smith 
12469d2c0f9SBarry Smith   PetscFunctionBegin;
12569d2c0f9SBarry Smith   ierr = PCReset_CP(pc);CHKERRQ(ierr);
12669d2c0f9SBarry Smith   ierr = PetscFree(cp->d);CHKERRQ(ierr);
12769d2c0f9SBarry Smith   ierr = PetscFree3(cp->a,cp->i,cp->j);CHKERRQ(ierr);
128c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
12924c02a0fSBarry Smith   PetscFunctionReturn(0);
13024c02a0fSBarry Smith }
13124c02a0fSBarry Smith 
13224c02a0fSBarry Smith #undef __FUNCT__
13324c02a0fSBarry Smith #define __FUNCT__ "PCSetFromOptions_CP"
13424c02a0fSBarry Smith static PetscErrorCode PCSetFromOptions_CP(PC pc)
13524c02a0fSBarry Smith {
13624c02a0fSBarry Smith   PetscFunctionBegin;
13724c02a0fSBarry Smith   PetscFunctionReturn(0);
13824c02a0fSBarry Smith }
13924c02a0fSBarry Smith 
14024c02a0fSBarry Smith 
14124c02a0fSBarry Smith /*MC
14224c02a0fSBarry Smith      PCCP - a "column-projection" preconditioner
14324c02a0fSBarry Smith 
14424c02a0fSBarry Smith      This is a terrible preconditioner and is not recommended, ever!
14524c02a0fSBarry Smith 
14624c02a0fSBarry Smith      Loops over the entries of x computing dx_i to
14724c02a0fSBarry Smith $
14824c02a0fSBarry Smith $        min || b - A(x + dx_i e_i ||_2
14924c02a0fSBarry Smith $        dx_i
15024c02a0fSBarry Smith $
15124c02a0fSBarry Smith $    That is, it changes a single entry of x to minimize the new residual.
15224c02a0fSBarry Smith $   Let A_i represent the ith column of A, then the minimization can be written as
15324c02a0fSBarry Smith $
15424c02a0fSBarry Smith $       min || r - (dx_i) A e_i ||_2
15524c02a0fSBarry Smith $       dx_i
15624c02a0fSBarry Smith $   or   min || r - (dx_i) A_i ||_2
15724c02a0fSBarry Smith $        dx_i
15824c02a0fSBarry Smith $
15924c02a0fSBarry Smith $    take the derivative with respect to dx_i to obtain
16024c02a0fSBarry Smith $        dx_i = (A_i^T A_i)^(-1) A_i^T r
16124c02a0fSBarry Smith $
16224c02a0fSBarry Smith $    This algorithm can be thought of as Gauss-Seidel on the normal equations
16324c02a0fSBarry Smith 
16424c02a0fSBarry Smith     Notes: This proceedure can also be done with block columns or any groups of columns
16524c02a0fSBarry Smith         but this is not coded.
16624c02a0fSBarry Smith 
16724c02a0fSBarry Smith       These "projections" can be done simultaneously for all columns (similar to Jacobi)
16824c02a0fSBarry Smith          or sequentially (similar to Gauss-Seidel/SOR). This is only coded for SOR type.
16924c02a0fSBarry Smith 
17024c02a0fSBarry Smith       This is related to, but not the same as "row projection" methods.
17124c02a0fSBarry Smith 
17224c02a0fSBarry Smith       This is currently coded only for SeqAIJ matrices in sequential (SOR) form.
17328529972SSatish Balay 
17428529972SSatish Balay   Level: intermediate
17524c02a0fSBarry Smith 
17624c02a0fSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PCJACOBI, PCSOR
17724c02a0fSBarry Smith 
17824c02a0fSBarry Smith M*/
17924c02a0fSBarry Smith 
18024c02a0fSBarry Smith EXTERN_C_BEGIN
18124c02a0fSBarry Smith #undef __FUNCT__
18224c02a0fSBarry Smith #define __FUNCT__ "PCCreate_CP"
1837087cfbeSBarry Smith PetscErrorCode  PCCreate_CP(PC pc)
18424c02a0fSBarry Smith {
18524c02a0fSBarry Smith   PC_CP          *cp;
18624c02a0fSBarry Smith   PetscErrorCode ierr;
18724c02a0fSBarry Smith 
18824c02a0fSBarry Smith   PetscFunctionBegin;
18924c02a0fSBarry Smith   ierr     = PetscNewLog(pc,PC_CP,&cp);CHKERRQ(ierr);
19024c02a0fSBarry Smith   pc->data = (void*)cp;
19124c02a0fSBarry Smith 
19224c02a0fSBarry Smith   pc->ops->apply           = PCApply_CP;
19324c02a0fSBarry Smith   pc->ops->applytranspose  = PCApply_CP;
19424c02a0fSBarry Smith   pc->ops->setup           = PCSetUp_CP;
19569d2c0f9SBarry Smith   pc->ops->reset           = PCReset_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