[petsc-users] KSP Solver giving inf's

Kaushik Kulkarni kaushikggg at gmail.com
Sat Apr 8 10:20:28 CDT 2017


Hello Matthew,
The same happened for:

PetscReal   loc[]  = { 0.0,  2.0,  1.0,
                           1.0, -2.0, -3.0,
                          -1.0,  1.0,  2.0};
PetscReal   b_array[] = { -8.0,
                                    0.0,
                                    3.0};

And this time I checked the matrix to be non-singular(determinant=1). I am
again getting (inf, inf, inf) as the solution.

On running with -ksp_converged_reason the following message comes:
Linear solve did not converge due to DIVERGED_PCSETUP_FAILED iterations 0
PCSETUP_FAILED due to FACTOR_NUMERIC_ZEROPIVOT

I am once again appending the whole code for reference.

Thanks.
---------

static char help[] = "KSP try.\n\n";


#include "petsc.h"
#include "petscksp.h"

int main(int argc, char *argv[])
{
    /// matrix creation variables.
    PetscInt    *idxm = new PetscInt[3];
    PetscInt    *idxn = new PetscInt[3];
    PetscReal   loc[]  = { 0.0,  2.0,  1.0,
                                 1.0, -2.0, -3.0,
                                -1.0,  1.0,  2.0};
    PetscReal   b_array[] = { -8.0,
                               0.0,
                               3.0};
    PetscInt i;
    KSP      ksp;

    /// Declaring the vectors
    Vec     x, b;

    // Declaring matrices
    Mat A;

    PetscInitialize(&argc,&argv,(char*)0,help);
    // Creating vectors
    VecCreateSeq(PETSC_COMM_SELF, 3, &x);
    VecCreateSeq(PETSC_COMM_SELF, 3, &b);
    // Creating matrix
    MatCreateSeqAIJ(PETSC_COMM_SELF, 3, 3, 3, NULL, &A);
    // Creating the indices
    for(i=0; i<3; i++) {
        idxm[i] = i;
        idxn[i] = i;
    }
    // Assembling the vector b and x
    VecSetValues(b, 3, idxm, b_array, INSERT_VALUES);
    VecAssemblyBegin(b);
    VecAssemblyEnd(b);

    //Assembling the Matrix
    MatSetValues(A, 3, idxm, 3, idxn, loc, INSERT_VALUES);
    MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
    MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);

    // KSP related operations
    KSPCreate(PETSC_COMM_SELF, &ksp);
    KSPSetType(ksp, KSPGMRES);
    KSPSetOperators(ksp, A, A);
    KSPSetFromOptions(ksp);
    KSPSolve(ksp,b,x);
    KSPDestroy(&ksp);

    VecView(x, PETSC_VIEWER_STDOUT_SELF);

    PetscFinalize();
    return 0;
}//End of file
--------


On Sat, Apr 8, 2017 at 7:01 PM, Matthew Knepley <knepley at gmail.com> wrote:

> On Sat, Apr 8, 2017 at 8:26 AM, Kaushik Kulkarni <kaushikggg at gmail.com>
> wrote:
>
>> I guess there is a miscalculation:
>> 2*row1 + row2 = (4, 0, 0)
>>
>
> You are right here, but I think it is still rank deficient.
>
>   Thanks,
>
>     Matt
>
>
>> Thanks,
>> Kaushik
>>
>> On Apr 8, 2017 6:52 PM, "Matthew Knepley" <knepley at gmail.com> wrote:
>>
>> On Sat, Apr 8, 2017 at 6:54 AM, Kaushik Kulkarni <kaushikggg at gmail.com>
>> wrote:
>>
>>> Hello,
>>> I just started with PETSc and was trying the KSP Solvers. I was trying
>>> for the problem -
>>>
>>> int main(int argc, char *argv[])
>>> {
>>>     /// matrix creation variables.
>>>     PetscInt    *idxm = new PetscInt[3];
>>>     PetscInt    *idxn = new PetscInt[3];
>>>
>>
>> Note that the problem below is rank deficient (2 * row 1 + row 2 = 0).
>> Its not clear to me whether b is
>> in the range space of the operator.
>>
>>   Thanks,
>>
>>     Matt
>>
>>
>>>     PetscReal   loc[]  = { 1.0, -2.0, -6.0,
>>>                                    2.0,  4.0, 12.0,
>>>                                    1.0, -4.0,-12.0};
>>>     PetscReal   b_array[] = {  12.0,
>>>                               -17.0,
>>>                                22.0};
>>>     PetscInt i;
>>>     KSP      ksp;
>>>
>>>     /// Declaring the vectors
>>>     Vec     x, b;
>>>
>>>     // Declaring matrices
>>>     Mat A;
>>>
>>>     PetscInitialize(&argc,&argv,(char*)0,help);
>>>     // Creating vectors
>>>     VecCreateSeq(PETSC_COMM_SELF, 3, &x);
>>>     VecCreateSeq(PETSC_COMM_SELF, 3, &b);
>>>     // Creating matrix
>>>     MatCreateSeqAIJ(PETSC_COMM_SELF, 3, 3, 3, NULL, &A);
>>>     // Creating the indices
>>>     for(i=0; i<3; i++) {
>>>         idxm[i] = i;
>>>         idxn[i] = i;
>>>     }
>>>     // Assembling the vector b and x
>>>     VecSetValues(b, 3, idxm, b_array, INSERT_VALUES);
>>>     VecAssemblyBegin(b);
>>>     VecAssemblyEnd(b);
>>>
>>>     //Assembling the Matrix
>>>     MatSetValues(A, 3, idxm, 3, idxn, loc, INSERT_VALUES);
>>>     MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
>>>     MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
>>>
>>>     // KSP related operations
>>>     KSPCreate(PETSC_COMM_SELF, &ksp);
>>>     KSPSetType(ksp, KSPGMRES);
>>>     KSPSetOperators(ksp, A, A);
>>>     KSPSetFromOptions(ksp);
>>>     KSPSolve(ksp,b,x);
>>>     KSPDestroy(&ksp);
>>>
>>>     VecView(x, PETSC_VIEWER_STDOUT_SELF);
>>>
>>>     PetscFinalize();
>>>     return 0;
>>> }
>>>
>>> But the obtained solution is found out to be- (inf, inf, inf).
>>>
>>> I wanted to know whether I am doing something wrong or is the problem
>>> inherently not solvable using GMRES. Currently I am running the code in a
>>> sequential manner(no parallelism intended).
>>>
>>> Thank you.
>>>
>>
>>
>>
>> --
>> 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
>



-- 
Kaushik Kulkarni
Fourth Year Undergraduate
Department of Mechanical Engineering
Indian Institute of Technology Bombay
Mumbai, India
https://kaushikcfd.github.io/About/
+91-9967687150
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20170408/9028fddc/attachment-0001.html>


More information about the petsc-users mailing list