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