[petsc-users] error related to nested vector

Matthew Knepley knepley at gmail.com
Tue Jan 14 06:08:41 CST 2020


On Tue, Jan 14, 2020 at 6:58 AM Y. Shidi <ys453 at cam.ac.uk> wrote:

> Thanks Matt,
>
> Instead of constructing the nested vector by VecCreateNest(), I use
> VecRestoreSubVector() to solve this issue.
>
> But I got problems for field splitting method,
> I use the following options:
> #PETSc Option Table entries:
> -fieldsplit_0_ksp_converged_reason true
> -fieldsplit_0_ksp_rtol 1e-12
> -fieldsplit_0_ksp_type cg
> -fieldsplit_0_pc_type ilu
> -fieldsplit_1_ksp_converged_reason true
> -fieldsplit_1_ksp_rtol 1e-12
> -fieldsplit_1_ksp_type cg
> -fieldsplit_1_pc_type hypre
> -ksp_converged_reason
> -ksp_rtol 1e-12
> -ksp_type fgmres
> -ksp_view
> -pc_fieldsplit_schur_fact_type upper
> -pc_fieldsplit_schur_precondition seflp
> -pc_type fieldsplit
> #End of PETSc Option Table entries
>
> There are some options that were not used.
>
> WARNING! There are options you set that were not used!
> WARNING! could be spelling mistake, etc!
> Option left: name:-fieldsplit_0_ksp_converged_reason value: true
> Option left: name:-fieldsplit_1_ksp_converged_reason value: true
> Option left: name:-pc_fieldsplit_schur_fact_type value: upper
> Option left: name:-pc_fieldsplit_schur_precondition value: seflp
>

PCFIELDSPLIT has to know how to divide your vector into fields. VecNest
provides that information
since it is divided into the nested pieces. Other ways you can provide this:

  1) Use a DM to describe the data layout. This is the best in my opinion.

  2) Fix the VecNest. If an overwrite is happening, its best to track it
down.

  3) Provide an IS to PCFIELDSPLIT using PCFieldSplitSetIS().

  4) Split by block component

  Thanks,

    Matt


> Kind Regards,
> Shidi
>
> On 2020-01-14 11:53, Matthew Knepley wrote:
> > It says that memory is being overwritten somewhere. You can track this
> > down using valgrind,
> > as it suggests in the error message.
> >
> >   Thanks,
> >
> >      Matt
> >
> > On Tue, Jan 14, 2020 at 5:44 AM Y. Shidi <ys453 at cam.ac.uk> wrote:
> >
> >> Dear developers,
> >>
> >> I have a 2x2 nested matrix and the corresponding nested vector.
> >> When I running the code with field splitting, it gets the following
> >> errors:
> >>
> >> [0]PETSC ERROR: PetscTrFreeDefault() called from
> >> VecRestoreArray_Nest()
> >> line 678 in
> >> /home/ys453/Sources/petsc/src/vec/vec/impls/nest/vecnest.c
> >> [0]PETSC ERROR: Block at address 0x3f95f60 is corrupted; cannot
> >> free;
> >> may be block not allocated with PetscMalloc()
> >> [0]PETSC ERROR: --------------------- Error Message
> >> --------------------------------------------------------------
> >> [0]PETSC ERROR: Memory corruption:
> >>
> > http://www.mcs.anl.gov/petsc/documentation/installation.html#valgrind
> >> [0]PETSC ERROR: Bad location or corrupted memory
> >> [0]PETSC ERROR: See
> >> http://www.mcs.anl.gov/petsc/documentation/faq.html
> >> for trouble shooting.
> >> [0]PETSC ERROR: Petsc Release Version 3.9.3, unknown
> >> [0]PETSC ERROR: 2DPetscSpuriousTest on a arch-linux2-c-debug named
> >> merlin by ys453 Tue Jan 14 10:36:53 2020
> >> [0]PETSC ERROR: Configure options --download-scalapack
> >> --download-mumps
> >> --download-parmetis --download-metis --download-ptscotch
> >> --download-superlu_dist --download-hypre
> >> [0]PETSC ERROR: #1 PetscTrFreeDefault() line 269 in
> >> /home/ys453/Sources/petsc/src/sys/memory/mtr.c
> >> [0]PETSC ERROR: #2 VecRestoreArray_Nest() line 678 in
> >> /home/ys453/Sources/petsc/src/vec/vec/impls/nest/vecnest.c
> >> [0]PETSC ERROR: #3 VecRestoreArrayRead() line 1835 in
> >> /home/ys453/Sources/petsc/src/vec/vec/interface/rvector.c
> >> [0]PETSC ERROR: #4 VecRestoreArrayPair() line 511 in
> >> /home/ys453/Sources/petsc/include/petscvec.h
> >> [0]PETSC ERROR: #5 VecScatterBegin_SSToSS() line 671 in
> >> /home/ys453/Sources/petsc/src/vec/vscat/impls/vscat.c
> >> [0]PETSC ERROR: #6 VecScatterBegin() line 1779 in
> >> /home/ys453/Sources/petsc/src/vec/vscat/impls/vscat.c
> >> [0]PETSC ERROR: #7 PCApply_FieldSplit() line 1010 in
> >> /home/ys453/Sources/petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c
> >> [0]PETSC ERROR: #8 PCApply() line 457 in
> >> /home/ys453/Sources/petsc/src/ksp/pc/interface/precon.c
> >> [0]PETSC ERROR: #9 KSP_PCApply() line 276 in
> >> /home/ys453/Sources/petsc/include/petsc/private/kspimpl.h
> >> [0]PETSC ERROR: #10 KSPFGMRESCycle() line 166 in
> >> /home/ys453/Sources/petsc/src/ksp/ksp/impls/gmres/fgmres/fgmres.c
> >> [0]PETSC ERROR: #11 KSPSolve_FGMRES() line 291 in
> >> /home/ys453/Sources/petsc/src/ksp/ksp/impls/gmres/fgmres/fgmres.c
> >> [0]PETSC ERROR: #12 KSPSolve() line 669 in
> >> /home/ys453/Sources/petsc/src/ksp/ksp/interface/itfunc.c
> >>
> >> I am not sure why it happens.
> >>
> >> Thank you for your time.
> >>
> >> Kind Regards,
> >> Shidi
> >
> > --
> >
> > 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
> >
> > https://www.cse.buffalo.edu/~knepley/ [1]
> >
> >
> > Links:
> > ------
> > [1] http://www.cse.buffalo.edu/~knepley/
>


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

https://www.cse.buffalo.edu/~knepley/ <http://www.cse.buffalo.edu/~knepley/>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20200114/f1cd08ca/attachment-0001.html>


More information about the petsc-users mailing list