PetscMalloc petsc-3.0.0-p4
Barry Smith
bsmith at mcs.anl.gov
Mon Mar 23 19:47:50 CDT 2009
[anlextwls097-161:~/Documents/RandD-100] barrysmith% python
>>> npixel = 256*256*5
>>> print npixel
327,680
>>> print npixel*npixel
107,374,182,400
The problem is that your n is too big to fit into a PetscInt hence bad
stuff happens.
This is a fatal flaw of C, it doesn't generation exceptions for
integer overflow.
If you want to have sizes this big you need to config/configure.py
PETSc with
--with-64-bit-indices then PETSc uses 64 bit integers for PETSc and
things will fit.
Of course, likely you still won't have enough memory to allocate such
a big dense matrix.
Barry
On Mar 23, 2009, at 6:41 PM, David Fuentes wrote:
>
>
> running on ubuntu 8.04.
>
> the following code should reproduce the problem.
>
>
>
> static char help[] = "MatCreateSeqDense\n\n";
>
> #include "petscpc.h"
>
> #undef __FUNCT__
> #define __FUNCT__ "main"
> int main(int argc,char **args)
> {
> Mat mat;
> PetscErrorCode ierr;
> PetscInt npixel = 256*256*5;
> PetscInt n= npixel * npixel;
>
> PetscInitialize(&argc,&args,(char *)0,help);
>
>
> /* Create and assemble matrix */
> ierr =
> MatCreateSeqDense(PETSC_COMM_SELF,n,n,PETSC_NULL,&mat);CHKERRQ(ierr);
> ierr = MatSetValue(mat,1,1,1.0,INSERT_VALUES);CHKERRQ(ierr);
> ierr = MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
> ierr = MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
>
> /* Free data structures */
> ierr = MatDestroy(mat);CHKERRQ(ierr);
> ierr = PetscFinalize();CHKERRQ(ierr);
> return 0;
> }
>
>
>
> below is my error message.
>
>
>
>
>
> [0]PETSC ERROR: --------------------- Error Message
> ------------------------------------
> [0]PETSC ERROR: Argument out of range!
> [0]PETSC ERROR: Row too large: row 1 max -1!
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
> [0]PETSC ERROR: Petsc Release Version 3.0.0, Patch 4, Fri Mar 6
> 14:46:08 CST 2009
> [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: ./ex4 on a intel-10. named setebos.ices.utexas.edu by
> fuentes Mon Mar 23 18:39:32 2009
> [0]PETSC ERROR: Libraries linked from
> /org/groups/oden/LIBRARIES/PETSC/petsc-3.0.0-p4/intel-10.1-
> mpich2-1.0.7-cxx-dbg/lib
> [0]PETSC ERROR: Configure run at Fri Mar 20 20:39:04 2009
> [0]PETSC ERROR: Configure options
> --with-mpi-dir=/org/groups/oden/LIBRARIES/MPI/mpich2-1.0.7-intel-10.1
> --with-clanguage=C++ --with-shared=1
> --with-blas-lapack-dir=/opt/intel/mkl/10.0.3.020 --download-mumps=1
> --download-parmetis=1 --download-scalapack=1 --download-blacs=1
> --download-plapack=1 --download-superlu_dist=1
> [0]PETSC ERROR:
> ------------------------------------------------------------------------
> [0]PETSC ERROR: MatSetValues_SeqDense() line 672 in
> src/mat/impls/dense/seq/dense.c
> [0]PETSC ERROR: MatSetValues() line 921 in src/mat/interface/matrix.c
> [0]PETSC ERROR: main() line 21 in src/ksp/pc/examples/tests/ex4.c
> application called MPI_Abort(MPI_COMM_WORLD, 63) - process 0[unset]:
> aborting job:
> application called MPI_Abort(MPI_COMM_WORLD, 63) - process 0
>
>
>
>
>
>
> df
>
>
>
>
> On Mon, 23 Mar 2009, Lisandro Dalcin wrote:
>
>> 1) Which system are you running on?
>>
>> 2) are you completely sure BOTH things happens, I mean... that malloc
>> returns null AND the error code is zero?
>>
>>
>> On Mon, Mar 23, 2009 at 7:49 PM, David Fuentes
>> <fuentesdt at gmail.com> wrote:
>>> In 3.0.0-p4
>>>
>>> PetscMalloc doesn't seem to return an error code if the requested
>>> memory
>>> wasn't allocated ? Is this correct ?
>>>
>>> I'm trying to allocate a dense matrix for which there was not
>>> enough memory
>>> and PetscMalloc is returning a null pointer but no errorcode to
>>> catch. which is causing seg faults later.
>>>
>>>
>>>
>>> thanks,
>>> df
>>>
>>>
>>
>>
>>
>> --
>> 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