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

libMesh::Problem_Interface::Problem_Interface ( NoxNonlinearSolver< Number > *  solver  )  [explicit]

Definition at line 81 of file trilinos_nox_nonlinear_solver.C.

00081                                                                         :
00082   _solver(solver)
00083 { }

libMesh::Problem_Interface::~Problem_Interface (  ) 

Definition at line 85 of file trilinos_nox_nonlinear_solver.C.

00086 { }


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::EpetraVector< T >::close(), libMesh::System::current_local_solution, libMesh::DofMap::enforce_constraints_exactly(), libMesh::err, libMesh::AutoPtr< Tp >::get(), 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::EpetraVector< T >::swap(), libMesh::NonlinearSolver< T >::system(), libMesh::System::update(), and libMesh::EpetraVector< T >::zero().

00090 {
00091   START_LOG("residual()", "TrilinosNoxNonlinearSolver");
00092 
00093   NonlinearImplicitSystem &sys = _solver->system();
00094 
00095   EpetraVector<Number> X_global(*const_cast<Epetra_Vector *>(&x)), R(r);
00096   EpetraVector<Number>& X_sys = *libmesh_cast_ptr<EpetraVector<Number>*>(sys.solution.get());
00097   EpetraVector<Number>& R_sys = *libmesh_cast_ptr<EpetraVector<Number>*>(sys.rhs);
00098 
00099   // Use the systems update() to get a good local version of the parallel solution
00100   X_global.swap(X_sys);
00101   R.swap(R_sys);
00102 
00103   sys.get_dof_map().enforce_constraints_exactly(sys);
00104   sys.update();
00105 
00106   // Swap back
00107   X_global.swap(X_sys);
00108   R.swap(R_sys);
00109 
00110   R.zero();
00111 
00112   //-----------------------------------------------------------------------------
00113   // if the user has provided both function pointers and objects only the pointer
00114   // will be used, so catch that as an error
00115 
00116   if (_solver->residual && _solver->residual_object)
00117   {
00118     libMesh::err << "ERROR: cannot specifiy both a function and object to compute the Residual!" << std::endl;
00119     libmesh_error();
00120   }
00121 
00122   if (_solver->matvec && _solver->residual_and_jacobian_object)
00123   {
00124     libMesh::err << "ERROR: cannot specifiy both a function and object to compute the combined Residual & Jacobian!" << std::endl;
00125     libmesh_error();
00126   }
00127   //-----------------------------------------------------------------------------
00128 
00129   if      (_solver->residual != NULL)                     _solver->residual                                            (*sys.current_local_solution.get(), R, sys);
00130   else if (_solver->residual_object != NULL)              _solver->residual_object->residual                           (*sys.current_local_solution.get(), R, sys);
00131   else if (_solver->matvec   != NULL)                     _solver->matvec                                              (*sys.current_local_solution.get(), &R, NULL, sys);
00132   else if (_solver->residual_and_jacobian_object != NULL) _solver->residual_and_jacobian_object->residual_and_jacobian (*sys.current_local_solution.get(), &R, NULL, sys);
00133   else return false;
00134 
00135   R.close();
00136   X_global.close();
00137 
00138   STOP_LOG("residual()", "TrilinosNoxNonlinearSolver");
00139 
00140   return true;
00141 }

bool libMesh::Problem_Interface::computeJacobian ( const Epetra_Vector &  x,
Epetra_Operator &  Jac 
)

Compute an explicit Jacobian.

Definition at line 143 of file trilinos_nox_nonlinear_solver.C.

References _solver, libMesh::SparseMatrix< T >::attach_dof_map(), libMesh::System::current_local_solution, libMesh::DofMap::enforce_constraints_exactly(), libMesh::err, libMesh::AutoPtr< Tp >::get(), libMesh::System::get_dof_map(), 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::NonlinearSolver< T >::system(), and libMesh::System::update().

00145 {
00146   START_LOG("jacobian()", "TrilinosNoxNonlinearSolver");
00147 
00148   NonlinearImplicitSystem &sys = _solver->system();
00149 
00150   EpetraMatrix<Number> Jac(&dynamic_cast<Epetra_FECrsMatrix &>(jac));
00151   EpetraVector<Number>& X_sys = *libmesh_cast_ptr<EpetraVector<Number>*>(sys.solution.get());
00152   EpetraVector<Number> X_global(*const_cast<Epetra_Vector *>(&x));
00153 
00154   // Set the dof maps
00155   Jac.attach_dof_map(sys.get_dof_map());
00156 
00157   // Use the systems update() to get a good local version of the parallel solution
00158   X_global.swap(X_sys);
00159 
00160   sys.get_dof_map().enforce_constraints_exactly(sys);
00161   sys.update();
00162 
00163   X_global.swap(X_sys);
00164 
00165   //-----------------------------------------------------------------------------
00166   // if the user has provided both function pointers and objects only the pointer
00167   // will be used, so catch that as an error
00168   if (_solver->jacobian && _solver->jacobian_object)
00169   {
00170     libMesh::err << "ERROR: cannot specify both a function and object to compute the Jacobian!" << std::endl;
00171     libmesh_error();
00172   }
00173 
00174   if (_solver->matvec && _solver->residual_and_jacobian_object)
00175   {
00176     libMesh::err << "ERROR: cannot specify both a function and object to compute the combined Residual & Jacobian!" << std::endl;
00177     libmesh_error();
00178   }
00179   //-----------------------------------------------------------------------------
00180 
00181   if      (_solver->jacobian != NULL)                     _solver->jacobian                                            (*sys.current_local_solution.get(), Jac, sys);
00182   else if (_solver->jacobian_object != NULL)              _solver->jacobian_object->jacobian                           (*sys.current_local_solution.get(), Jac, sys);
00183   else if (_solver->matvec   != NULL)                     _solver->matvec                                              (*sys.current_local_solution.get(), NULL, &Jac, sys);
00184   else if (_solver->residual_and_jacobian_object != NULL) _solver->residual_and_jacobian_object->residual_and_jacobian (*sys.current_local_solution.get(), NULL, &Jac, sys);
00185   else libmesh_error();
00186 
00187   Jac.close();
00188   X_global.close();
00189 
00190   STOP_LOG("jacobian()", "TrilinosNoxNonlinearSolver");
00191 
00192   return true;
00193 }

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 195 of file trilinos_nox_nonlinear_solver.C.

00196 {
00197 //   cout << "ERROR: Problem_Interface::preconditionVector() - Use Explicit Jacobian only for this test problem!" << endl;
00198    throw 1;
00199 }

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 201 of file trilinos_nox_nonlinear_solver.C.

References _solver, libMesh::TrilinosPreconditioner< T >::compute(), libMesh::System::current_local_solution, libMesh::DofMap::enforce_constraints_exactly(), libMesh::err, libMesh::AutoPtr< Tp >::get(), libMesh::System::get_dof_map(), 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::NonlinearSolver< T >::system(), and libMesh::System::update().

00204 {
00205   START_LOG("preconditioner()", "TrilinosNoxNonlinearSolver");
00206 
00207   NonlinearImplicitSystem &sys = _solver->system();
00208   TrilinosPreconditioner<Number> & tpc = dynamic_cast<TrilinosPreconditioner<Number> &>(prec);
00209 
00210   EpetraMatrix<Number> Jac(dynamic_cast<Epetra_FECrsMatrix *>(tpc.mat()));
00211   EpetraVector<Number>& X_sys = *libmesh_cast_ptr<EpetraVector<Number>*>(sys.solution.get());
00212   EpetraVector<Number> X_global(*const_cast<Epetra_Vector *>(&x));
00213 
00214   // Set the dof maps
00215   Jac.attach_dof_map(sys.get_dof_map());
00216 
00217   // Use the systems update() to get a good local version of the parallel solution
00218   X_global.swap(X_sys);
00219 
00220   sys.get_dof_map().enforce_constraints_exactly(sys);
00221   sys.update();
00222 
00223   X_global.swap(X_sys);
00224 
00225   //-----------------------------------------------------------------------------
00226   // if the user has provided both function pointers and objects only the pointer
00227   // will be used, so catch that as an error
00228   if (_solver->jacobian && _solver->jacobian_object)
00229   {
00230     libMesh::err << "ERROR: cannot specify both a function and object to compute the Jacobian!" << std::endl;
00231     libmesh_error();
00232   }
00233 
00234   if (_solver->matvec && _solver->residual_and_jacobian_object)
00235   {
00236     libMesh::err << "ERROR: cannot specify both a function and object to compute the combined Residual & Jacobian!" << std::endl;
00237     libmesh_error();
00238   }
00239   //-----------------------------------------------------------------------------
00240 
00241   if      (_solver->jacobian != NULL)                     _solver->jacobian                                            (*sys.current_local_solution.get(), Jac, sys);
00242   else if (_solver->jacobian_object != NULL)              _solver->jacobian_object->jacobian                           (*sys.current_local_solution.get(), Jac, sys);
00243   else if (_solver->matvec   != NULL)                     _solver->matvec                                              (*sys.current_local_solution.get(), NULL, &Jac, sys);
00244   else if (_solver->residual_and_jacobian_object != NULL) _solver->residual_and_jacobian_object->residual_and_jacobian (*sys.current_local_solution.get(), NULL, &Jac, sys);
00245   else libmesh_error();
00246 
00247   Jac.close();
00248   X_global.close();
00249 
00250   tpc.compute();
00251 
00252   STOP_LOG("preconditioner()", "TrilinosNoxNonlinearSolver");
00253 
00254   return true;
00255 }


Member Data Documentation


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

Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:06 UTC

Hosted By:
SourceForge.net Logo