[petsc-users] Some clarification on PCMG

Aulisa, Eugenio eugenio.aulisa at ttu.edu
Fri Feb 19 11:38:13 CST 2016


Thanks Berry,

In Amat we have the weak form of a Jacobian with all Neumann Boundary conditions, 
while in Pmat we have the same Jacobian  with the correct  Dirichlet nodes imposed by penalty, 
So Pmat is the real problem and the fact that the preconditioned residual converges 
tell us that we are solving  the problem correctly. We know that for many other reasons.

The  residual update with a Neumann Amat on the Dirichlet nodes 
is used only to generate the correct residual restriction on the coarse nodes
whose support overlaps with the support of the fine boundary nodes.

The fact that the true residual does not converge is then correct,
because it involves the residual of those fake Neumann boundary nodes ,
that are never solved. 

Now, what puzzles us is that using Pmat or Amat does 
not change the convergence of the Multigrid not even of a single digit.
We would expect a better convergence using Amat, but it is not the case.

Thanks again,
Eugenio

Eugenio Aulisa

Department of Mathematics and Statistics,
Texas Tech University
Lubbock TX, 79409-1042
room: 226
http://www.math.ttu.edu/~eaulisa/
phone: (806) 834-6684
fax: (806) 742-1112




________________________________________
From: Barry Smith [bsmith at mcs.anl.gov]
Sent: Friday, February 19, 2016 10:42 AM
To: Aulisa, Eugenio
Cc: PETSc; simone.bna at gmail.com; Bornia, Giorgio
Subject: Re: [petsc-users] Some clarification on PCMG

> On Feb 19, 2016, at 8:55 AM, Aulisa, Eugenio <eugenio.aulisa at ttu.edu> wrote:
>
> Hi,
>
> I would like a clarification how PCMG updates
> the residual on the level-l before restricting it to the level l-1
> in the down-cycle.
>
> In particular I would expect that what is projected down is the true-residual,
> the one updated using Amat_l, however from what I see in my test it seams
> to be the residual updated using Pmat_l.
>
> My test is very simple:
> When setting up  the ksp object on the l-level, in one case I use
>
> KSPSetOperators(_ksp_l, _Amat_l, _Pmat_l);
>
> in the other case I use
>
> KSPSetOperators(_ksp_l, _Pmat_l, _Pmat_l);
>
> when I run my application I get exactly the same results for the preconditioned residual:
>
> these are the very last 2 outputs:
>
> with Amat,Pmat
> 9 KSP preconditioned resid norm 2.939296481331e-16 true resid norm 1.460300485277e-02 ||r(i)||/||b|| 9.999999999999e-01

 It uses the Amat to update the residual. Note that this has not worked essentially at all the true residual norm 1.46e-2 is huge while using the Pmat below it is 5e-15.  The fact that the "residual norm" is the same is a red herring.

  I suspect that the Amat, and the Pmat are "too far apart" hence using Amat to define the operator and Pmat to define the smoother produces a multigrid that really doesn't work.

  As a test try to make the Pmat be just a small perturbation of Amat (or if it is easier the other way around Amat a small perturbation of Pmat) then run and look at the true residual history.

  Barry


>
> with Pmat, Pmat:
> 9 KSP preconditioned resid norm 2.939296481331e-16 true resid norm 5.165834461813e-15 ||r(i)||/||b|| 3.537514719672e-13
>
> Since Amat and Pmat are 2 different matrices, and PCMG does not seams to see any difference,
> the only explanation I gave  myself is that the residual projected down to the coarser level is the
> the one updated using Pmat.
>
> Now, if this is the case, is there any way I can use Amat_l to update the residual?
>
> If this is not the case, and Amat_l is already used,
> do you have any other explanation why the results are identical?
>
> Thanks,
> Eugenio
>
>
>
>
>
>
>
>
>
> Eugenio Aulisa
>
> Department of Mathematics and Statistics,
> Texas Tech University
> Lubbock TX, 79409-1042
> room: 226
> http://www.math.ttu.edu/~eaulisa/
> phone: (806) 834-6684
> fax: (806) 742-1112
>
>
>



More information about the petsc-users mailing list