2: #include <petsc/private/pcimpl.h> 3: #include <../src/mat/impls/aij/seq/aij.h> 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;
17: static PetscErrorCode PCSetUp_CP(PC pc) 18: {
19: PC_CP *cp = (PC_CP*)pc->data;
20: PetscInt i,j,*colcnt;
22: PetscBool flg;
23: Mat_SeqAIJ *aij = (Mat_SeqAIJ*)pc->pmat->data;
26: PetscObjectTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg);
27: if (!flg) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Currently only handles SeqAIJ matrices");
29: MatGetLocalSize(pc->pmat,&cp->m,&cp->n);
30: if (cp->m != cp->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Currently only for square matrices");
32: if (!cp->work) {MatCreateVecs(pc->pmat,&cp->work,NULL);}
33: if (!cp->d) {PetscMalloc1(cp->n,&cp->d);}
34: if (cp->a && pc->flag != SAME_NONZERO_PATTERN) {
35: PetscFree3(cp->a,cp->i,cp->j);
36: cp->a = 0;
37: }
39: /* convert to column format */
40: if (!cp->a) {
41: PetscMalloc3(aij->nz,&cp->a,cp->n+1,&cp->i,aij->nz,&cp->j);
42: }
43: PetscCalloc1(cp->n,&colcnt);
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: PetscMemzero(colcnt,cp->n*sizeof(PetscInt));
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: PetscFree(colcnt);
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: return(0);
64: }
65: /* -------------------------------------------------------------------------- */
66: static PetscErrorCode PCApply_CP(PC pc,Vec bb,Vec xx) 67: {
68: PC_CP *cp = (PC_CP*)pc->data;
70: PetscScalar *b,*x,xt;
71: PetscInt i,j;
74: VecCopy(bb,cp->work);
75: VecGetArray(cp->work,&b);
76: VecGetArray(xx,&x);
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: }
93: VecRestoreArray(cp->work,&b);
94: VecRestoreArray(xx,&x);
95: return(0);
96: }
97: /* -------------------------------------------------------------------------- */
98: static PetscErrorCode PCReset_CP(PC pc) 99: {
100: PC_CP *cp = (PC_CP*)pc->data;
104: PetscFree(cp->d);
105: VecDestroy(&cp->work);
106: PetscFree3(cp->a,cp->i,cp->j);
107: return(0);
108: }
110: static PetscErrorCode PCDestroy_CP(PC pc)111: {
112: PC_CP *cp = (PC_CP*)pc->data;
116: PCReset_CP(pc);
117: PetscFree(cp->d);
118: PetscFree3(cp->a,cp->i,cp->j);
119: PetscFree(pc->data);
120: return(0);
121: }
123: static PetscErrorCode PCSetFromOptions_CP(PetscOptionItems *PetscOptionsObject,PC pc)124: {
126: return(0);
127: }
129: /*MC
130: PCCP - a "column-projection" preconditioner
132: This is a terrible preconditioner and is not recommended, ever!
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
152: Notes:
153: This proceedure can also be done with block columns or any groups of columns
154: but this is not coded.
156: These "projections" can be done simultaneously for all columns (similar to Jacobi)
157: or sequentially (similar to Gauss-Seidel/SOR). This is only coded for SOR type.
159: This is related to, but not the same as "row projection" methods.
161: This is currently coded only for SeqAIJ matrices in sequential (SOR) form.
163: Level: intermediate
165: .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCJACOBI, PCSOR167: M*/
169: PETSC_EXTERN PetscErrorCode PCCreate_CP(PC pc)170: {
171: PC_CP *cp;
175: PetscNewLog(pc,&cp);
176: pc->data = (void*)cp;
178: pc->ops->apply = PCApply_CP;
179: pc->ops->applytranspose = PCApply_CP;
180: pc->ops->setup = PCSetUp_CP;
181: pc->ops->reset = PCReset_CP;
182: pc->ops->destroy = PCDestroy_CP;
183: pc->ops->setfromoptions = PCSetFromOptions_CP;
184: pc->ops->view = 0;
185: pc->ops->applyrichardson = 0;
186: return(0);
187: }