libMesh::AztecLinearSolver< T > Class Template Reference

#include <trilinos_aztec_linear_solver.h>

Inheritance diagram for libMesh::AztecLinearSolver< T >:

Public Member Functions

 AztecLinearSolver (const libMesh::Parallel::Communicator &comm)
 
 ~AztecLinearSolver ()
 
virtual void clear () override
 
virtual void init (const char *name=nullptr) override
 
virtual std::pair< unsigned int, Realsolve (SparseMatrix< T > &matrix_in, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const double tol, const unsigned int m_its) override
 
virtual std::pair< unsigned int, Realsolve (SparseMatrix< T > &matrix, SparseMatrix< T > &preconditioner, NumericVector< T > &solution, NumericVector< T > &rhs, const double tol, const unsigned int m_its) override
 
virtual std::pair< unsigned int, Realsolve (const ShellMatrix< T > &shell_matrix, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const double tol, const unsigned int m_its) override
 
virtual std::pair< unsigned int, Realsolve (const ShellMatrix< T > &shell_matrix, const SparseMatrix< T > &precond_matrix, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const double tol, const unsigned int m_its) override
 
void get_residual_history (std::vector< double > &hist)
 
Real get_initial_residual ()
 
virtual void print_converged_reason () const override
 
virtual LinearConvergenceReason get_converged_reason () const override
 
bool initialized () const
 
virtual void init_names (const System &)
 
SolverType solver_type () const
 
void set_solver_type (const SolverType st)
 
PreconditionerType preconditioner_type () const
 
void set_preconditioner_type (const PreconditionerType pct)
 
void attach_preconditioner (Preconditioner< T > *preconditioner)
 
virtual void reuse_preconditioner (bool)
 
bool get_same_preconditioner ()
 
virtual void restrict_solve_to (const std::vector< unsigned int > *const dofs, const SubsetSolveMode subset_solve_mode=SUBSET_ZERO)
 
std::pair< unsigned int, Realsolve (SparseMatrix< T > &matrix, SparseMatrix< T > *precond_matrix, NumericVector< T > &, NumericVector< T > &, const double, const unsigned int)
 
std::pair< unsigned int, Realsolve (const ShellMatrix< T > &matrix, const SparseMatrix< T > *precond_matrix, NumericVector< T > &, NumericVector< T > &, const double, const unsigned int)
 
virtual std::pair< unsigned int, Realadjoint_solve (SparseMatrix< T > &, NumericVector< T > &, NumericVector< T > &, const double, const unsigned int)
 
void set_solver_configuration (SolverConfiguration &solver_configuration)
 
const Parallel::Communicatorcomm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Static Public Member Functions

static std::unique_ptr< LinearSolver< T > > build (const libMesh::Parallel::Communicator &comm_in, const SolverPackage solver_package=libMesh::default_solver_package())
 
static std::string get_info ()
 
static void print_info (std::ostream &out=libMesh::out)
 
static unsigned int n_objects ()
 
static void enable_print_counter_info ()
 
static void disable_print_counter_info ()
 

Protected Types

typedef std::map< std::string, std::pair< unsigned int, unsigned int > > Counts
 

Protected Member Functions

void increment_constructor_count (const std::string &name)
 
void increment_destructor_count (const std::string &name)
 

Protected Attributes

SolverType _solver_type
 
PreconditionerType _preconditioner_type
 
bool _is_initialized
 
Preconditioner< T > * _preconditioner
 
bool same_preconditioner
 
SolverConfiguration_solver_configuration
 
const Parallel::Communicator_communicator
 

Static Protected Attributes

static Counts _counts
 
static Threads::atomic< unsigned int > _n_objects
 
static Threads::spin_mutex _mutex
 
static bool _enable_print_counter = true
 

Private Member Functions

void set_solver_type ()
 

Private Attributes

Epetra_LinearProblem * _linear_problem
 
AztecOO * _linear_solver
 

Detailed Description

template<typename T>
class libMesh::AztecLinearSolver< T >

This class provides an interface to AztecOO iterative solvers that is compatible with the libMesh LinearSolver<>

Author
Benjamin Kirk
Date
2002-2008

Definition at line 53 of file trilinos_aztec_linear_solver.h.

Member Typedef Documentation

◆ Counts

typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts
protectedinherited

Data structure to log the information. The log is identified by the class name.

Definition at line 117 of file reference_counter.h.

Constructor & Destructor Documentation

◆ AztecLinearSolver()

template<typename T >
libMesh::AztecLinearSolver< T >::AztecLinearSolver ( const libMesh::Parallel::Communicator comm)

Constructor. Initializes Aztec data structures

Definition at line 41 of file trilinos_aztec_linear_solver.C.

References libMesh::LinearSolver< T >::_preconditioner_type, libMesh::BLOCK_JACOBI_PRECOND, libMesh::ILU_PRECOND, and libMesh::ParallelObject::n_processors().

41  :
42  LinearSolver<T>(comm)
43 {
44  if (this->n_processors() == 1)
46  else
48 }
const Parallel::Communicator & comm() const
processor_id_type n_processors() const
PreconditionerType _preconditioner_type

◆ ~AztecLinearSolver()

template<typename T >
libMesh::AztecLinearSolver< T >::~AztecLinearSolver ( )
inline

Destructor.

Definition at line 177 of file trilinos_aztec_linear_solver.h.

178 {
179  this->clear ();
180 }

Member Function Documentation

◆ adjoint_solve()

template<typename T>
std::pair< unsigned int, Real > libMesh::LinearSolver< T >::adjoint_solve ( SparseMatrix< T > &  mat,
NumericVector< T > &  sol,
NumericVector< T > &  rhs,
const double  tol,
const unsigned int  n_iter 
)
virtualinherited

Function to solve the adjoint system.

Note
This method will compute the preconditioner from the system matrix. This is not a pure virtual function and is defined linear_solver.C

Reimplemented in libMesh::PetscLinearSolver< T >, libMesh::LaspackLinearSolver< T >, and libMesh::EigenSparseLinearSolver< T >.

Definition at line 145 of file linear_solver.C.

Referenced by libMesh::ImplicitSystem::adjoint_solve().

150 {
151  // Log how long the linear solve takes.
152  LOG_SCOPE("adjoint_solve()", "LinearSolver");
153 
154  // Take the discrete adjoint
155  mat.close();
156  mat.get_transpose(mat);
157 
158  // Call the solve function for the relevant linear algebra library and
159  // solve the transpose matrix
160  const std::pair<unsigned int, Real> totalrval = this->solve (mat, sol, rhs, tol, n_iter);
161 
162  // Now transpose back and restore the original matrix
163  // by taking the discrete adjoint
164  mat.get_transpose(mat);
165 
166  return totalrval;
167 }
virtual std::pair< unsigned int, Real > solve(SparseMatrix< T > &, NumericVector< T > &, NumericVector< T > &, const double, const unsigned int)=0

◆ attach_preconditioner()

template<typename T>
void libMesh::LinearSolver< T >::attach_preconditioner ( Preconditioner< T > *  preconditioner)
inherited

Attaches a Preconditioner object to be used

Definition at line 118 of file linear_solver.C.

119 {
120  if (this->_is_initialized)
121  libmesh_error_msg("Preconditioner must be attached before the solver is initialized!");
122 
124  _preconditioner = preconditioner;
125 }
Preconditioner< T > * _preconditioner
PreconditionerType _preconditioner_type

◆ build()

template<typename T >
std::unique_ptr< LinearSolver< T > > libMesh::LinearSolver< T >::build ( const libMesh::Parallel::Communicator comm_in,
const SolverPackage  solver_package = libMesh::default_solver_package() 
)
staticinherited

Builds a LinearSolver using the linear solver package specified by solver_package

Definition at line 57 of file linear_solver.C.

Referenced by libMesh::ImplicitSystem::get_linear_solver(), and libMesh::TimeSolver::init().

59 {
60  // Avoid unused parameter warnings when no solver packages are enabled.
62 
63  // Build the appropriate solver
64  switch (solver_package)
65  {
66 #ifdef LIBMESH_HAVE_LASPACK
67  case LASPACK_SOLVERS:
68  return libmesh_make_unique<LaspackLinearSolver<T>>(comm);
69 #endif
70 
71 
72 #ifdef LIBMESH_HAVE_PETSC
73  case PETSC_SOLVERS:
74  return libmesh_make_unique<PetscLinearSolver<T>>(comm);
75 #endif
76 
77 
78 #ifdef LIBMESH_TRILINOS_HAVE_AZTECOO
79  case TRILINOS_SOLVERS:
80  return libmesh_make_unique<AztecLinearSolver<T>>(comm);
81 #endif
82 
83 
84 #ifdef LIBMESH_HAVE_EIGEN
85  case EIGEN_SOLVERS:
86  return libmesh_make_unique<EigenSparseLinearSolver<T>>(comm);
87 #endif
88 
89  default:
90  libmesh_error_msg("ERROR: Unrecognized solver package: " << solver_package);
91  }
92 
93  return std::unique_ptr<LinearSolver<T>>();
94 }
EIGEN_SOLVERS
Definition: libmesh.C:246
TRILINOS_SOLVERS
Definition: libmesh.C:244
const Parallel::Communicator & comm() const
void libmesh_ignore(const Args &...)
LASPACK_SOLVERS
Definition: libmesh.C:248

◆ clear()

template<typename T >
void libMesh::AztecLinearSolver< T >::clear ( )
overridevirtual

Release all memory and clear data structures.

Reimplemented from libMesh::LinearSolver< T >.

Definition at line 53 of file trilinos_aztec_linear_solver.C.

References libMesh::libMeshPrivateData::_is_initialized, libMesh::BLOCK_JACOBI_PRECOND, libMesh::GMRES, libMesh::ILU_PRECOND, and libMesh::initialized().

54 {
55  if (this->initialized())
56  {
57  this->_is_initialized = false;
58 
59  // Mimic PETSc default solver and preconditioner
60  this->_solver_type = GMRES;
61 
62  if (this->n_processors() == 1)
64  else
66  }
67 }
processor_id_type n_processors() const
bool initialized() const
Definition: linear_solver.h:96
PreconditionerType _preconditioner_type

◆ comm()

const Parallel::Communicator& libMesh::ParallelObject::comm ( ) const
inlineinherited
Returns
A reference to the Parallel::Communicator object used by this mesh.

Definition at line 89 of file parallel_object.h.

References libMesh::ParallelObject::_communicator.

Referenced by libMesh::__libmesh_petsc_diff_solver_jacobian(), libMesh::__libmesh_petsc_diff_solver_monitor(), libMesh::__libmesh_petsc_diff_solver_residual(), libMesh::__libmesh_tao_equality_constraints(), libMesh::__libmesh_tao_equality_constraints_jacobian(), libMesh::__libmesh_tao_gradient(), libMesh::__libmesh_tao_hessian(), libMesh::__libmesh_tao_inequality_constraints(), libMesh::__libmesh_tao_inequality_constraints_jacobian(), libMesh::__libmesh_tao_objective(), libMesh::MeshRefinement::_coarsen_elements(), libMesh::ExactSolution::_compute_error(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::BoundaryInfo::_find_id_maps(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult_add(), libMesh::EquationSystems::_read_impl(), libMesh::MeshRefinement::_refine_elements(), libMesh::MeshRefinement::_smooth_flags(), libMesh::PetscDMWrapper::add_dofs_helper(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::ImplicitSystem::add_matrix(), libMesh::System::add_vector(), libMesh::UnstructuredMesh::all_second_order(), libMesh::MeshTools::Modification::all_tri(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assemble_qoi(), libMesh::MeshCommunication::assign_global_indices(), libMesh::DofMap::attach_matrix(), libMesh::MeshTools::Generation::build_extrusion(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::PetscDMWrapper::build_section(), libMesh::PetscDMWrapper::build_sf(), libMesh::MeshBase::cache_elem_dims(), libMesh::System::calculate_norm(), libMesh::DofMap::check_dirichlet_bcid_consistency(), libMesh::PetscDMWrapper::check_section_n_dofs(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Problem_Interface::computeF(), libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), libMesh::ExodusII_IO::copy_elemental_solution(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::MeshTools::create_bounding_box(), libMesh::MeshTools::create_nodal_bounding_box(), libMesh::MeshRefinement::create_parent_error_vector(), libMesh::MeshTools::create_processor_bounding_box(), libMesh::MeshTools::create_subdomain_bounding_box(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), DMlibMeshFunction(), DMlibMeshJacobian(), DMlibMeshSetSystem_libMesh(), DMVariableBounds_libMesh(), libMesh::MeshRefinement::eliminate_unrefined_patches(), libMesh::EpetraVector< T >::EpetraVector(), libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::ExactErrorEstimator::estimate_error(), libMesh::MeshRefinement::flag_elements_by_elem_fraction(), libMesh::MeshRefinement::flag_elements_by_error_fraction(), libMesh::MeshRefinement::flag_elements_by_nelem_target(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::DofMap::get_info(), libMesh::ImplicitSystem::get_linear_solver(), libMesh::LocationMap< T >::init(), libMesh::TimeSolver::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_and_attach_petscdm(), libMesh::EigenSystem::init_data(), libMesh::EigenSystem::init_matrices(), libMesh::OptimizationSystem::initialize_equality_constraints_storage(), libMesh::OptimizationSystem::initialize_inequality_constraints_storage(), libMesh::MeshTools::libmesh_assert_consistent_distributed(), libMesh::MeshTools::libmesh_assert_consistent_distributed_nodes(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_new_node_procids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_flags(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_p_levels(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::MeshTools::libmesh_assert_valid_unique_ids(), libMesh::libmesh_petsc_snes_fd_residual(), libMesh::libmesh_petsc_snes_jacobian(), libMesh::libmesh_petsc_snes_mffd_residual(), libMesh::libmesh_petsc_snes_postcheck(), libMesh::libmesh_petsc_snes_residual(), libMesh::libmesh_petsc_snes_residual_helper(), libMesh::MeshRefinement::limit_level_mismatch_at_edge(), libMesh::MeshRefinement::limit_level_mismatch_at_node(), libMesh::MeshRefinement::limit_overrefined_boundary(), libMesh::MeshRefinement::limit_underrefined_boundary(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshCommunication::make_elems_parallel_consistent(), libMesh::MeshRefinement::make_flags_parallel_consistent(), libMesh::MeshCommunication::make_new_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_new_nodes_parallel_consistent(), libMesh::MeshCommunication::make_node_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_unique_ids_parallel_consistent(), libMesh::MeshCommunication::make_nodes_parallel_consistent(), libMesh::MeshCommunication::make_p_levels_parallel_consistent(), libMesh::MeshRefinement::make_refinement_compatible(), libMesh::FEMSystem::mesh_position_set(), libMesh::DistributedMesh::n_active_elem(), libMesh::MeshTools::n_active_levels(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::CondensedEigenSystem::n_global_non_condensed_dofs(), libMesh::MeshTools::n_levels(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::MeshTools::n_p_levels(), libMesh::BoundaryInfo::n_shellface_conds(), libMesh::DistributedMesh::parallel_max_elem_id(), libMesh::DistributedMesh::parallel_max_node_id(), libMesh::ReplicatedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_n_elem(), libMesh::DistributedMesh::parallel_n_nodes(), libMesh::SparsityPattern::Build::parallel_sync(), libMesh::MeshTools::paranoid_n_levels(), libMesh::petsc_auto_fieldsplit(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshBase::prepare_for_use(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::MeshBase::recalculate_n_partitions(), libMesh::MeshRefinement::refine_and_coarsen_elements(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::CheckpointIO::select_split_config(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::PetscDiffSolver::setup_petsc_data(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::split_mesh(), libMesh::MeshBase::subdomain_ids(), libMesh::BoundaryInfo::sync(), libMesh::MeshRefinement::test_level_one(), libMesh::MeshRefinement::test_unflagged(), libMesh::MeshTools::total_weight(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().

90  { return _communicator; }
const Parallel::Communicator & _communicator

◆ disable_print_counter_info()

void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 106 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

Referenced by libMesh::LibMeshInit::LibMeshInit().

107 {
108  _enable_print_counter = false;
109  return;
110 }

◆ enable_print_counter_info()

void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

Methods to enable/disable the reference counter output from print_info()

Definition at line 100 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

101 {
102  _enable_print_counter = true;
103  return;
104 }

◆ get_converged_reason()

template<typename T >
LinearConvergenceReason libMesh::AztecLinearSolver< T >::get_converged_reason ( ) const
overridevirtual
Returns
The solver's convergence flag

Implements libMesh::LinearSolver< T >.

Definition at line 232 of file trilinos_aztec_linear_solver.C.

References libMesh::CONVERGED_RTOL_NORMAL, libMesh::DIVERGED_ITS, and libMesh::DIVERGED_NULL.

233 {
234  const double *status = _linear_solver->GetAztecStatus();
235 
236  switch (static_cast<int>(status[AZ_why]))
237  {
238  case AZ_normal :
239  return CONVERGED_RTOL_NORMAL;
240  case AZ_maxits :
241  return DIVERGED_ITS;
242  default :
243  break;
244  }
245  return DIVERGED_NULL;
246 }
MPI_Status status
Definition: status.h:41

◆ get_info()

std::string libMesh::ReferenceCounter::get_info ( )
staticinherited

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().

Referenced by libMesh::ReferenceCounter::print_info().

48 {
49 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
50 
51  std::ostringstream oss;
52 
53  oss << '\n'
54  << " ---------------------------------------------------------------------------- \n"
55  << "| Reference count information |\n"
56  << " ---------------------------------------------------------------------------- \n";
57 
58  for (const auto & pr : _counts)
59  {
60  const std::string name(pr.first);
61  const unsigned int creations = pr.second.first;
62  const unsigned int destructions = pr.second.second;
63 
64  oss << "| " << name << " reference count information:\n"
65  << "| Creations: " << creations << '\n'
66  << "| Destructions: " << destructions << '\n';
67  }
68 
69  oss << " ---------------------------------------------------------------------------- \n";
70 
71  return oss.str();
72 
73 #else
74 
75  return "";
76 
77 #endif
78 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42

◆ get_initial_residual()

template<typename T >
Real libMesh::AztecLinearSolver< T >::get_initial_residual ( )
Returns
The initial residual for the solve just completed with this interface.

Use this method instead of the one above if you just want the starting residual and not the entire history.

Definition at line 191 of file trilinos_aztec_linear_solver.C.

192 {
193  return _linear_solver->TrueResidual();
194 }

◆ get_residual_history()

template<typename T >
void libMesh::AztecLinearSolver< T >::get_residual_history ( std::vector< double > &  hist)

Fills the input vector with the sequence of residual norms from the latest iterative solve.

Definition at line 182 of file trilinos_aztec_linear_solver.C.

183 {
184  libmesh_not_implemented();
185 }

◆ get_same_preconditioner()

template<typename T >
bool libMesh::LinearSolver< T >::get_same_preconditioner ( )
inlineinherited
Returns
same_preconditioner, which indicates if we reuse the same preconditioner for subsequent solves.

Definition at line 322 of file linear_solver.h.

323 {
324  return same_preconditioner;
325 }

◆ increment_constructor_count()

void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name)
inlineprotectedinherited

Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.

Definition at line 181 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().

182 {
183  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
184  std::pair<unsigned int, unsigned int> & p = _counts[name];
185 
186  p.first++;
187 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
spin_mutex spin_mtx
Definition: threads.C:29

◆ increment_destructor_count()

void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name)
inlineprotectedinherited

Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.

Definition at line 194 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().

195 {
196  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
197  std::pair<unsigned int, unsigned int> & p = _counts[name];
198 
199  p.second++;
200 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
spin_mutex spin_mtx
Definition: threads.C:29

◆ init()

template<typename T >
void libMesh::AztecLinearSolver< T >::init ( const char *  name = nullptr)
overridevirtual

Initialize data structures if not done so already.

Implements libMesh::LinearSolver< T >.

Definition at line 72 of file trilinos_aztec_linear_solver.C.

References libMesh::libMeshPrivateData::_is_initialized, libMesh::BLOCK_JACOBI_PRECOND, libMesh::ICC_PRECOND, libMesh::ILU_PRECOND, libMesh::initialized(), and libMesh::LU_PRECOND.

73 {
74  // Initialize the data structures if not done so already.
75  if (!this->initialized())
76  {
77  this->_is_initialized = true;
78 
79  _linear_solver = new AztecOO();
80 
82 
83  switch(this->_preconditioner_type)
84  {
85  case ILU_PRECOND:
86  _linear_solver->SetAztecOption(AZ_precond,AZ_dom_decomp);
87  _linear_solver->SetAztecOption(AZ_subdomain_solve,AZ_ilu);
88  break;
89 
91  _linear_solver->SetAztecOption(AZ_precond,AZ_Jacobi);
92  break;
93 
94  case ICC_PRECOND:
95  _linear_solver->SetAztecOption(AZ_precond,AZ_dom_decomp);
96  _linear_solver->SetAztecOption(AZ_subdomain_solve,AZ_icc);
97  break;
98 
99  case LU_PRECOND:
100  _linear_solver->SetAztecOption(AZ_precond,AZ_dom_decomp);
101  _linear_solver->SetAztecOption(AZ_subdomain_solve,AZ_lu);
102  break;
103 
104  default:
105  _linear_solver->SetAztecOption(AZ_precond,AZ_dom_decomp);
106  _linear_solver->SetAztecOption(AZ_subdomain_solve,AZ_ilu);
107  }
108  }
109 }
bool initialized() const
Definition: linear_solver.h:96
PreconditionerType _preconditioner_type

◆ init_names()

template<typename T>
virtual void libMesh::LinearSolver< T >::init_names ( const System )
inlinevirtualinherited

Apply names to the system to be solved. For most packages this is a no-op; for PETSc this sets an option prefix from the system name and sets field names from the system's variable names.

Since field names are applied to DoF numberings, this method must be called again after any System reinit.

Reimplemented in libMesh::PetscLinearSolver< T >.

Definition at line 117 of file linear_solver.h.

117 {}

◆ initialized()

template<typename T>
bool libMesh::LinearSolver< T >::initialized ( ) const
inlineinherited
Returns
true if the data structures are initialized, false otherwise.

Definition at line 96 of file linear_solver.h.

96 { return _is_initialized; }

◆ n_objects()

static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited

Prints the number of outstanding (created, but not yet destroyed) objects.

Definition at line 83 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

84  { return _n_objects; }
static Threads::atomic< unsigned int > _n_objects

◆ n_processors()

processor_id_type libMesh::ParallelObject::n_processors ( ) const
inlineinherited
Returns
The number of processors in the group.

Definition at line 95 of file parallel_object.h.

References libMesh::ParallelObject::_communicator, and libMesh::Parallel::Communicator::size().

Referenced by libMesh::BoundaryInfo::_find_id_maps(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::DistributedMesh::add_elem(), libMesh::DistributedMesh::add_node(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assembly(), libMesh::AztecLinearSolver< T >::AztecLinearSolver(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::DistributedMesh::clear(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::UnstructuredMesh::create_pid_mesh(), libMesh::MeshTools::create_processor_bounding_box(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::EnsightIO::EnsightIO(), libMesh::MeshBase::get_info(), libMesh::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_and_attach_petscdm(), libMesh::Nemesis_IO_Helper::initialize(), libMesh::DistributedMesh::insert_elem(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_new_node_procids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_elem_on_proc(), libMesh::MeshBase::n_elem_on_proc(), libMesh::MeshBase::n_nodes_on_proc(), libMesh::MeshBase::partition(), libMesh::PetscLinearSolver< T >::PetscLinearSolver(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::CheckpointIO::read_connectivity(), libMesh::XdrIO::read_header(), libMesh::CheckpointIO::read_nodes(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::System::read_serialized_vector(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().

96  { return cast_int<processor_id_type>(_communicator.size()); }
processor_id_type size() const
Definition: communicator.h:175
const Parallel::Communicator & _communicator

◆ preconditioner_type()

template<typename T >
PreconditionerType libMesh::LinearSolver< T >::preconditioner_type ( ) const
inherited
Returns
The type of preconditioner to use.

Definition at line 98 of file linear_solver.C.

99 {
100  if (_preconditioner)
101  return _preconditioner->type();
102 
103  return _preconditioner_type;
104 }
Preconditioner< T > * _preconditioner
PreconditionerType _preconditioner_type

◆ print_converged_reason()

template<typename T >
void libMesh::AztecLinearSolver< T >::print_converged_reason ( ) const
overridevirtual

Prints a useful message about why the latest linear solve con(di)verged.

Reimplemented from libMesh::LinearSolver< T >.

Definition at line 199 of file trilinos_aztec_linear_solver.C.

References libMesh::out.

200 {
201  const double *status = _linear_solver->GetAztecStatus();
202 
203  switch (static_cast<int>(status[AZ_why]))
204  {
205  case AZ_normal :
206  libMesh::out << "AztecOO converged.\n";
207  break;
208  case AZ_maxits :
209  libMesh::out << "AztecOO failed to converge within maximum iterations.\n";
210  break;
211  case AZ_param :
212  libMesh::out << "AztecOO failed to support a user-requested parameter.\n";
213  break;
214  case AZ_breakdown :
215  libMesh::out << "AztecOO encountered numerical breakdown.\n";
216  break;
217  case AZ_loss :
218  libMesh::out << "AztecOO encountered numerical loss of precision.\n";
219  break;
220  case AZ_ill_cond :
221  libMesh::out << "AztecOO encountered an ill-conditioned GMRES Hessian.\n";
222  break;
223  default:
224  libMesh::out << "AztecOO reported an unrecognized condition.\n";
225  break;
226  }
227 }
MPI_Status status
Definition: status.h:41
OStreamProxy out(std::cout)

◆ print_info()

void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out)
staticinherited

Prints the reference information, by default to libMesh::out.

Definition at line 87 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().

88 {
90  out_stream << ReferenceCounter::get_info();
91 }
static std::string get_info()

◆ processor_id()

processor_id_type libMesh::ParallelObject::processor_id ( ) const
inlineinherited
Returns
The rank of this processor in the group.

Definition at line 101 of file parallel_object.h.

References libMesh::ParallelObject::_communicator, and libMesh::Parallel::Communicator::rank().

Referenced by libMesh::BoundaryInfo::_find_id_maps(), libMesh::EquationSystems::_read_impl(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::DistributedMesh::add_elem(), libMesh::BoundaryInfo::add_elements(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::DistributedMesh::add_node(), libMesh::UnstructuredMesh::all_second_order(), libMesh::MeshTools::Modification::all_tri(), libMesh::FEMSystem::assembly(), libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::InfElemBuilder::build_inf_elem(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::DistributedMesh::clear(), libMesh::ExodusII_IO_Helper::close(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::compute_communication_map_parameters(), libMesh::Nemesis_IO_Helper::compute_internal_and_border_elems_and_internal_nodes(), libMesh::Nemesis_IO_Helper::compute_node_communication_maps(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::ExodusII_IO_Helper::create(), libMesh::DistributedMesh::delete_elem(), libMesh::DistributedMesh::delete_node(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::DistributedMesh::DistributedMesh(), libMesh::DofMap::end_dof(), libMesh::DofMap::end_old_dof(), libMesh::EnsightIO::EnsightIO(), libMesh::MeshFunction::find_element(), libMesh::MeshFunction::find_elements(), libMesh::UnstructuredMesh::find_neighbors(), libMesh::DofMap::first_dof(), libMesh::DofMap::first_old_dof(), libMesh::Nemesis_IO_Helper::get_cmap_params(), libMesh::Nemesis_IO_Helper::get_eb_info_global(), libMesh::Nemesis_IO_Helper::get_elem_cmap(), libMesh::Nemesis_IO_Helper::get_elem_map(), libMesh::MeshBase::get_info(), libMesh::DofMap::get_info(), libMesh::Nemesis_IO_Helper::get_init_global(), libMesh::Nemesis_IO_Helper::get_init_info(), libMesh::Nemesis_IO_Helper::get_loadbal_param(), libMesh::Nemesis_IO_Helper::get_node_cmap(), libMesh::Nemesis_IO_Helper::get_node_map(), libMesh::Nemesis_IO_Helper::get_ns_param_global(), libMesh::Nemesis_IO_Helper::get_ss_param_global(), libMesh::SparsityPattern::Build::handle_vi_vj(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::ExodusII_IO_Helper::initialize_element_variables(), libMesh::ExodusII_IO_Helper::initialize_global_variables(), libMesh::ExodusII_IO_Helper::initialize_nodal_variables(), libMesh::DistributedMesh::insert_elem(), libMesh::DofMap::is_evaluable(), libMesh::SparsityPattern::Build::join(), libMesh::DofMap::last_dof(), libMesh::MeshTools::libmesh_assert_consistent_distributed(), libMesh::MeshTools::libmesh_assert_consistent_distributed_nodes(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_local_elem(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::DofMap::n_local_dofs(), libMesh::System::n_local_dofs(), libMesh::MeshBase::n_local_elem(), libMesh::MeshBase::n_local_nodes(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::BoundaryInfo::n_shellface_conds(), libMesh::SparsityPattern::Build::operator()(), libMesh::DistributedMesh::own_node(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::Nemesis_IO_Helper::put_cmap_params(), libMesh::Nemesis_IO_Helper::put_elem_cmap(), libMesh::Nemesis_IO_Helper::put_elem_map(), libMesh::Nemesis_IO_Helper::put_loadbal_param(), libMesh::Nemesis_IO_Helper::put_node_cmap(), libMesh::Nemesis_IO_Helper::put_node_map(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read(), libMesh::ExodusII_IO_Helper::read_elem_num_map(), libMesh::ExodusII_IO_Helper::read_global_values(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::ExodusII_IO_Helper::read_node_num_map(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::System::read_serialized_data(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::System::read_serialized_vectors(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::CheckpointIO::select_split_config(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::MeshTools::total_weight(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::Parallel::Packing< T >::unpack(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::CheckpointIO::write(), libMesh::EquationSystems::write(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::ExodusII_IO::write_element_data(), libMesh::ExodusII_IO_Helper::write_element_values(), libMesh::ExodusII_IO_Helper::write_elements(), libMesh::ExodusII_IO::write_global_data(), libMesh::ExodusII_IO_Helper::write_global_values(), libMesh::System::write_header(), libMesh::ExodusII_IO::write_information_records(), libMesh::ExodusII_IO_Helper::write_information_records(), libMesh::ExodusII_IO_Helper::write_nodal_coordinates(), libMesh::UCDIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), libMesh::ExodusII_IO_Helper::write_nodal_values(), libMesh::Nemesis_IO_Helper::write_nodesets(), libMesh::ExodusII_IO_Helper::write_nodesets(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bc_names(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::System::write_serialized_data(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), libMesh::System::write_serialized_vector(), libMesh::System::write_serialized_vectors(), libMesh::Nemesis_IO_Helper::write_sidesets(), libMesh::ExodusII_IO_Helper::write_sidesets(), libMesh::ExodusII_IO::write_timestep(), libMesh::ExodusII_IO_Helper::write_timestep(), and libMesh::ExodusII_IO::write_timestep_discontinuous().

102  { return cast_int<processor_id_type>(_communicator.rank()); }
const Parallel::Communicator & _communicator
processor_id_type rank() const
Definition: communicator.h:173

◆ restrict_solve_to()

template<typename T >
void libMesh::LinearSolver< T >::restrict_solve_to ( const std::vector< unsigned int > *const  dofs,
const SubsetSolveMode  subset_solve_mode = SUBSET_ZERO 
)
virtualinherited

After calling this method, all successive solves will be restricted to the given set of dofs, which must contain local dofs on each processor only and not contain any duplicates. This mode can be disabled by calling this method with dofs being a nullptr.

Reimplemented in libMesh::PetscLinearSolver< T >.

Definition at line 136 of file linear_solver.C.

138 {
139  if (dofs != nullptr)
140  libmesh_not_implemented();
141 }

◆ reuse_preconditioner()

template<typename T >
void libMesh::LinearSolver< T >::reuse_preconditioner ( bool  reuse_flag)
virtualinherited

Set the same_preconditioner flag, which indicates if we reuse the same preconditioner for subsequent solves.

Definition at line 129 of file linear_solver.C.

Referenced by libMesh::ImplicitSystem::disable_cache().

130 {
131  same_preconditioner = reuse_flag;
132 }

◆ set_preconditioner_type()

template<typename T >
void libMesh::LinearSolver< T >::set_preconditioner_type ( const PreconditionerType  pct)
inherited

Sets the type of preconditioner to use.

Definition at line 108 of file linear_solver.C.

109 {
110  if (_preconditioner)
111  _preconditioner->set_type(pct);
112  else
113  _preconditioner_type = pct;
114 }
Preconditioner< T > * _preconditioner
PreconditionerType _preconditioner_type

◆ set_solver_configuration()

template<typename T >
void libMesh::LinearSolver< T >::set_solver_configuration ( SolverConfiguration solver_configuration)
inherited

Set the solver configuration object.

Definition at line 177 of file linear_solver.C.

178 {
179  _solver_configuration = &solver_configuration;
180 }
SolverConfiguration * _solver_configuration

◆ set_solver_type() [1/2]

template<typename T>
void libMesh::LinearSolver< T >::set_solver_type ( const SolverType  st)
inlineinherited

Sets the type of solver to use.

Definition at line 127 of file linear_solver.h.

128  { _solver_type = st; }

◆ set_solver_type() [2/2]

template<typename T >
void libMesh::AztecLinearSolver< T >::set_solver_type ( )
private

Tells AztecOO to use the user-specified solver stored in _solver_type

Definition at line 251 of file trilinos_aztec_linear_solver.C.

References libMesh::BICGSTAB, libMesh::CG, libMesh::CGS, libMesh::Utility::enum_to_string(), libMesh::err, libMesh::GMRES, and libMesh::TFQMR.

252 {
253  switch (this->_solver_type)
254  {
255  case CG:
256  _linear_solver->SetAztecOption(AZ_solver, AZ_cg); return;
257 
258  case CGS:
259  _linear_solver->SetAztecOption(AZ_solver, AZ_cgs); return;
260 
261  case TFQMR:
262  _linear_solver->SetAztecOption(AZ_solver, AZ_tfqmr); return;
263 
264  case BICGSTAB:
265  _linear_solver->SetAztecOption(AZ_solver, AZ_bicgstab); return;
266 
267  case GMRES:
268  _linear_solver->SetAztecOption(AZ_solver, AZ_gmres); return;
269 
270  default:
271  libMesh::err << "ERROR: Unsupported AztecOO Solver: "
272  << Utility::enum_to_string(this->_solver_type) << std::endl
273  << "Continuing with AztecOO defaults" << std::endl;
274  }
275 }
OStreamProxy err(std::cerr)
std::string enum_to_string(const T e)

◆ solve() [1/6]

template<typename T >
virtual std::pair<unsigned int, Real> libMesh::AztecLinearSolver< T >::solve ( SparseMatrix< T > &  matrix_in,
NumericVector< T > &  solution_in,
NumericVector< T > &  rhs_in,
const double  tol,
const unsigned int  m_its 
)
inlineoverridevirtual

Call the Aztec solver. It calls the method below, using the same matrix for the system and preconditioner matrices.

Implements libMesh::LinearSolver< T >.

Definition at line 81 of file trilinos_aztec_linear_solver.h.

86  {
87  return this->solve(matrix_in, matrix_in, solution_in, rhs_in, tol, m_its);
88  }
virtual std::pair< unsigned int, Real > solve(SparseMatrix< T > &matrix_in, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const double tol, const unsigned int m_its) override

◆ solve() [2/6]

template<typename T >
std::pair< unsigned int, Real > libMesh::AztecLinearSolver< T >::solve ( SparseMatrix< T > &  matrix,
SparseMatrix< T > &  preconditioner,
NumericVector< T > &  solution,
NumericVector< T > &  rhs,
const double  tol,
const unsigned int  m_its 
)
overridevirtual

This method allows you to call a linear solver while specifying the matrix to use as the (left) preconditioning matrix.

Note
The linear solver will not compute a preconditioner in this case, and will instead premultiply by the matrix you provide.

Implements libMesh::LinearSolver< T >.

Definition at line 116 of file trilinos_aztec_linear_solver.C.

References libMesh::EpetraMatrix< T >::close(), libMesh::TriangleWrapper::init(), and libMesh::EpetraMatrix< T >::mat().

122 {
123  LOG_SCOPE("solve()", "AztecLinearSolver");
124 
125  // Make sure the data passed in are really of Epetra types
126  EpetraMatrix<T> * matrix = cast_ptr<EpetraMatrix<T> *>(&matrix_in);
127  EpetraMatrix<T> * precond = cast_ptr<EpetraMatrix<T> *>(&precond_in);
128  EpetraVector<T> * solution = cast_ptr<EpetraVector<T> *>(&solution_in);
129  EpetraVector<T> * rhs = cast_ptr<EpetraVector<T> *>(&rhs_in);
130 
131  this->init();
132 
133  // Close the matrices and vectors in case this wasn't already done.
134  matrix->close ();
135  precond->close ();
136  solution->close ();
137  rhs->close ();
138 
139  _linear_solver->SetAztecOption(AZ_max_iter,m_its);
140  _linear_solver->SetAztecParam(AZ_tol,tol);
141 
142  Epetra_FECrsMatrix * emat = matrix->mat();
143  Epetra_Vector * esol = solution->vec();
144  Epetra_Vector * erhs = rhs->vec();
145 
146  _linear_solver->Iterate(emat, esol, erhs, m_its, tol);
147 
148  // return the # of its. and the final residual norm.
149  return std::make_pair(_linear_solver->NumIters(), _linear_solver->TrueResidual());
150 }
virtual void init(const char *name=nullptr) override

◆ solve() [3/6]

template<typename T >
std::pair< unsigned int, Real > libMesh::AztecLinearSolver< T >::solve ( const ShellMatrix< T > &  shell_matrix,
NumericVector< T > &  solution_in,
NumericVector< T > &  rhs_in,
const double  tol,
const unsigned int  m_its 
)
overridevirtual

This function solves a system whose matrix is a shell matrix.

Implements libMesh::LinearSolver< T >.

Definition at line 156 of file trilinos_aztec_linear_solver.C.

161 {
162  libmesh_not_implemented();
163 }

◆ solve() [4/6]

template<typename T >
std::pair< unsigned int, Real > libMesh::AztecLinearSolver< T >::solve ( const ShellMatrix< T > &  shell_matrix,
const SparseMatrix< T > &  precond_matrix,
NumericVector< T > &  solution_in,
NumericVector< T > &  rhs_in,
const double  tol,
const unsigned int  m_its 
)
overridevirtual

This function solves a system whose matrix is a shell matrix, but a sparse matrix is used as preconditioning matrix, this allowing other preconditioners than JACOBI.

Implements libMesh::LinearSolver< T >.

Definition at line 169 of file trilinos_aztec_linear_solver.C.

175 {
176  libmesh_not_implemented();
177 }

◆ solve() [5/6]

template<typename T>
std::pair< unsigned int, Real > libMesh::LinearSolver< T >::solve ( SparseMatrix< T > &  matrix,
SparseMatrix< T > *  precond_matrix,
NumericVector< T > &  sol,
NumericVector< T > &  rhs,
const double  tol,
const unsigned int  n_iter 
)
inlineinherited

This function calls the solver "_solver_type" preconditioned with the "_preconditioner_type" preconditioner. The preconditioning matrix is used if it is provided, or the system matrix is used if precond_matrix is null

Definition at line 330 of file linear_solver.h.

336 {
337  if (pc_mat)
338  return this->solve(mat, *pc_mat, sol, rhs, tol, n_iter);
339  else
340  return this->solve(mat, sol, rhs, tol, n_iter);
341 }
virtual std::pair< unsigned int, Real > solve(SparseMatrix< T > &, NumericVector< T > &, NumericVector< T > &, const double, const unsigned int)=0

◆ solve() [6/6]

template<typename T>
std::pair< unsigned int, Real > libMesh::LinearSolver< T >::solve ( const ShellMatrix< T > &  matrix,
const SparseMatrix< T > *  precond_matrix,
NumericVector< T > &  sol,
NumericVector< T > &  rhs,
const double  tol,
const unsigned int  n_iter 
)
inlineinherited

This function solves a system whose matrix is a shell matrix, but an optional sparse matrix may be used as preconditioning matrix.

Definition at line 347 of file linear_solver.h.

353 {
354  if (pc_mat)
355  return this->solve(mat, *pc_mat, sol, rhs, tol, n_iter);
356  else
357  return this->solve(mat, sol, rhs, tol, n_iter);
358 }
virtual std::pair< unsigned int, Real > solve(SparseMatrix< T > &, NumericVector< T > &, NumericVector< T > &, const double, const unsigned int)=0

◆ solver_type()

template<typename T>
SolverType libMesh::LinearSolver< T >::solver_type ( ) const
inlineinherited
Returns
The type of solver to use.

Definition at line 122 of file linear_solver.h.

122 { return _solver_type; }

Member Data Documentation

◆ _communicator

◆ _counts

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited

◆ _enable_print_counter

bool libMesh::ReferenceCounter::_enable_print_counter = true
staticprotectedinherited

Flag to control whether reference count information is printed when print_info is called.

Definition at line 141 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().

◆ _is_initialized

template<typename T>
bool libMesh::LinearSolver< T >::_is_initialized
protectedinherited

Flag indicating if the data structures have been initialized.

Definition at line 287 of file linear_solver.h.

Referenced by libMesh::LinearSolver< Number >::initialized().

◆ _linear_problem

template<typename T >
Epetra_LinearProblem* libMesh::AztecLinearSolver< T >::_linear_problem
private

The Epetra linear problem object.

Definition at line 165 of file trilinos_aztec_linear_solver.h.

◆ _linear_solver

template<typename T >
AztecOO* libMesh::AztecLinearSolver< T >::_linear_solver
private

The AztecOO solver object

Definition at line 170 of file trilinos_aztec_linear_solver.h.

◆ _mutex

Threads::spin_mutex libMesh::ReferenceCounter::_mutex
staticprotectedinherited

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 135 of file reference_counter.h.

◆ _n_objects

Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects
staticprotectedinherited

The number of objects. Print the reference count information when the number returns to 0.

Definition at line 130 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().

◆ _preconditioner

template<typename T>
Preconditioner<T>* libMesh::LinearSolver< T >::_preconditioner
protectedinherited

Holds the Preconditioner object to be used for the linear solves.

Definition at line 292 of file linear_solver.h.

◆ _preconditioner_type

template<typename T>
PreconditionerType libMesh::LinearSolver< T >::_preconditioner_type
protectedinherited

Enum stating with type of preconditioner to use.

Definition at line 282 of file linear_solver.h.

Referenced by libMesh::AztecLinearSolver< T >::AztecLinearSolver(), and libMesh::PetscLinearSolver< T >::PetscLinearSolver().

◆ _solver_configuration

template<typename T>
SolverConfiguration* libMesh::LinearSolver< T >::_solver_configuration
protectedinherited

Optionally store a SolverOptions object that can be used to set parameters like solver type, tolerances and iteration limits.

Definition at line 306 of file linear_solver.h.

◆ _solver_type

template<typename T>
SolverType libMesh::LinearSolver< T >::_solver_type
protectedinherited

◆ same_preconditioner

template<typename T>
bool libMesh::LinearSolver< T >::same_preconditioner
protectedinherited

Boolean flag to indicate whether we want to use an identical preconditioner to the previous solve. This can save substantial work in the cases where the system matrix is the same for successive solves.

Definition at line 300 of file linear_solver.h.


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