xref: /petsc/src/ksp/pc/impls/cp/cp.c (revision 24676f2a9731969d4e2872adbcfa46b47a966276)
1af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/
2c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>
324c02a0fSBarry Smith 
424c02a0fSBarry Smith /*
524c02a0fSBarry Smith    Private context (data structure) for the CP preconditioner.
624c02a0fSBarry Smith */
724c02a0fSBarry Smith typedef struct {
824c02a0fSBarry Smith   PetscInt     n, m;
924c02a0fSBarry Smith   Vec          work;
1024c02a0fSBarry Smith   PetscScalar *d;     /* sum of squares of each column */
1124c02a0fSBarry Smith   PetscScalar *a;     /* non-zeros by column */
1224c02a0fSBarry Smith   PetscInt    *i, *j; /* offsets of nonzeros by column, non-zero indices by column */
1324c02a0fSBarry Smith } PC_CP;
1424c02a0fSBarry Smith 
PCSetUp_CP(PC pc)15d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_CP(PC pc)
16d71ae5a4SJacob Faibussowitsch {
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   }
583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5924c02a0fSBarry Smith }
60f1580f4eSBarry Smith 
PCApply_CP(PC pc,Vec bb,Vec xx)61d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_CP(PC pc, Vec bb, Vec xx)
62d71ae5a4SJacob Faibussowitsch {
6324c02a0fSBarry Smith   PC_CP       *cp = (PC_CP *)pc->data;
6424c02a0fSBarry Smith   PetscScalar *b, *x, xt;
6524c02a0fSBarry Smith   PetscInt     i, j;
6624c02a0fSBarry Smith 
6724c02a0fSBarry Smith   PetscFunctionBegin;
689566063dSJacob Faibussowitsch   PetscCall(VecCopy(bb, cp->work));
699566063dSJacob Faibussowitsch   PetscCall(VecGetArray(cp->work, &b));
709566063dSJacob Faibussowitsch   PetscCall(VecGetArray(xx, &x));
7124c02a0fSBarry Smith 
7224c02a0fSBarry Smith   for (i = 0; i < cp->n; i++) { /* over columns */
7375567043SBarry Smith     xt = 0.;
742fa5cd67SKarl Rupp     for (j = cp->i[i]; j < cp->i[i + 1]; j++) xt += cp->a[j] * b[cp->j[j]]; /* over rows in column */
7524c02a0fSBarry Smith     xt *= cp->d[i];
7624c02a0fSBarry Smith     x[i] = xt;
772fa5cd67SKarl 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*/
7824c02a0fSBarry Smith   }
7924c02a0fSBarry Smith   for (i = cp->n - 1; i > -1; i--) { /* over columns */
8075567043SBarry Smith     xt = 0.;
812fa5cd67SKarl Rupp     for (j = cp->i[i]; j < cp->i[i + 1]; j++) xt += cp->a[j] * b[cp->j[j]]; /* over rows in column */
8224c02a0fSBarry Smith     xt *= cp->d[i];
8324c02a0fSBarry Smith     x[i] = xt;
842fa5cd67SKarl 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*/
8524c02a0fSBarry Smith   }
8624c02a0fSBarry Smith 
879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(cp->work, &b));
889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(xx, &x));
893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9024c02a0fSBarry Smith }
91f1580f4eSBarry Smith 
PCReset_CP(PC pc)92d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_CP(PC pc)
93d71ae5a4SJacob Faibussowitsch {
9424c02a0fSBarry Smith   PC_CP *cp = (PC_CP *)pc->data;
9524c02a0fSBarry Smith 
9624c02a0fSBarry Smith   PetscFunctionBegin;
979566063dSJacob Faibussowitsch   PetscCall(PetscFree(cp->d));
989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&cp->work));
999566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cp->a, cp->i, cp->j));
1003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10169d2c0f9SBarry Smith }
10269d2c0f9SBarry Smith 
PCDestroy_CP(PC pc)103d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_CP(PC pc)
104d71ae5a4SJacob Faibussowitsch {
10569d2c0f9SBarry Smith   PC_CP *cp = (PC_CP *)pc->data;
10669d2c0f9SBarry Smith 
10769d2c0f9SBarry Smith   PetscFunctionBegin;
1089566063dSJacob Faibussowitsch   PetscCall(PCReset_CP(pc));
1099566063dSJacob Faibussowitsch   PetscCall(PetscFree(cp->d));
1109566063dSJacob Faibussowitsch   PetscCall(PetscFree3(cp->a, cp->i, cp->j));
1119566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
1123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
11324c02a0fSBarry Smith }
11424c02a0fSBarry Smith 
11524c02a0fSBarry Smith /*MC
116*0b4b7b1cSBarry Smith      PCCP - a "column-projection" preconditioner. Iteratively projects the current residual onto the one dimensional spaces
117*0b4b7b1cSBarry Smith             spanned by each of the columns of the matrix.
11824c02a0fSBarry Smith 
11924c02a0fSBarry Smith      This is a terrible preconditioner and is not recommended, ever!
12024c02a0fSBarry Smith 
12173a58da7SBarry Smith      Loops over the entries of x computing dx_i (e_i is the unit vector in the ith direction) to
122f1580f4eSBarry Smith .vb
123f1580f4eSBarry Smith 
124f1580f4eSBarry Smith         min || b - A(x + dx_i e_i ||_2
125f1580f4eSBarry Smith         dx_i
126f1580f4eSBarry Smith 
127f1580f4eSBarry Smith     That is, it changes a single entry of x to minimize the new residual norm.
128f1580f4eSBarry Smith    Let A_i represent the ith column of A, then the minimization can be written as
129f1580f4eSBarry Smith 
130f1580f4eSBarry Smith        min || r - (dx_i) A e_i ||_2
131f1580f4eSBarry Smith        dx_i
132f1580f4eSBarry Smith    or   min || r - (dx_i) A_i ||_2
133f1580f4eSBarry Smith         dx_i
134f1580f4eSBarry Smith 
135f1580f4eSBarry Smith     take the derivative with respect to dx_i to obtain
136f1580f4eSBarry Smith         dx_i = (A_i^T A_i)^(-1) A_i^T r
137f1580f4eSBarry Smith 
138*0b4b7b1cSBarry Smith     This is equivalent to using Gauss-Seidel on the normal equations
139f1580f4eSBarry Smith .ve
14024c02a0fSBarry Smith 
14195452b02SPatrick Sanan     Notes:
142da81f932SPierre Jolivet     This procedure can also be done with block columns or any groups of columns
14324c02a0fSBarry Smith     but this is not coded.
14424c02a0fSBarry Smith 
145*0b4b7b1cSBarry Smith     These "projections" can be done simultaneously for all columns (similar to the Jacobi method)
14624c02a0fSBarry Smith     or sequentially (similar to Gauss-Seidel/SOR). This is only coded for SOR type.
14724c02a0fSBarry Smith 
14824c02a0fSBarry Smith     This is related to, but not the same as "row projection" methods.
14924c02a0fSBarry Smith 
150f1580f4eSBarry Smith     This is currently coded only for `MATSEQAIJ` matrices
15128529972SSatish Balay 
15228529972SSatish Balay   Level: intermediate
15324c02a0fSBarry Smith 
154562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PCJACOBI`, `PCSOR`
15524c02a0fSBarry Smith M*/
15624c02a0fSBarry Smith 
PCCreate_CP(PC pc)157d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_CP(PC pc)
158d71ae5a4SJacob Faibussowitsch {
15924c02a0fSBarry Smith   PC_CP *cp;
16024c02a0fSBarry Smith 
16124c02a0fSBarry Smith   PetscFunctionBegin;
1624dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&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;
1700a545947SLisandro Dalcin   pc->ops->view            = NULL;
1710a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
1723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17324c02a0fSBarry Smith }
174