1 2 #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 3 #include <../src/mat/impls/aij/seq/aij.h> 4 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; 15 16 static PetscErrorCode PCSetUp_CP(PC pc) { 17 PC_CP *cp = (PC_CP *)pc->data; 18 PetscInt i, j, *colcnt; 19 PetscBool flg; 20 Mat_SeqAIJ *aij = (Mat_SeqAIJ *)pc->pmat->data; 21 22 PetscFunctionBegin; 23 PetscCall(PetscObjectTypeCompare((PetscObject)pc->pmat, MATSEQAIJ, &flg)); 24 PetscCheck(flg, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Currently only handles SeqAIJ matrices"); 25 26 PetscCall(MatGetLocalSize(pc->pmat, &cp->m, &cp->n)); 27 PetscCheck(cp->m == cp->n, PETSC_COMM_SELF, PETSC_ERR_SUP, "Currently only for square matrices"); 28 29 if (!cp->work) PetscCall(MatCreateVecs(pc->pmat, &cp->work, NULL)); 30 if (!cp->d) PetscCall(PetscMalloc1(cp->n, &cp->d)); 31 if (cp->a && pc->flag != SAME_NONZERO_PATTERN) { 32 PetscCall(PetscFree3(cp->a, cp->i, cp->j)); 33 cp->a = NULL; 34 } 35 36 /* convert to column format */ 37 if (!cp->a) { PetscCall(PetscMalloc3(aij->nz, &cp->a, cp->n + 1, &cp->i, aij->nz, &cp->j)); } 38 PetscCall(PetscCalloc1(cp->n, &colcnt)); 39 40 for (i = 0; i < aij->nz; i++) colcnt[aij->j[i]]++; 41 cp->i[0] = 0; 42 for (i = 0; i < cp->n; i++) cp->i[i + 1] = cp->i[i] + colcnt[i]; 43 PetscCall(PetscArrayzero(colcnt, cp->n)); 44 for (i = 0; i < cp->m; i++) { /* over rows */ 45 for (j = aij->i[i]; j < aij->i[i + 1]; j++) { /* over columns in row */ 46 cp->j[cp->i[aij->j[j]] + colcnt[aij->j[j]]] = i; 47 cp->a[cp->i[aij->j[j]] + colcnt[aij->j[j]]++] = aij->a[j]; 48 } 49 } 50 PetscCall(PetscFree(colcnt)); 51 52 /* compute sum of squares of each column d[] */ 53 for (i = 0; i < cp->n; i++) { /* over columns */ 54 cp->d[i] = 0.; 55 for (j = cp->i[i]; j < cp->i[i + 1]; j++) cp->d[i] += cp->a[j] * cp->a[j]; /* over rows in column */ 56 cp->d[i] = 1.0 / cp->d[i]; 57 } 58 PetscFunctionReturn(0); 59 } 60 /* -------------------------------------------------------------------------- */ 61 static PetscErrorCode PCApply_CP(PC pc, Vec bb, Vec xx) { 62 PC_CP *cp = (PC_CP *)pc->data; 63 PetscScalar *b, *x, xt; 64 PetscInt i, j; 65 66 PetscFunctionBegin; 67 PetscCall(VecCopy(bb, cp->work)); 68 PetscCall(VecGetArray(cp->work, &b)); 69 PetscCall(VecGetArray(xx, &x)); 70 71 for (i = 0; i < cp->n; i++) { /* over columns */ 72 xt = 0.; 73 for (j = cp->i[i]; j < cp->i[i + 1]; j++) xt += cp->a[j] * b[cp->j[j]]; /* over rows in column */ 74 xt *= cp->d[i]; 75 x[i] = xt; 76 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*/ 77 } 78 for (i = cp->n - 1; i > -1; 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 86 PetscCall(VecRestoreArray(cp->work, &b)); 87 PetscCall(VecRestoreArray(xx, &x)); 88 PetscFunctionReturn(0); 89 } 90 /* -------------------------------------------------------------------------- */ 91 static PetscErrorCode PCReset_CP(PC pc) { 92 PC_CP *cp = (PC_CP *)pc->data; 93 94 PetscFunctionBegin; 95 PetscCall(PetscFree(cp->d)); 96 PetscCall(VecDestroy(&cp->work)); 97 PetscCall(PetscFree3(cp->a, cp->i, cp->j)); 98 PetscFunctionReturn(0); 99 } 100 101 static PetscErrorCode PCDestroy_CP(PC pc) { 102 PC_CP *cp = (PC_CP *)pc->data; 103 104 PetscFunctionBegin; 105 PetscCall(PCReset_CP(pc)); 106 PetscCall(PetscFree(cp->d)); 107 PetscCall(PetscFree3(cp->a, cp->i, cp->j)); 108 PetscCall(PetscFree(pc->data)); 109 PetscFunctionReturn(0); 110 } 111 112 static PetscErrorCode PCSetFromOptions_CP(PC pc, PetscOptionItems *PetscOptionsObject) { 113 PetscFunctionBegin; 114 PetscFunctionReturn(0); 115 } 116 117 /*MC 118 PCCP - a "column-projection" preconditioner 119 120 This is a terrible preconditioner and is not recommended, ever! 121 122 Loops over the entries of x computing dx_i (e_i is the unit vector in the ith direction) to 123 $ 124 $ min || b - A(x + dx_i e_i ||_2 125 $ dx_i 126 $ 127 $ That is, it changes a single entry of x to minimize the new residual norm. 128 $ Let A_i represent the ith column of A, then the minimization can be written as 129 $ 130 $ min || r - (dx_i) A e_i ||_2 131 $ dx_i 132 $ or min || r - (dx_i) A_i ||_2 133 $ dx_i 134 $ 135 $ take the derivative with respect to dx_i to obtain 136 $ dx_i = (A_i^T A_i)^(-1) A_i^T r 137 $ 138 $ This algorithm can be thought of as Gauss-Seidel on the normal equations 139 140 Notes: 141 This proceedure can also be done with block columns or any groups of columns 142 but this is not coded. 143 144 These "projections" can be done simultaneously for all columns (similar to Jacobi) 145 or sequentially (similar to Gauss-Seidel/SOR). This is only coded for SOR type. 146 147 This is related to, but not the same as "row projection" methods. 148 149 This is currently coded only for SeqAIJ matrices in sequential (SOR) form. 150 151 Level: intermediate 152 153 .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PCJACOBI`, `PCSOR` 154 155 M*/ 156 157 PETSC_EXTERN PetscErrorCode PCCreate_CP(PC pc) { 158 PC_CP *cp; 159 160 PetscFunctionBegin; 161 PetscCall(PetscNewLog(pc, &cp)); 162 pc->data = (void *)cp; 163 164 pc->ops->apply = PCApply_CP; 165 pc->ops->applytranspose = PCApply_CP; 166 pc->ops->setup = PCSetUp_CP; 167 pc->ops->reset = PCReset_CP; 168 pc->ops->destroy = PCDestroy_CP; 169 pc->ops->setfromoptions = PCSetFromOptions_CP; 170 pc->ops->view = NULL; 171 pc->ops->applyrichardson = NULL; 172 PetscFunctionReturn(0); 173 } 174