bool Problem_Interface::computePreconditioner(const Epetra_Vector & x,
                                              Epetra_Operator & prec,
                                              Teuchos::ParameterList * /*p*/)
{
  START_LOG("preconditioner()", "TrilinosNoxNonlinearSolver");

  NonlinearImplicitSystem & sys = _solver->system();
  TrilinosPreconditioner<Number> & tpc = dynamic_cast<TrilinosPreconditioner<Number> &>(prec);

  EpetraMatrix<Number> Jac(dynamic_cast<Epetra_FECrsMatrix *>(tpc.mat()),sys.comm());
  EpetraVector<Number> & X_sys = *cast_ptr<EpetraVector<Number> *>(sys.solution.get());
  EpetraVector<Number> X_global(*const_cast<Epetra_Vector *>(&x), sys.comm());

  // Set the dof maps
  Jac.attach_dof_map(sys.get_dof_map());

  // Use the systems update() to get a good local version of the parallel solution
  X_global.swap(X_sys);

  sys.get_dof_map().enforce_constraints_exactly(sys);
  sys.update();

  X_global.swap(X_sys);

  //-----------------------------------------------------------------------------
  // if the user has provided both function pointers and objects only the pointer
  // will be used, so catch that as an error
  if (_solver->jacobian && _solver->jacobian_object)
    libmesh_error_msg("ERROR: cannot specify both a function and object to compute the Jacobian!");

  if (_solver->matvec && _solver->residual_and_jacobian_object)
    libmesh_error_msg("ERROR: cannot specify both a function and object to compute the combined Residual & Jacobian!");

  if      (_solver->jacobian != NULL)                     _solver->jacobian                                            (*sys.current_local_solution.get(), Jac, sys);
  else if (_solver->jacobian_object != NULL)              _solver->jacobian_object->jacobian                           (*sys.current_local_solution.get(), Jac, sys);
  else if (_solver->matvec   != NULL)                     _solver->matvec                                              (*sys.current_local_solution.get(), NULL, &Jac, sys);
  else if (_solver->residual_and_jacobian_object != NULL) _solver->residual_and_jacobian_object->residual_and_jacobian (*sys.current_local_solution.get(), NULL, &Jac, sys);
  else
    libmesh_error_msg("Error! Unable to compute residual and/or Jacobian!");

  Jac.close();
  X_global.close();

  tpc.compute();

  STOP_LOG("preconditioner()", "TrilinosNoxNonlinearSolver");

  return true;
}
bool Problem_Interface::computeF(const Epetra_Vector & x,
                                 Epetra_Vector & r,
                                 NOX::Epetra::Interface::Required::FillType /*fillType*/)
{
  START_LOG("residual()", "TrilinosNoxNonlinearSolver");

  NonlinearImplicitSystem & sys = _solver->system();

  EpetraVector<Number> X_global(*const_cast<Epetra_Vector *>(&x), sys.comm()), R(r, sys.comm());
  EpetraVector<Number> & X_sys = *cast_ptr<EpetraVector<Number> *>(sys.solution.get());
  EpetraVector<Number> & R_sys = *cast_ptr<EpetraVector<Number> *>(sys.rhs);

  // Use the systems update() to get a good local version of the parallel solution
  X_global.swap(X_sys);
  R.swap(R_sys);

  sys.get_dof_map().enforce_constraints_exactly(sys);
  sys.update();

  // Swap back
  X_global.swap(X_sys);
  R.swap(R_sys);

  R.zero();

  //-----------------------------------------------------------------------------
  // if the user has provided both function pointers and objects only the pointer
  // will be used, so catch that as an error

  if (_solver->residual && _solver->residual_object)
    libmesh_error_msg("ERROR: cannot specifiy both a function and object to compute the Residual!");

  if (_solver->matvec && _solver->residual_and_jacobian_object)
    libmesh_error_msg("ERROR: cannot specifiy both a function and object to compute the combined Residual & Jacobian!");

  if      (_solver->residual != NULL)                     _solver->residual                                            (*sys.current_local_solution.get(), R, sys);
  else if (_solver->residual_object != NULL)              _solver->residual_object->residual                           (*sys.current_local_solution.get(), R, sys);
  else if (_solver->matvec   != NULL)                     _solver->matvec                                              (*sys.current_local_solution.get(), &R, NULL, sys);
  else if (_solver->residual_and_jacobian_object != NULL) _solver->residual_and_jacobian_object->residual_and_jacobian (*sys.current_local_solution.get(), &R, NULL, sys);
  else return false;

  R.close();
  X_global.close();

  STOP_LOG("residual()", "TrilinosNoxNonlinearSolver");

  return true;
}
Beispiel #3
0
  //---------------------------------------------------------------
  // this function is called by PETSc to evaluate the residual at X
  PetscErrorCode
  __libmesh_petsc_snes_residual (SNES snes, Vec x, Vec r, void *ctx)
  {
    START_LOG("residual()", "PetscNonlinearSolver");

    PetscErrorCode ierr=0;

    libmesh_assert(x);
    libmesh_assert(r);
    libmesh_assert(ctx);

    PetscNonlinearSolver<Number>* solver =
      static_cast<PetscNonlinearSolver<Number>*> (ctx);

    // Get the current iteration number from the snes object,
    // store it in the PetscNonlinearSolver object for possible use
    // by the user's residual function.
    {
      PetscInt n_iterations = 0;
      ierr = SNESGetIterationNumber(snes, &n_iterations);
      CHKERRABORT(solver->comm().get(),ierr);
      solver->set_current_nonlinear_iteration_number( static_cast<unsigned>(n_iterations) );
    }

    NonlinearImplicitSystem &sys = solver->system();

    PetscVector<Number>& X_sys = *libmesh_cast_ptr<PetscVector<Number>*>(sys.solution.get());
    PetscVector<Number>& R_sys = *libmesh_cast_ptr<PetscVector<Number>*>(sys.rhs);
    PetscVector<Number> X_global(x, sys.comm()), R(r, sys.comm());

    // Use the systems update() to get a good local version of the parallel solution
    X_global.swap(X_sys);
    R.swap(R_sys);

    sys.get_dof_map().enforce_constraints_exactly(sys);

    sys.update();

    //Swap back
    X_global.swap(X_sys);
    R.swap(R_sys);

    R.zero();

    //-----------------------------------------------------------------------------
    // if the user has provided both function pointers and objects only the pointer
    // will be used, so catch that as an error
    if (solver->residual && solver->residual_object)
      {
	libMesh::err << "ERROR: cannot specifiy both a function and object to compute the Residual!" << std::endl;
	libmesh_error();
      }

    if (solver->matvec && solver->residual_and_jacobian_object)
      {
	libMesh::err << "ERROR: cannot specifiy both a function and object to compute the combined Residual & Jacobian!" << std::endl;
	libmesh_error();
      }
    //-----------------------------------------------------------------------------

    if      (solver->residual != NULL)                     solver->residual                                            (*sys.current_local_solution.get(), R, sys);
    else if (solver->residual_object != NULL)              solver->residual_object->residual                           (*sys.current_local_solution.get(), R, sys);
    else if (solver->matvec   != NULL)                     solver->matvec                                              (*sys.current_local_solution.get(), &R, NULL, sys);
    else if (solver->residual_and_jacobian_object != NULL) solver->residual_and_jacobian_object->residual_and_jacobian (*sys.current_local_solution.get(), &R, NULL, sys);
    else libmesh_error();

    R.close();
    X_global.close();

    STOP_LOG("residual()", "PetscNonlinearSolver");

    return ierr;
  }
Beispiel #4
0
  //---------------------------------------------------------------
  // this function is called by PETSc to evaluate the Jacobian at X
  PetscErrorCode
  __libmesh_petsc_snes_jacobian (SNES snes, Vec x, Mat *jac, Mat *pc, MatStructure *msflag, void *ctx)
  {
    START_LOG("jacobian()", "PetscNonlinearSolver");

    PetscErrorCode ierr=0;

    libmesh_assert(ctx);

    PetscNonlinearSolver<Number>* solver =
      static_cast<PetscNonlinearSolver<Number>*> (ctx);

    // Get the current iteration number from the snes object,
    // store it in the PetscNonlinearSolver object for possible use
    // by the user's Jacobian function.
    {
      PetscInt n_iterations = 0;
      ierr = SNESGetIterationNumber(snes, &n_iterations);
      CHKERRABORT(solver->comm().get(),ierr);
      solver->set_current_nonlinear_iteration_number( static_cast<unsigned>(n_iterations) );
    }

    NonlinearImplicitSystem &sys = solver->system();

    PetscMatrix<Number> PC(*pc, sys.comm());
    PetscMatrix<Number> Jac(*jac, sys.comm());
    PetscVector<Number>& X_sys = *libmesh_cast_ptr<PetscVector<Number>*>(sys.solution.get());
    PetscMatrix<Number>& Jac_sys = *libmesh_cast_ptr<PetscMatrix<Number>*>(sys.matrix);
    PetscVector<Number> X_global(x, sys.comm());

    // Set the dof maps
    PC.attach_dof_map(sys.get_dof_map());
    Jac.attach_dof_map(sys.get_dof_map());

    // Use the systems update() to get a good local version of the parallel solution
    X_global.swap(X_sys);
    Jac.swap(Jac_sys);

    sys.get_dof_map().enforce_constraints_exactly(sys);
    sys.update();

    X_global.swap(X_sys);
    Jac.swap(Jac_sys);

    PC.zero();

    //-----------------------------------------------------------------------------
    // if the user has provided both function pointers and objects only the pointer
    // will be used, so catch that as an error
    if (solver->jacobian && solver->jacobian_object)
      {
	libMesh::err << "ERROR: cannot specifiy both a function and object to compute the Jacobian!" << std::endl;
	libmesh_error();
      }

    if (solver->matvec && solver->residual_and_jacobian_object)
      {
	libMesh::err << "ERROR: cannot specifiy both a function and object to compute the combined Residual & Jacobian!" << std::endl;
	libmesh_error();
      }
    //-----------------------------------------------------------------------------

    if      (solver->jacobian != NULL)                     solver->jacobian                                            (*sys.current_local_solution.get(), PC, sys);
    else if (solver->jacobian_object != NULL)              solver->jacobian_object->jacobian                           (*sys.current_local_solution.get(), PC, sys);
    else if (solver->matvec   != NULL)                     solver->matvec                                              (*sys.current_local_solution.get(), NULL, &PC, sys);
    else if (solver->residual_and_jacobian_object != NULL) solver->residual_and_jacobian_object->residual_and_jacobian (*sys.current_local_solution.get(), NULL, &PC, sys);
    else libmesh_error();

    PC.close();
    Jac.close();
    X_global.close();

    *msflag = SAME_NONZERO_PATTERN;

    STOP_LOG("jacobian()", "PetscNonlinearSolver");

    return ierr;
  }
Beispiel #5
0
  //---------------------------------------------------------------
  // this function is called by PETSc to evaluate the residual at X
  PetscErrorCode
  __libmesh_petsc_snes_residual (SNES snes, Vec x, Vec r, void * ctx)
  {
    LOG_SCOPE("residual()", "PetscNonlinearSolver");

    PetscErrorCode ierr=0;

    libmesh_assert(x);
    libmesh_assert(r);
    libmesh_assert(ctx);

    // No way to safety-check this cast, since we got a void *...
    PetscNonlinearSolver<Number> * solver =
      static_cast<PetscNonlinearSolver<Number> *> (ctx);

    // Get the current iteration number from the snes object,
    // store it in the PetscNonlinearSolver object for possible use
    // by the user's residual function.
    {
      PetscInt n_iterations = 0;
      ierr = SNESGetIterationNumber(snes, &n_iterations);
      CHKERRABORT(solver->comm().get(),ierr);
      solver->_current_nonlinear_iteration_number = cast_int<unsigned>(n_iterations);
    }

    NonlinearImplicitSystem & sys = solver->system();

    PetscVector<Number> & X_sys = *cast_ptr<PetscVector<Number> *>(sys.solution.get());
    PetscVector<Number> X_global(x, sys.comm()), R(r, sys.comm());

    // Use the system's update() to get a good local version of the
    // parallel solution.  This operation does not modify the incoming
    // "x" vector, it only localizes information from "x" into
    // sys.current_local_solution.
    X_global.swap(X_sys);
    sys.update();
    X_global.swap(X_sys);

    // Enforce constraints (if any) exactly on the
    // current_local_solution.  This is the solution vector that is
    // actually used in the computation of the residual below, and is
    // not locked by debug-enabled PETSc the way that "x" is.
    sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());

    if (solver->_zero_out_residual)
      R.zero();

    //-----------------------------------------------------------------------------
    // if the user has provided both function pointers and objects only the pointer
    // will be used, so catch that as an error
    if (solver->residual && solver->residual_object)
      libmesh_error_msg("ERROR: cannot specifiy both a function and object to compute the Residual!");

    if (solver->matvec && solver->residual_and_jacobian_object)
      libmesh_error_msg("ERROR: cannot specifiy both a function and object to compute the combined Residual & Jacobian!");

    if (solver->residual != libmesh_nullptr)
      solver->residual(*sys.current_local_solution.get(), R, sys);

    else if (solver->residual_object != libmesh_nullptr)
      solver->residual_object->residual(*sys.current_local_solution.get(), R, sys);

    else if (solver->matvec != libmesh_nullptr)
      solver->matvec (*sys.current_local_solution.get(), &R, libmesh_nullptr, sys);

    else if (solver->residual_and_jacobian_object != libmesh_nullptr)
      solver->residual_and_jacobian_object->residual_and_jacobian (*sys.current_local_solution.get(), &R, libmesh_nullptr, sys);

    else
      libmesh_error_msg("Error! Unable to compute residual and/or Jacobian!");

    R.close();

    return ierr;
  }