[petsc-users] Problem with least squares commutators(LSC)

keguoyi coyigg at hotmail.com
Thu Feb 4 09:45:50 CST 2016


Dear Users,

Thanks a lot for the reply. It really helps a lot.

The code for -fieldsplit-pressure is here, and we just let PETSc compute it from the blocks of the full matrix we provide.
if( _preconditioner == LSC_PRECOND ) {
      _rtol = 1.e-3;
      _abstol = 1.e-20;
      _dtol = 1.e+50;
      _maxits = 1;

      SetPetscSolverType(ksp);
      PC pc;
      KSPGetPC( ksp, &pc );
      KSPSetTolerances( ksp, _rtol, _abstol, _dtol, _maxits );
      KSPSetFromOptions( ksp );
      PCSetType( pc, PCLSC );
    }

We think what our problem is

This is a NavierStokes problem with Dirchhlet Boundary conditions everywhere but one point (to fix the pressure).  We impose BC using penalty method, thus changing the diagonal of A00, but we leave A01 unchanged. So in evaluating A10*A01 it is like we are having a Neumann Problem everywhere.   

We have to rethink how to impose the BC.

If you have other suggestions, please let us know. Thank you so much!!

Best,
Guoyi





> Subject: Re: [petsc-users] Problem with least squares commutators(LSC)
> From: bsmith at mcs.anl.gov
> Date: Thu, 4 Feb 2016 00:09:57 -0600
> CC: petsc-users at mcs.anl.gov
> To: coyigg at hotmail.com
> 
> 
>   It is failing due to zero diagonal entry in factoring the matrix Lp which is A10*A01 which I think in your case is B * B^T
> 
>   Do you use code like? 
> 
>    PetscObjectCompose((PetscObject)Sp,"LSC_L",(PetscObject)L);
>    PetscObjectCompose((PetscObject)Sp,"LSC_Lp",(PetscObject)Lp);
> 
>    And then your Jacobian assembly would look like
>    PetscObjectQuery((PetscObject)Sp,"LSC_L",(PetscObject*)&L);
>    PetscObjectQuery((PetscObject)Sp,"LSC_Lp",(PetscObject*)&Lp);
>    if (L) { assembly L }
>    if (Lp) { assemble Lp }
> 
>   That is are you providing a matrix that contains B * B^T already computed or are you just letting PETSc compute it from the blocks of the 
> full matrix you provide? (I am guessing you let PETSc compute it). 
> 
>    PETSc computes it in the routine PCSetUp_LSC() 
> 
>     ierr = MatSchurComplementGetSubMatrices(pc->mat,NULL,NULL,&B,&C,NULL);CHKERRQ(ierr);
>     if (!lsc->L) {
>       ierr = MatMatMult(C,B,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&lsc->L);CHKERRQ(ierr);
>     } else {
>       ierr = MatMatMult(C,B,MAT_REUSE_MATRIX,PETSC_DEFAULT,&lsc->L);CHKERRQ(ierr);
>     }
>     Lp = L = lsc->L;
> 
> It would be useful to see what Lp is. You can run with -start_in_debugger put a breakpoint in PCSetUp_LSC then at line Lp = L = lsc->L;
> and call MatView(Lp,0) and MatView(B,0) and MatView(C,0) at that line to see what the matrices look like. Of course run on 1 process with a small size for the matrices.
> 
>   Do you maybe have a periodic domain or Neumann boundary conditions so that B * B^T is singular?
> 
>   Barry
> 
> 
> 
> 
> > On Feb 3, 2016, at 11:30 PM, keguoyi <coyigg at hotmail.com> wrote:
> > 
> > Dear Users,
> > 
> > The following is the full error message. Thank you so much!
> > 
> >  *************** Linear iteration 1 ***********
> > [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
> > [0]PETSC ERROR: Object is in wrong state
> > [0]PETSC ERROR: Matrix is missing diagonal entry 0
> > [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
> > [0]PETSC ERROR: Petsc Release Version 3.6.1, unknown
> > [0]PETSC ERROR: ./PCFieldSplit_ex13 on a arch-linux2-cxx-opt named linux-394h by coyigg Wed Feb  3 14:57:21 2016
> > [0]PETSC ERROR: Configure options --with-clanguage=cxx --with-debugging=0 --download-fblaslapack=1 --download-hdf5=1 --download-openmpi=1 --download-metis=1 --download-parmetis=1 --with-shared-libraries=1 --with-cc=gcc --with-fc=gfortran --with-cxx=g++ --download-blacs=1 --download-scalapack=1 --download-mumps=1
> > [0]PETSC ERROR: #1 MatILUFactorSymbolic_SeqAIJ() line 1729 in /home/eaulisa/software/petsc/src/mat/impls/aij/seq/aijfact.c
> > [0]PETSC ERROR: #2 MatILUFactorSymbolic() line 6457 in /home/eaulisa/software/petsc/src/mat/interface/matrix.c
> > [0]PETSC ERROR: #3 PCSetUp_ILU() line 204 in /home/eaulisa/software/petsc/src/ksp/pc/impls/factor/ilu/ilu.c
> > [0]PETSC ERROR: #4 PCSetUp() line 983 in /home/eaulisa/software/petsc/src/ksp/pc/interface/precon.c
> > [0]PETSC ERROR: #5 KSPSetUp() line 332 in /home/eaulisa/software/petsc/src/ksp/ksp/interface/itfunc.c
> > [0]PETSC ERROR: #6 KSPSolve() line 546 in /home/eaulisa/software/petsc/src/ksp/ksp/interface/itfunc.c
> > [0]PETSC ERROR: #7 PCApply_LSC() line 83 in /home/eaulisa/software/petsc/src/ksp/pc/impls/lsc/lsc.c
> > [0]PETSC ERROR: #8 PCApply() line 483 in /home/eaulisa/software/petsc/src/ksp/pc/interface/precon.c
> > [0]PETSC ERROR: #9 KSP_PCApply() line 242 in /home/eaulisa/software/petsc/include/petsc/private/kspimpl.h
> > [0]PETSC ERROR: #10 KSPSolve_PREONLY() line 26 in /home/eaulisa/software/petsc/src/ksp/ksp/impls/preonly/preonly.c
> > [0]PETSC ERROR: #11 KSPSolve() line 604 in /home/eaulisa/software/petsc/src/ksp/ksp/interface/itfunc.c
> > [0]PETSC ERROR: #12 PCApply_FieldSplit_Schur() line 874 in /home/eaulisa/software/petsc/src/ksp/pc/impls/fieldsplit/fieldsplit.c
> > [0]PETSC ERROR: #13 PCApply() line 483 in /home/eaulisa/software/petsc/src/ksp/pc/interface/precon.c
> > [0]PETSC ERROR: #14 KSP_PCApply() line 242 in /home/eaulisa/software/petsc/include/petsc/private/kspimpl.h
> > [0]PETSC ERROR: #15 KSPInitialResidual() line 63 in /home/eaulisa/software/petsc/src/ksp/ksp/interface/itres.c
> > [0]PETSC ERROR: #16 KSPSolve_GMRES() line 235 in /home/eaulisa/software/petsc/src/ksp/ksp/impls/gmres/gmres.c
> > [0]PETSC ERROR: #17 KSPSolve() line 604 in /home/eaulisa/software/petsc/src/ksp/ksp/interface/itfunc.c
> > [0]PETSC ERROR: #18 PCMGMCycle_Private() line 19 in /home/eaulisa/software/petsc/src/ksp/pc/impls/mg/mg.c
> > [0]PETSC ERROR: #19 PCApply_MG() line 340 in /home/eaulisa/software/petsc/src/ksp/pc/impls/mg/mg.c
> > [0]PETSC ERROR: #20 PCApply() line 483 in /home/eaulisa/software/petsc/src/ksp/pc/interface/precon.c
> > [0]PETSC ERROR: #21 KSPSolve() line 582 in /home/eaulisa/software/petsc/src/ksp/ksp/interface/itfunc.c
> >  *************** MG linear solver time:    0.002
> >  *************** Number of outer ksp solver iterations = 0
> >  *************** Convergence reason = 0
> >  *************** Residual norm =          0
> >  *************** Level Max 2  Linear Res  L2norm U = 4.033131e-02
> >  *************** Level Max 2  Linear Res  L2norm V = 5.330749e-02
> >  *************** Level Max 2  Linear Res  L2norm P = 1.396813e-01
> > 
> > 
> > > Subject: Re: [petsc-users] Problem with least squares commutators(LSC)
> > > From: bsmith at mcs.anl.gov
> > > Date: Wed, 3 Feb 2016 22:26:42 -0600
> > > CC: petsc-users at mcs.anl.gov
> > > To: coyigg at hotmail.com
> > > 
> > > 
> > > We'd need the full error message to determine what is happening.
> > > 
> > > But note that since the lower right block is all 0 any attempt to use a factorization of that inside a preconditioner will fail with zero pivot. 
> > > 
> > > 
> > > > On Feb 2, 2016, at 10:40 PM, keguoyi <coyigg at hotmail.com> wrote:
> > > > 
> > > > Dear PETSc users,
> > > > 
> > > > This is Guoyi Ke, a graduate student of Texas Tech University. I have a 2D Navier Stokes problem that has block matrices: J=[F B^T; B 0]. I use Schur complement preconditioner for block J. The code was set as:
> > > > 
> > > > PCFieldSplitSetType( pc, PC_COMPOSITE_SCHUR );
> > > > PCFieldSplitSetSchurFactType(pc, PC_FIELDSPLIT_SCHUR_FACT_LOWER);
> > > > PCFieldSplitSetSchurPre(pc,PC_FIELDSPLIT_SCHUR_PRE_SELFP,NULL);
> > > > 
> > > > For both velocity block (fieldsplit: 0) and pressure block (fieldsplit: 1) , we set the KSP as KSPPREONLY and PC as PCILU. For this case, it works fine for us. 
> > > > However, in order to use LSC we change the code as:
> > > > 
> > > > PCFieldSplitSetType( pc, PC_COMPOSITE_SCHUR );
> > > > PCFieldSplitSetSchurFactType(pc, PC_FIELDSPLIT_SCHUR_FACT_LOWER);
> > > > PCFieldSplitSetSchurPre(pc,PC_FIELDSPLIT_SCHUR_PRE_SELF,NULL);
> > > > 
> > > > We keep the same setup for velocity block (fieldsplit: 0) as above. For pressure block (fieldsplit: 1), we set the KSP as KSPGMRES and PC as PCLSU. I got the following errors:
> > > > 
> > > > [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
> > > > [0]PETSC ERROR: Object is in wrong state
> > > > [0]PETSC ERROR: Matrix is missing diagonal entry 0
> > > > 
> > > > 
> > > > Any suggestion and help will be highly appreciated. Thank you so much!
> > > > 
> > > > Best,
> > > > Guoyi 
> > > 
> 
 		 	   		  
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20160204/8afc230b/attachment.html>


More information about the petsc-users mailing list