SNES problem
JoaoFlavio Vasconcellos
joaoflavio.vasconcellos at jcu.edu.au
Wed Oct 3 16:32:14 CDT 2007
Dear Lisandro
FormFunction1 is called only once. I'm quite sure that F is a non zero vector. I'd printed it and its values are 100% rigth.
I've been used snes functions as written above:
...
...
....
ierr = SNESCreate(PETSC_COMM_WORLD, &snes);CHKERRQ(ierr);
ierr = SNESSetFunction(snes, resid, FormFunction1, this);CHKERRQ(ierr);
....
.... // Calculating c_ctx.A
....
ierr = FormJacobian1(snes, resid, &c_ctx.A, &c_ctx.A, &mat_str, this); CHKERRQ(ierr);
ierr = SNESSetJacobian(snes, c_ctx.A, c_ctx.A, FormJacobian1, this); CHKERRQ(ierr);
ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
cout << "Calling SNESSolve" << endl;
PetscScalar zero = 0;
ierr = VecSet(C, zero);CHKERRQ(ierr);
ierr = SNESSolve(snes,PETSC_NULL, C);CHKERRQ(ierr);
ierr = SNESGetIterationNumber(snes, &its);CHKERRQ(ierr);
ierr = SNESGetSolutionUpdate(snes, &C);CHKERRQ(ierr);
VecView (C, 0);
...
...
...
All outputs of VevView in FormFunction1 are ok. The output of VecView (C, 0); after call SNESSolve is an unexpected zero vector.
Thanks
JF
Joao Flavio Vieira de Vasconcellos
School of Engineering
James Cook University
Phone: +61 7 4781 4340
Fax: +61 7 4775 1184
Email: jflavio at iprj.uerj.br
joaoflavio.vasconcellos at jcu.edu.au
skype: jflavio_vasconcellos
msn: jf_vasconcellos at hotmail.com
---- Original message ----
>Date: Wed, 3 Oct 2007 12:40:28 -0300
>From: "Lisandro Dalcin" <dalcinl at gmail.com>
>Subject: Re: SNES problem
>To: petsc-users at mcs.anl.gov
>
>Flavio, are you completelly sure your FormFunction1() is actually
>being called? I can see the SNES monitor at iter 0, but I did not see
>in the the outputs of VecView() or your 'cout << "Residual " << endl'.
>Did you called SNESSetFunction() to set your FormFunction1()??
>
>On 10/3/07, JoaoFlavio Vasconcellos <joaoflavio.vasconcellos at jcu.edu.au> wrote:
>> Hello
>>
>> I'm trying to use SNES to solve a problem. I wrote the code above to evaluate the residual of the linear system. In this routine norm_2 of F is different from zero but outside of this routine norm_2 of F is zero.
>>
>> #undef __FUNCT__
>> #define __FUNCT__ "FormFunction1"
>> PetscErrorCode FormFunction1(SNES snes, Vec X, Vec F, void *_simulator) {
>> Simulator *simulator = (Simulator*) _simulator;
>> PetscErrorCode ierr;
>>
>> VecView (X, 0);
>> ierr = simulator->CalculateMatrixCoefficient(X); CHKERRQ(ierr);
>> ierr = simulator->CalculateSourceTerm(X); CHKERRQ(ierr);
>> ierr = simulator->CalculateResidual (X, &F); CHKERRQ(ierr);
>>
>> cout << "Residual " << endl;
>> VecView (F, 0); // norm_2 of F is != zero
>> return 0;
>>
>> }
>>
>> // End of routine
>>
>>
>>
>>
>> Messages from -snes_view:
>>
>> 0 SNES Function norm 0.000000000000e+00
>> SNES Object:
>> type: ls
>> line search variant: SNESLineSearchCubic
>> alpha=0.0001, maxstep=1e+08, steptol=1e-12
>> maximum iterations=50, maximum function evaluations=10000
>> tolerances: relative=1e-08, absolute=1e-50, solution=1e-08
>> total number of linear solver iterations=0
>> total number of function evaluations=1
>> KSP Object:
>> type: gmres
>> GMRES: restart=30, using Classical (unmodified) Gram-Schmidt Orthogonalization with no iterative refinement
>> GMRES: happy breakdown tolerance 1e-30
>> maximum iterations=10000, initial guess is zero
>> tolerances: relative=1e-05, absolute=1e-50, divergence=10000
>> left preconditioning
>> PC Object:
>> type: ilu
>> ILU: 0 levels of fill
>> ILU: factor fill ratio allocated 1
>> ILU: tolerance for zero pivot 1e-12
>> out-of-place factorization
>> matrix ordering: natural
>>
>> // end of message
>>
>> I dont know why "SNES Function norm 0.000000000000e+00" if
>> F's norm_2 inside my routine is not. Has anyone that could sugest how to fix this problem ?
>>
>>
>> Thank you in advance.
>>
>>
>> Joao Flavio Vieira de Vasconcellos
>> School of Engineering
>> James Cook University
>>
>> Phone: +61 7 4781 4340
>> Fax: +61 7 4775 1184
>> Email: jflavio at iprj.uerj.br
>> joaoflavio.vasconcellos at jcu.edu.au
>> skype: jflavio_vasconcellos
>> msn: jf_vasconcellos at hotmail.com
>>
>>
>
>
>--
>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
>
More information about the petsc-users
mailing list