1 #define PETSCKSP_DLL 2 3 #include "private/pcimpl.h" /*I "petscpc.h" I*/ 4 #include "src/mat/impls/aij/seq/aij.h" 5 6 /* 7 Private context (data structure) for the CP preconditioner. 8 */ 9 typedef struct { 10 PetscInt n,m; 11 Vec work; 12 PetscScalar *d; /* sum of squares of each column */ 13 PetscScalar *a; /* non-zeros by column */ 14 PetscInt *i,*j; /* offsets of nonzeros by column, non-zero indices by column */ 15 } PC_CP; 16 17 18 #undef __FUNCT__ 19 #define __FUNCT__ "PCSetUp_CP" 20 static PetscErrorCode PCSetUp_CP(PC pc) 21 { 22 PC_CP *cp = (PC_CP*)pc->data; 23 PetscInt i,j,*colcnt; 24 PetscErrorCode ierr; 25 PetscTruth flg; 26 Mat_SeqAIJ *aij = (Mat_SeqAIJ*)pc->pmat->data; 27 28 PetscFunctionBegin; 29 ierr = PetscTypeCompare((PetscObject)pc->pmat,MATSEQAIJ,&flg);CHKERRQ(ierr); 30 if (!flg) SETERRQ(PETSC_ERR_SUP,"Currently only handles SeqAIJ matrices"); 31 32 ierr = MatGetLocalSize(pc->pmat,&cp->m,&cp->n);CHKERRQ(ierr); 33 if (cp->m != cp->n) SETERRQ(PETSC_ERR_SUP,"Currently only for square matrices"); 34 35 if (!cp->work) {ierr = MatGetVecs(pc->pmat,&cp->work,PETSC_NULL);CHKERRQ(ierr);} 36 if (!cp->d) {ierr = PetscMalloc(cp->n*sizeof(PetscScalar),&cp->d);CHKERRQ(ierr);} 37 if (cp->a && pc->flag != SAME_NONZERO_PATTERN) { 38 ierr = PetscFree3(cp->a,cp->i,cp->j);CHKERRQ(ierr); 39 cp->a = 0; 40 } 41 42 /* convert to column format */ 43 if (!cp->a) { 44 ierr = PetscMalloc3(aij->nz,PetscScalar,&cp->a,cp->n+1,PetscInt,&cp->i,aij->nz,PetscInt,&cp->j);CHKERRQ(ierr); 45 } 46 ierr = PetscMalloc(cp->n*sizeof(PetscInt),&colcnt);CHKERRQ(ierr); 47 ierr = PetscMemzero(colcnt,cp->n*sizeof(PetscInt));CHKERRQ(ierr); 48 49 for (i=0; i<aij->nz; i++) { 50 colcnt[aij->j[i]]++; 51 } 52 cp->i[0] = 0; 53 for (i=0; i<cp->n; i++) { 54 cp->i[i+1] = cp->i[i] + colcnt[i]; 55 } 56 ierr = PetscMemzero(colcnt,cp->n*sizeof(PetscInt));CHKERRQ(ierr); 57 for (i=0; i<cp->m; i++) { /* over rows */ 58 for (j=aij->i[i]; j<aij->i[i+1]; j++) { /* over columns in row */ 59 cp->j[cp->i[aij->j[j]]+colcnt[aij->j[j]]] = i; 60 cp->a[cp->i[aij->j[j]]+colcnt[aij->j[j]]++] = aij->a[j]; 61 } 62 } 63 ierr = PetscFree(colcnt);CHKERRQ(ierr); 64 65 /* compute sum of squares of each column d[] */ 66 for (i=0; i<cp->n; i++) { /* over columns */ 67 cp->d[i] = 0; 68 for (j=cp->i[i]; j<cp->i[i+1]; j++) { /* over rows in column */ 69 cp->d[i] += cp->a[j]*cp->a[j]; 70 } 71 cp->d[i] = 1.0/cp->d[i]; 72 } 73 PetscFunctionReturn(0); 74 } 75 /* -------------------------------------------------------------------------- */ 76 #undef __FUNCT__ 77 #define __FUNCT__ "PCApply_CP" 78 static PetscErrorCode PCApply_CP(PC pc,Vec bb,Vec xx) 79 { 80 PC_CP *cp = (PC_CP*)pc->data; 81 PetscErrorCode ierr; 82 PetscScalar *b,*x,xt; 83 PetscInt i,j; 84 85 PetscFunctionBegin; 86 ierr = VecCopy(bb,cp->work);CHKERRQ(ierr); 87 ierr = VecGetArray(cp->work,&b);CHKERRQ(ierr); 88 ierr = VecGetArray(xx,&x);CHKERRQ(ierr); 89 90 for (i=0; i<cp->n; i++) { /* over columns */ 91 xt = 0; 92 for (j=cp->i[i]; j<cp->i[i+1]; j++) { /* over rows in column */ 93 xt += cp->a[j]*b[cp->j[j]]; 94 } 95 xt *= cp->d[i]; 96 x[i] = xt; 97 for (j=cp->i[i]; j<cp->i[i+1]; j++) { /* over rows in column updating b*/ 98 b[cp->j[j]] -= xt*cp->a[j]; 99 } 100 } 101 for (i=cp->n-1; i>-1; i--) { /* over columns */ 102 xt = 0; 103 for (j=cp->i[i]; j<cp->i[i+1]; j++) { /* over rows in column */ 104 xt += cp->a[j]*b[cp->j[j]]; 105 } 106 xt *= cp->d[i]; 107 x[i] = xt; 108 for (j=cp->i[i]; j<cp->i[i+1]; j++) { /* over rows in column updating b*/ 109 b[cp->j[j]] -= xt*cp->a[j]; 110 } 111 } 112 113 ierr = VecRestoreArray(cp->work,&b);CHKERRQ(ierr); 114 ierr = VecRestoreArray(xx,&x);CHKERRQ(ierr); 115 PetscFunctionReturn(0); 116 } 117 /* -------------------------------------------------------------------------- */ 118 #undef __FUNCT__ 119 #define __FUNCT__ "PCDestroy_CP" 120 static PetscErrorCode PCDestroy_CP(PC pc) 121 { 122 PC_CP *cp = (PC_CP*)pc->data; 123 PetscErrorCode ierr; 124 125 PetscFunctionBegin; 126 ierr = PetscFree(cp->d);CHKERRQ(ierr); 127 if (cp->work) {ierr = VecDestroy(cp->work);CHKERRQ(ierr);} 128 ierr = PetscFree3(cp->a,cp->i,cp->j);CHKERRQ(ierr); 129 ierr = PetscFree(cp);CHKERRQ(ierr); 130 PetscFunctionReturn(0); 131 } 132 133 #undef __FUNCT__ 134 #define __FUNCT__ "PCSetFromOptions_CP" 135 static PetscErrorCode PCSetFromOptions_CP(PC pc) 136 { 137 PetscFunctionBegin; 138 PetscFunctionReturn(0); 139 } 140 141 142 /*MC 143 PCCP - a "column-projection" preconditioner 144 145 This is a terrible preconditioner and is not recommended, ever! 146 147 Loops over the entries of x computing dx_i to 148 $ 149 $ min || b - A(x + dx_i e_i ||_2 150 $ dx_i 151 $ 152 $ That is, it changes a single entry of x to minimize the new residual. 153 $ Let A_i represent the ith column of A, then the minimization can be written as 154 $ 155 $ min || r - (dx_i) A e_i ||_2 156 $ dx_i 157 $ or min || r - (dx_i) A_i ||_2 158 $ dx_i 159 $ 160 $ take the derivative with respect to dx_i to obtain 161 $ dx_i = (A_i^T A_i)^(-1) A_i^T r 162 $ 163 $ This algorithm can be thought of as Gauss-Seidel on the normal equations 164 165 Notes: This proceedure can also be done with block columns or any groups of columns 166 but this is not coded. 167 168 These "projections" can be done simultaneously for all columns (similar to Jacobi) 169 or sequentially (similar to Gauss-Seidel/SOR). This is only coded for SOR type. 170 171 This is related to, but not the same as "row projection" methods. 172 173 This is currently coded only for SeqAIJ matrices in sequential (SOR) form. 174 175 Level: intermediate 176 177 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCJACOBI, PCSOR 178 179 M*/ 180 181 EXTERN_C_BEGIN 182 #undef __FUNCT__ 183 #define __FUNCT__ "PCCreate_CP" 184 PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_CP(PC pc) 185 { 186 PC_CP *cp; 187 PetscErrorCode ierr; 188 189 PetscFunctionBegin; 190 ierr = PetscNewLog(pc,PC_CP,&cp);CHKERRQ(ierr); 191 pc->data = (void*)cp; 192 193 pc->ops->apply = PCApply_CP; 194 pc->ops->applytranspose = PCApply_CP; 195 pc->ops->setup = PCSetUp_CP; 196 pc->ops->destroy = PCDestroy_CP; 197 pc->ops->setfromoptions = PCSetFromOptions_CP; 198 pc->ops->view = 0; 199 pc->ops->applyrichardson = 0; 200 PetscFunctionReturn(0); 201 } 202 EXTERN_C_END 203 204 205