[petsc-users] Segmentation violation issue when using MatSetVals in Fortran
Matthew Knepley
knepley at gmail.com
Thu Feb 2 13:32:15 CST 2017
On Thu, Feb 2, 2017 at 1:25 PM, Austin Herrema <aherrema at iastate.edu> wrote:
> Hello all,
>
> I'm having a strange issue when using the MatSetVals function in Fortran.
> The code is quite large, but relevant snippets of what I'm trying to do are
> below.
>
> ! PETSc variable declaration
> Mat LHS
> PetscErrorCode ier
> PetscInt idx1(nsd*nsd*nshl*nshl), idx2(nsd*nsd*nshl*nshl)
> PetscScalar vals(nsd*nsd*nshl*nshl)
>
vals should be of size nidx1 * nidx2.
Matt
>
> ! Code for filling in index and value arrays... I've verified that
> there are 2304 elements in each array.
>
> ! Fill in values
> call MatSetValues(LHS, nvals, idx1, nvals, idx2, vals, ADD_VALUES, ier)
>
> Upon executing this, I get a general segmentation violation error. Running
> in Valgrind I get:
>
> ==87676== Invalid read of size 8
> ==87676== at 0x1008D1F45: MatSetValues_SeqAIJ (aij.c:461)
> ==87676== by 0x1004C523C: MatSetValues (matrix.c:1190)
> ==87676== by 0x10053F1A5: matsetvalues_ (matrixf.c:696)
> ==87676== by 0x10005F47D: petsc_fillmat_ (petsc_fillMat.F90:78)
> ==87676== by 0x10005E5C9: intelmass_shell_ (shell_IntElmAss.F90:330)
> ==87676== by 0x100057259: solflowgmres_ (solflow.F90:125)
> ==87676== by 0x100051051: MAIN__ (driver.F90:140)
> ==87676== by 0x1000545C4: main (driver.F90:7)
> ==87676== Address 0x10a4dac60 is 0 bytes after a block of size 18,432
> alloc'd
> ==87676== at 0x100088681: malloc (in /usr/local/Cellar/valgrind/3.
> 12.0/lib/valgrind/vgpreload_memcheck-amd64-darwin.so)
> ==87676== by 0x10005EF47: petsc_fillmat_ (petsc_fillMat.F90:23)
> ==87676== by 0x10005E5C9: intelmass_shell_ (shell_IntElmAss.F90:330)
> ==87676== by 0x100057259: solflowgmres_ (solflow.F90:125)
> ==87676== by 0x100051051: MAIN__ (driver.F90:140)
> ==87676== by 0x1000545C4: main (driver.F90:7)
> ==87676==
> ==87676== Conditional jump or move depends on uninitialised value(s)
> ==87676== at 0x1008D1F57: MatSetValues_SeqAIJ (aij.c:463)
> ==87676== by 0x1004C523C: MatSetValues (matrix.c:1190)
> ==87676== by 0x10053F1A5: matsetvalues_ (matrixf.c:696)
> ==87676== by 0x10005F47D: petsc_fillmat_ (petsc_fillMat.F90:78)
> ==87676== by 0x10005E5C9: intelmass_shell_ (shell_IntElmAss.F90:330)
> ==87676== by 0x100057259: solflowgmres_ (solflow.F90:125)
> ==87676== by 0x100051051: MAIN__ (driver.F90:140)
> ==87676== by 0x1000545C4: main (driver.F90:7)
> ==87676==
> ==87676== Conditional jump or move depends on uninitialised value(s)
> ==87676== at 0x1008D1F62: MatSetValues_SeqAIJ (aij.c:463)
> ==87676== by 0x1004C523C: MatSetValues (matrix.c:1190)
> ==87676== by 0x10053F1A5: matsetvalues_ (matrixf.c:696)
> ==87676== by 0x10005F47D: petsc_fillmat_ (petsc_fillMat.F90:78)
> ==87676== by 0x10005E5C9: intelmass_shell_ (shell_IntElmAss.F90:330)
> ==87676== by 0x100057259: solflowgmres_ (solflow.F90:125)
> ==87676== by 0x100051051: MAIN__ (driver.F90:140)
> ==87676== by 0x1000545C4: main (driver.F90:7)
> ==87676==
>
> I'm not terribly experienced with Valgrind but I can't find anything too
> helpful in that output. It shows the traceback but I still cannot
> understand *why *the error is being produced. Strangely, if I try to set
> the matrix values in a few smaller chunks, i.e.
>
> call MatSetValues(LHS, 500, idx1(1:500), 500, idx2(1:500), vals(1:500),
> ADD_VALUES, ier)
> call MatSetValues(LHS, 500, idx1(501:1000), 500, idx2(501:1000),
> vals(501:1000), ADD_VALUES, ier)
> call MatSetValues(LHS, 500, idx1(1001:1500), 500, idx2(1001:1500),
> vals(1001:1500), ADD_VALUES, ier)
> call MatSetValues(LHS, 500, idx1(1501:2000), 500, idx2(1501:2000),
> vals(1501:2000), ADD_VALUES, ier)
> call MatSetValues(LHS, 304, idx1(2001:2304), 304, idx2(2001:2304),
> vals(2001:2304), ADD_VALUES, ier)
>
> Then everything works with no errors. Is there some kind of limitation to
> how many values you can insert into a matrix at once? Would strike me as
> strange. I can assign elements one at a time but reducing the number of
> calls to MatSetVals by a factor of 2000 would be fantastic. I would be
> happy to provide more details of the matrix setup if necessary.
>
> Thank you for any help!
> Austin
>
>
> --
> *Austin Herrema*
> PhD Student | Graduate Research Assistant | Iowa State University
> Wind Energy Science, Engineering, and Policy | Mechanical Engineering
>
--
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/20170202/118c4d94/attachment-0001.html>
More information about the petsc-users
mailing list