[petsc-users] VecValidValues() reports NaN found

linjing bo francium87 at hotmail.com
Thu May 8 05:48:05 CDT 2014


Sorry, forgot to attatch the matrix file
 
Date: Thu, 8 May 2014 05:20:52 -0500
Subject: Re: [petsc-users] VecValidValues() reports NaN found
From: knepley at gmail.com
To: francium87 at hotmail.com
CC: petsc-users at mcs.anl.gov

On Thu, May 8, 2014 at 2:49 AM, linjing bo <francium87 at hotmail.com> wrote:







Yes, I called KSPDestroy(). I have reproduce the problem using a small C code, this code with default ilu preconditioner will show an error
Can you also send your matrix so I can run it?

  Thanks,
     Matt  
[0]PETSC ERROR: --------------------- Error Message ------------------------------------
[0]PETSC ERROR: Argument out of range!
[0]PETSC ERROR: Cannot log negative flops!
[0]PETSC ERROR: ------------------------------------------------------------------------

[0]PETSC ERROR: Petsc Release Version 3.4.4, Mar, 13, 2014
[0]PETSC ERROR: See docs/changes/index.html for recent updates.
[0]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
[0]PETSC ERROR: See docs/index.html for manual pages.

[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: ./bptest_c on a arch-linux2-c-debug named node2.indac.info by jlin Thu May  8 15:43:13 2014

[0]PETSC ERROR: Libraries linked from /opt/sfw/petsc/3.4.4/intel/openmpi/lib
[0]PETSC ERROR: Configure run at Sat Apr 26 20:19:41 2014
[0]PETSC ERROR: Configure options --prefix=/opt/sfw/petsc/3.4.4/intel/openmpi --with-mpi-dir=/opt/sfw/openmpi/1.6.3/intel --with-blas-lapack-dir=/opt/sfw/intel/composer_xe_2011_sp1.7.256/mkl/lib/intel64 --with-mpiexec=mpiexec

[0]PETSC ERROR: ------------------------------------------------------------------------
[0]PETSC ERROR: PetscLogFlops() line 204 in /tmp/petsc-3.4.4/include/petsclog.h
[0]PETSC ERROR: MatLUFactorNumeric_SeqAIJ() line 552 in /tmp/petsc-3.4.4/src/mat/impls/aij/seq/aijfact.c

[0]PETSC ERROR: MatLUFactorNumeric() line 2889 in /tmp/petsc-3.4.4/src/mat/interface/matrix.c
[0]PETSC ERROR: PCSetUp_ILU() line 232 in /tmp/petsc-3.4.4/src/ksp/pc/impls/factor/ilu/ilu.c
[0]PETSC ERROR: PCSetUp() line 890 in /tmp/petsc-3.4.4/src/ksp/pc/interface/precon.c

[0]PETSC ERROR: KSPSetUp() line 278 in /tmp/petsc-3.4.4/src/ksp/ksp/interface/itfunc.c
[0]PETSC ERROR: KSPSolve() line 399 in /tmp/petsc-3.4.4/src/ksp/ksp/interface/itfunc.c
[0]PETSC ERROR: main() line 46 in "unknowndirectory/"test.c

===================================================Below is the code-------------------------------------------------static char help[] = "Solve";

#include <petsc.h>int main(int argc, char **args){
        Vec x,b,u;
        Mat A;
        KSP ksp;
        PC  pc;
        PetscViewer fd;
        PetscErrorCode ierr;
        PetscReal tol=1.e-4;

        PetscScalar one = 1.0;
        PetscInt n=1023;        PetscInitialize(&argc,&args,(char*)0,help);        ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);

        ierr = PetscObjectSetName((PetscObject) x, "Solution");CHKERRQ(ierr);
        ierr = VecSetSizes(x,PETSC_DECIDE,n);CHKERRQ(ierr);
        ierr = VecSetFromOptions(x);CHKERRQ(ierr);
        ierr = VecDuplicate(x,&b);CHKERRQ(ierr);
        PetscViewerBinaryOpen(PETSC_COMM_WORLD, "tor0bp.bin", FILE_MODE_READ, &fd);
        ierr = MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n, \
                        11,PETSC_NULL,11,PETSC_NULL,&A);

        ierr = MatLoad(A, fd);
        PetscViewerDestroy(&fd);        VecSet( b, one);
        VecSet( x, one);
        VecAssemblyBegin(b);
        VecAssemblyEnd(b);
        VecAssemblyBegin(x);

        VecAssemblyEnd(x);        ierr = KSPCreate(PETSC_COMM_WORLD,&ksp);CHKERRQ(ierr);        ierr = KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr);
        ierr = KSPSetTolerances(ksp,1.e-5,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr);        ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);        ierr = KSPSetInitialGuessNonzero(ksp,PETSC_TRUE);CHKERRQ(ierr);
        ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr);        ierr = KSPView(ksp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);

        KSPDestroy(&ksp);
        VecDestroy(&x);

        VecDestroy(&b);
        PetscFinalize();
        return 0;
}-----------------------------------------------------------

Date: Wed, 7 May 2014 05:52:13 -0500

Subject: Re: [petsc-users] VecValidValues() reports NaN found
From: knepley at gmail.com
To: francium87 at hotmail.com

CC: petsc-users at mcs.anl.gov

On Tue, May 6, 2014 at 10:53 PM, linjing bo <francium87 at hotmail.com> wrote:





The Valgrind shows memory leak in memalign() called by KSPSetup and PCSetup. Is that normal?

Did you call KSPDestroy()?


  Matt  
---------------------------------------------------------------------------------------------


==31551== 136 bytes in 1 blocks are definitely lost in loss record 2,636 of 3,327
==31551==    at 0x4A05458: memalign (vg_replace_malloc.c:727)
==31551==    by 0x5498E79: PetscMallocAlign (mal.c:27)
==31551==    by 0x5DE9D5E: KSPSetUp_GMRES (gmres.c:75)


==31551==    by 0x5E1C41F: KSPSetUp (itfunc.c:239)
==31551==    by 0x5E1821D: KSPSolve (itfunc.c:399)
==31551==    by 0x5CAD7F1: kspsolve_ (itfuncf.c:219)
==31551==    by 0x8D9842: petsc_solver_defi_ (petsc_defi.F90:434)


==31551==    by 0x831618: field_solver_defi_ (field_defi.F90:57)
==31551==    by 0x46CB4E: MAIN__ (main.F90:96)
==31551==    by 0x41852B: main (in /home/jlin/defi_field/ftest/gtc)
 
...
 
==31551== 4,092 bytes in 1 blocks are definitely lost in loss record 3,113 of 3,327


==31551==    at 0x4A05458: memalign (vg_replace_malloc.c:727)
==31551==    by 0x5498E79: PetscMallocAlign (mal.c:27)
==31551==    by 0x5A085AE: MatDuplicateNoCreate_SeqAIJ (aij.c:4011)
==31551==    by 0x5A2D24F: MatILUFactorSymbolic_SeqAIJ_ilu0 (aijfact.c:1655)


==31551==    by 0x5A2B7CB: MatILUFactorSymbolic_SeqAIJ (aijfact.c:1756)
==31551==    by 0x573D152: MatILUFactorSymbolic (matrix.c:6240)
==31551==    by 0x5CFB843: PCSetUp_ILU (ilu.c:204)
==31551==    by 0x5D8BAA7: PCSetUp (precon.c:890)


==31551==    by 0x5E1C639: KSPSetUp (itfunc.c:278)
==31551==    by 0x5E1821D: KSPSolve (itfunc.c:399)
==31551==    by 0x5CAD7F1: kspsolve_ (itfuncf.c:219)
==31551==    by 0x8D9842: petsc_solver_defi_ (petsc_defi.F90:434)


==31551==    by 0x831618: field_solver_defi_ (field_defi.F90:57)
==31551==    by 0x46CB4E: MAIN__ (main.F90:96)
==31551==    by 0x41852B: main (in /home/jlin/defi_field/ftest/gtc)

-------------------------------------------------------------------------------------------------



 
From: francium87 at hotmail.com
To: knepley at gmail.com
CC: petsc-users at mcs.anl.gov


Subject: RE: [petsc-users] VecValidValues() reports NaN found
Date: Mon, 5 May 2014 13:15:38 +0000




Ok, I will try it . Thanks for your advise. 

Date: Mon, 5 May 2014 08:12:05 -0500
Subject: Re: [petsc-users] VecValidValues() reports NaN found
From: knepley at gmail.com


To: francium87 at hotmail.com
CC: petsc-users at mcs.anl.gov

On Mon, May 5, 2014 at 7:56 AM, linjing bo <francium87 at hotmail.com> wrote:






I use JACOBI. The message showed is with JACOBI. 


Wired situation is that the backtrack information shows the location
 is before actually apply PC, so I guess the rhs vec is not changed at 
this point.

Another wired thing is : Because the original code is
 to complex. I write out the A matrix in Ax=b, and write a small test 
code to read in this matrix and solve it, no error showed. The KSP, PC 
are all set to be the same.

When I try to using ILU, more wired error happens, the backtrack info shows it died in a Flops logging function:

1) Run in serial until it works



2) It looks like you have memory overwriting problems. Run with valgrind
   Matt 


[2]PETSC ERROR: --------------------- Error Message ------------------------------------
[2]PETSC ERROR: Argument out of range!                                                  
[2]PETSC ERROR: Cannot log negative flops!                                              



[2]PETSC ERROR: ------------------------------------------------------------------------
[2]PETSC ERROR: Petsc Release Version 3.4.4, Mar, 13, 2014                              
[2]PETSC ERROR: See docs/changes/index.html for recent updates.                         



[2]PETSC ERROR: See docs/faq.html for hints about trouble shooting.                     
[2]PETSC ERROR: See docs/index.html for manual pages.                                   
[2]PETSC ERROR: ------------------------------------------------------------------------



[2]PETSC ERROR: ./gtc on a arch-linux2-c-debug named node2.indac.info by jlin Mon May  5 20:51:27 2014                                                                                                



[2]PETSC ERROR: Libraries linked from /opt/sfw/petsc/3.4.4/intel/openmpi/lib                       
[2]PETSC ERROR: Configure run at Sat Apr 26 20:19:41 2014                                          
[2]PETSC ERROR: Configure options --prefix=/opt/sfw/petsc/3.4.4/intel/openmpi --with-mpi-dir=/opt/sfw/openmpi/1.6.3/intel --with-blas-lapack-dir=/opt/sfw/intel/composer_xe_2011_sp1.7.256/mkl/lib/intel64 --with-mpiexec=mpiexec                                                                        



[2]PETSC ERROR: ------------------------------------------------------------------------           
[2]PETSC ERROR: PetscLogFlops() line 204 in /tmp/petsc-3.4.4/include/petsclog.h                    
[2]PETSC ERROR: MatLUFactorNumeric_SeqAIJ() line 552 in /tmp/petsc-3.4.4/src/mat/impls/aij/seq/aijfact.c                                                                                              



[2]PETSC ERROR: MatLUFactorNumeric() line 2889 in /tmp/petsc-3.4.4/src/mat/interface/matrix.c      
[2]PETSC ERROR: PCSetUp_ILU() line 232 in /tmp/petsc-3.4.4/src/ksp/pc/impls/factor/ilu/ilu.c       
[2]PETSC ERROR: PCSetUp() line 890 in /tmp/petsc-3.4.4/src/ksp/pc/interface/precon.c               



[2]PETSC ERROR: KSPSetUp() line 278 in /tmp/petsc-3.4.4/src/ksp/ksp/interface/itfunc.c             
[2]PETSC ERROR: KSPSolve() line 399 in /tmp/petsc-3.4.4/src/ksp/ksp/interface/itfunc.c             







Date: Mon, 5 May 2014 07:27:52 -0500
Subject: Re: [petsc-users] VecValidValues() reports NaN found
From: knepley at gmail.com
To: francium87 at hotmail.com



CC: petsc-users at mcs.anl.gov

On Mon, May 5, 2014 at 7:25 AM, linjing bo <francium87 at hotmail.com> wrote:










Hi, I'm trying to use PETSc's ksp method to solve a linear system. When running, Error is reported by VecValidValues() that NaN or Inf is found with error message listed below 


[3]PETSC ERROR: --------------------- Error Message ------------------------------------           




[3]PETSC ERROR: Floating point exception!                                                          
[3]PETSC ERROR: Vec entry at local location 0 is not-a-number or infinite at beginning of function: Parameter number 2!                                                                               




[3]PETSC ERROR: ------------------------------------------------------------------------           
[3]PETSC ERROR: Petsc Release Version 3.4.4, Mar, 13, 2014                                         
[3]PETSC ERROR: See docs/changes/index.html for recent updates.                                    




[3]PETSC ERROR: See docs/faq.html for hints about trouble shooting.                                
[3]PETSC ERROR: See docs/index.html for manual pages.                                              
[3]PETSC ERROR: ------------------------------------------------------------------------           




[3]PETSC ERROR: ./gtc on a arch-linux2-c-debug named node2.indac.info by jlin Mon May  5 20:03:20 2014                                                                                                




[3]PETSC ERROR: Libraries linked from /opt/sfw/petsc/3.4.4/intel/openmpi/lib                       
[3]PETSC ERROR: Configure run at Sat Apr 26 20:19:41 2014                                          
[3]PETSC ERROR: Configure options --prefix=/opt/sfw/petsc/3.4.4/intel/openmpi --with-mpi-dir=/opt/sfw/openmpi/1.6.3/intel --with-blas-lapack-dir=/opt/sfw/intel/composer_xe_2011_sp1.7.256/mkl/lib/intel64 --with-mpiexec=mpiexec                                                                        




[3]PETSC ERROR: ------------------------------------------------------------------------           
[3]PETSC ERROR: VecValidValues() line 28 in /tmp/petsc-3.4.4/src/vec/vec/interface/rvector.c  




It looks like the vector after preconditioner application is bad. What is the preconditioner?
  Matt 



     
[3]PETSC ERROR: PCApply() line 436 in /tmp/petsc-3.4.4/src/ksp/pc/interface/precon.c               
[3]PETSC ERROR: KSP_PCApply() line 227 in /tmp/petsc-3.4.4/include/petsc-private/kspimpl.h         




[3]PETSC ERROR: KSPInitialResidual() line 64 in /tmp/petsc-3.4.4/src/ksp/ksp/interface/itres.c     
[3]PETSC ERROR: KSPSolve_GMRES() line 239 in /tmp/petsc-3.4.4/src/ksp/ksp/impls/gmres/gmres.c      
[3]PETSC ERROR: KSPSolve() line 441 in /tmp/petsc-3.4.4/src/ksp/ksp/interface/itfunc.c             






After read the source code shown by backtrack informations, I realize the problem is in the right hand side vector. So I make a trial of set right hand side vector to ONE by VecSet, But the program still shows error message above, and using VecView or VecGetValue to investigate the first value of rhs vec shows the value is 1.0 as I set it to. Hope I clearly describe the problem. The code related is listed below





---------------------------Solver section--------------------------

  call VecSet( pet_bp_b, one, ierr)

  vecidx=[0,1]
  call VecGetValues( pet_bp_b, 2, vecidx, first, ierr)
  write(*,*) ' first two values ', first(1), first(2)





  call KSPSetInitialGuessNonzero(solver_bp,Petsc_True,ierr)
  call KSPSolve(solver_bp,pet_bp_b,pet_bp_x,ierr)
  call KSPView(solver_bp, PETSC_VIEWER_STDOUT_SELF,ierr)
  CHKERRQ(ierr)



 		 	   		  


-- 
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
 		 	   		   		 	   		  


-- 
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
 		 	   		  
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140508/bd01922a/attachment-0001.html>
-------------- next part --------------
A non-text attachment was scrubbed...
Name: matrix.zip
Type: application/zip
Size: 50130 bytes
Desc: not available
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140508/bd01922a/attachment-0001.zip>


More information about the petsc-users mailing list