[petsc-users] PCFIELDSPLIT question

Matthew Knepley knepley at gmail.com
Tue Feb 2 14:42:16 CST 2016


On Tue, Feb 2, 2016 at 2:35 PM, Hom Nath Gharti <hng.email at gmail.com> wrote:

> Thank you so much Matt.
>
> I now tried the following:
>
> ======================================================
> call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)
> call KSPGetPC(ksp,pc,ierr)
>
> call PCFieldSplitSetIS(pc,"0",gdofu_is,ierr);
> call ISDestroy(gdofu_is,ierr)
> call PCFieldSplitSetIS(pc,"1",gdofchi_is,ierr);
> call ISDestroy(gdofchi_is,ierr)
> call PCFieldSplitSetIS(pc,"2",gdofp_is,ierr);
> call ISDestroy(gdofp_is,ierr)
> call PCFieldSplitSetIS(pc,"3",gdofphi_is,ierr);
> call ISDestroy(gdofphi_is,ierr)
>
> ! Amat changes in each time steps
> call KSPSetOperators(ksp,Amat,Amat,ierr) !version >= 3.5.0
>
> call KSPSolve(ksp,bvec,xvec,ierr)
> ======================================================
>

I am guessing that "0" is not a valid string for your Fortran compiler. Are
you sure
its not single quotes '0'?

  Matt


> But I get the following error:
>
> [0]PETSC ERROR: --------------------- Error Message
> --------------------------------------------------------------
> [0]PETSC ERROR: Petsc has generated inconsistent data
> [0]PETSC ERROR: Unhandled case, must have at least two fields, not 1
> [0]PETSC ERROR: See
> http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble
> shooting.
> [0]PETSC ERROR: Petsc Release Version 3.6.3, Dec, 03, 2015
> [0]PETSC ERROR:
> /tigress/hgharti/gitwork/SPECFEM3D_GLOBEVSPP/./bin/xspecfem3D on a
> arch-linux2-c-debug named tiger-r3c1n7 by hgharti Tue Feb  2 15:
> 29:30 2016
> [0]PETSC ERROR: Configure options
> --with-blas-lapack-dir=/opt/intel/composer_xe_2015.2.164/mkl/lib/intel64/
> --with-cc=mpicc --with-cxx=mpicxx --wit
> h-fc=mpif90 --with-mpiexec=mpiexec --with-debugging=1
> --download-scalapack --download-mumps --download-pastix
> --download-superlu --download-superlu_dist --download-metis
> --download-parmetis --download-ptscotch --download-hypre
> [0]PETSC ERROR: #1 PCFieldSplitSetDefaults() line 469 in
>
> /home/hgharti/lsoft/petsc-3.6.3-intel16-mpi5/src/ksp/pc/impls/fieldsplit/fieldsplit.c
> [0]PETSC ERROR: #2 PCSetUp_FieldSplit() line 486 in
>
> /home/hgharti/lsoft/petsc-3.6.3-intel16-mpi5/src/ksp/pc/impls/fieldsplit/fieldsplit.c
> [0]PETSC ERROR: #3 PCSetUp() line 983 in
> /home/hgharti/lsoft/petsc-3.6.3-intel16-mpi5/src/ksp/pc/interface/precon.c
> [0]PETSC ERROR: #4 KSPSetUp() line 332 in
> /home/hgharti/lsoft/petsc-3.6.3-intel16-mpi5/src/ksp/ksp/interface/itfunc.c
> [0]PETSC ERROR: #5 KSPSolve() line 546 in
> /home/hgharti/lsoft/petsc-3.6.3-intel16-mpi5/src/ksp/ksp/interface/itfunc.c
> forrtl: error (76): Abort trap signal
>
> Am I missing something?
>
> Thanks,
> Hom Nath
>
> On Tue, Feb 2, 2016 at 3:24 PM, Matthew Knepley <knepley at gmail.com> wrote:
> > On Tue, Feb 2, 2016 at 2:20 PM, Hom Nath Gharti <hng.email at gmail.com>
> wrote:
> >>
> >> Hi Matt, hi all,
> >>
> >> I am trying to use PCFIELDSPLIT for my solver which consists of 4
> >> different variables, namely, u (displacement vector), \chi
> >> (displacement potential), p(pressure), and \phi (gravity potential).
> >>
> >> My code segment looks like the following:
> >> ======================================================
> >> call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)
> >> call KSPGetPC(ksp,pc,ierr)
> >> nsplit=4
> >> call PCFieldSplitSetBlockSize(pc,nsplit,ierr);
> >
> >
> > You do not need the statement above.
> >
> >>
> >> call PCFieldSplitSetIS(pc,"0",gdofu_is,ierr);
> >> call ISDestroy(gdofu_is,ierr)
> >> call PCFieldSplitSetIS(pc,"1",gdofchi_is,ierr);
> >> call ISDestroy(gdofchi_is,ierr)
> >> call PCFieldSplitSetIS(pc,"2",gdofp_is,ierr);
> >> call ISDestroy(gdofp_is,ierr)
> >> call PCFieldSplitSetIS(pc,"3",gdofphi_is,ierr);
> >> call ISDestroy(gdofphi_is,ierr)
> >>
> >> call PCFieldSplitGetSubKSP(pc,nsplit,subksp,ierr);
> >
> >
> > These SetOperators() calls are wrong. I am not sure why you want them
> here.
> > Also, that means you do not need the call above.
> >
> >   Thanks,
> >
> >      Matt
> >
> >>
> >> call KSPSetOperators(subksp(1),Amat,Amat,ierr);
> >> call KSPSetOperators(subksp(2),Amat,Amat,ierr);
> >> call KSPSetOperators(subksp(3),Amat,Amat,ierr);
> >> call KSPSetOperators(subksp(4),Amat,Amat,ierr);
> >>
> >> call KSPSolve(ksp,bvec,xvec,ierr)
> >> ======================================================
> >>
> >> But I am getting the following error:
> >> [79]PETSC ERROR: Null argument, when expecting valid pointer
> >> [79]PETSC ERROR: Null Object: Parameter # 1
> >> [79]PETSC ERROR: #1 KSPSetOperators() line 536 in
> >> /home/hgharti/lsoft/petsc-3.6.3-intel16-mpi5/src/ksp/ksp/interf
> >>
> >> It looks like I am doing something wrong in "call KSPSetOperators" but
> >> I could not figure that out.
> >>
> >> Could anybody help me to fix this problem? I looked into almost all
> >> related examples but I could not really figure out the correct steps
> >> after "call PCFieldSplitSetIS".
> >>
> >> Any help will be greatly appreciated.
> >>
> >> Best,
> >> Hom nath
> >>
> >> On Sun, Jan 24, 2016 at 7:14 PM, Hom Nath Gharti <hng.email at gmail.com>
> >> wrote:
> >> > Thank you so much Matt! I will try.
> >> >
> >> > Hom Nath
> >> >
> >> > On Sun, Jan 24, 2016 at 6:26 AM, Matthew Knepley <knepley at gmail.com>
> >> > wrote:
> >> >> On Fri, Jan 22, 2016 at 2:19 PM, Hom Nath Gharti <
> hng.email at gmail.com>
> >> >> wrote:
> >> >>>
> >> >>> Dear all,
> >> >>>
> >> >>> I am new to PcFieldSplit.
> >> >>>
> >> >>> I have a matrix formed using MATMPIAIJ. Is it possible to use
> >> >>> PCFIELDSPLIT operations in this type of matrix? Or does it have to
> be
> >> >>> MATMPIBIJ or MATNEST format?
> >> >>
> >> >>
> >> >> Yes, you can split AIJ.
> >> >>
> >> >>>
> >> >>> If possible for MATMPIAIJ, could anybody provide me a simple example
> >> >>> or few steps? Variables in the equations are displacement vector,
> >> >>> scalar potential and pressure.
> >> >>
> >> >>
> >> >> If you do not have a collocated discretization, then you have to use
> >> >>
> >> >>
> >> >>
> >> >>
> http://www.mcs.anl.gov/petsc/petsc-current/docs/manualpages/PC/PCFieldSplitSetIS.html
> >> >>
> >> >>   Thanks,
> >> >>
> >> >>      Matt
> >> >>
> >> >>>
> >> >>> Thanks for help.
> >> >>>
> >> >>> Hom Nath
> >> >>
> >> >>
> >> >>
> >> >>
> >> >> --
> >> >> 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/20160202/eb9b0710/attachment.html>


More information about the petsc-users mailing list