[petsc-users] SNES convergence issue

Tian(ICT) rongtian at ncic.ac.cn
Mon May 16 09:50:13 CDT 2011


./pgfem -pc_type none does not work either, leading to negative jacobian of element.
The SNES solve is soling a meshfree method.
I attached the segment of code.

What is the output for -mat_view, are all the output the jacobian, or do they inlcude preconditioning matrix?
I found the output matrices do not have the same nonzero structure.

PetscErrorCode FormFunction(SNES snes,Vec u, Vec r, void *ctx)
{
  Log::Send(Log::Info, " Form Function ...");
  ApplicationCtx *user = (ApplicationCtx*) ctx;
  //NOTE: r is the nonlinea function, NOT the r.h.s term
  //the r.h.s is '-r'
  //u is the solution during the current Newton-Raphson iteration
  //not an increment 

  //VecSet(u,0.001);
  //user->pDomain->pPhyDomain->GetMesh(0)->GetEBoundary()->Treatment(NULL,&u);

  //printf("solutoin vector.................\n");
  //VecView(u,PETSC_VIEWER_STDOUT_WORLD);

  //dof = du
  VecAYPX(user->oldu, -1, u);
  user->pDomain->pPhyDomain->UpdateDofs(user->oldu);
  VecCopy(u,user->oldu);

  //x = x + du
  user->pDomain->pPhyDomain->UpdateConfiguration(user->oldu);

  //QPoint's gradients in current config, new jacobian, new locations etc.
  user->pDomain->pPhyDomain->UpdateQPoints();

  //mfree connectivity in the new configuration
  user->pDomain->pPhyDomain->MFreeConnectivity();

  //F = (dX/dx)^-1
  user->pDomain->pPhyDomain->DeformGradients();

  //stress
  user->pDomain->pPhyDomain->UpdateStress();

  //f(x) of J(x) dx = f(x)
  FormInitialGuess(r);
  user->pDomain->pPhyDomain->InternalForce(r);
  VecAssemblyBegin(r);
  VecAssemblyEnd(r);

  user->pDomain->pPhyDomain->Loads(user->pDomain->contrl.loadContlType, -1.0*user->lambda, &u, &r);

  //printf("residual force.................\n");
  //VecView(r,PETSC_VIEWER_STDOUT_WORLD);

  PetscFunctionReturn(0);
}


PetscErrorCode FormJacobian(SNES snes, Vec x, Mat *jac/*jacobian*/, Mat *B/*precon*/,MatStructure*flag,void *ctx)
{
  Log::Send(Log::Info, " Form Jacobian ...");
  ApplicationCtx *user = (ApplicationCtx*) ctx;
  
  //initialize jacobian matrix
  MatZeroEntries(*B);

  //assemble
  user->pDomain->pPhyDomain->TangentMatrix(*B);
  MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
  MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);

  //boundary treatment, mainly MatZeroRows()
  user->pDomain->pPhyDomain->Constraints(user->pDomain->contrl.loadContlType, B);

  //printf("Global matrix........\n");
  //MatView(*B,PETSC_VIEWER_STDOUT_WORLD);

  if (*jac != *B){
    MatAssemblyBegin(*jac,MAT_FINAL_ASSEMBLY);
    MatAssemblyEnd(*jac,MAT_FINAL_ASSEMBLY);
  }

  ::DumpMatrixToDisk(*B);

  *flag = SAME_NONZERO_PATTERN;

  PetscFunctionReturn(0);
}




  ----- Original Message ----- 
  From: Jed Brown 
  To: PETSc users list 
  Sent: Monday, May 16, 2011 9:04 PM
  Subject: Re: [petsc-users] SNES convergence issue


  On Mon, May 16, 2011 at 14:54, Tian(ICT) <rongtian at ncic.ac.cn> wrote:

     12 KSP preconditioned resid norm 5.794625326821e+000 true resid norm 5.490686672936e+006 ||Ae||/||Ax|| 7.765003561961e+004
     13 KSP preconditioned resid norm 7.239823936251e-001 true resid norm 4.262009749173e+005 ||Ae||/||Ax|| 6.027391992073e+003
     14 KSP preconditioned resid norm 6.718858761053e-001 true resid norm 7.089916382344e+005 ||Ae||/||Ax|| 1.002665590704e+004
     15 KSP preconditioned resid norm 6.284315560280e-001 true resid norm 4.441924320726e+005 ||Ae||/||Ax|| 6.281829619309e+003
     16 KSP preconditioned resid norm 6.861179511523e-011 true resid norm 1.032548577293e+006 ||Ae||/||Ax|| 1.460244202260e+004


  Notice how the true residual is blowing up while the GMRES estimate converges. (The language here is not correct, "preconditioned resid norm" in this case is actually the unpreconditioned norm as estimated by FGMRES.) This indicates that the preconditioner is singular which could be due to lack of pivoting. If you run with -pc_type none, it should converge correctly. Do that and see if SNES is converging properly. Also try -pc_type lu -pc_factor_mat_solver_package superlu (or mumps or umfpack) if you have any of these packages installed.


  What equations are you solving and with what discretization?
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20110516/60000c65/attachment-0001.htm>


More information about the petsc-users mailing list