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: 167 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: 197 Either this or PCGalerkinSetRestriction() or both must be called 198 199 Level: Intermediate 200 201 .keywords: PC, set, Galerkin preconditioner 202 203 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN, 204 PCGalerkinSetRestriction(), PCGalerkinGetKSP() 205 206 @*/ 207 PetscErrorCode PCGalerkinSetInterpolation(PC pc,Mat P) 208 { 209 PetscErrorCode ierr; 210 211 PetscFunctionBegin; 212 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 213 ierr = PetscTryMethod(pc,"PCGalerkinSetInterpolation_C",(PC,Mat),(pc,P));CHKERRQ(ierr); 214 PetscFunctionReturn(0); 215 } 216 217 /*@ 218 PCGalerkinSetComputeSubmatrix - Provide a routine that will be called to compute the Galerkin submatrix 219 220 Logically Collective 221 222 Input Parameter: 223 + pc - the preconditioner context 224 . computeAsub - routine that computes the submatrix from the global matrix 225 - ctx - context used by the routine, or NULL 226 227 Calling sequence of computeAsub: 228 $ computeAsub(PC pc,Mat A, Mat Ap, Mat *cAP,void *ctx); 229 230 + PC - the Galerkin PC 231 . A - the matrix in the Galerkin PC 232 . Ap - the computed submatrix from any previous computation, if NULL it has not previously been computed 233 . cAp - the submatrix computed by this routine 234 - ctx - optional user-defined function context 235 236 Level: Intermediate 237 238 Notes: 239 Instead of providing this routine you can call PCGalerkinGetKSP() and then KSPSetOperators() to provide the submatrix, 240 but that will not work for multiple KSPSolves with different matrices unless you call it for each solve. 241 242 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 243 matrix and computes its values in cAp. On each subsequent call the routine should up the Ap matrix. 244 245 Developer Notes: 246 If the user does not call this routine nor call PCGalerkinGetKSP() and KSPSetOperators() then PCGalerkin could 247 could automatically compute the submatrix via calls to MatGalerkin() or MatRARt() 248 249 .keywords: PC, get, Galerkin preconditioner, sub preconditioner 250 251 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN, 252 PCGalerkinSetRestriction(), PCGalerkinSetInterpolation(), PCGalerkinGetKSP() 253 254 @*/ 255 PetscErrorCode PCGalerkinSetComputeSubmatrix(PC pc,PetscErrorCode (*computeAsub)(PC,Mat,Mat,Mat*,void*),void *ctx) 256 { 257 PetscErrorCode ierr; 258 259 PetscFunctionBegin; 260 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 261 ierr = PetscTryMethod(pc,"PCGalerkinSetComputeSubmatrix_C",(PC,PetscErrorCode (*)(PC,Mat,Mat,Mat*,void*),void*),(pc,computeAsub,ctx));CHKERRQ(ierr); 262 PetscFunctionReturn(0); 263 } 264 265 /*@ 266 PCGalerkinGetKSP - Gets the KSP object in the Galerkin PC. 267 268 Not Collective 269 270 Input Parameter: 271 . pc - the preconditioner context 272 273 Output Parameters: 274 . ksp - the KSP object 275 276 Level: Intermediate 277 278 Notes: 279 Once you have called this routine you can call KSPSetOperators() on the resulting ksp to provide the operator for the Galerkin problem, 280 an alternative is to use PCGalerkinSetComputeSubmatrix() to provide a routine that computes the submatrix as needed. 281 282 .keywords: PC, get, Galerkin preconditioner, sub preconditioner 283 284 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCGALERKIN, 285 PCGalerkinSetRestriction(), PCGalerkinSetInterpolation(), PCGalerkinSetComputeSubmatrix() 286 287 @*/ 288 PetscErrorCode PCGalerkinGetKSP(PC pc,KSP *ksp) 289 { 290 PetscErrorCode ierr; 291 292 PetscFunctionBegin; 293 PetscValidHeaderSpecific(pc,PC_CLASSID,1); 294 PetscValidPointer(ksp,2); 295 ierr = PetscUseMethod(pc,"PCGalerkinGetKSP_C",(PC,KSP*),(pc,ksp));CHKERRQ(ierr); 296 PetscFunctionReturn(0); 297 } 298 299 static PetscErrorCode PCSetFromOptions_Galerkin(PetscOptionItems *PetscOptionsObject,PC pc) 300 { 301 PC_Galerkin *jac = (PC_Galerkin*)pc->data; 302 PetscErrorCode ierr; 303 const char *prefix; 304 PetscBool flg; 305 306 PetscFunctionBegin; 307 ierr = KSPGetOptionsPrefix(jac->ksp,&prefix);CHKERRQ(ierr); 308 ierr = PetscStrendswith(prefix,"galerkin_",&flg);CHKERRQ(ierr); 309 if (!flg) { 310 ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 311 ierr = KSPSetOptionsPrefix(jac->ksp,prefix);CHKERRQ(ierr); 312 ierr = KSPAppendOptionsPrefix(jac->ksp,"galerkin_");CHKERRQ(ierr); 313 } 314 315 ierr = PetscOptionsHead(PetscOptionsObject,"Galerkin options");CHKERRQ(ierr); 316 if (jac->ksp) { 317 ierr = KSPSetFromOptions(jac->ksp);CHKERRQ(ierr); 318 } 319 ierr = PetscOptionsTail();CHKERRQ(ierr); 320 PetscFunctionReturn(0); 321 } 322 323 /* -------------------------------------------------------------------------------------------*/ 324 325 /*MC 326 PCGALERKIN - Build (part of) a preconditioner by P S R (where P is often R^T) 327 328 $ Use PCGalerkinSetRestriction(pc,R) and/or PCGalerkinSetInterpolation(pc,P) followed by 329 $ PCGalerkinGetKSP(pc,&ksp); KSPSetOperators(ksp,A,....) 330 331 Level: intermediate 332 333 Developer Note: If KSPSetOperators() has not been called on the inner KSP then PCGALERKIN could use MatRARt() or MatPtAP() to compute 334 the operators automatically. 335 Should there be a prefix for the inner KSP. 336 There is no KSPSetFromOptions_Galerkin() that calls KSPSetFromOptions() on the inner KSP 337 338 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 339 PCSHELL, PCKSP, PCGalerkinSetRestriction(), PCGalerkinSetInterpolation(), PCGalerkinGetKSP() 340 341 M*/ 342 343 PETSC_EXTERN PetscErrorCode PCCreate_Galerkin(PC pc) 344 { 345 PetscErrorCode ierr; 346 PC_Galerkin *jac; 347 348 PetscFunctionBegin; 349 ierr = PetscNewLog(pc,&jac);CHKERRQ(ierr); 350 351 pc->ops->apply = PCApply_Galerkin; 352 pc->ops->setup = PCSetUp_Galerkin; 353 pc->ops->reset = PCReset_Galerkin; 354 pc->ops->destroy = PCDestroy_Galerkin; 355 pc->ops->view = PCView_Galerkin; 356 pc->ops->setfromoptions = PCSetFromOptions_Galerkin; 357 pc->ops->applyrichardson = NULL; 358 359 ierr = KSPCreate(PetscObjectComm((PetscObject)pc),&jac->ksp);CHKERRQ(ierr); 360 ierr = KSPSetErrorIfNotConverged(jac->ksp,pc->erroriffailure);CHKERRQ(ierr); 361 ierr = PetscObjectIncrementTabLevel((PetscObject)jac->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 362 363 pc->data = (void*)jac; 364 365 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGalerkinSetRestriction_C",PCGalerkinSetRestriction_Galerkin);CHKERRQ(ierr); 366 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGalerkinSetInterpolation_C",PCGalerkinSetInterpolation_Galerkin);CHKERRQ(ierr); 367 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGalerkinGetKSP_C",PCGalerkinGetKSP_Galerkin);CHKERRQ(ierr); 368 ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGalerkinSetComputeSubmatrix_C",PCGalerkinSetComputeSubmatrix_Galerkin);CHKERRQ(ierr); 369 PetscFunctionReturn(0); 370 } 371 372