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 16d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_CP(PC pc) 17d71ae5a4SJacob Faibussowitsch { 1824c02a0fSBarry Smith PC_CP *cp = (PC_CP *)pc->data; 1924c02a0fSBarry Smith PetscInt i, j, *colcnt; 20ace3abfcSBarry Smith PetscBool flg; 2124c02a0fSBarry Smith Mat_SeqAIJ *aij = (Mat_SeqAIJ *)pc->pmat->data; 2224c02a0fSBarry Smith 2324c02a0fSBarry Smith PetscFunctionBegin; 249566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATSEQAIJ, &flg)); 2528b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Currently only handles SeqAIJ matrices"); 2624c02a0fSBarry Smith 279566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(pc->pmat, &cp->m, &cp->n)); 2808401ef6SPierre Jolivet PetscCheck(cp->m == cp->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Currently only for square matrices"); 2924c02a0fSBarry Smith 309566063dSJacob Faibussowitsch if (!cp->work) PetscCall(MatCreateVecs(pc->pmat, &cp->work, NULL)); 319566063dSJacob Faibussowitsch if (!cp->d) PetscCall(PetscMalloc1(cp->n, &cp->d)); 3224c02a0fSBarry Smith if (cp->a && pc->flag != SAME_NONZERO_PATTERN) { 339566063dSJacob Faibussowitsch PetscCall(PetscFree3(cp->a, cp->i, cp->j)); 340a545947SLisandro Dalcin cp->a = NULL; 3524c02a0fSBarry Smith } 3624c02a0fSBarry Smith 3724c02a0fSBarry Smith /* convert to column format */ 3848a46eb9SPierre Jolivet if (!cp->a) PetscCall(PetscMalloc3(aij->nz, &cp->a, cp->n + 1, &cp->i, aij->nz, &cp->j)); 399566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(cp->n, &colcnt)); 4024c02a0fSBarry Smith 412fa5cd67SKarl Rupp for (i = 0; i < aij->nz; i++) colcnt[aij->j[i]]++; 42e60cf9a0SBarry Smith cp->i[0] = 0; 432fa5cd67SKarl Rupp for (i = 0; i < cp->n; i++) cp->i[i + 1] = cp->i[i] + colcnt[i]; 449566063dSJacob Faibussowitsch PetscCall(PetscArrayzero(colcnt, cp->n)); 4524c02a0fSBarry Smith for (i = 0; i < cp->m; i++) { /* over rows */ 4624c02a0fSBarry Smith for (j = aij->i[i]; j < aij->i[i + 1]; j++) { /* over columns in row */ 4724c02a0fSBarry Smith cp->j[cp->i[aij->j[j]] + colcnt[aij->j[j]]] = i; 4824c02a0fSBarry Smith cp->a[cp->i[aij->j[j]] + colcnt[aij->j[j]]++] = aij->a[j]; 4924c02a0fSBarry Smith } 5024c02a0fSBarry Smith } 519566063dSJacob Faibussowitsch PetscCall(PetscFree(colcnt)); 5224c02a0fSBarry Smith 5324c02a0fSBarry Smith /* compute sum of squares of each column d[] */ 5424c02a0fSBarry Smith for (i = 0; i < cp->n; i++) { /* over columns */ 5575567043SBarry Smith cp->d[i] = 0.; 562fa5cd67SKarl 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 */ 5724c02a0fSBarry Smith cp->d[i] = 1.0 / cp->d[i]; 5824c02a0fSBarry Smith } 59*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6024c02a0fSBarry Smith } 61f1580f4eSBarry Smith 62d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_CP(PC pc, Vec bb, Vec xx) 63d71ae5a4SJacob Faibussowitsch { 6424c02a0fSBarry Smith PC_CP *cp = (PC_CP *)pc->data; 6524c02a0fSBarry Smith PetscScalar *b, *x, xt; 6624c02a0fSBarry Smith PetscInt i, j; 6724c02a0fSBarry Smith 6824c02a0fSBarry Smith PetscFunctionBegin; 699566063dSJacob Faibussowitsch PetscCall(VecCopy(bb, cp->work)); 709566063dSJacob Faibussowitsch PetscCall(VecGetArray(cp->work, &b)); 719566063dSJacob Faibussowitsch PetscCall(VecGetArray(xx, &x)); 7224c02a0fSBarry Smith 7324c02a0fSBarry Smith for (i = 0; i < cp->n; i++) { /* over columns */ 7475567043SBarry Smith xt = 0.; 752fa5cd67SKarl Rupp for (j = cp->i[i]; j < cp->i[i + 1]; j++) xt += cp->a[j] * b[cp->j[j]]; /* over rows in column */ 7624c02a0fSBarry Smith xt *= cp->d[i]; 7724c02a0fSBarry Smith x[i] = xt; 782fa5cd67SKarl 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*/ 7924c02a0fSBarry Smith } 8024c02a0fSBarry Smith for (i = cp->n - 1; i > -1; i--) { /* over columns */ 8175567043SBarry Smith xt = 0.; 822fa5cd67SKarl Rupp for (j = cp->i[i]; j < cp->i[i + 1]; j++) xt += cp->a[j] * b[cp->j[j]]; /* over rows in column */ 8324c02a0fSBarry Smith xt *= cp->d[i]; 8424c02a0fSBarry Smith x[i] = xt; 852fa5cd67SKarl 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*/ 8624c02a0fSBarry Smith } 8724c02a0fSBarry Smith 889566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(cp->work, &b)); 899566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(xx, &x)); 90*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9124c02a0fSBarry Smith } 92f1580f4eSBarry Smith 93d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_CP(PC pc) 94d71ae5a4SJacob Faibussowitsch { 9524c02a0fSBarry Smith PC_CP *cp = (PC_CP *)pc->data; 9624c02a0fSBarry Smith 9724c02a0fSBarry Smith PetscFunctionBegin; 989566063dSJacob Faibussowitsch PetscCall(PetscFree(cp->d)); 999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&cp->work)); 1009566063dSJacob Faibussowitsch PetscCall(PetscFree3(cp->a, cp->i, cp->j)); 101*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10269d2c0f9SBarry Smith } 10369d2c0f9SBarry Smith 104d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_CP(PC pc) 105d71ae5a4SJacob Faibussowitsch { 10669d2c0f9SBarry Smith PC_CP *cp = (PC_CP *)pc->data; 10769d2c0f9SBarry Smith 10869d2c0f9SBarry Smith PetscFunctionBegin; 1099566063dSJacob Faibussowitsch PetscCall(PCReset_CP(pc)); 1109566063dSJacob Faibussowitsch PetscCall(PetscFree(cp->d)); 1119566063dSJacob Faibussowitsch PetscCall(PetscFree3(cp->a, cp->i, cp->j)); 1129566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 113*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11424c02a0fSBarry Smith } 11524c02a0fSBarry Smith 116d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_CP(PC pc, PetscOptionItems *PetscOptionsObject) 117d71ae5a4SJacob Faibussowitsch { 11824c02a0fSBarry Smith PetscFunctionBegin; 119*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12024c02a0fSBarry Smith } 12124c02a0fSBarry Smith 12224c02a0fSBarry Smith /*MC 12324c02a0fSBarry Smith PCCP - a "column-projection" preconditioner 12424c02a0fSBarry Smith 12524c02a0fSBarry Smith This is a terrible preconditioner and is not recommended, ever! 12624c02a0fSBarry Smith 12773a58da7SBarry Smith Loops over the entries of x computing dx_i (e_i is the unit vector in the ith direction) to 128f1580f4eSBarry Smith .vb 129f1580f4eSBarry Smith 130f1580f4eSBarry Smith min || b - A(x + dx_i e_i ||_2 131f1580f4eSBarry Smith dx_i 132f1580f4eSBarry Smith 133f1580f4eSBarry Smith That is, it changes a single entry of x to minimize the new residual norm. 134f1580f4eSBarry Smith Let A_i represent the ith column of A, then the minimization can be written as 135f1580f4eSBarry Smith 136f1580f4eSBarry Smith min || r - (dx_i) A e_i ||_2 137f1580f4eSBarry Smith dx_i 138f1580f4eSBarry Smith or min || r - (dx_i) A_i ||_2 139f1580f4eSBarry Smith dx_i 140f1580f4eSBarry Smith 141f1580f4eSBarry Smith take the derivative with respect to dx_i to obtain 142f1580f4eSBarry Smith dx_i = (A_i^T A_i)^(-1) A_i^T r 143f1580f4eSBarry Smith 144f1580f4eSBarry Smith This algorithm can be thought of as Gauss-Seidel on the normal equations 145f1580f4eSBarry Smith .ve 14624c02a0fSBarry Smith 14795452b02SPatrick Sanan Notes: 14895452b02SPatrick Sanan This proceedure can also be done with block columns or any groups of columns 14924c02a0fSBarry Smith but this is not coded. 15024c02a0fSBarry Smith 15124c02a0fSBarry Smith These "projections" can be done simultaneously for all columns (similar to Jacobi) 15224c02a0fSBarry Smith or sequentially (similar to Gauss-Seidel/SOR). This is only coded for SOR type. 15324c02a0fSBarry Smith 15424c02a0fSBarry Smith This is related to, but not the same as "row projection" methods. 15524c02a0fSBarry Smith 156f1580f4eSBarry Smith This is currently coded only for `MATSEQAIJ` matrices 15728529972SSatish Balay 15828529972SSatish Balay Level: intermediate 15924c02a0fSBarry Smith 160db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PCJACOBI`, `PCSOR` 16124c02a0fSBarry Smith M*/ 16224c02a0fSBarry Smith 163d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_CP(PC pc) 164d71ae5a4SJacob Faibussowitsch { 16524c02a0fSBarry Smith PC_CP *cp; 16624c02a0fSBarry Smith 16724c02a0fSBarry Smith PetscFunctionBegin; 1684dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&cp)); 16924c02a0fSBarry Smith pc->data = (void *)cp; 17024c02a0fSBarry Smith 17124c02a0fSBarry Smith pc->ops->apply = PCApply_CP; 17224c02a0fSBarry Smith pc->ops->applytranspose = PCApply_CP; 17324c02a0fSBarry Smith pc->ops->setup = PCSetUp_CP; 17469d2c0f9SBarry Smith pc->ops->reset = PCReset_CP; 17524c02a0fSBarry Smith pc->ops->destroy = PCDestroy_CP; 17624c02a0fSBarry Smith pc->ops->setfromoptions = PCSetFromOptions_CP; 1770a545947SLisandro Dalcin pc->ops->view = NULL; 1780a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 179*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18024c02a0fSBarry Smith } 181