$extrastylesheet
libMesh::Problem_Interface Class Reference

List of all members.

Public Member Functions

 Problem_Interface (NoxNonlinearSolver< Number > *solver)
 ~Problem_Interface ()
bool computeF (const Epetra_Vector &x, Epetra_Vector &FVec, NOX::Epetra::Interface::Required::FillType fillType)
 Compute and return F.
bool computeJacobian (const Epetra_Vector &x, Epetra_Operator &Jac)
 Compute an explicit Jacobian.
bool computePrecMatrix (const Epetra_Vector &x, Epetra_RowMatrix &M)
 Compute the Epetra_RowMatrix M, that will be used by the Aztec preconditioner instead of the Jacobian. This is used when there is no explicit Jacobian present (i.e. Matrix-Free Newton-Krylov). This MUST BE an Epetra_RowMatrix since the Aztec preconditioners need to know the sparsity pattern of the matrix. Returns true if computation was successful.
bool computePreconditioner (const Epetra_Vector &x, Epetra_Operator &Prec, Teuchos::ParameterList *p)
 Computes a user supplied preconditioner based on input vector x. Returns true if computation was successful.

Public Attributes

NoxNonlinearSolver< Number > * _solver

Detailed Description

Definition at line 53 of file trilinos_nox_nonlinear_solver.C.


Constructor & Destructor Documentation


Member Function Documentation

bool libMesh::Problem_Interface::computeF ( const Epetra_Vector &  x,
Epetra_Vector &  FVec,
NOX::Epetra::Interface::Required::FillType  fillType 
)

Compute and return F.

Definition at line 88 of file trilinos_nox_nonlinear_solver.C.

References _solver, libMesh::ParallelObject::comm(), libMesh::System::current_local_solution, libMesh::DofMap::enforce_constraints_exactly(), libMesh::System::get_dof_map(), libMesh::NonlinearSolver< T >::matvec, libMesh::NonlinearImplicitSystem::ComputeResidual::residual(), libMesh::NonlinearSolver< T >::residual, libMesh::NonlinearImplicitSystem::ComputeResidualandJacobian::residual_and_jacobian(), libMesh::NonlinearSolver< T >::residual_and_jacobian_object, libMesh::NonlinearSolver< T >::residual_object, libMesh::ExplicitSystem::rhs, libMesh::System::solution, libMesh::START_LOG(), libMesh::sys, libMesh::NonlinearSolver< T >::system(), libMesh::System::update(), libMesh::X_global(), and libMesh::X_sys.

{
  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;
}
bool libMesh::Problem_Interface::computeJacobian ( const Epetra_Vector &  x,
Epetra_Operator &  Jac 
)

Compute an explicit Jacobian.

Definition at line 136 of file trilinos_nox_nonlinear_solver.C.

References _solver, libMesh::ParallelObject::comm(), libMesh::System::current_local_solution, libMesh::DofMap::enforce_constraints_exactly(), libMesh::System::get_dof_map(), libMesh::Jac(), libMesh::NonlinearImplicitSystem::ComputeJacobian::jacobian(), libMesh::NonlinearSolver< T >::jacobian, libMesh::NonlinearSolver< T >::jacobian_object, libMesh::NonlinearSolver< T >::matvec, libMesh::NonlinearImplicitSystem::ComputeResidualandJacobian::residual_and_jacobian(), libMesh::NonlinearSolver< T >::residual_and_jacobian_object, libMesh::System::solution, libMesh::START_LOG(), libMesh::sys, libMesh::NonlinearSolver< T >::system(), libMesh::System::update(), libMesh::X_global(), and libMesh::X_sys.

{
  START_LOG("jacobian()", "TrilinosNoxNonlinearSolver");

  NonlinearImplicitSystem &sys = _solver->system();

  EpetraMatrix<Number> Jac(&dynamic_cast<Epetra_FECrsMatrix &>(jac), 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();

  STOP_LOG("jacobian()", "TrilinosNoxNonlinearSolver");

  return true;
}
bool libMesh::Problem_Interface::computePrecMatrix ( const Epetra_Vector &  x,
Epetra_RowMatrix &  M 
)

Compute the Epetra_RowMatrix M, that will be used by the Aztec preconditioner instead of the Jacobian. This is used when there is no explicit Jacobian present (i.e. Matrix-Free Newton-Krylov). This MUST BE an Epetra_RowMatrix since the Aztec preconditioners need to know the sparsity pattern of the matrix. Returns true if computation was successful.

Definition at line 182 of file trilinos_nox_nonlinear_solver.C.

{
  //   libMesh::out << "ERROR: Problem_Interface::preconditionVector() - Use Explicit Jacobian only for this test problem!" << endl;
  throw 1;
}
bool libMesh::Problem_Interface::computePreconditioner ( const Epetra_Vector &  x,
Epetra_Operator &  Prec,
Teuchos::ParameterList *  p 
)

Computes a user supplied preconditioner based on input vector x. Returns true if computation was successful.

Definition at line 188 of file trilinos_nox_nonlinear_solver.C.

References _solver, libMesh::ParallelObject::comm(), libMesh::TrilinosPreconditioner< T >::compute(), libMesh::System::current_local_solution, libMesh::DofMap::enforce_constraints_exactly(), libMesh::System::get_dof_map(), libMesh::Jac(), libMesh::NonlinearImplicitSystem::ComputeJacobian::jacobian(), libMesh::NonlinearSolver< T >::jacobian, libMesh::NonlinearSolver< T >::jacobian_object, libMesh::TrilinosPreconditioner< T >::mat(), libMesh::NonlinearSolver< T >::matvec, libMesh::NonlinearImplicitSystem::ComputeResidualandJacobian::residual_and_jacobian(), libMesh::NonlinearSolver< T >::residual_and_jacobian_object, libMesh::System::solution, libMesh::START_LOG(), libMesh::sys, libMesh::NonlinearSolver< T >::system(), libMesh::System::update(), libMesh::X_global(), and libMesh::X_sys.

{
  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;
}

Member Data Documentation


The documentation for this class was generated from the following file: