[petsc-users] SNESSetConvergenceTest
Matthew Knepley
knepley at gmail.com
Wed Jan 22 17:45:03 CST 2014
On Wed, Jan 22, 2014 at 5:37 PM, Dharmendar Reddy
<dharmareddy84 at gmail.com>wrote:
> On Wed, Jan 22, 2014 at 5:33 PM, Matthew Knepley <knepley at gmail.com>
> wrote:
> >
> > On Wed, Jan 22, 2014 at 4:50 PM, Dharmendar Reddy <
> dharmareddy84 at gmail.com>
> > wrote:
> >>
> >> I take back my earlier email...
> >
> >
> > Do you have the threadcomm stuff turned on?
>
> I do not know what this means. Is related to compile time flags for my
> fortran code ? I link to mkl which is threaded. But otherwise, i
> guess, the answer is no.
>
The weird business looks like its in pthreads, so could you reconfigure
without MKL?
$PETSC_ARCH/conf/reconfigure-$PETSC_ARCH.py --download-f-blas-lapack
--PETSC_ARCH=arch-testing
It would not surprise me if this was the problem.
Matt
>
> >
> > Matt
> >
> >>
> >> I see the error irrespective of valgrind options. I did a clean
> >> compile to test it few times now.
> >> See the attached valgrind.log
> >>
> >> Now the difference between, ex5f90.F and my code is that i use
> >> DMSNESSetFunctionLocal in the case where the code fails.
> >>
> >> I have also attached the solver.F90 which wraps the calls to snes, Now
> >> the solver has to interface (1) for equations which are handled with a
> >> DM and (2) for equation with out DM , here i use SNESSetFunction, and
> >> this part of the code runs fine with the test case in the code which
> >> can be compiled with -DTestSolver and standard petsc makefile
> >>
> >> The case with DM is linked deeply into my code, i am trying to get to
> >> a test case for that.
> >>
> >> thanks
> >> Reddy
> >>
> >> On Wed, Jan 22, 2014 at 4:30 PM, Barry Smith <bsmith at mcs.anl.gov>
> wrote:
> >> >
> >> > The thing to use in valgrind is -q --tool=memcheck
> >> >
> >> > This checks for corrupted and illegally accessed memory,
> >> > uninitialized memory etc It shows the exact lines that memory is used
> >> > incorrectly so you can fix your code.
> >> >
> >> > Barry
> >> >
> >> > On Jan 22, 2014, at 4:02 PM, Dharmendar Reddy <
> dharmareddy84 at gmail.com>
> >> > wrote:
> >> >
> >> >> Hello,
> >> >> The code works fine if i run with valgrind --leak-check=yes
> >> >> exename args
> >> >> it fails if i run valgrind exename args
> >> >>
> >> >> Any suggestions on what to look for my code?
> >> >>
> >> >> Thanks
> >> >> Reddy
> >> >>
> >> >> On Wed, Jan 22, 2014 at 3:31 PM, Dharmendar Reddy
> >> >> <dharmareddy84 at gmail.com> wrote:
> >> >>> Hello,
> >> >>> I tried using the SNESSetConvergenceTest in ex5f90.F using my
> >> >>> subroutine. Works as expected. I will try to run my code through
> >> >>> valgrind.
> >> >>>
> >> >>> Thanks
> >> >>> Reddy
> >> >>>
> >> >>> On Wed, Jan 22, 2014 at 2:33 PM, Dharmendar Reddy
> >> >>> <dharmareddy84 at gmail.com> wrote:
> >> >>>> Hello,
> >> >>>> I will also try to create a stand alone test case.
> >> >>>>
> >> >>>> Thanks
> >> >>>> Reddy
> >> >>>>
> >> >>>> On Wed, Jan 22, 2014 at 1:53 PM, Matthew Knepley <
> knepley at gmail.com>
> >> >>>> wrote:
> >> >>>>> On Wed, Jan 22, 2014 at 1:52 PM, Barry Smith <bsmith at mcs.anl.gov>
> >> >>>>> wrote:
> >> >>>>>>
> >> >>>>>>
> >> >>>>>> We need more information to track this down. Could you make a
> >> >>>>>> simple
> >> >>>>>> standalone code that reproduces the problem? Just creates a this,
> >> >>>>>> creates a
> >> >>>>>> SNES, sets the function, sets the test code and then crashes when
> >> >>>>>> runs.
> >> >>>>>>
> >> >>>>>> But since it is not crashing in the SNESConvergenceTest_snes
> >> >>>>>> routine I
> >> >>>>>> suspect it may be memory corruption. Have you run your code under
> >> >>>>>> valgrind?
> >> >>>>>> http://www.mcs.anl.gov/petsc/documentation/faq.html#valgrind
> >> >>>>>
> >> >>>>>
> >> >>>>> I am putting his convergence test into ex5f right now. This can be
> >> >>>>> valgrinded. I will send an update when it runs.
> >> >>>>>
> >> >>>>> Matt
> >> >>>>>
> >> >>>>>>
> >> >>>>>>
> >> >>>>>> Barry
> >> >>>>>>
> >> >>>>>> On Jan 20, 2014, at 7:14 PM, Dharmendar Reddy
> >> >>>>>> <dharmareddy84 at gmail.com>
> >> >>>>>> wrote:
> >> >>>>>>
> >> >>>>>>> Hello,
> >> >>>>>>> I am getting a segmentation fault when i use
> >> >>>>>>> SNESSetConvergenceTest from Fortran.
> >> >>>>>>>
> >> >>>>>>> I use the following:
> >> >>>>>>>
> >> >>>>>>> call SNESSetConvergenceTest(this%snes, SNESConvergenceTest_snes,
> >> >>>>>>> this,PETSC_NULL_FUNCTION, ierr) (Line 1)
> >> >>>>>>>
> >> >>>>>>> where the soubtoutine is
> >> >>>>>>>
> >> >>>>>>> subroutine
> >> >>>>>>> SNESConvergenceTest_snes(snes,it,xnorm,gnorm,fnorm,reason,this,
> >> >>>>>>> ierr)
> >> >>>>>>> implicit none
> >> >>>>>>> #include "finclude/petsc.h"
> >> >>>>>>> ! IO Variables
> >> >>>>>>> SNES :: snes
> >> >>>>>>> integer :: it
> >> >>>>>>> PetscReal :: xnorm
> >> >>>>>>> PetscReal :: gnorm
> >> >>>>>>> PetscReal :: fnorm
> >> >>>>>>> SNESConvergedReason :: reason
> >> >>>>>>> type(Solver_t) :: this
> >> >>>>>>> PetscErrorCode :: ierr
> >> >>>>>>> ! Local Variables
> >> >>>>>>> Vec :: X, dX
> >> >>>>>>> Vec :: F, W, G
> >> >>>>>>> real(WP) :: lambda
> >> >>>>>>> print*, 'Calling custom test',it,xnorm,fnorm,gnorm
> >> >>>>>>> !call SNESLineSearchGetLambda(this%linesearch, lambda, ierr)
> >> >>>>>>> !call SNESLineSearchGetVecs(this%linesearch, X, F, dX, W, G,
> ierr)
> >> >>>>>>> !call this%SNESConvergenceTest(it, fnorm, lambda, X, dX, reason,
> >> >>>>>>> ierr)
> >> >>>>>>> reason = SNES_DIVERGED_FUNCTION_COUNT
> >> >>>>>>> end subroutine SNESConvergenceTest_snes
> >> >>>>>>>
> >> >>>>>>> Code works fine if the call at Line 1 is commented i.e, if i
> >> >>>>>>> donot
> >> >>>>>>> set the convergencetest
> >> >>>>>>>
> >> >>>>>>> I have attached the stack trace when the code fails
> >> >>>>>>> (traceFail.log)
> >> >>>>>>> and stack trace of the code which runs well (tracePass.log).
> >> >>>>>>>
> >> >>>>>>>
> >> >>>>>>>
> >> >>>>>>> Thanks
> >> >>>>>>> Reddy
> >> >>>>>>> <traceFail.log><tracePass.log>
> >> >>>>>>
> >> >>>>>
> >> >>>>>
> >> >>>>>
> >> >>>>> --
> >> >>>>> 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
> >> >>>>
> >> >>>>
> >> >>>>
> >> >>>> --
> >> >>>> -----------------------------------------------------
> >> >>>> Dharmendar Reddy Palle
> >> >>>> Graduate Student
> >> >>>> Microelectronics Research center,
> >> >>>> University of Texas at Austin,
> >> >>>> 10100 Burnet Road, Bldg. 160
> >> >>>> MER 2.608F, TX 78758-4445
> >> >>>> e-mail: dharmareddy84 at gmail.com
> >> >>>> Phone: +1-512-350-9082
> >> >>>> United States of America.
> >> >>>> Homepage: https://webspace.utexas.edu/~dpr342
> >> >>>
> >> >>>
> >> >>>
> >> >>> --
> >> >>> -----------------------------------------------------
> >> >>> Dharmendar Reddy Palle
> >> >>> Graduate Student
> >> >>> Microelectronics Research center,
> >> >>> University of Texas at Austin,
> >> >>> 10100 Burnet Road, Bldg. 160
> >> >>> MER 2.608F, TX 78758-4445
> >> >>> e-mail: dharmareddy84 at gmail.com
> >> >>> Phone: +1-512-350-9082
> >> >>> United States of America.
> >> >>> Homepage: https://webspace.utexas.edu/~dpr342
> >> >>
> >> >>
> >> >>
> >> >> --
> >> >> -----------------------------------------------------
> >> >> Dharmendar Reddy Palle
> >> >> Graduate Student
> >> >> Microelectronics Research center,
> >> >> University of Texas at Austin,
> >> >> 10100 Burnet Road, Bldg. 160
> >> >> MER 2.608F, TX 78758-4445
> >> >> e-mail: dharmareddy84 at gmail.com
> >> >> Phone: +1-512-350-9082
> >> >> United States of America.
> >> >> Homepage: https://webspace.utexas.edu/~dpr342
> >> >
> >>
> >>
> >>
> >> --
> >> -----------------------------------------------------
> >> Dharmendar Reddy Palle
> >> Graduate Student
> >> Microelectronics Research center,
> >> University of Texas at Austin,
> >> 10100 Burnet Road, Bldg. 160
> >> MER 2.608F, TX 78758-4445
> >> e-mail: dharmareddy84 at gmail.com
> >> Phone: +1-512-350-9082
> >> United States of America.
> >> Homepage: https://webspace.utexas.edu/~dpr342
> >
> >
> >
> >
> > --
> > 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
>
>
>
> --
> -----------------------------------------------------
> Dharmendar Reddy Palle
> Graduate Student
> Microelectronics Research center,
> University of Texas at Austin,
> 10100 Burnet Road, Bldg. 160
> MER 2.608F, TX 78758-4445
> e-mail: dharmareddy84 at gmail.com
> Phone: +1-512-350-9082
> United States of America.
> Homepage: https://webspace.utexas.edu/~dpr342
>
--
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
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140122/39d46b17/attachment-0001.html>
More information about the petsc-users
mailing list