xref: /petsc/src/ksp/pc/impls/cp/cp.c (revision f1580f4e3ce5d5b2393648fd039d0d41b440385d)
124c02a0fSBarry Smith 
2af0996ceSBarry Smith #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 
169371c9d4SSatish Balay static PetscErrorCode PCSetUp_CP(PC pc) {
1724c02a0fSBarry Smith   PC_CP      *cp = (PC_CP *)pc->data;
1824c02a0fSBarry Smith   PetscInt    i, j, *colcnt;
19ace3abfcSBarry Smith   PetscBool   flg;
2024c02a0fSBarry Smith   Mat_SeqAIJ *aij = (Mat_SeqAIJ *)pc->pmat->data;
2124c02a0fSBarry Smith 
2224c02a0fSBarry Smith   PetscFunctionBegin;
239566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATSEQAIJ, &flg));
2428b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Currently only handles SeqAIJ matrices");
2524c02a0fSBarry Smith 
269566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(pc->pmat, &cp->m, &cp->n));
2708401ef6SPierre Jolivet   PetscCheck(cp->m == cp->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Currently only for square matrices");
2824c02a0fSBarry Smith 
299566063dSJacob Faibussowitsch   if (!cp->work) PetscCall(MatCreateVecs(pc->pmat, &cp->work, NULL));
309566063dSJacob Faibussowitsch   if (!cp->d) PetscCall(PetscMalloc1(cp->n, &cp->d));
3124c02a0fSBarry Smith   if (cp->a && pc->flag != SAME_NONZERO_PATTERN) {
329566063dSJacob Faibussowitsch     PetscCall(PetscFree3(cp->a, cp->i, cp->j));
330a545947SLisandro Dalcin     cp->a = NULL;
3424c02a0fSBarry Smith   }
3524c02a0fSBarry Smith 
3624c02a0fSBarry Smith   /* convert to column format */
3748a46eb9SPierre Jolivet   if (!cp->a) PetscCall(PetscMalloc3(aij->nz, &cp->a, cp->n + 1, &cp->i, aij->nz, &cp->j));
389566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(cp->n, &colcnt));
3924c02a0fSBarry Smith 
402fa5cd67SKarl Rupp   for (i = 0; i < aij->nz; i++) colcnt[aij->j[i]]++;
41e60cf9a0SBarry Smith   cp->i[0] = 0;
422fa5cd67SKarl Rupp   for (i = 0; i < cp->n; i++) cp->i[i + 1] = cp->i[i] + colcnt[i];
439566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(colcnt, cp->n));
4424c02a0fSBarry Smith   for (i = 0; i < cp->m; i++) {                   /* over rows */
4524c02a0fSBarry Smith     for (j = aij->i[i]; j < aij->i[i + 1]; j++) { /* over columns in row */
4624c02a0fSBarry Smith       cp->j[cp->i[aij->j[j]] + colcnt[aij->j[j]]]   = i;
4724c02a0fSBarry Smith       cp->a[cp->i[aij->j[j]] + colcnt[aij->j[j]]++] = aij->a[j];
4824c02a0fSBarry Smith     }
4924c02a0fSBarry Smith   }
509566063dSJacob Faibussowitsch   PetscCall(PetscFree(colcnt));
5124c02a0fSBarry Smith 
5224c02a0fSBarry Smith   /* compute sum of squares of each column d[] */
5324c02a0fSBarry Smith   for (i = 0; i < cp->n; i++) { /* over columns */
5475567043SBarry Smith     cp->d[i] = 0.;
552fa5cd67SKarl 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 */
5624c02a0fSBarry Smith     cp->d[i] = 1.0 / cp->d[i];
5724c02a0fSBarry Smith   }
5824c02a0fSBarry Smith   PetscFunctionReturn(0);
5924c02a0fSBarry Smith }
60*f1580f4eSBarry Smith 
619371c9d4SSatish Balay static PetscErrorCode PCApply_CP(PC pc, Vec bb, Vec xx) {
6224c02a0fSBarry Smith   PC_CP       *cp = (PC_CP *)pc->data;
6324c02a0fSBarry Smith   PetscScalar *b, *x, xt;
6424c02a0fSBarry Smith   PetscInt     i, j;
6524c02a0fSBarry Smith 
6624c02a0fSBarry Smith   PetscFunctionBegin;
679566063dSJacob Faibussowitsch   PetscCall(VecCopy(bb, cp->work));
689566063dSJacob Faibussowitsch   PetscCall(VecGetArray(cp->work, &b));
699566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
7024c02a0fSBarry Smith 
7124c02a0fSBarry Smith   for (i = 0; i < cp->n; i++) { /* over columns */
7275567043SBarry Smith     xt = 0.;
732fa5cd67SKarl Rupp     for (j = cp->i[i]; j < cp->i[i + 1]; j++) xt += cp->a[j] * b[cp->j[j]]; /* over rows in column */
7424c02a0fSBarry Smith     xt *= cp->d[i];
7524c02a0fSBarry Smith     x[i] = xt;
762fa5cd67SKarl 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*/
7724c02a0fSBarry Smith   }
7824c02a0fSBarry Smith   for (i = cp->n - 1; i > -1; i--) { /* over columns */
7975567043SBarry Smith     xt = 0.;
802fa5cd67SKarl Rupp     for (j = cp->i[i]; j < cp->i[i + 1]; j++) xt += cp->a[j] * b[cp->j[j]]; /* over rows in column */
8124c02a0fSBarry Smith     xt *= cp->d[i];
8224c02a0fSBarry Smith     x[i] = xt;
832fa5cd67SKarl 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*/
8424c02a0fSBarry Smith   }
8524c02a0fSBarry Smith 
869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(cp->work, &b));
879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
8824c02a0fSBarry Smith   PetscFunctionReturn(0);
8924c02a0fSBarry Smith }
90*f1580f4eSBarry Smith 
919371c9d4SSatish Balay static PetscErrorCode PCReset_CP(PC pc) {
9224c02a0fSBarry Smith   PC_CP *cp = (PC_CP *)pc->data;
9324c02a0fSBarry Smith 
9424c02a0fSBarry Smith   PetscFunctionBegin;
959566063dSJacob Faibussowitsch   PetscCall(PetscFree(cp->d));
969566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&cp->work));
979566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cp->a, cp->i, cp->j));
9869d2c0f9SBarry Smith   PetscFunctionReturn(0);
9969d2c0f9SBarry Smith }
10069d2c0f9SBarry Smith 
1019371c9d4SSatish Balay static PetscErrorCode PCDestroy_CP(PC pc) {
10269d2c0f9SBarry Smith   PC_CP *cp = (PC_CP *)pc->data;
10369d2c0f9SBarry Smith 
10469d2c0f9SBarry Smith   PetscFunctionBegin;
1059566063dSJacob Faibussowitsch   PetscCall(PCReset_CP(pc));
1069566063dSJacob Faibussowitsch   PetscCall(PetscFree(cp->d));
1079566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cp->a, cp->i, cp->j));
1089566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
10924c02a0fSBarry Smith   PetscFunctionReturn(0);
11024c02a0fSBarry Smith }
11124c02a0fSBarry Smith 
1129371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_CP(PC pc, PetscOptionItems *PetscOptionsObject) {
11324c02a0fSBarry Smith   PetscFunctionBegin;
11424c02a0fSBarry Smith   PetscFunctionReturn(0);
11524c02a0fSBarry Smith }
11624c02a0fSBarry Smith 
11724c02a0fSBarry Smith /*MC
11824c02a0fSBarry Smith      PCCP - a "column-projection" preconditioner
11924c02a0fSBarry Smith 
12024c02a0fSBarry Smith      This is a terrible preconditioner and is not recommended, ever!
12124c02a0fSBarry Smith 
12273a58da7SBarry Smith      Loops over the entries of x computing dx_i (e_i is the unit vector in the ith direction) to
123*f1580f4eSBarry Smith .vb
124*f1580f4eSBarry Smith 
125*f1580f4eSBarry Smith         min || b - A(x + dx_i e_i ||_2
126*f1580f4eSBarry Smith         dx_i
127*f1580f4eSBarry Smith 
128*f1580f4eSBarry Smith     That is, it changes a single entry of x to minimize the new residual norm.
129*f1580f4eSBarry Smith    Let A_i represent the ith column of A, then the minimization can be written as
130*f1580f4eSBarry Smith 
131*f1580f4eSBarry Smith        min || r - (dx_i) A e_i ||_2
132*f1580f4eSBarry Smith        dx_i
133*f1580f4eSBarry Smith    or   min || r - (dx_i) A_i ||_2
134*f1580f4eSBarry Smith         dx_i
135*f1580f4eSBarry Smith 
136*f1580f4eSBarry Smith     take the derivative with respect to dx_i to obtain
137*f1580f4eSBarry Smith         dx_i = (A_i^T A_i)^(-1) A_i^T r
138*f1580f4eSBarry Smith 
139*f1580f4eSBarry Smith     This algorithm can be thought of as Gauss-Seidel on the normal equations
140*f1580f4eSBarry Smith .ve
14124c02a0fSBarry Smith 
14295452b02SPatrick Sanan     Notes:
14395452b02SPatrick Sanan     This proceedure can also be done with block columns or any groups of columns
14424c02a0fSBarry Smith     but this is not coded.
14524c02a0fSBarry Smith 
14624c02a0fSBarry Smith     These "projections" can be done simultaneously for all columns (similar to Jacobi)
14724c02a0fSBarry Smith     or sequentially (similar to Gauss-Seidel/SOR). This is only coded for SOR type.
14824c02a0fSBarry Smith 
14924c02a0fSBarry Smith     This is related to, but not the same as "row projection" methods.
15024c02a0fSBarry Smith 
151*f1580f4eSBarry Smith     This is currently coded only for `MATSEQAIJ` matrices
15228529972SSatish Balay 
15328529972SSatish Balay   Level: intermediate
15424c02a0fSBarry Smith 
155db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PCJACOBI`, `PCSOR`
15624c02a0fSBarry Smith M*/
15724c02a0fSBarry Smith 
1589371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_CP(PC pc) {
15924c02a0fSBarry Smith   PC_CP *cp;
16024c02a0fSBarry Smith 
16124c02a0fSBarry Smith   PetscFunctionBegin;
1629566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc, &cp));
16324c02a0fSBarry Smith   pc->data = (void *)cp;
16424c02a0fSBarry Smith 
16524c02a0fSBarry Smith   pc->ops->apply           = PCApply_CP;
16624c02a0fSBarry Smith   pc->ops->applytranspose  = PCApply_CP;
16724c02a0fSBarry Smith   pc->ops->setup           = PCSetUp_CP;
16869d2c0f9SBarry Smith   pc->ops->reset           = PCReset_CP;
16924c02a0fSBarry Smith   pc->ops->destroy         = PCDestroy_CP;
17024c02a0fSBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_CP;
1710a545947SLisandro Dalcin   pc->ops->view            = NULL;
1720a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
17324c02a0fSBarry Smith   PetscFunctionReturn(0);
17424c02a0fSBarry Smith }
175