[petsc-users] [petsc-3.2] LSQR new convergence test

Satish Balay balay at mcs.anl.gov
Thu Sep 22 07:54:02 CDT 2011


Loks like Barry sent you the wrong file. zitfuncf90.c belongs to
src/ksp/ksp/interface/f90-custom [and is unchanged]

Try the attached patch

cd petsc-3.2-p2
patch -Np1 < ksp-conv.patch

Satish

On Thu, 22 Sep 2011, Alexander Grayver wrote:

> Barry,
> 
> Unfortunately, the problem is still there. Even though I put zitfuncf90.c into
> src/ksp/ksp/interface/ftn-custom, modified makefile there and made
> configure/make from scratch.
> When are you going to release new patch?
> 
> Regards,
> Alexander
> 
> On 21.09.2011 21:58, Barry Smith wrote:
> >     Alexander,
> > 
> >      Yes, a problem with trying to handle C pointers on the Fortran side. It
> > will be fixed in the next patch release. You can drop the attached file into
> > src/ksp/ksp/interface/ftn-custom run make in that directory and then relink
> > your program and it should work.
> > 
> >      Barry
> > [see attached file: zitfuncf90.c]
> > On Sep 21, 2011, at 11:00 AM, Alexander Grayver wrote:
> > 
> > > Barry,
> > > 
> > > I think the error is related to this
> > > 
> > > PetscFortranAddr cctx
> > > 
> > > definition. Speaking in terms of C++ notation, it seems that PETSc
> > > takes&cctx as KSPDefaultConvergedCtx, however real KSPDefaultConvergedCtx
> > > is *cctx.
> > > Not sure if I'm clear now, but it's easy to see if you trace code like
> > > that:
> > > 
> > >   PetscFortranAddr      cctx
> > >   KSP  ksp
> > > 
> > >   call KSPCreate(MPI_COMM_WORLD,ksp,ierr)
> > >   call KSPSetType(ksp,KSPLSQR,ierr)
> > >   call KSPDefaultConvergedCreate(cctx,ierr)
> > >   call
> > > KSPSetConvergenceTest(ksp,KSPDefaultConverged,cctx,PETSC_NULL_FUNCTION,ierr)
> > >   call KSPDestroy(ksp,ierr)
> > > 
> > > I guess the error should be reproducible from the code above.
> > > 
> > > Regards,
> > > Alexander
> > > 
> > > On 21.09.2011 17:17, Alexander Grayver wrote:
> > > > Thanks Barry, I added your code and solver works well now, but when I
> > > > call KSPDestroy afterward it crashes:
> > > > 
> > > > [8]PETSC ERROR: --------------------- Error Message
> > > > ------------------------------------
> > > > [8]PETSC ERROR: Invalid argument!
> > > > [8]PETSC ERROR: Wrong type of object: Parameter # 1!
> > > > [8]PETSC ERROR:
> > > > ------------------------------------------------------------------------
> > > > [8]PETSC ERROR: Petsc Release Version 3.2.0, Patch 2, Fri Sep 16
> > > > 10:10:45 CDT 2011
> > > > [8]PETSC ERROR: See docs/changes/index.html for recent updates.
> > > > [8]PETSC ERROR: See docs/faq.html for hints about trouble shooting.
> > > > [8]PETSC ERROR: See docs/index.html for manual pages.
> > > > [8]PETSC ERROR:
> > > > ------------------------------------------------------------------------
> > > > [8]PETSC ERROR: /home/test on a openmpi-i named node226 by user Wed Sep
> > > > 21 17:01:23 2011
> > > > [8]PETSC ERROR: Libraries linked from
> > > > /home/lib/petsc-3.2-p2/openmpi-intel-complex-debug-f/lib
> > > > [8]PETSC ERROR: Configure run at Tue Sep 20 09:28:29 2011
> > > > [8]PETSC ERROR: Configure options
> > > > --with-petsc-arch=openmpi-intel-complex-debug-f
> > > > --with-fortran-interfaces=1 --with-mpi-dir=/opt/mpi/intel/openmpi-1.4.2
> > > > --with-scalar-type=complex
> > > > --with-blas-lapack-dir=/opt/intel/Compiler/11.1/072/mkl/lib/em64t
> > > > --with-precision=double --with-x=0
> > > > [8]PETSC ERROR:
> > > > ------------------------------------------------------------------------
> > > > [8]PETSC ERROR: VecDestroy() line 566 in
> > > > /home/lib/petsc-3.2-p2/src/vec/vec/interface/vector.c
> > > > [8]PETSC ERROR: KSPDefaultConvergedDestroy() line 669 in
> > > > /home/lib/petsc-3.2-p2/src/ksp/ksp/interface/iterativ.c
> > > > [8]PETSC ERROR: KSPDestroy() line 758 in
> > > > /home/lib/petsc-3.2-p2/src/ksp/ksp/interface/itfunc.c
> > > > 
> > > > Is it my fault?
> > > > 
> > > > Regards,
> > > > Alexander
> > > > 
> > > > On 21.09.2011 15:36, Barry Smith wrote:
> > > > > On Sep 21, 2011, at 2:39 AM, Alexander Grayver wrote:
> > > > > 
> > > > > > Hong, Barry,
> > > > > > 
> > > > > > Thanks! Sorry for not to be clear, but Barry right, my question was
> > > > > > how to get back to DefaultConvergenceTest in case of LSQR solver.
> > > > > > 
> > > > > > Barry, yes, both versions give the same two norms.
> > > > > > 
> > > > > > Can you clarify please how to implement that in Fortran:
> > > > > > void *cctx;
> > > > > > KSPDefaultConvergedCreate(&cctx);
> > > > > > KSPSetConvergenceTest(ksp,KSPDefaultConverged,cctx);
> > > > >     It should be almost the same
> > > > > 
> > > > >      PetscFortranAddr cctx
> > > > > 
> > > > >      call KSPDefaultConvergedCreate(cctx,ierr)
> > > > >      call
> > > > > KSPSetConvergenceTest(ksp,KSPDefaultConverged,cctx,PETSC_NULL_FUNCTION,ierr)
> > > > > 
> > > > >     Barry
> > > > > 
> > > > > 
> > > > > > I'm a little bit confused about cctx parameter.
> > > > > > 
> > > > > > Regards,
> > > > > > Alexander
> > > > > > 
> > > > > > 
> > > > > > On 21.09.2011 01:10, Barry Smith wrote:
> > > > > > > It has its own monitor that provides additional information
> > > > > > > -ksp_monitor_lsqr
> > > > > > > 
> > > > > > > You can also remove the new convergence test and get back the old
> > > > > > > one with code like
> > > > > > > 
> > > > > > > void *cctx;
> > > > > > > KSPDefaultConvergedCreate(&cctx);
> > > > > > > KSPSetConvergenceTest(ksp,KSPDefaultConverged,cctx);
> > > > > > > 
> > > > > > > after the KSPType is set to LSQR. So if you are happy with the old
> > > > > > > test.
> > > > > > > 
> > > > > > > 
> > > > > > >    Do both versions give the same first two norms?
> > > > > > > 
> > > > > > > > 0 KSP Residual norm 9.386670021557e-17
> > > > > > > >   1 KSP Residual norm 8.258308650175e-17
> > > > > > >    Barry
> > > > > > > 
> > > > > > > 
> > > > > > > On Sep 20, 2011, at 4:40 AM, Alexander Grayver wrote:
> > > > > > > 
> > > > > > > > Hello,
> > > > > > > > 
> > > > > > > > In comparison with petsc-3.1-p7 in the latest petsc LSQR solver
> > > > > > > > behaves differently.
> > > > > > > > 
> > > > > > > > In petsc31 I had:
> > > > > > > > 
> > > > > > > >   0 KSP Residual norm 9.386670021557e-17
> > > > > > > > ...
> > > > > > > > 95 KSP Residual norm 9.341367317075e-18
> > > > > > > > Linear solve converged due to CONVERGED_RTOL iterations 95
> > > > > > > > KSP Object:
> > > > > > > >   type: lsqr
> > > > > > > >   maximum iterations=200, initial guess is zero
> > > > > > > >   tolerances:  relative=0.1, absolute=1e-50, divergence=10000
> > > > > > > >   left preconditioning
> > > > > > > >   using PRECONDITIONED norm type for convergence test
> > > > > > > > PC Object:
> > > > > > > >   type: none
> > > > > > > >   linear system matrix = precond matrix:
> > > > > > > >   Matrix Object:
> > > > > > > >     type=shell, rows=19584, cols=19584
> > > > > > > > 
> > > > > > > > relative residual norm =    0.9951737192872134E-01
> > > > > > > > 
> > > > > > > > Now, with petsc32 I have:
> > > > > > > > 
> > > > > > > >   0 KSP Residual norm 9.386670021557e-17
> > > > > > > >   1 KSP Residual norm 8.258308650175e-17
> > > > > > > > Linear solve converged due to CONVERGED_RTOL_NORMAL iterations 1
> > > > > > > > KSP Object: 12 MPI processes
> > > > > > > >   type: lsqr
> > > > > > > >   maximum iterations=200, initial guess is zero
> > > > > > > >   tolerances:  relative=0.1, absolute=1e-50, divergence=10000
> > > > > > > >   left preconditioning
> > > > > > > >   using UNPRECONDITIONED norm type for convergence test
> > > > > > > > PC Object: 12 MPI processes
> > > > > > > >   type: none
> > > > > > > >   linear system matrix = precond matrix:
> > > > > > > >   Matrix Object:   12 MPI processes
> > > > > > > >     type: shell
> > > > > > > >     rows=19584, cols=19584
> > > > > > > > relative residual norm =    0.8797910900468064E+00
> > > > > > > > 
> > > > > > > > So I found new convergence test here which sets
> > > > > > > > CONVERGED_RTOL_NORMAL:
> > > > > > > > http://www.mcs.anl.gov/petsc/petsc-as/snapshots/petsc-current/docs/manualpages/KSP/KSPLSQRDefaultConverged.html
> > > > > > > > 
> > > > > > > > The question is, how to interpret this new test and make it
> > > > > > > > works properly for me?
> > > > > > > > Thanks in advance.
> > > > > > > > 
> > > > > > > > Regards,
> > > > > > > > Alexander
> 
> 
-------------- next part --------------
tree 80c3b1f3210d
parent 98f24f0f10f1
author Barry Smith <bsmith at mcs.anl.gov> 1316634856 18000
committer Barry Smith <bsmith at mcs.anl.gov> 1316634856 18000
revision 20311
branch default

fix for using KSPSetConvergenceTest() with context from KSPDefaultConvergedCreate() in Fortran
diff --git a/src/ksp/ksp/interface/ftn-custom/zitfuncf.c b/src/ksp/ksp/interface/ftn-custom/zitfuncf.c
--- a/src/ksp/ksp/interface/ftn-custom/zitfuncf.c
+++ b/src/ksp/ksp/interface/ftn-custom/zitfuncf.c
@@ -8,6 +8,7 @@
 
 #define kspdefaultconverged_       KSPDEFAULTCONVERGED
 #define kspdefaultconvergedcreate_  KSPDEFAULTCONVERGEDCREATE
+#define kspdefaultconvergeddestroy_  KSPDEFAULTCONVERGEDDESTROY
 #define kspskipconverged_          KSPSKIPCONVERGED
 #define kspgmresmonitorkrylov_     KSPGMRESMONITORKRYLOV
 #define kspmonitordefault_         KSPMONITORDEFAULT
@@ -24,6 +25,7 @@
 #define kspgetresidualhistory_         kspgetresidualhistory
 #define kspdefaultconverged_           kspdefaultconverged
 #define kspdefaultconvergedcreate_     kspdefaultconvergedcreate
+#define kspdefaultconvergeddestroy_    kspdefaultconvergeddestroy
 #define kspskipconverged_              kspskipconverged
 #define kspmonitorsingularvalue_       kspmonitorsingularvalue
 #define kspgmresmonitorkrylov_         kspgmresmonitorkrylov
@@ -163,7 +165,7 @@
 }
 
 void PETSC_STDCALL kspsetconvergencetest_(KSP *ksp,
-      void (PETSC_STDCALL *converge)(KSP*,PetscInt*,PetscReal*,KSPConvergedReason*,void*,PetscErrorCode*),void *cctx,
+      void (PETSC_STDCALL *converge)(KSP*,PetscInt*,PetscReal*,KSPConvergedReason*,void*,PetscErrorCode*),void **cctx,
       void (PETSC_STDCALL *destroy)(void*,PetscErrorCode*),PetscErrorCode *ierr)
 {
   CHKFORTRANNULLOBJECT(cctx);
@@ -171,7 +173,7 @@
 
   PetscObjectAllocateFortranPointers(*ksp,6);
   if ((PetscVoidFunction)converge == (PetscVoidFunction)kspdefaultconverged_) {
-    *ierr = KSPSetConvergenceTest(*ksp,KSPDefaultConverged,cctx,KSPDefaultConvergedDestroy);
+    *ierr = KSPSetConvergenceTest(*ksp,KSPDefaultConverged,*cctx,KSPDefaultConvergedDestroy);
   } else if ((PetscVoidFunction)converge == (PetscVoidFunction)kspskipconverged_) {
     *ierr = KSPSetConvergenceTest(*ksp,KSPSkipConverged,0,0);
   } else {
@@ -191,6 +193,11 @@
   *ierr = KSPDefaultConvergedCreate((void**)ctx);
 }
 
+void PETSC_STDCALL kspdefaultconvergeddestroy_(PetscFortranAddr *ctx,PetscErrorCode *ierr)
+{
+  *ierr = KSPDefaultConvergedDestroy(*(void**)ctx);
+}
+
 void PETSC_STDCALL kspgetresidualhistory_(KSP *ksp,PetscInt *na,PetscErrorCode *ierr)
 {
   *ierr = KSPGetResidualHistory(*ksp,PETSC_NULL,na);


More information about the petsc-users mailing list