[petsc-dev] [petsc-maint #78191] Re: [petsc-users] Reusing ML preconditioner

Barry Smith bsmith at mcs.anl.gov
Fri Jul 1 22:48:17 CDT 2011


  John,

    I have determined the issue and it is NOT an error in our PCML. It is doing what it is suppose to do, it is using the entire ML preconditioner that was built during the first build of the preconditioner. 

    The reason for the inconsistent behavior and and the difference between the two runs is due to the fact that the "built" ML preconditioner contains as part of itself the original matrix; it is used for smoothing on the finest level and for computing the residual on the finest level.  What your code is doing is changing the original matrix to contain the new matrix values so you are actually changing the ML preconditioner; it is exactly the same as before on the coarser levels but on the finest level it is now using the new matrix for smoothing and residual computations. Naturally it works much better as a preconditioner for the new problem.

    The currently "correct" way for you to get the effect you want is to instead manually change the smoother and residual computation on the finest level before using the preconditioner again: For example I added the following to your code and it results in the same behavior for both "versions" of your code.

if (i)
      {
          ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
          PC pc;
          ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
          PetscInt nlevels;
          ierr = PCMGGetLevels(pc,&nlevels);CHKERRQ(ierr);
          KSP sksp;
          ierr = PCMGGetSmoother(pc,nlevels-1,&sksp);CHKERRQ(ierr);
          ierr = KSPSetOperators(sksp,A[i],A[i],SAME_NONZERO_PATTERN);CHKERRQ(ierr);
          ierr = PCMGSetResidual(pc,nlevels-1,PCMGDefaultResidual,A[i]);CHKERRQ(ierr);

Note that this is not technically reusing the exact old preconditioner, it is actually partially modifying the old preconditioner with new matrix information.


But now thinking about it, it would be nice if this was the default behavior and the user doesn't need to do the trick above.  The question is how to achieve this I could

1) provide a PCMGSetNewFinest() that does what the code above does or 
2) somehow have it happen automatically. This is not possible because when SAME_PRECONDITIONER is used the PC_MG doesn't even know that the operator has changed and so cannot change 
       the smoother on the finest level, plus it really isn't a SAME_PRECONDITIONER so is philosophically wrong or
3) have a new flag for KSPSetOperators  ALMOST_SAME_PRECONDITIONER that causes 1 to happen.

  Note that in the normal model of using KSP where one changes values in the matrix and calls KSPSetOperators() repeatedly with that same matrix the new smoother will be automatically used with SOR smoothing (though not with ILU) since the PC_MG points to the same matrix that has the new values in it.

  Remind me why are you passing entirely new matrices to KSPSetOperators instead of the same matrix with new values in it?

  Thanks for bringing this to our attention, we need to think about how to proceed. 


   Barry



On Jul 1, 2011, at 11:08 AM, John Fettig via RT wrote:

> Just in case the attachment didn't make it to you (I should have
> looked at the filesize), I've put it on dropbox:
> 
> http://dl.dropbox.com/u/29777859/reuse_ml.tar.gz
> 
> John
> 
> On Fri, Jul 1, 2011 at 12:04 PM, John Fettig <john.fettig at gmail.com> wrote:
>> Here is a small code that reproduces the problem I'm seeing.  I've
>> written out two matrices, one from the first timestep and one from the
>> second timestep.  This code reads them in and solves them in two
>> different ways.  If you run with
>> 
>> mpirun -np 1 ./reuse_ml -ksp_type cg -pc_type ml
>> 
>> you will see that the first time through (which corresponds to
>> "version 2" of my code in the petsc-users correspondence), the second
>> system takes significantly more iterations than the second time
>> through (which corresponds to "version 1" of my code).  Another thing
>> that really bothers me is that the residual from Solve 0 is different
>> between the two versions....even though both the matrix and the right
>> hand side passed to the KSP are identical and everything has been
>> destroyed and re-created.  I would love to know how to reuse this
>> preconditioner effectively.
>> 
>> John
>> 
>> 
>> On Thu, Jun 30, 2011 at 5:20 PM, John Fettig <john.fettig at gmail.com> wrote:
>>> On Thu, Jun 30, 2011 at 5:11 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:
>>>> 
>>>> On Jun 29, 2011, at 10:35 PM, John Fettig wrote:
>>>> 
>>>>> On Wed, Jun 29, 2011 at 11:21 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:
>>>>>>  If you pass in SAME_PRECONDITIONER to the KSPSetOperators() it is simply impossible that it will rebuild the preconditioner. There must be another explanation as to what is going on. Are you calling anything else on the KSP?
>>>>> 
>>>>> You are right.  I was playing around with options and had set it to
>>>>> SAME_NONZERO_PATTERN for the first iteration at each timestep to make
>>>>> sure that would work (it does), and forgot to change it back.
>>>>> 
>>>>> When it is actually set to SAME_PRECONDITIONER and I duplicate Amat to
>>>>> Pmat and pass them to KSPSetOperators() before every KSPSolve, it
>>>>> definitely uses the same preconditioner (as reported by -info), but
>>>>> unfortunately it reverts back to the behavior I saw previously.
>>>> 
>>>>   What do you mean? Reverts back to? After some iterations? Copying the matrix and passing it into KSPSetOperators() as the Pmat cannot make any difference.
>>> 
>>> I mean that it doesn't make any difference, as you expect.
>>> 
>>>>   Can you send us a code that reproduces the problem you think is happening (petsc-maint at mcs.anl.gov), I cannot get a consistent handle on the issue.
>>> 
>>> I will see what I can do.
>>> 
>>> John
>>> 
>> 
> 




More information about the petsc-dev mailing list