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