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