[petsc-users] Preconditioner in multigrid solver

Smith, Barry F. bsmith at mcs.anl.gov
Mon Mar 11 13:05:39 CDT 2019



> On Mar 11, 2019, at 9:42 AM, Pietro Benedusi via petsc-users <petsc-users at mcs.anl.gov> wrote:
> 
> Dear Petsc team,
> 
> I have a question about the setting up of a multigrid solver.
> 
> I would like yo use a PCG smoother, preconditioned with a mass matrix, just on the fine level.
> But when add the line for preconditioning the CG with the mass matrix my MG diverges.
> 
> I have implemented the same solver in MATLAB and it converges fine. Also the operators in PETSc are the same and the PCG applied directly on the problem (without MG) works the same in both PETSC and MATLAB. 
> 
> This is what I do in PETSC for 2 levels:
> 
> 
>     KSP space_solver;       
> 
>     ierr = KSPCreate(PETSC_COMM_WORLD,&space_solver);CHKERRQ(ierr);
>     ierr = KSPSetTolerances(space_solver,1e-8,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);          
>     ierr = KSPSetOperators(space_solver, K, K);CHKERRQ(ierr); 
>     ierr = KSPSetNormType(space_solver,  KSP_NORM_UNPRECONDITIONED );CHKERRQ(ierr);  
> 
>     ierr = KSPSetType(space_solver,KSPRICHARDSON);CHKERRQ(ierr);    
>     ierr = KSPSetFromOptions(space_solver);CHKERRQ(ierr);         
>     ierr = KSPSetUp(space_solver);CHKERRQ(ierr);        
> 
>     PC pcmg;
>     ierr = KSPGetPC(space_solver, &pcmg);
>     ierr = PCSetType(pcmg, PCMG);    
>     ierr = PCMGSetLevels(pcmg,levels, NULL);CHKERRQ(ierr);
>     ierr = PCMGSetGalerkin(pcmg,PC_MG_GALERKIN_BOTH);CHKERRQ(ierr);
> 
>     // smoothers
>     for (int i = 1; i < levels; ++i)  
>     {
>         KSP smoother;
>         ierr = PCMGGetSmoother(pcmg, i, &smoother);CHKERRQ(ierr);
> 
>             ierr = KSPSetType(smoother, KSPCG);CHKERRQ(ierr); 
>             ierr = KSPSetOperators(smoother, K, M);CHKERRQ(ierr);   


   I'm not sure what you mean by "preconditioned with a mass matrix" but I don't think this line makes sense. The last argument to this call is the matrix from which the preconditioner is constructed.  I don't think it ever makes sense to precondition a Laplacian with a mass matrix; they are very different beasts.


   Barry


>          
> 
>             // ierr = KSPSetUp(smoother);CHKERRQ(ierr);              
> 
>             ierr = KSPSetTolerances(smoother,1e-12,PETSC_DEFAULT,PETSC_DEFAULT,s_p);CHKERRQ(ierr);               
>             ierr = KSPSetNormType(smoother, KSP_NORM_NONE);CHKERRQ(ierr);  
> 
>             PC sm;
>             ierr = KSPGetPC(smoother, &sm);CHKERRQ(ierr);            
>             ierr = PCSetType(sm,   PCLU);CHKERRQ(ierr);       
>                         
>         ierr = PCMGSetInterpolation(pcmg, i, interpolation_operators[i-1]);CHKERRQ(ierr);                   
>     }
> 
> 
> 
> I think there is a problem with the PETSc syntax, because I checked everything else and it is fine.
> 
> Do you any ideas?
> 
> Thank you very much!
> 
> Best,
> Pietro 
> 
> ~~~~~~~~~~~~
> Pietro Benedusi
> 
> Numerical Simulation in Science,
> Medicine and Engineering research group
> ICS, Institute of Computational Science
> USI, Università della Svizzera Italiana
> Via Giuseppe Buffi, 13
> CH - 6900 Lugano
> benedp at usi.ch
> 
> 



More information about the petsc-users mailing list