xref: /petsc/src/ksp/pc/impls/cp/cp.c (revision 5b6bfdb9644f185dbf5e5a09b808ec241507e1e7)
1 
2 #include <petsc/private/pcimpl.h>   /*I "petscpc.h" I*/
3 #include <../src/mat/impls/aij/seq/aij.h>
4 
5 /*
6    Private context (data structure) for the CP preconditioner.
7 */
8 typedef struct {
9   PetscInt    n,m;
10   Vec         work;
11   PetscScalar *d;       /* sum of squares of each column */
12   PetscScalar *a;       /* non-zeros by column */
13   PetscInt    *i,*j;    /* offsets of nonzeros by column, non-zero indices by column */
14 } PC_CP;
15 
16 
17 static PetscErrorCode PCSetUp_CP(PC pc)
18 {
19   PC_CP          *cp = (PC_CP*)pc->data;
20   PetscInt       i,j,*colcnt;
21   PetscErrorCode ierr;
22   PetscBool      flg;
23   Mat_SeqAIJ     *aij = (Mat_SeqAIJ*)pc->pmat->data;
24 
25   PetscFunctionBegin;
26   ierr = PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg);CHKERRQ(ierr);
27   if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Currently only handles SeqAIJ matrices");
28 
29   ierr = MatGetLocalSize(pc->pmat,&cp->m,&cp->n);CHKERRQ(ierr);
30   if (cp->m != cp->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently only for square matrices");
31 
32   if (!cp->work) {ierr = MatCreateVecs(pc->pmat,&cp->work,NULL);CHKERRQ(ierr);}
33   if (!cp->d) {ierr = PetscMalloc1(cp->n,&cp->d);CHKERRQ(ierr);}
34   if (cp->a && pc->flag != SAME_NONZERO_PATTERN) {
35     ierr  = PetscFree3(cp->a,cp->i,cp->j);CHKERRQ(ierr);
36     cp->a = 0;
37   }
38 
39   /* convert to column format */
40   if (!cp->a) {
41     ierr = PetscMalloc3(aij->nz,&cp->a,cp->n+1,&cp->i,aij->nz,&cp->j);CHKERRQ(ierr);
42   }
43   ierr = PetscCalloc1(cp->n,&colcnt);CHKERRQ(ierr);
44 
45   for (i=0; i<aij->nz; i++) colcnt[aij->j[i]]++;
46   cp->i[0] = 0;
47   for (i=0; i<cp->n; i++) cp->i[i+1] = cp->i[i] + colcnt[i];
48   ierr = PetscMemzero(colcnt,cp->n*sizeof(PetscInt));CHKERRQ(ierr);
49   for (i=0; i<cp->m; i++) {  /* over rows */
50     for (j=aij->i[i]; j<aij->i[i+1]; j++) {  /* over columns in row */
51       cp->j[cp->i[aij->j[j]]+colcnt[aij->j[j]]]   = i;
52       cp->a[cp->i[aij->j[j]]+colcnt[aij->j[j]]++] = aij->a[j];
53     }
54   }
55   ierr = PetscFree(colcnt);CHKERRQ(ierr);
56 
57   /* compute sum of squares of each column d[] */
58   for (i=0; i<cp->n; i++) {  /* over columns */
59     cp->d[i] = 0.;
60     for (j=cp->i[i]; j<cp->i[i+1]; j++) cp->d[i] += cp->a[j]*cp->a[j]; /* over rows in column */
61     cp->d[i] = 1.0/cp->d[i];
62   }
63   PetscFunctionReturn(0);
64 }
65 /* -------------------------------------------------------------------------- */
66 static PetscErrorCode PCApply_CP(PC pc,Vec bb,Vec xx)
67 {
68   PC_CP          *cp = (PC_CP*)pc->data;
69   PetscErrorCode ierr;
70   PetscScalar    *b,*x,xt;
71   PetscInt       i,j;
72 
73   PetscFunctionBegin;
74   ierr = VecCopy(bb,cp->work);CHKERRQ(ierr);
75   ierr = VecGetArray(cp->work,&b);CHKERRQ(ierr);
76   ierr = VecGetArray(xx,&x);CHKERRQ(ierr);
77 
78   for (i=0; i<cp->n; i++) {  /* over columns */
79     xt = 0.;
80     for (j=cp->i[i]; j<cp->i[i+1]; j++) xt += cp->a[j]*b[cp->j[j]]; /* over rows in column */
81     xt  *= cp->d[i];
82     x[i] = xt;
83     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*/
84   }
85   for (i=cp->n-1; i>-1; i--) {  /* over columns */
86     xt = 0.;
87     for (j=cp->i[i]; j<cp->i[i+1]; j++) xt += cp->a[j]*b[cp->j[j]]; /* over rows in column */
88     xt  *= cp->d[i];
89     x[i] = xt;
90     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*/
91   }
92 
93   ierr = VecRestoreArray(cp->work,&b);CHKERRQ(ierr);
94   ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr);
95   PetscFunctionReturn(0);
96 }
97 /* -------------------------------------------------------------------------- */
98 static PetscErrorCode PCReset_CP(PC pc)
99 {
100   PC_CP          *cp = (PC_CP*)pc->data;
101   PetscErrorCode ierr;
102 
103   PetscFunctionBegin;
104   ierr = PetscFree(cp->d);CHKERRQ(ierr);
105   ierr = VecDestroy(&cp->work);CHKERRQ(ierr);
106   ierr = PetscFree3(cp->a,cp->i,cp->j);CHKERRQ(ierr);
107   PetscFunctionReturn(0);
108 }
109 
110 static PetscErrorCode PCDestroy_CP(PC pc)
111 {
112   PC_CP          *cp = (PC_CP*)pc->data;
113   PetscErrorCode ierr;
114 
115   PetscFunctionBegin;
116   ierr = PCReset_CP(pc);CHKERRQ(ierr);
117   ierr = PetscFree(cp->d);CHKERRQ(ierr);
118   ierr = PetscFree3(cp->a,cp->i,cp->j);CHKERRQ(ierr);
119   ierr = PetscFree(pc->data);CHKERRQ(ierr);
120   PetscFunctionReturn(0);
121 }
122 
123 static PetscErrorCode PCSetFromOptions_CP(PetscOptionItems *PetscOptionsObject,PC pc)
124 {
125   PetscFunctionBegin;
126   PetscFunctionReturn(0);
127 }
128 
129 /*MC
130      PCCP - a "column-projection" preconditioner
131 
132      This is a terrible preconditioner and is not recommended, ever!
133 
134      Loops over the entries of x computing dx_i (e_i is the unit vector in the ith direction) to
135 $
136 $        min || b - A(x + dx_i e_i ||_2
137 $        dx_i
138 $
139 $    That is, it changes a single entry of x to minimize the new residual norm.
140 $   Let A_i represent the ith column of A, then the minimization can be written as
141 $
142 $       min || r - (dx_i) A e_i ||_2
143 $       dx_i
144 $   or   min || r - (dx_i) A_i ||_2
145 $        dx_i
146 $
147 $    take the derivative with respect to dx_i to obtain
148 $        dx_i = (A_i^T A_i)^(-1) A_i^T r
149 $
150 $    This algorithm can be thought of as Gauss-Seidel on the normal equations
151 
152     Notes: This proceedure can also be done with block columns or any groups of columns
153         but this is not coded.
154 
155       These "projections" can be done simultaneously for all columns (similar to Jacobi)
156          or sequentially (similar to Gauss-Seidel/SOR). This is only coded for SOR type.
157 
158       This is related to, but not the same as "row projection" methods.
159 
160       This is currently coded only for SeqAIJ matrices in sequential (SOR) form.
161 
162   Level: intermediate
163 
164 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PCJACOBI, PCSOR
165 
166 M*/
167 
168 PETSC_EXTERN PetscErrorCode PCCreate_CP(PC pc)
169 {
170   PC_CP          *cp;
171   PetscErrorCode ierr;
172 
173   PetscFunctionBegin;
174   ierr     = PetscNewLog(pc,&cp);CHKERRQ(ierr);
175   pc->data = (void*)cp;
176 
177   pc->ops->apply           = PCApply_CP;
178   pc->ops->applytranspose  = PCApply_CP;
179   pc->ops->setup           = PCSetUp_CP;
180   pc->ops->reset           = PCReset_CP;
181   pc->ops->destroy         = PCDestroy_CP;
182   pc->ops->setfromoptions  = PCSetFromOptions_CP;
183   pc->ops->view            = 0;
184   pc->ops->applyrichardson = 0;
185   PetscFunctionReturn(0);
186 }
187 
188 
189