1 2 /* 3 Defines a preconditioner defined by R^T S R 4 */ 5 #include <petsc/private/pcimpl.h> 6 #include <petscksp.h> /*I "petscksp.h" I*/ 7 8 typedef struct { 9 KSP ksp; 10 Mat R,P; 11 Vec b,x; 12 PetscErrorCode (*computeasub)(PC,Mat,Mat,Mat*,void*); 13 void *computeasub_ctx; 14 } PC_Galerkin; 15 16 static PetscErrorCode PCApply_Galerkin(PC pc,Vec x,Vec y) 17 { 18 PetscErrorCode ierr; 19 PC_Galerkin *jac = (PC_Galerkin*)pc->data; 20 21 PetscFunctionBegin; 22 if (jac->R) { 23 ierr = MatRestrict(jac->R,x,jac->b);CHKERRQ(ierr); 24 } else { 25 ierr = MatRestrict(jac->P,x,jac->b);CHKERRQ(ierr); 26 } 27 ierr = KSPSolve(jac->ksp,jac->b,jac->x);CHKERRQ(ierr); 28 if (jac->P) { 29 ierr = MatInterpolate(jac->P,jac->x,y);CHKERRQ(ierr); 30 } else { 31 ierr = MatInterpolate(jac->R,jac->x,y);CHKERRQ(ierr); 32 } 33 PetscFunctionReturn(0); 34 } 35 36 static PetscErrorCode PCSetUp_Galerkin(PC pc) 37 { 38 PetscErrorCode ierr; 39 PC_Galerkin *jac = (PC_Galerkin*)pc->data; 40 PetscBool a; 41 Vec *xx,*yy; 42 43 PetscFunctionBegin; 44 if (jac->computeasub) { 45 Mat Ap; 46 if (!pc->setupcalled) { 47 ierr = (*jac->computeasub)(pc,pc->pmat,NULL,&Ap,jac->computeasub_ctx);CHKERRQ(ierr); 48 ierr = KSPSetOperators(jac->ksp,Ap,Ap);CHKERRQ(ierr); 49 ierr = MatDestroy(&Ap);CHKERRQ(ierr); 50 } else { 51 ierr = KSPGetOperators(jac->ksp,NULL,&Ap);CHKERRQ(ierr); 52 ierr = (*jac->computeasub)(pc,pc->pmat,Ap,NULL,jac->computeasub_ctx);CHKERRQ(ierr); 53 } 54 } 55 56 if (!jac->x) { 57 ierr = KSPGetOperatorsSet(jac->ksp,&a,NULL);CHKERRQ(ierr); 58 if (!a) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set operator of PCGALERKIN KSP with PCGalerkinGetKSP()/KSPSetOperators()"); 59 ierr = KSPCreateVecs(jac->ksp,1,&xx,1,&yy);CHKERRQ(ierr); 60 jac->x = *xx; 61 jac->b = *yy; 62 ierr = PetscFree(xx);CHKERRQ(ierr); 63 ierr = PetscFree(yy);CHKERRQ(ierr); 64 } 65 if (!jac->R && !jac->P) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or interpolation of PCGALERKIN with PCGalerkinSetRestriction()/Interpolation()"); 66 /* should check here that sizes of R/P match size of a */ 67 68 PetscFunctionReturn(0); 69 } 70 71 static PetscErrorCode PCReset_Galerkin(PC pc) 72 { 73 PC_Galerkin *jac = (PC_Galerkin*)pc->data; 74 PetscErrorCode ierr; 75 76 PetscFunctionBegin; 77 ierr = MatDestroy(&jac->R);CHKERRQ(ierr); 78 ierr = MatDestroy(&jac->P);CHKERRQ(ierr); 79 ierr = VecDestroy(&jac->x);CHKERRQ(ierr); 80 ierr = VecDestroy(&jac->b);CHKERRQ(ierr); 81 ierr = KSPReset(jac->ksp);CHKERRQ(ierr); 82 PetscFunctionReturn(0); 83 } 84 85 static PetscErrorCode PCDestroy_Galerkin(PC pc) 86 { 87 PC_Galerkin *jac = (PC_Galerkin*)pc->data; 88 PetscErrorCode ierr; 89 90 PetscFunctionBegin; 91 ierr = PCReset_Galerkin(pc);CHKERRQ(ierr); 92 ierr = KSPDestroy(&jac->ksp);CHKERRQ(ierr); 93 ierr = PetscFree(pc->data);CHKERRQ(ierr); 94 PetscFunctionReturn(0); 95 } 96 97 static PetscErrorCode PCView_Galerkin(PC pc,PetscViewer viewer) 98 { 99 PC_Galerkin *jac = (PC_Galerkin*)pc->data; 100 PetscErrorCode ierr; 101 PetscBool iascii; 102 103 PetscFunctionBegin; 104 ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 105 if (iascii) { 106 ierr = PetscViewerASCIIPrintf(viewer,"Galerkin PC\n");CHKERRQ(ierr); 107 ierr = PetscViewerASCIIPrintf(viewer,"KSP on Galerkin follow\n");CHKERRQ(ierr); 108 ierr = PetscViewerASCIIPrintf(viewer,"---------------------------------\n");CHKERRQ(ierr); 109 } 110 ierr = KSPView(jac->ksp,viewer);CHKERRQ(ierr); 111 PetscFunctionReturn(0); 112 } 113 114 static PetscErrorCode PCGalerkinGetKSP_Galerkin(PC pc,KSP *ksp) 115 { 116 PC_Galerkin *jac = (PC_Galerkin*)pc->data; 117 118 PetscFunctionBegin; 119 *ksp = jac->ksp; 120 PetscFunctionReturn(0); 121 } 122 123 static PetscErrorCode PCGalerkinSetRestriction_Galerkin(PC pc,Mat R) 124 { 125 PC_Galerkin *jac = (PC_Galerkin*)pc->data; 126 PetscErrorCode ierr; 127 128 PetscFunctionBegin; 129 ierr = PetscObjectReference((PetscObject)R);CHKERRQ(ierr); 130 ierr = MatDestroy(&jac->R);CHKERRQ(ierr); 131 jac->R = R; 132 PetscFunctionReturn(0); 133 } 134 135 static PetscErrorCode PCGalerkinSetInterpolation_Galerkin(PC pc,Mat P) 136 { 137 PC_Galerkin *jac = (PC_Galerkin*)pc->data; 138 PetscErrorCode ierr; 139 140 PetscFunctionBegin; 141 ierr = PetscObjectReference((PetscObject)P);CHKERRQ(ierr); 142 ierr = MatDestroy(&jac->P);CHKERRQ(ierr); 143 jac->P = P; 144 PetscFunctionReturn(0); 145 } 146 147 static PetscErrorCode PCGalerkinSetComputeSubmatrix_Galerkin(PC pc,PetscErrorCode (*computeAsub)(PC,Mat,Mat,Mat*,void*),void *ctx) 148 { 149 PC_Galerkin *jac = (PC_Galerkin*)pc->data; 150 151 PetscFunctionBegin; 152 jac->computeasub = computeAsub; 153 jac->computeasub_ctx = ctx; 154 PetscFunctionReturn(0); 155 } 156 157 /* -------------------------------------------------------------------------------- */ 158 /*@ 159 PCGalerkinSetRestriction - Sets the restriction operator for the "Galerkin-type" preconditioner 160 161 Logically Collective on PC 162 163 Input Parameter: 164 + pc - the preconditioner context 165 - R - the restriction operator 166 167 Notes: Either this or PCGalerkinSetInterpolation() or both must be called 168 169 Level: Intermediate 170 171 .keywords: PC, set, Galerkin preconditioner 172 173 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN, 174 PCGalerkinSetInterpolation(), PCGalerkinGetKSP() 175 176 @*/ 177 PetscErrorCode PCGalerkinSetRestriction(PC pc,Mat R) 178 { 179 PetscErrorCode ierr; 180 181 PetscFunctionBegin; 182 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 183 ierr = PetscTryMethod(pc,"PCGalerkinSetRestriction_C",(PC,Mat),(pc,R));CHKERRQ(ierr); 184 PetscFunctionReturn(0); 185 } 186 187 /*@ 188 PCGalerkinSetInterpolation - Sets the interpolation operator for the "Galerkin-type" preconditioner 189 190 Logically Collective on PC 191 192 Input Parameter: 193 + pc - the preconditioner context 194 - R - the interpolation operator 195 196 Notes: Either this or PCGalerkinSetRestriction() or both must be called 197 198 Level: Intermediate 199 200 .keywords: PC, set, Galerkin preconditioner 201 202 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN, 203 PCGalerkinSetRestriction(), PCGalerkinGetKSP() 204 205 @*/ 206 PetscErrorCode PCGalerkinSetInterpolation(PC pc,Mat P) 207 { 208 PetscErrorCode ierr; 209 210 PetscFunctionBegin; 211 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 212 ierr = PetscTryMethod(pc,"PCGalerkinSetInterpolation_C",(PC,Mat),(pc,P));CHKERRQ(ierr); 213 PetscFunctionReturn(0); 214 } 215 216 /*@ 217 PCGalerkinSetComputeSubmatrix - Provide a routine that will be called to compute the Galerkin submatrix 218 219 Logically Collective 220 221 Input Parameter: 222 + pc - the preconditioner context 223 . computeAsub - routine that computes the submatrix from the global matrix 224 - ctx - context used by the routine, or NULL 225 226 Calling sequence of computeAsub: 227 $ computeAsub(PC pc,Mat A, Mat Ap, Mat *cAP,void *ctx); 228 229 + PC - the Galerkin PC 230 . A - the matrix in the Galerkin PC 231 . Ap - the computed submatrix from any previous computation, if NULL it has not previously been computed 232 . cAp - the submatrix computed by this routine 233 - ctx - optional user-defined function context 234 235 Level: Intermediate 236 237 Notes: Instead of providing this routine you can call PCGalerkinGetKSP() and then KSPSetOperators() to provide the submatrix, 238 but that will not work for multiple KSPSolves with different matrices unless you call it for each solve. 239 240 This routine is called each time the outter matrix is changed. In the first call the Ap argument is NULL and the routine should create the 241 matrix and computes its values in cAp. On each subsequent call the routine should up the Ap matrix. 242 243 Developer Notes: If the user does not call this routine nor call PCGalerkinGetKSP() and KSPSetOperators() then PCGalerkin could 244 could automatically compute the submatrix via calls to MatGalerkin() or MatRARt() 245 246 .keywords: PC, get, Galerkin preconditioner, sub preconditioner 247 248 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN, 249 PCGalerkinSetRestriction(), PCGalerkinSetInterpolation(), PCGalerkinGetKSP() 250 251 @*/ 252 PetscErrorCode PCGalerkinSetComputeSubmatrix(PC pc,PetscErrorCode (*computeAsub)(PC,Mat,Mat,Mat*,void*),void *ctx) 253 { 254 PetscErrorCode ierr; 255 256 PetscFunctionBegin; 257 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 258 ierr = PetscTryMethod(pc,"PCGalerkinSetComputeSubmatrix_C",(PC,PetscErrorCode (*)(PC,Mat,Mat,Mat*,void*),void*),(pc,computeAsub,ctx));CHKERRQ(ierr); 259 PetscFunctionReturn(0); 260 } 261 262 /*@ 263 PCGalerkinGetKSP - Gets the KSP object in the Galerkin PC. 264 265 Not Collective 266 267 Input Parameter: 268 . pc - the preconditioner context 269 270 Output Parameters: 271 . ksp - the KSP object 272 273 Level: Intermediate 274 275 Notes: Once you have called this routine you can call KSPSetOperators() on the resulting ksp to provide the operator for the Galerkin problem, 276 an alternative is to use PCGalerkinSetComputeSubmatrix() to provide a routine that computes the submatrix as needed. 277 278 .keywords: PC, get, Galerkin preconditioner, sub preconditioner 279 280 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN, 281 PCGalerkinSetRestriction(), PCGalerkinSetInterpolation(), PCGalerkinSetComputeSubmatrix() 282 283 @*/ 284 PetscErrorCode PCGalerkinGetKSP(PC pc,KSP *ksp) 285 { 286 PetscErrorCode ierr; 287 288 PetscFunctionBegin; 289 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 290 PetscValidPointer(ksp,2); 291 ierr = PetscUseMethod(pc,"PCGalerkinGetKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr); 292 PetscFunctionReturn(0); 293 } 294 295 static PetscErrorCode PCSetFromOptions_Galerkin(PetscOptionItems *PetscOptionsObject,PC pc) 296 { 297 PC_Galerkin *jac = (PC_Galerkin*)pc->data; 298 PetscErrorCode ierr; 299 const char *prefix; 300 PetscBool flg; 301 302 PetscFunctionBegin; 303 ierr = KSPGetOptionsPrefix(jac->ksp,&prefix);CHKERRQ(ierr); 304 ierr = PetscStrendswith(prefix,"galerkin_",&flg);CHKERRQ(ierr); 305 if (!flg) { 306 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 307 ierr = KSPSetOptionsPrefix(jac->ksp,prefix);CHKERRQ(ierr); 308 ierr = KSPAppendOptionsPrefix(jac->ksp,"galerkin_");CHKERRQ(ierr); 309 } 310 311 ierr = PetscOptionsHead(PetscOptionsObject,"Galerkin options");CHKERRQ(ierr); 312 if (jac->ksp) { 313 ierr = KSPSetFromOptions(jac->ksp);CHKERRQ(ierr); 314 } 315 ierr = PetscOptionsTail();CHKERRQ(ierr); 316 PetscFunctionReturn(0); 317 } 318 319 /* -------------------------------------------------------------------------------------------*/ 320 321 /*MC 322 PCGALERKIN - Build (part of) a preconditioner by P S R (where P is often R^T) 323 324 $ Use PCGalerkinSetRestriction(pc,R) and/or PCGalerkinSetInterpolation(pc,P) followed by 325 $ PCGalerkinGetKSP(pc,&ksp); KSPSetOperators(ksp,A,....) 326 327 Level: intermediate 328 329 Developer Note: If KSPSetOperators() has not been called on the inner KSP then PCGALERKIN could use MatRARt() or MatPtAP() to compute 330 the operators automatically. 331 Should there be a prefix for the inner KSP. 332 There is no KSPSetFromOptions_Galerkin() that calls KSPSetFromOptions() on the inner KSP 333 334 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 335 PCSHELL, PCKSP, PCGalerkinSetRestriction(), PCGalerkinSetInterpolation(), PCGalerkinGetKSP() 336 337 M*/ 338 339 PETSC_EXTERN PetscErrorCode PCCreate_Galerkin(PC pc) 340 { 341 PetscErrorCode ierr; 342 PC_Galerkin *jac; 343 344 PetscFunctionBegin; 345 ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr); 346 347 pc->ops->apply = PCApply_Galerkin; 348 pc->ops->setup = PCSetUp_Galerkin; 349 pc->ops->reset = PCReset_Galerkin; 350 pc->ops->destroy = PCDestroy_Galerkin; 351 pc->ops->view = PCView_Galerkin; 352 pc->ops->setfromoptions = PCSetFromOptions_Galerkin; 353 pc->ops->applyrichardson = NULL; 354 355 ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&jac->ksp);CHKERRQ(ierr); 356 ierr = KSPSetErrorIfNotConverged(jac->ksp,pc->erroriffailure);CHKERRQ(ierr); 357 ierr = PetscObjectIncrementTabLevel((PetscObject)jac->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 358 359 pc->data = (void*)jac; 360 361 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGalerkinSetRestriction_C",PCGalerkinSetRestriction_Galerkin);CHKERRQ(ierr); 362 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGalerkinSetInterpolation_C",PCGalerkinSetInterpolation_Galerkin);CHKERRQ(ierr); 363 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGalerkinGetKSP_C",PCGalerkinGetKSP_Galerkin);CHKERRQ(ierr); 364 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGalerkinSetComputeSubmatrix_C",PCGalerkinSetComputeSubmatrix_Galerkin);CHKERRQ(ierr); 365 PetscFunctionReturn(0); 366 } 367 368