[petsc-users] eps_target malloc issue

jifeng zhao jifengzhao2015 at u.northwestern.edu
Mon Jul 28 16:48:12 CDT 2014


Hi Barry,

Thanks for the reply. I see. I am still testing it and looking for the bug.


It is pretty weird that the error appears, since the code I am running is
ex7 and I didn't change anything.

My matrices are assembled in binary files separately. As I explained
earlier, I have a series of matrices to solve, each of them have the same
patterns but with different values. The way I assembly the matrix has
proven to be correct for all the other cases. Only for this one, when the
matrix has an zero eigenvalue, ex7 failed to solve it.

In fact, I don't understand why this "Out of range" error could possibly
appear? 1. my matrices have the correct parrellel lay up. 2. In ex7, I am
not explicitly accessing to any elements at all! All I did is EPSSolve();

Any ideas?

Best regards,
Jifeng Zhao




On Mon, Jul 28, 2014 at 12:07 PM, Barry Smith <bsmith at mcs.anl.gov> wrote:

>
> On Jul 28, 2014, at 10:27 AM, jifeng zhao <
> jifengzhao2015 at u.northwestern.edu> wrote:
>
> > Hello everyone,
> >
> > Are there more ideas why when using MatSetOption(mat,
> MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE), the processor seems trying to
> access elements that are out of local range?
>
>    Add the option has NOTHING to do with the “processor seems trying to
> access elements out of local range”. It is just that adding that option
> allows the code to run further until it hits this other bug.  So they
> second bug needs to be tracked down.
>
>   Barry
>
> >
> > Anybody can help? Thank you so much.
> >
> > Best regards,
> > Jifeng Zhao
> >
> >
> > On Sun, Jul 27, 2014 at 3:37 PM, jifeng zhao <
> jifengzhao2015 at u.northwestern.edu> wrote:
> > Thank you all for the insights.
> >
> > I tried to add MatSetOption(mat,
> MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE). It does solve the malloc
> error. However there appears to be another error:
> >
> > In addition, 66324 appears to be the number of index that is assigned to
> each processor.
> >
> > Any more insights? Thanks a lot!
> >
> > PETSC ERROR: Argument out of range!
> > PETSC ERROR: Column too large: col 935748 max 66324!
> > PETSC ERROR:
> ------------------------------------------------------------------------
> > PETSC ERROR: Petsc Release Version 3.4.4, Mar, 13, 2014
> > PETSC ERROR: See docs/changes/index.html for recent updates.
> > PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> > PETSC ERROR: See docs/index.html for manual pages.
> > PETSC ERROR:
> ------------------------------------------------------------------------
> > PETSC ERROR: eigen_solver20 on a arch-linux2-c-debug named qnode4091 by
> jzh953 Sun Jul 27 15:17:28 2014
> > PETSC ERROR: Libraries linked from /software/petsc/3.4.4/lib
> > PETSC ERROR: Configure run at Tue Jul  8 14:20:27 2014
> > ]PETSC ERROR: Configure options --prefix=/software/petsc/3.4.4
> --with-blas-lib=/software/lapack/3.5.0/lib/libblas.so
> --with-lapack-lib=/software/lapack/3.5.0/liblapack.so
> --with-hdf5-dir=/software/hdf5/1.8.12 with-boost=1
> --with-boost-dir=/software/hdf5/1.8.12 --with-fftw=1
> --with-fftw-dir=/software/FFTW/3.3.3-gcc-threads
> --with-valgrind-dir=/software/valgrind/3.8.1 --with-netcdf=1
> --with-netcdf-dir=/software/netCDF/4.2 --with-hdf5=1 --with-debugging=1
> --with-mpi-dir=/software/mpi/openmpi-1.6.3-gcc-4.6.3-trq4
> --download-scalapack --download-mumps --download-hypre --download-parmetis
> --download-metis
> >
> > PETSC ERROR:
> ------------------------------------------------------------------------
> > PETSC ERROR: MatSetValues_SeqSBAIJ() line 986 in
> /hpc/software/sources/builds/petsc-3.4.4/src/mat/impls/sbaij/seq/sbaij.c
> > PETSC ERROR: MatDisAssemble_MPISBAIJ() line 356 in
> /hpc/software/sources/builds/petsc-3.4.4/src/mat/impls/sbaij/mpi/mmsbaij.c
> > PETSC ERROR: MatSetValues_MPISBAIJ() line 213 in
> /hpc/software/sources/builds/petsc-3.4.4/src/mat/impls/sbaij/mpi/mpisbaij.c
> > PETSC ERROR: MatSetValues() line 1106 in
> /hpc/software/sources/builds/petsc-3.4.4/src/mat/interface/matrix.c
> > PETSC ERROR: MatAXPY_Basic() line 78 in
> /hpc/software/sources/builds/petsc-3.4.4/src/mat/utils/axpy.c
> > PETSC ERROR: MatAXPY_MPISBAIJ() line 1398 in
> /hpc/software/sources/builds/petsc-3.4.4/src/mat/impls/sbaij/mpi/mpisbaij.c
> > PETSC ERROR: MatAXPY() line 39 in
> /hpc/software/sources/builds/petsc-3.4.4/src/mat/utils/axpy.c
> > PETSC ERROR: STMatGAXPY_Private() line 375 in src/st/interface/stsolve.c
> > PETSC ERROR: STSetUp_Sinvert() line 138 in src/st/impls/sinvert/sinvert.c
> > PETSC ERROR: STSetUp() line 294 in src/st/interface/stsolve.c
> > PETSC ERROR: EPSSetUp() line 215 in src/eps/interface/setup.c
> > PETSC ERROR: EPSSolve() line 90 in src/eps/interface/solve.c
> > PETSC ERROR:
> --------------------------------------------------------------------------
> >
> >
> > On Sat, Jul 26, 2014 at 12:35 PM, Jose E. Roman <jroman at dsic.upv.es>
> wrote:
> > The problem is that your stiffness and mass matrices have different
> nonzero pattern.
> >
> > The -eps_target configuration computes a MatAXPY on the first matrix
> (stiffness), so MatSetOption(mat,
> MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE) could be called on this matrix,
> as Barry suggests.
> >
> > Another workaround would be to force both of them having the same
> nonzero pattern, maybe using
> MatSetOption(mat,MAT_IGNORE_ZERO_ENTRIES,PETSC_FALSE).
> >
> > Also, this problem would not appear with target=0.0.
> >
> > See details on section 3.4.2 of SLEPc's users guide.
> >
> > Jose
> >
> >
> >
> > El 26/07/2014, a las 19:04, Barry Smith escribió:
> >
> > >
> > >  We need the entire error message. A matrix is being formed that
> wasn’t fully preallocated for so you need to determine which matrix and
> make sure either 1) it gets fully preallocated or 2) MatSetOption(mat,
> MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE) is called on the matrix
> > >
> > >  Barry
> > >
> > > On Jul 26, 2014, at 10:34 AM, jifeng zhao <
> jifengzhao2015 at u.northwestern.edu> wrote:
> > >
> > >> Hello all,
> > >>
> > >> I am trying to use ex7 to solve series of generalized problem. The
> running commands are:
> > >>
> > >> ./ex7 -f1 petsc_stiff19.dat -f2 petsc_mass19.dat -mat_type sbaij
> -eps_gen_hermitian -eps_type krylovschur -eps_smallest_magnitude.
> -eps_monitor_conv -st_k    sp_type bcgs -st_pc_type bjacobi -st_sub_pc_type
> icc -st_ksp_rtol 1.e-4 -eps_tol 1.e-4 -eps_nev 40 -st_type sinvert
> > >>
> > >> It works fine, but when my problem contains zero eigenvalue, it fails
> to converge.
> > >>
> > >> Supposedly, I should use eps_target -0.1 instead of
> eps_smallest_magnitude. However if eps_target is used, I run into a memory
> issue:
> > >>
> > >> PETSC ERROR: --------------------- Error Message
> ------------------------------------
> > >> PETSC ERROR: Argument out of range!
> > >> PETSC ERROR: New nonzero at (0,6) caused a malloc!
> > >> PETSC ERROR:
> ------------------------------------------------------------------------
> > >>
> > >> Does anybody know how to solve this? Thanks!
> > >>
> > >> Jifeng Zhao
> > >> PhD candidate at Northwestern University, US
> > >> Theoretical and Applied Mechanics Program
> > >
> >
> >
> >
> >
> > --
> > Jifeng Zhao
> > PhD candidate at Northwestern University, US
> > Theoretical and Applied Mechanics Program
> >
> >
> >
> > --
> > Jifeng Zhao
> > PhD candidate at Northwestern University, US
> > Theoretical and Applied Mechanics Program
>
>


-- 
Jifeng Zhao
PhD candidate at Northwestern University, US
Theoretical and Applied Mechanics Program
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20140728/e3a96c94/attachment.html>


More information about the petsc-users mailing list