docs for MINRES, only for positive definite?

Barry Smith bsmith at mcs.anl.gov
Thu Feb 28 13:11:08 CST 2008


On Feb 28, 2008, at 11:44 AM, Matthew Knepley wrote:

> On Thu, Feb 28, 2008 at 11:35 AM, Barry Smith <bsmith at mcs.anl.gov>  
> wrote:
>>
>>   symmetry has nothing to do with it. Yes the matrix and
>> preconditioner must be
>> symmetric. The point is that the preconditioner has to also be
>> positive definite.
>> Because B is used in the algorithm to define a norm.
>
> I disagree. Here is a weaker condition:
>
>  1) A is symmetric indefinite and so is B
>
>  2) BA is also symmetric
>
> Since A and B are full rank, so is BA. Thus BA is symmetric indefinite
> and MINRES will work.


    But that is a different algorithm. That is running minres on the  
operator
BA with the preconditioner I which is very different than running
minres on A with the preconditioner B.


>
> My comment was that if B is SPD, then BA is equivalent to B^1/2 A
> B^1/2 and thus symmetric.
>
>   Matt
>
>>    Barry
>>
>>
>>
>> On Feb 28, 2008, at 11:28 AM, Matthew Knepley wrote:
>>
>>> On Thu, Feb 28, 2008 at 11:25 AM, Lisandro Dalcin
>>> <dalcinl at gmail.com> wrote:
>>>> Good point, the current code seems to require that...
>>>
>>> Its not the code, its the algorithm. It requires symmetry.
>>>
>>> Matt
>>>
>>>> ierr = KSP_PCApply(ksp,R,Z);CHKERRQ(ierr); /*     z  <- B*r        
>>>> */
>>>> ierr = VecDot(R,Z,&dp);CHKERRQ(ierr);
>>>> /*...*/
>>>> if (dp < 0.0) {
>>>>   ksp->reason = KSP_DIVERGED_INDEFINITE_PC;
>>>>   PetscFunctionReturn(0);
>>>> }
>>>>
>>>> Indeed, the following (simple minded, diagonal matrix) test fails
>>>> with
>>>> -pc_type jacobi, but success with -pc_type none
>>>>
>>>> import sys, petsc4py
>>>> petsc4py.init(sys.argv)
>>>> from petsc4py import PETSc
>>>> import numpy as N
>>>> A = PETSc.Mat().createAIJ([10,10])
>>>> for i in range(0,5):
>>>>   A[i,i] = -(i + 1)
>>>> for i in range(5,10):
>>>>   A[i,i] = +(i + 1)
>>>> A.assemble()
>>>> A.view()
>>>> x, b= A.getVecs()
>>>> b.set(1)
>>>> ksp = PETSc.KSP().create()
>>>> ksp.type = 'minres'
>>>> ksp.setOperators(A)
>>>> ksp.setFromOptions()
>>>> ksp.solve(b,x)
>>>>
>>>>
>>>>
>>>>
>>>> On 2/28/08, Barry Smith <bsmith at mcs.anl.gov> wrote:
>>>>>
>>>>>   But does it require a positive definite preconditioner?
>>>>>
>>>>>
>>>>
>>>>
>>>>>   Barry
>>>>>
>>>>>
>>>>> On Feb 28, 2008, at 9:32 AM, Matthew Knepley wrote:
>>>>>
>>>>>> Docs are wrong.
>>>>>>
>>>>>> Matt
>>>>>>
>>>>>> 2008/2/28 Lisandro Dalcin <dalcinl at gmail.com>:
>>>>>>> I've noticed that the docs for MINRES say that the operator and
>>>>>>> the
>>>>>>> preconditioner must be POSITIVE DEFINITE. But I understand
>>>>>>> MINRES is
>>>>>>> tailored for the symmetric/hermitian-indefinite case.
>>>>>>>
>>>>>>> Are the docs wrong? Or the actual code is a (very peculiar)  
>>>>>>> MINRES
>>>>>>> variant?
>>>>>>>
>>>>>>>
>>>>>>> --
>>>>>>> Lisandro Dalcín
>>>>>>> ---------------
>>>>>>> Centro Internacional de Métodos Computacionales en Ingeniería
>>>>>>> (CIMEC)
>>>>>>> Instituto de Desarrollo Tecnológico para la Industria Química
>>>>>>> (INTEC)
>>>>>>> Consejo Nacional de Investigaciones Científicas y Técnicas
>>>>>>> (CONICET)
>>>>>>> PTLC - Güemes 3450, (3000) Santa Fe, Argentina
>>>>>>> Tel/Fax: +54-(0)342-451.1594
>>>>>>>
>>>>>>>
>>>>>>
>>>>>>
>>>>>>
>>>>>> --
>>>>>> What most experimenters take for granted before they begin their
>>>>>> experiments is infinitely more interesting than any results to
>>>>>> which
>>>>>> their experiments lead.
>>>>>> -- Norbert Wiener
>>>>>>
>>>>>
>>>>>
>>>>
>>>>
>>>> --
>>>> Lisandro Dalcín
>>>> ---------------
>>>> Centro Internacional de Métodos Computacionales en Ingeniería  
>>>> (CIMEC)
>>>> Instituto de Desarrollo Tecnológico para la Industria Química  
>>>> (INTEC)
>>>> Consejo Nacional de Investigaciones Científicas y Técnicas  
>>>> (CONICET)
>>>> PTLC - Güemes 3450, (3000) Santa Fe, Argentina
>>>> Tel/Fax: +54-(0)342-451.1594
>>>>
>>>>
>>>
>>>
>>>
>>> --
>>> What most experimenters take for granted before they begin their
>>> experiments is infinitely more interesting than any results to which
>>> their experiments lead.
>>> -- Norbert Wiener
>>>
>>
>>
>
>
>
> -- 
> What most experimenters take for granted before they begin their
> experiments is infinitely more interesting than any results to which
> their experiments lead.
> -- Norbert Wiener
>




More information about the petsc-dev mailing list