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," KSP on Galerkin follow\n");CHKERRQ(ierr); 107 ierr = PetscViewerASCIIPrintf(viewer," ---------------------------------\n");CHKERRQ(ierr); 108 } 109 ierr = KSPView(jac->ksp,viewer);CHKERRQ(ierr); 110 PetscFunctionReturn(0); 111 } 112 113 static PetscErrorCode PCGalerkinGetKSP_Galerkin(PC pc,KSP *ksp) 114 { 115 PC_Galerkin *jac = (PC_Galerkin*)pc->data; 116 117 PetscFunctionBegin; 118 *ksp = jac->ksp; 119 PetscFunctionReturn(0); 120 } 121 122 static PetscErrorCode PCGalerkinSetRestriction_Galerkin(PC pc,Mat R) 123 { 124 PC_Galerkin *jac = (PC_Galerkin*)pc->data; 125 PetscErrorCode ierr; 126 127 PetscFunctionBegin; 128 ierr = PetscObjectReference((PetscObject)R);CHKERRQ(ierr); 129 ierr = MatDestroy(&jac->R);CHKERRQ(ierr); 130 jac->R = R; 131 PetscFunctionReturn(0); 132 } 133 134 static PetscErrorCode PCGalerkinSetInterpolation_Galerkin(PC pc,Mat P) 135 { 136 PC_Galerkin *jac = (PC_Galerkin*)pc->data; 137 PetscErrorCode ierr; 138 139 PetscFunctionBegin; 140 ierr = PetscObjectReference((PetscObject)P);CHKERRQ(ierr); 141 ierr = MatDestroy(&jac->P);CHKERRQ(ierr); 142 jac->P = P; 143 PetscFunctionReturn(0); 144 } 145 146 static PetscErrorCode PCGalerkinSetComputeSubmatrix_Galerkin(PC pc,PetscErrorCode (*computeAsub)(PC,Mat,Mat,Mat*,void*),void *ctx) 147 { 148 PC_Galerkin *jac = (PC_Galerkin*)pc->data; 149 150 PetscFunctionBegin; 151 jac->computeasub = computeAsub; 152 jac->computeasub_ctx = ctx; 153 PetscFunctionReturn(0); 154 } 155 156 /* -------------------------------------------------------------------------------- */ 157 /*@ 158 PCGalerkinSetRestriction - Sets the restriction operator for the "Galerkin-type" preconditioner 159 160 Logically Collective on PC 161 162 Input Parameter: 163 + pc - the preconditioner context 164 - R - the restriction operator 165 166 Notes: Either this or PCGalerkinSetInterpolation() or both must be called 167 168 Level: Intermediate 169 170 .keywords: PC, set, Galerkin preconditioner 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 Parameter: 192 + pc - the preconditioner context 193 - R - the interpolation operator 194 195 Notes: Either this or PCGalerkinSetRestriction() or both must be called 196 197 Level: Intermediate 198 199 .keywords: PC, set, Galerkin preconditioner 200 201 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN, 202 PCGalerkinSetRestriction(), PCGalerkinGetKSP() 203 204 @*/ 205 PetscErrorCode PCGalerkinSetInterpolation(PC pc,Mat P) 206 { 207 PetscErrorCode ierr; 208 209 PetscFunctionBegin; 210 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 211 ierr = PetscTryMethod(pc,"PCGalerkinSetInterpolation_C",(PC,Mat),(pc,P));CHKERRQ(ierr); 212 PetscFunctionReturn(0); 213 } 214 215 /*@ 216 PCGalerkinSetComputeSubmatrix - Provide a routine that will be called to compute the Galerkin submatrix 217 218 Logically Collective 219 220 Input Parameter: 221 + pc - the preconditioner context 222 . computeAsub - routine that computes the submatrix from the global matrix 223 - ctx - context used by the routine, or NULL 224 225 Calling sequence of computeAsub: 226 $ computeAsub(PC pc,Mat A, Mat Ap, Mat *cAP,void *ctx); 227 228 + PC - the Galerkin PC 229 . A - the matrix in the Galerkin PC 230 . Ap - the computed submatrix from any previous computation, if NULL it has not previously been computed 231 . cAp - the submatrix computed by this routine 232 - ctx - optional user-defined function context 233 234 Level: Intermediate 235 236 Notes: 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 outter 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: If the user does not call this routine nor call PCGalerkinGetKSP() and KSPSetOperators() then PCGalerkin could 243 could automatically compute the submatrix via calls to MatGalerkin() or MatRARt() 244 245 .keywords: PC, get, Galerkin preconditioner, sub preconditioner 246 247 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN, 248 PCGalerkinSetRestriction(), PCGalerkinSetInterpolation(), PCGalerkinGetKSP() 249 250 @*/ 251 PetscErrorCode PCGalerkinSetComputeSubmatrix(PC pc,PetscErrorCode (*computeAsub)(PC,Mat,Mat,Mat*,void*),void *ctx) 252 { 253 PetscErrorCode ierr; 254 255 PetscFunctionBegin; 256 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 257 ierr = PetscTryMethod(pc,"PCGalerkinSetComputeSubmatrix_C",(PC,PetscErrorCode (*)(PC,Mat,Mat,Mat*,void*),void*),(pc,computeAsub,ctx));CHKERRQ(ierr); 258 PetscFunctionReturn(0); 259 } 260 261 /*@ 262 PCGalerkinGetKSP - Gets the KSP object in the Galerkin PC. 263 264 Not Collective 265 266 Input Parameter: 267 . pc - the preconditioner context 268 269 Output Parameters: 270 . ksp - the KSP object 271 272 Level: Intermediate 273 274 Notes: 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 .keywords: PC, get, Galerkin preconditioner, sub preconditioner 278 279 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN, 280 PCGalerkinSetRestriction(), PCGalerkinSetInterpolation(), PCGalerkinSetComputeSubmatrix() 281 282 @*/ 283 PetscErrorCode PCGalerkinGetKSP(PC pc,KSP *ksp) 284 { 285 PetscErrorCode ierr; 286 287 PetscFunctionBegin; 288 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 289 PetscValidPointer(ksp,2); 290 ierr = PetscUseMethod(pc,"PCGalerkinGetKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr); 291 PetscFunctionReturn(0); 292 } 293 294 static PetscErrorCode PCSetFromOptions_Galerkin(PetscOptionItems *PetscOptionsObject,PC pc) 295 { 296 PC_Galerkin *jac = (PC_Galerkin*)pc->data; 297 PetscErrorCode ierr; 298 const char *prefix; 299 PetscBool flg; 300 301 PetscFunctionBegin; 302 ierr = KSPGetOptionsPrefix(jac->ksp,&prefix);CHKERRQ(ierr); 303 ierr = PetscStrendswith(prefix,"galerkin_",&flg);CHKERRQ(ierr); 304 if (!flg) { 305 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 306 ierr = KSPSetOptionsPrefix(jac->ksp,prefix);CHKERRQ(ierr); 307 ierr = KSPAppendOptionsPrefix(jac->ksp,"galerkin_");CHKERRQ(ierr); 308 } 309 310 ierr = PetscOptionsHead(PetscOptionsObject,"Galerkin options");CHKERRQ(ierr); 311 if (jac->ksp) { 312 ierr = KSPSetFromOptions(jac->ksp);CHKERRQ(ierr); 313 } 314 ierr = PetscOptionsTail();CHKERRQ(ierr); 315 PetscFunctionReturn(0); 316 } 317 318 /* -------------------------------------------------------------------------------------------*/ 319 320 /*MC 321 PCGALERKIN - Build (part of) a preconditioner by P S R (where P is often R^T) 322 323 $ Use PCGalerkinSetRestriction(pc,R) and/or PCGalerkinSetInterpolation(pc,P) followed by 324 $ PCGalerkinGetKSP(pc,&ksp); KSPSetOperators(ksp,A,....) 325 326 Level: intermediate 327 328 Developer Note: If KSPSetOperators() has not been called on the inner KSP then PCGALERKIN could use MatRARt() or MatPtAP() to compute 329 the operators automatically. 330 Should there be a prefix for the inner KSP. 331 There is no KSPSetFromOptions_Galerkin() that calls KSPSetFromOptions() on the inner KSP 332 333 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 334 PCSHELL, PCKSP, PCGalerkinSetRestriction(), PCGalerkinSetInterpolation(), PCGalerkinGetKSP() 335 336 M*/ 337 338 PETSC_EXTERN PetscErrorCode PCCreate_Galerkin(PC pc) 339 { 340 PetscErrorCode ierr; 341 PC_Galerkin *jac; 342 343 PetscFunctionBegin; 344 ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr); 345 346 pc->ops->apply = PCApply_Galerkin; 347 pc->ops->setup = PCSetUp_Galerkin; 348 pc->ops->reset = PCReset_Galerkin; 349 pc->ops->destroy = PCDestroy_Galerkin; 350 pc->ops->view = PCView_Galerkin; 351 pc->ops->setfromoptions = PCSetFromOptions_Galerkin; 352 pc->ops->applyrichardson = NULL; 353 354 ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&jac->ksp);CHKERRQ(ierr); 355 ierr = KSPSetErrorIfNotConverged(jac->ksp,pc->erroriffailure);CHKERRQ(ierr); 356 ierr = PetscObjectIncrementTabLevel((PetscObject)jac->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 357 358 pc->data = (void*)jac; 359 360 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGalerkinSetRestriction_C",PCGalerkinSetRestriction_Galerkin);CHKERRQ(ierr); 361 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGalerkinSetInterpolation_C",PCGalerkinSetInterpolation_Galerkin);CHKERRQ(ierr); 362 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGalerkinGetKSP_C",PCGalerkinGetKSP_Galerkin);CHKERRQ(ierr); 363 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGalerkinSetComputeSubmatrix_C",PCGalerkinSetComputeSubmatrix_Galerkin);CHKERRQ(ierr); 364 PetscFunctionReturn(0); 365 } 366 367