xref: /petsc/src/ksp/pc/impls/cp/cp.c (revision b00a91154f763f12aa55f3d53a3f2776f15f49e3)
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);
29ce94432eSBarry Smith   if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),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 
340298fd71SBarry Smith   if (!cp->work) {ierr = MatGetVecs(pc->pmat,&cp->work,NULL);CHKERRQ(ierr);}
35785e854fSJed Brown   if (!cp->d) {ierr = PetscMalloc1(cp->n,&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) {
43dcca6d9dSJed Brown     ierr = PetscMalloc3(aij->nz,&cp->a,cp->n+1,&cp->i,aij->nz,&cp->j);CHKERRQ(ierr);
4424c02a0fSBarry Smith   }
451795a4d1SJed Brown   ierr = PetscCalloc1(cp->n,&colcnt);CHKERRQ(ierr);
4624c02a0fSBarry Smith 
472fa5cd67SKarl Rupp   for (i=0; i<aij->nz; i++) colcnt[aij->j[i]]++;
48e60cf9a0SBarry Smith   cp->i[0] = 0;
492fa5cd67SKarl Rupp   for (i=0; i<cp->n; i++) cp->i[i+1] = cp->i[i] + colcnt[i];
5024c02a0fSBarry Smith   ierr = PetscMemzero(colcnt,cp->n*sizeof(PetscInt));CHKERRQ(ierr);
5124c02a0fSBarry Smith   for (i=0; i<cp->m; i++) {  /* over rows */
5224c02a0fSBarry Smith     for (j=aij->i[i]; j<aij->i[i+1]; j++) {  /* over columns in row */
5324c02a0fSBarry Smith       cp->j[cp->i[aij->j[j]]+colcnt[aij->j[j]]]   = i;
5424c02a0fSBarry Smith       cp->a[cp->i[aij->j[j]]+colcnt[aij->j[j]]++] = aij->a[j];
5524c02a0fSBarry Smith     }
5624c02a0fSBarry Smith   }
5724c02a0fSBarry Smith   ierr = PetscFree(colcnt);CHKERRQ(ierr);
5824c02a0fSBarry Smith 
5924c02a0fSBarry Smith   /* compute sum of squares of each column d[] */
6024c02a0fSBarry Smith   for (i=0; i<cp->n; i++) {  /* over columns */
6175567043SBarry Smith     cp->d[i] = 0.;
622fa5cd67SKarl 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 */
6324c02a0fSBarry Smith     cp->d[i] = 1.0/cp->d[i];
6424c02a0fSBarry Smith   }
6524c02a0fSBarry Smith   PetscFunctionReturn(0);
6624c02a0fSBarry Smith }
6724c02a0fSBarry Smith /* -------------------------------------------------------------------------- */
6824c02a0fSBarry Smith #undef __FUNCT__
6924c02a0fSBarry Smith #define __FUNCT__ "PCApply_CP"
7024c02a0fSBarry Smith static PetscErrorCode PCApply_CP(PC pc,Vec bb,Vec xx)
7124c02a0fSBarry Smith {
7224c02a0fSBarry Smith   PC_CP          *cp = (PC_CP*)pc->data;
7324c02a0fSBarry Smith   PetscErrorCode ierr;
7424c02a0fSBarry Smith   PetscScalar    *b,*x,xt;
7524c02a0fSBarry Smith   PetscInt       i,j;
7624c02a0fSBarry Smith 
7724c02a0fSBarry Smith   PetscFunctionBegin;
7824c02a0fSBarry Smith   ierr = VecCopy(bb,cp->work);CHKERRQ(ierr);
7924c02a0fSBarry Smith   ierr = VecGetArray(cp->work,&b);CHKERRQ(ierr);
8024c02a0fSBarry Smith   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
8124c02a0fSBarry Smith 
8224c02a0fSBarry Smith   for (i=0; i<cp->n; 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   for (i=cp->n-1; i>-1; i--) {  /* over columns */
9075567043SBarry Smith     xt = 0.;
912fa5cd67SKarl Rupp     for (j=cp->i[i]; j<cp->i[i+1]; j++) xt += cp->a[j]*b[cp->j[j]]; /* over rows in column */
9224c02a0fSBarry Smith     xt  *= cp->d[i];
9324c02a0fSBarry Smith     x[i] = xt;
942fa5cd67SKarl 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*/
9524c02a0fSBarry Smith   }
9624c02a0fSBarry Smith 
9724c02a0fSBarry Smith   ierr = VecRestoreArray(cp->work,&b);CHKERRQ(ierr);
9824c02a0fSBarry Smith   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
9924c02a0fSBarry Smith   PetscFunctionReturn(0);
10024c02a0fSBarry Smith }
10124c02a0fSBarry Smith /* -------------------------------------------------------------------------- */
10224c02a0fSBarry Smith #undef __FUNCT__
10369d2c0f9SBarry Smith #define __FUNCT__ "PCReset_CP"
10469d2c0f9SBarry Smith static PetscErrorCode PCReset_CP(PC pc)
10524c02a0fSBarry Smith {
10624c02a0fSBarry Smith   PC_CP          *cp = (PC_CP*)pc->data;
10724c02a0fSBarry Smith   PetscErrorCode ierr;
10824c02a0fSBarry Smith 
10924c02a0fSBarry Smith   PetscFunctionBegin;
11024c02a0fSBarry Smith   ierr = PetscFree(cp->d);CHKERRQ(ierr);
1116bf464f9SBarry Smith   ierr = VecDestroy(&cp->work);CHKERRQ(ierr);
11224c02a0fSBarry Smith   ierr = PetscFree3(cp->a,cp->i,cp->j);CHKERRQ(ierr);
11369d2c0f9SBarry Smith   PetscFunctionReturn(0);
11469d2c0f9SBarry Smith }
11569d2c0f9SBarry Smith 
11669d2c0f9SBarry Smith #undef __FUNCT__
11769d2c0f9SBarry Smith #define __FUNCT__ "PCDestroy_CP"
11869d2c0f9SBarry Smith static PetscErrorCode PCDestroy_CP(PC pc)
11969d2c0f9SBarry Smith {
12069d2c0f9SBarry Smith   PC_CP          *cp = (PC_CP*)pc->data;
12169d2c0f9SBarry Smith   PetscErrorCode ierr;
12269d2c0f9SBarry Smith 
12369d2c0f9SBarry Smith   PetscFunctionBegin;
12469d2c0f9SBarry Smith   ierr = PCReset_CP(pc);CHKERRQ(ierr);
12569d2c0f9SBarry Smith   ierr = PetscFree(cp->d);CHKERRQ(ierr);
12669d2c0f9SBarry Smith   ierr = PetscFree3(cp->a,cp->i,cp->j);CHKERRQ(ierr);
127c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
12824c02a0fSBarry Smith   PetscFunctionReturn(0);
12924c02a0fSBarry Smith }
13024c02a0fSBarry Smith 
13124c02a0fSBarry Smith #undef __FUNCT__
13224c02a0fSBarry Smith #define __FUNCT__ "PCSetFromOptions_CP"
13324c02a0fSBarry Smith static PetscErrorCode PCSetFromOptions_CP(PC pc)
13424c02a0fSBarry Smith {
13524c02a0fSBarry Smith   PetscFunctionBegin;
13624c02a0fSBarry Smith   PetscFunctionReturn(0);
13724c02a0fSBarry Smith }
13824c02a0fSBarry Smith 
13924c02a0fSBarry Smith 
14024c02a0fSBarry Smith /*MC
14124c02a0fSBarry Smith      PCCP - a "column-projection" preconditioner
14224c02a0fSBarry Smith 
14324c02a0fSBarry Smith      This is a terrible preconditioner and is not recommended, ever!
14424c02a0fSBarry Smith 
14524c02a0fSBarry Smith      Loops over the entries of x computing dx_i to
14624c02a0fSBarry Smith $
14724c02a0fSBarry Smith $        min || b - A(x + dx_i e_i ||_2
14824c02a0fSBarry Smith $        dx_i
14924c02a0fSBarry Smith $
15024c02a0fSBarry Smith $    That is, it changes a single entry of x to minimize the new residual.
15124c02a0fSBarry Smith $   Let A_i represent the ith column of A, then the minimization can be written as
15224c02a0fSBarry Smith $
15324c02a0fSBarry Smith $       min || r - (dx_i) A e_i ||_2
15424c02a0fSBarry Smith $       dx_i
15524c02a0fSBarry Smith $   or   min || r - (dx_i) A_i ||_2
15624c02a0fSBarry Smith $        dx_i
15724c02a0fSBarry Smith $
15824c02a0fSBarry Smith $    take the derivative with respect to dx_i to obtain
15924c02a0fSBarry Smith $        dx_i = (A_i^T A_i)^(-1) A_i^T r
16024c02a0fSBarry Smith $
16124c02a0fSBarry Smith $    This algorithm can be thought of as Gauss-Seidel on the normal equations
16224c02a0fSBarry Smith 
16324c02a0fSBarry Smith     Notes: This proceedure can also be done with block columns or any groups of columns
16424c02a0fSBarry Smith         but this is not coded.
16524c02a0fSBarry Smith 
16624c02a0fSBarry Smith       These "projections" can be done simultaneously for all columns (similar to Jacobi)
16724c02a0fSBarry Smith          or sequentially (similar to Gauss-Seidel/SOR). This is only coded for SOR type.
16824c02a0fSBarry Smith 
16924c02a0fSBarry Smith       This is related to, but not the same as "row projection" methods.
17024c02a0fSBarry Smith 
17124c02a0fSBarry Smith       This is currently coded only for SeqAIJ matrices in sequential (SOR) form.
17228529972SSatish Balay 
17328529972SSatish Balay   Level: intermediate
17424c02a0fSBarry Smith 
17524c02a0fSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PCJACOBI, PCSOR
17624c02a0fSBarry Smith 
17724c02a0fSBarry Smith M*/
17824c02a0fSBarry Smith 
17924c02a0fSBarry Smith #undef __FUNCT__
18024c02a0fSBarry Smith #define __FUNCT__ "PCCreate_CP"
1818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_CP(PC pc)
18224c02a0fSBarry Smith {
18324c02a0fSBarry Smith   PC_CP          *cp;
18424c02a0fSBarry Smith   PetscErrorCode ierr;
18524c02a0fSBarry Smith 
18624c02a0fSBarry Smith   PetscFunctionBegin;
187*b00a9115SJed Brown   ierr     = PetscNewLog(pc,&cp);CHKERRQ(ierr);
18824c02a0fSBarry Smith   pc->data = (void*)cp;
18924c02a0fSBarry Smith 
19024c02a0fSBarry Smith   pc->ops->apply           = PCApply_CP;
19124c02a0fSBarry Smith   pc->ops->applytranspose  = PCApply_CP;
19224c02a0fSBarry Smith   pc->ops->setup           = PCSetUp_CP;
19369d2c0f9SBarry Smith   pc->ops->reset           = PCReset_CP;
19424c02a0fSBarry Smith   pc->ops->destroy         = PCDestroy_CP;
19524c02a0fSBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_CP;
19624c02a0fSBarry Smith   pc->ops->view            = 0;
19724c02a0fSBarry Smith   pc->ops->applyrichardson = 0;
19824c02a0fSBarry Smith   PetscFunctionReturn(0);
19924c02a0fSBarry Smith }
20024c02a0fSBarry Smith 
20124c02a0fSBarry Smith 
202