[petsc-users] user provided local preconditioner with additive schwarz preconditioner

Barry Smith bsmith at mcs.anl.gov
Fri Jul 1 12:25:49 CDT 2016


> On Jun 30, 2016, at 11:48 PM, Duan Zhaowen <dzw.pku at gmail.com> wrote:
> 
> Thank you Barry. I'll try and follow the code. In my code the global matrix A was partitioned into CSR format. The local preconditioner I want to use only effect on local part of matrix A (diagonal part).

   What do you mean? In additive Schwarz the subproblems contain overlapping sets of variables that are solved for an then updated. 

> So in MyApplyFunc(PC, Vec x, Vec y) of the shell preconditioner, should I only take care of local part of vector y and leave alone its non-local part (or overlap)? 

   The local solve updates the entire y (which has overlap with other y from other subdomains). 

    There are several variants of overlapping Schwarz 

/*E
    PCASMType - Type of additive Schwarz method to use

$  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
$                        and computed values in ghost regions are added together.
$                        Classical standard additive Schwarz.
$  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
$                        region are discarded.
$                        Default.
$  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
$                        region are added back in.
$  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
$                        discarded.
$                        Not very good.

   Level: beginner

.seealso: PCASMSetType()
E*/
typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;

but this is all handled inside the PCASM code. Your local function doesn't know or care which of the variants is used. It is the job of your local function to solve for all the values in the y output based on all the values in the x input.

  Barry

> 
> Thank you again. If I have more problem I will let you know.
> 
> Zhaowen
> 
> On Thu, Jun 30, 2016 at 6:09 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:
> 
>    I don't think we have an example that does exactly that.
> 
>    If you are working with KSP directly and not SNES here is how to proceed
> 
>    KSPGetPC(ksp,&pc);
>    PCSetType(pc,PCASM);
>    KSPSetOperators()
>    KSPSetUp()                        <--- this must be called before the code below otherwise the subksps don't exist yet
> 
>    PetscInt n_local;
>    KSP     *subksps;
> 
>    PCASMGetSubKSP(pc,&n_local,NULL,&subksps);
>    for (i=0; i<n_local; i++) {
>       PC subpc;
> 
>       KSPGetPC(subksps[i],&subpc);
>       PCSetType(subpc,PCSHELL);
>       PCShellSetApply(subpc,yourapplyfunction);
>       anything else you need to set for your shell preconditioner here
>    }
>    KSPSetUpOnBlocks(ksp);
> 
>    KSPSolve();
> 
>    Now if you want to solve with a different right hand side or different entries in you matrix just call
>    KSPSolve() again  you don't need to repeat the code above.
> 
>   Barry
> 
>  Note that any of PETSc's preconditioners can be used on the subdomains so normally you can just -sub_pc_type typeyouwant and you don't need to mess with shell preconditioners.
> 
> 
> 
> 
> 
> 
> 
> > On Jun 30, 2016, at 5:49 PM, Duan Zhaowen <dzw.pku at gmail.com> wrote:
> >
> > Hi,
> >
> > I was trying to define a shell preconditioner for local partition, and let it work with global additive schwarz preconditioner for parallel computing. Is any one can give an example on this kind of preconditioners combination. Thanks!
> >
> > ZW
> 
> 



More information about the petsc-users mailing list