libMesh / libmesh

libMesh github repository
http://libmesh.github.io
GNU Lesser General Public License v2.1
651 stars 285 forks source link

DM type 'libmesh' did not attach the DM to the matrix #3802

Open jltuhnu opened 6 months ago

jltuhnu commented 6 months ago
[1]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
[1]PETSC ERROR: Petsc has generated inconsistent data
[1]PETSC ERROR: DM type 'libmesh' did not attach the DM to the matrix
[1]PETSC ERROR: WARNING! There are option(s) set that were not used! Could be the program crashed before they were used or a spelling mistake, etc!
[1]PETSC ERROR: Option left: name:-ksp_converged_reason (no value)
[1]PETSC ERROR: Option left: name:-snes_converged_reason (no value)

I want to use petsc_auto_fieldsplit(pc, system), The following code runs successfully on my mac (libmesh+petsc-3.19.0), but when I move it to the supercomputing (libmesh+petsc-3.18.0), it throws the petsc error.

void inbsolve(EquationSystems &es, int &n_iters, int &l_iters, SNESConvergedReason &reason)
{
  TransientNonlinearImplicitSystem &system =
      es.get_system<TransientNonlinearImplicitSystem>("twophase");
  NonlinearSolver<Number> *nonlinear_solver = system.nonlinear_solver.get();
  PetscNonlinearSolver<Number> *petsc_nonlinear_solver =
      cast_ptr<PetscNonlinearSolver<Number> *>(nonlinear_solver);

  petsc_nonlinear_solver->init();
  PetscVector<Number> &x =
      *(cast_ptr<PetscVector<Number> *>(system.solution.get()));
  PetscMatrix<Number> &jac =
      *(cast_ptr<PetscMatrix<Number> *>(system.matrix));
  PetscVector<Number> &r =
      *(cast_ptr<PetscVector<Number> *>(system.rhs));

  PetscErrorCode ierr = 0;
  // snes settings
  SNES snes = petsc_nonlinear_solver->snes();

  ierr = SNESSetFromOptions(snes);
  LIBMESH_CHKERR(ierr);
  // ksp settings
  KSP ksp;
  ierr = SNESGetKSP(snes, &ksp);
  LIBMESH_CHKERR(ierr);
  ierr = KSPSetFromOptions(ksp);
  LIBMESH_CHKERR(ierr);
  // pc settings
  PC pc;
  ierr = KSPGetPC(ksp, &pc);
  LIBMESH_CHKERR(ierr);
  petsc_auto_fieldsplit(pc, system);
  ierr = PCSetFromOptions(pc);
  LIBMESH_CHKERR(ierr);

  ierr = SNESSolve(snes, nullptr, x.vec());
  LIBMESH_CHKERR(ierr);
  ierr = SNESGetIterationNumber(snes, &n_iters);
  LIBMESH_CHKERR(ierr);

  ierr = SNESGetLinearSolveIterations(snes, &l_iters);
  LIBMESH_CHKERR(ierr);
  ierr = SNESGetConvergedReason(snes, &reason);
  LIBMESH_CHKERR(ierr);
  system.update();
}
jwpeterson commented 6 months ago

but when I move it to the supercomputing

Just to be clear, you are using the same version of libmesh in both places, and when you say "move it" you mean that you re-complied both your application and libmesh itself using the native compilers of each system?

Which version (or git hash) of libmesh are you using exactly?

jltuhnu commented 6 months ago

you mean that you re-complied both your application and libmesh itself using the native compilers of each system? yes, I use the same version of libmesh, and re-complied both my application and libmesh in different environments. There is a significant version gap in MPI and GCC between the two instances. Which version (or git hash) of libmesh are you using exactly? the git HEAD of libmesh :

commit cd2192759881961864dce40ecbc66af2b7c9c1fd (HEAD -> devel, origin/master, origin/devel, origin/HEAD)
Merge: 7729913f8 3905097b5
Author: John W. Peterson <jwpeterson@gmail.com>
Date:   Thu Feb 29 17:59:17 2024 -0600

I want to emphasize that when I use equation_systems.get_system("twophase").solve(), everything works fine. but i need much more infomation of the slover so I write the void inbsolve()function.

jwpeterson commented 6 months ago

OK, well I see some relevant looking commented out code in nonlinear_implicit_system.C:

  // FIXME - this is necessary for petsc_auto_fieldsplit
  // nonlinear_solver->init_names(*this);

The commented lines were added in 7a24259a92 by @roystgnr but there's no hint as to why they are commented out? This doesn't necessarily explain why you would see different behavior from different versions of PETSc, though. I don't know much about how the petsc_auto_fieldsplit code works, and I'm not even sure if we have tests or examples of using it with nonlinear solvers.

jltuhnu commented 6 months ago

@jwpeterson It seems that this issue is not caused by petsc_auto_fieldsplit. Even after removing petsc_auto_fieldsplit(pc, system);, the error persists for non-fieldsplit preconditioners. I am wondering if there are any missing elements related to the DM type 'libmesh' in my code that could be causing this issue.

jwpeterson commented 6 months ago

@jltuhnu are you saying that you are not using fieldsplit preconditioner (even via PETSc command line args) and you are still getting that same error message (DM type 'libmesh' did not attach the DM to the matrix)? If so, I'm not sure how that's possible since I don't think anything other than the fieldsplit preconditioner uses the DM...

jltuhnu commented 6 months ago

@jwpeterson are you saying that you are not using fieldsplit preconditioner (even via PETSc command line args) and you are still getting that same error message (DM type 'libmesh' did not attach the DM to the matrix) Yes, when I use asm preconditioner also get that error, Although I resolved my issue by modifying the libmesh nonlinear solve() function, I am still curious about the cause of the error with different Petsc versions.