libMesh::TaoOptimizationSolver< T > Class Template Reference

#include <tao_optimization_solver.h>

Inheritance diagram for libMesh::TaoOptimizationSolver< T >:

Public Types

typedef OptimizationSystem sys_type
 

Public Member Functions

 TaoOptimizationSolver (sys_type &system)
 
 ~TaoOptimizationSolver ()
 
virtual void clear () override
 
virtual void init () override
 
Tao tao ()
 
virtual void solve () override
 
virtual void get_dual_variables () override
 
virtual void print_converged_reason () override
 
virtual int get_converged_reason () override
 
bool initialized () const
 
const sys_typesystem () const
 
sys_typesystem ()
 
const Parallel::Communicatorcomm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Static Public Member Functions

static std::unique_ptr< OptimizationSolver< T > > build (sys_type &s, 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 ()
 

Public Attributes

OptimizationSystem::ComputeObjectiveobjective_object
 
OptimizationSystem::ComputeGradientgradient_object
 
OptimizationSystem::ComputeHessianhessian_object
 
OptimizationSystem::ComputeEqualityConstraintsequality_constraints_object
 
OptimizationSystem::ComputeEqualityConstraintsJacobianequality_constraints_jacobian_object
 
OptimizationSystem::ComputeInequalityConstraintsinequality_constraints_object
 
OptimizationSystem::ComputeInequalityConstraintsJacobianinequality_constraints_jacobian_object
 
OptimizationSystem::ComputeLowerAndUpperBoundslower_and_upper_bounds_object
 
unsigned int max_objective_function_evaluations
 
Real objective_function_relative_tolerance
 
bool verbose
 

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

Tao _tao
 
TaoConvergedReason _reason
 
sys_type_system
 
bool _is_initialized
 
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
 

Friends

PetscErrorCode __libmesh_tao_objective (Tao tao, Vec x, PetscReal *objective, void *ctx)
 
PetscErrorCode __libmesh_tao_gradient (Tao tao, Vec x, Vec g, void *ctx)
 
PetscErrorCode __libmesh_tao_hessian (Tao tao, Vec x, Mat h, Mat pc, void *ctx)
 
PetscErrorCode __libmesh_tao_equality_constraints (Tao tao, Vec x, Vec ce, void *ctx)
 
PetscErrorCode __libmesh_tao_equality_constraints_jacobian (Tao tao, Vec x, Mat J, Mat Jpre, void *ctx)
 
PetscErrorCode __libmesh_tao_inequality_constraints (Tao tao, Vec x, Vec cineq, void *ctx)
 
PetscErrorCode __libmesh_tao_inequality_constraints_jacobian (Tao tao, Vec x, Mat J, Mat Jpre, void *ctx)
 

Detailed Description

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

This class provides an interface to the Tao optimization solvers.

Author
David Knezevic
Date
2015

Definition at line 58 of file tao_optimization_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.

◆ sys_type

template<typename T>
typedef OptimizationSystem libMesh::TaoOptimizationSolver< T >::sys_type

The type of system that we use in conjunction with this solver.

Definition at line 65 of file tao_optimization_solver.h.

Constructor & Destructor Documentation

◆ TaoOptimizationSolver()

template<typename T >
libMesh::TaoOptimizationSolver< T >::TaoOptimizationSolver ( sys_type system)
explicit

Constructor. Initializes Tao data structures.

Definition at line 412 of file tao_optimization_solver.C.

412  :
413  OptimizationSolver<T>(system_in),
414  _reason(TAO_CONVERGED_USER) // Arbitrary initial value...
415 {
416 }

◆ ~TaoOptimizationSolver()

template<typename T >
libMesh::TaoOptimizationSolver< T >::~TaoOptimizationSolver ( )

Destructor.

Definition at line 421 of file tao_optimization_solver.C.

422 {
423  this->clear ();
424 }

Member Function Documentation

◆ build()

template<typename T >
std::unique_ptr< OptimizationSolver< T > > libMesh::OptimizationSolver< T >::build ( sys_type s,
const SolverPackage  solver_package = libMesh::default_solver_package() 
)
staticinherited

Builds an OptimizationSolver using the package specified by solver_package

Definition at line 62 of file optimization_solver.C.

References libMesh::libmesh_ignore(), libMesh::NLOPT_SOLVERS, and libMesh::PETSC_SOLVERS.

63 {
64  // Prevent unused variables warnings when Tao is not available
65  libmesh_ignore(s);
66 
67  // Build the appropriate solver
68  switch (solver_package)
69  {
70 
71 #if defined(LIBMESH_HAVE_PETSC_TAO) && !defined(LIBMESH_USE_COMPLEX_NUMBERS)
72  case PETSC_SOLVERS:
73  return libmesh_make_unique<TaoOptimizationSolver<T>>(s);
74 #endif // #if defined(LIBMESH_HAVE_PETSC_TAO) && !defined(LIBMESH_USE_COMPLEX_NUMBERS)
75 
76 #if defined(LIBMESH_HAVE_NLOPT) && !defined(LIBMESH_USE_COMPLEX_NUMBERS)
77  case NLOPT_SOLVERS:
78  return libmesh_make_unique<NloptOptimizationSolver<T>>(s);
79 #endif // #if defined(LIBMESH_HAVE_NLOPT) && !defined(LIBMESH_USE_COMPLEX_NUMBERS)
80 
81  default:
82  libmesh_error_msg("ERROR: Unrecognized solver package: " << solver_package);
83  }
84 }
void libmesh_ignore(const Args &...)

◆ clear()

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

Release all memory and clear data structures.

Reimplemented from libMesh::OptimizationSolver< T >.

Definition at line 429 of file tao_optimization_solver.C.

References libMesh::libMeshPrivateData::_is_initialized, ierr, and libMesh::initialized().

430 {
431  if (this->initialized())
432  {
433  this->_is_initialized = false;
434 
435  PetscErrorCode ierr=0;
436 
437  ierr = TaoDestroy(&_tao);
438  LIBMESH_CHKERR(ierr);
439  }
440 }
PetscErrorCode ierr

◆ 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 >
int libMesh::TaoOptimizationSolver< T >::get_converged_reason ( )
overridevirtual
Returns
The currently-available (or most recently obtained, if the Tao object has been destroyed) convergence reason.

Refer to Tao docs for the meaning of different TaoConvergedReason.

Reimplemented from libMesh::OptimizationSolver< T >.

Definition at line 640 of file tao_optimization_solver.C.

References ierr, and libMesh::initialized().

641 {
642  PetscErrorCode ierr=0;
643 
644  if (this->initialized())
645  {
646  ierr = TaoGetConvergedReason(_tao, &_reason);
647  LIBMESH_CHKERR(ierr);
648  }
649 
650  return static_cast<int>(_reason);
651 }
PetscErrorCode ierr

◆ get_dual_variables()

template<typename T >
void libMesh::TaoOptimizationSolver< T >::get_dual_variables ( )
overridevirtual

Get the current values of dual variables associated with inequality and equality constraints. The variables will be stored in _system.lambda_eq and _system.lambda_ineq.

Reimplemented from libMesh::OptimizationSolver< T >.

Definition at line 610 of file tao_optimization_solver.C.

References ierr, and libMesh::PetscVector< T >::vec().

611 {
612  LOG_SCOPE("get_dual_variables()", "TaoOptimizationSolver");
613 
614  PetscVector<T> * lambda_eq_petsc =
615  cast_ptr<PetscVector<T> *>(this->system().lambda_eq.get());
616  PetscVector<T> * lambda_ineq_petsc =
617  cast_ptr<PetscVector<T> *>(this->system().lambda_ineq.get());
618 
619  Vec lambda_eq_petsc_vec = lambda_eq_petsc->vec();
620  Vec lambda_ineq_petsc_vec = lambda_ineq_petsc->vec();
621 
622  PetscErrorCode ierr = 0;
623  ierr = TaoGetDualVariables(_tao,
624  &lambda_eq_petsc_vec,
625  &lambda_ineq_petsc_vec);
626  LIBMESH_CHKERR(ierr);
627 }
std::unique_ptr< NumericVector< Number > > lambda_ineq
const sys_type & system() const
std::unique_ptr< NumericVector< Number > > lambda_eq
PetscErrorCode ierr

◆ 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

◆ 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::TaoOptimizationSolver< T >::init ( )
overridevirtual

Initialize data structures if not done so already.

Implements libMesh::OptimizationSolver< T >.

Definition at line 445 of file tao_optimization_solver.C.

References libMesh::libMeshPrivateData::_is_initialized, ierr, and libMesh::initialized().

Referenced by libMesh::TaoOptimizationSolver< T >::tao().

446 {
447  // Initialize the data structures if not done so already.
448  if (!this->initialized())
449  {
450  this->_is_initialized = true;
451 
452  PetscErrorCode ierr=0;
453 
454  ierr = TaoCreate(this->comm().get(),&_tao);
455  LIBMESH_CHKERR(ierr);
456  }
457 }
const Parallel::Communicator & comm() const
PetscErrorCode ierr

◆ initialized()

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

Definition at line 92 of file optimization_solver.h.

References libMesh::OptimizationSolver< T >::_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

◆ print_converged_reason()

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

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

Reimplemented from libMesh::OptimizationSolver< T >.

Definition at line 631 of file tao_optimization_solver.C.

References libMesh::out.

632 {
633  libMesh::out << "Tao optimization solver convergence/divergence reason: "
634  << TaoConvergedReasons[this->get_converged_reason()] << std::endl;
635 }
virtual int get_converged_reason() override
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

◆ solve()

template<typename T >
void libMesh::TaoOptimizationSolver< T >::solve ( )
overridevirtual

Call the Tao solver.

Implements libMesh::OptimizationSolver< T >.

Definition at line 460 of file tao_optimization_solver.C.

References 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(), ierr, libMesh::TriangleWrapper::init(), libMesh::PetscMatrix< T >::mat(), and libMesh::PetscVector< T >::zero().

461 {
462  LOG_SCOPE("solve()", "TaoOptimizationSolver");
463 
464  this->init ();
465 
466  this->system().solution->zero();
467 
468  PetscMatrix<T> * hessian = cast_ptr<PetscMatrix<T> *>(this->system().matrix);
469  // PetscVector<T> * gradient = cast_ptr<PetscVector<T> *>(this->system().rhs);
470  PetscVector<T> * x = cast_ptr<PetscVector<T> *>(this->system().solution.get());
471  PetscVector<T> * ceq = cast_ptr<PetscVector<T> *>(this->system().C_eq.get());
472  PetscMatrix<T> * ceq_jac = cast_ptr<PetscMatrix<T> *>(this->system().C_eq_jac.get());
473  PetscVector<T> * cineq = cast_ptr<PetscVector<T> *>(this->system().C_ineq.get());
474  PetscMatrix<T> * cineq_jac = cast_ptr<PetscMatrix<T> *>(this->system().C_ineq_jac.get());
475  PetscVector<T> * lb = cast_ptr<PetscVector<T> *>(&this->system().get_vector("lower_bounds"));
476  PetscVector<T> * ub = cast_ptr<PetscVector<T> *>(&this->system().get_vector("upper_bounds"));
477 
478  // Set the starting guess to zero.
479  x->zero();
480 
481  PetscErrorCode ierr = 0;
482 
483  // Workaround for bug where TaoSetFromOptions *reset*
484  // programmatically set tolerance and max. function evaluation
485  // values when "-tao_type ipm" was specified on the command line: we
486  // call TaoSetFromOptions twice (both before and after setting
487  // custom options programmatically)
488  ierr = TaoSetFromOptions(_tao);
489  LIBMESH_CHKERR(ierr);
490 
491  // Set convergence tolerances
492  // f(X) - f(X*) (estimated) <= fatol
493  // |f(X) - f(X*)| (estimated) / |f(X)| <= frtol
494  // ||g(X)|| <= gatol
495  // ||g(X)|| / |f(X)| <= grtol
496  // ||g(X)|| / ||g(X0)|| <= gttol
497  // Command line equivalents: -tao_fatol, -tao_frtol, -tao_gatol, -tao_grtol, -tao_gttol
498  ierr = TaoSetTolerances(_tao,
499 #if PETSC_RELEASE_LESS_THAN(3,7,0)
500  // Releases up to 3.X.Y had fatol and frtol, after that they were removed.
501  // Hopefully we'll be able to know X and Y soon. Guessing at 3.7.0.
502  /*fatol=*/PETSC_DEFAULT,
503  /*frtol=*/PETSC_DEFAULT,
504 #endif
505  /*gatol=*/PETSC_DEFAULT,
507  /*gttol=*/PETSC_DEFAULT);
508  LIBMESH_CHKERR(ierr);
509 
510  // Set the max-allowed number of objective function evaluations
511  // Command line equivalent: -tao_max_funcs
512  ierr = TaoSetMaximumFunctionEvaluations(_tao, this->max_objective_function_evaluations);
513  LIBMESH_CHKERR(ierr);
514 
515  // Set the max-allowed number of optimization iterations.
516  // Command line equivalent: -tao_max_it
517  // Not implemented for now as it seems fairly similar to
518  // ierr = TaoSetMaximumIterations(_tao, 4);
519  // LIBMESH_CHKERR(ierr);
520 
521  // Set solution vec and an initial guess
522  ierr = TaoSetInitialVector(_tao, x->vec());
523  LIBMESH_CHKERR(ierr);
524 
525  // We have to have an objective function
526  libmesh_assert( this->objective_object );
527 
528  // Set routines for objective, gradient, hessian evaluation
529  ierr = TaoSetObjectiveRoutine(_tao, __libmesh_tao_objective, this);
530  LIBMESH_CHKERR(ierr);
531 
532  if (this->gradient_object)
533  {
534  ierr = TaoSetGradientRoutine(_tao, __libmesh_tao_gradient, this);
535  LIBMESH_CHKERR(ierr);
536  }
537 
538  if (this->hessian_object)
539  {
540  ierr = TaoSetHessianRoutine(_tao, hessian->mat(), hessian->mat(), __libmesh_tao_hessian, this);
541  LIBMESH_CHKERR(ierr);
542  }
543 
545  {
546  // Need to actually compute the bounds vectors first
548 
549  ierr = TaoSetVariableBounds(_tao,
550  lb->vec(),
551  ub->vec());
552  LIBMESH_CHKERR(ierr);
553  }
554 
555  if (this->equality_constraints_object)
556  {
557  ierr = TaoSetEqualityConstraintsRoutine(_tao, ceq->vec(), __libmesh_tao_equality_constraints, this);
558  LIBMESH_CHKERR(ierr);
559  }
560 
562  {
563  ierr = TaoSetJacobianEqualityRoutine(_tao,
564  ceq_jac->mat(),
565  ceq_jac->mat(),
567  this);
568  LIBMESH_CHKERR(ierr);
569  }
570 
571  // Optionally set inequality constraints
573  {
574  ierr = TaoSetInequalityConstraintsRoutine(_tao, cineq->vec(), __libmesh_tao_inequality_constraints, this);
575  LIBMESH_CHKERR(ierr);
576  }
577 
578  // Optionally set inequality constraints Jacobian
580  {
581  ierr = TaoSetJacobianInequalityRoutine(_tao,
582  cineq_jac->mat(),
583  cineq_jac->mat(),
585  this);
586  LIBMESH_CHKERR(ierr);
587  }
588 
589  // Check for Tao command line options
590  ierr = TaoSetFromOptions(_tao);
591  LIBMESH_CHKERR(ierr);
592 
593  // Perform the optimization
594  ierr = TaoSolve(_tao);
595  LIBMESH_CHKERR(ierr);
596 
597  // Enforce constraints exactly now that the solve is done. We have
598  // been enforcing them on the current_local_solution during the
599  // solve, but now need to be sure they are enforced on the parallel
600  // solution vector as well.
602 
603  // Store the convergence/divergence reason
604  ierr = TaoGetConvergedReason(_tao, &_reason);
605  LIBMESH_CHKERR(ierr);
606 }
OptimizationSystem::ComputeObjective * objective_object
friend PetscErrorCode __libmesh_tao_inequality_constraints(Tao tao, Vec x, Vec cineq, void *ctx)
friend PetscErrorCode __libmesh_tao_hessian(Tao tao, Vec x, Mat h, Mat pc, void *ctx)
std::unique_ptr< NumericVector< Number > > C_eq
friend PetscErrorCode __libmesh_tao_equality_constraints(Tao tao, Vec x, Vec ce, void *ctx)
OptimizationSystem::ComputeGradient * gradient_object
std::unique_ptr< SparseMatrix< Number > > C_ineq_jac
OptimizationSystem::ComputeHessian * hessian_object
friend PetscErrorCode __libmesh_tao_equality_constraints_jacobian(Tao tao, Vec x, Mat J, Mat Jpre, void *ctx)
OptimizationSystem::ComputeInequalityConstraints * inequality_constraints_object
OptimizationSystem::ComputeEqualityConstraintsJacobian * equality_constraints_jacobian_object
friend PetscErrorCode __libmesh_tao_inequality_constraints_jacobian(Tao tao, Vec x, Mat J, Mat Jpre, void *ctx)
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:774
std::unique_ptr< NumericVector< Number > > solution
Definition: system.h:1523
friend PetscErrorCode __libmesh_tao_objective(Tao tao, Vec x, PetscReal *objective, void *ctx)
unsigned int max_objective_function_evaluations
friend PetscErrorCode __libmesh_tao_gradient(Tao tao, Vec x, Vec g, void *ctx)
std::unique_ptr< SparseMatrix< Number > > C_eq_jac
const sys_type & system() const
OptimizationSystem::ComputeLowerAndUpperBounds * lower_and_upper_bounds_object
PetscErrorCode ierr
SparseMatrix< Number > * matrix
const DofMap & get_dof_map() const
Definition: system.h:2049
OptimizationSystem::ComputeInequalityConstraintsJacobian * inequality_constraints_jacobian_object
std::unique_ptr< NumericVector< Number > > C_ineq
OptimizationSystem::ComputeEqualityConstraints * equality_constraints_object
void enforce_constraints_exactly(const System &system, NumericVector< Number > *v=nullptr, bool homogeneous=false) const

◆ system() [1/2]

◆ system() [2/2]

template<typename T >
sys_type& libMesh::OptimizationSolver< T >::system ( )
inlineinherited
Returns
A writable reference to the system we are using to define the optimization problem.

Definition at line 189 of file optimization_solver.h.

References libMesh::OptimizationSolver< T >::_system.

189 { return _system; }

◆ tao()

template<typename T>
Tao libMesh::TaoOptimizationSolver< T >::tao ( )
inline
Returns
The raw PETSc Tao context pointer.

Definition at line 91 of file tao_optimization_solver.h.

References libMesh::TaoOptimizationSolver< T >::_tao, and libMesh::TaoOptimizationSolver< T >::init().

Friends And Related Function Documentation

◆ __libmesh_tao_equality_constraints

template<typename T>
PetscErrorCode __libmesh_tao_equality_constraints ( Tao  tao,
Vec  x,
Vec  ce,
void *  ctx 
)
friend

Definition at line 205 of file tao_optimization_solver.C.

206  {
207  LOG_SCOPE("equality_constraints()", "TaoOptimizationSolver");
208 
209  PetscErrorCode ierr = 0;
210 
211  libmesh_assert(x);
212  libmesh_assert(ce);
213  libmesh_assert(ctx);
214 
215  // ctx should be a pointer to the solver (it was passed in as void *)
216  TaoOptimizationSolver<Number> * solver =
217  static_cast<TaoOptimizationSolver<Number> *> (ctx);
218 
219  OptimizationSystem & sys = solver->system();
220 
221  // We'll use current_local_solution below, so let's ensure that it's consistent
222  // with the vector x that was passed in.
223  PetscVector<Number> & X_sys = *cast_ptr<PetscVector<Number> *>(sys.solution.get());
224  PetscVector<Number> X(x, sys.comm());
225 
226  // Perform a swap so that sys.solution points to the input vector
227  // "x", update sys.current_local_solution based on "x", then swap
228  // back.
229  X.swap(X_sys);
230  sys.update();
231  X.swap(X_sys);
232 
233  // We'll also pass the constraints vector ce into the assembly routine
234  // so let's make a PETSc vector for that too.
235  PetscVector<Number> eq_constraints(ce, sys.comm());
236 
237  // Clear the gradient prior to assembly
238  eq_constraints.zero();
239 
240  // Enforce constraints exactly on the current_local_solution.
241  sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
242 
243  if (solver->equality_constraints_object != nullptr)
244  solver->equality_constraints_object->equality_constraints(*(sys.current_local_solution), eq_constraints, sys);
245  else
246  libmesh_error_msg("Constraints function not defined in __libmesh_tao_equality_constraints");
247 
248  eq_constraints.close();
249 
250  return ierr;
251  }
PetscErrorCode ierr

◆ __libmesh_tao_equality_constraints_jacobian

template<typename T>
PetscErrorCode __libmesh_tao_equality_constraints_jacobian ( Tao  tao,
Vec  x,
Mat  J,
Mat  Jpre,
void *  ctx 
)
friend

Definition at line 257 of file tao_optimization_solver.C.

258  {
259  LOG_SCOPE("equality_constraints_jacobian()", "TaoOptimizationSolver");
260 
261  PetscErrorCode ierr = 0;
262 
263  libmesh_assert(x);
264  libmesh_assert(J);
265  libmesh_assert(Jpre);
266 
267  // ctx should be a pointer to the solver (it was passed in as void *)
268  TaoOptimizationSolver<Number> * solver =
269  static_cast<TaoOptimizationSolver<Number> *> (ctx);
270 
271  OptimizationSystem & sys = solver->system();
272 
273  // We'll use current_local_solution below, so let's ensure that it's consistent
274  // with the vector x that was passed in.
275  PetscVector<Number> & X_sys = *cast_ptr<PetscVector<Number> *>(sys.solution.get());
276  PetscVector<Number> X(x, sys.comm());
277 
278  // Perform a swap so that sys.solution points to the input vector
279  // "x", update sys.current_local_solution based on "x", then swap
280  // back.
281  X.swap(X_sys);
282  sys.update();
283  X.swap(X_sys);
284 
285  // Let's also wrap J and Jpre in PetscMatrix objects for convenience
286  PetscMatrix<Number> J_petsc(J, sys.comm());
287  PetscMatrix<Number> Jpre_petsc(Jpre, sys.comm());
288 
289  // Enforce constraints exactly on the current_local_solution.
290  sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
291 
292  if (solver->equality_constraints_jacobian_object != nullptr)
293  solver->equality_constraints_jacobian_object->equality_constraints_jacobian(*(sys.current_local_solution), J_petsc, sys);
294  else
295  libmesh_error_msg("Constraints function not defined in __libmesh_tao_equality_constraints_jacobian");
296 
297  J_petsc.close();
298  Jpre_petsc.close();
299 
300  return ierr;
301  }
PetscErrorCode ierr

◆ __libmesh_tao_gradient

template<typename T>
PetscErrorCode __libmesh_tao_gradient ( Tao  tao,
Vec  x,
Vec  g,
void *  ctx 
)
friend

Definition at line 97 of file tao_optimization_solver.C.

98  {
99  LOG_SCOPE("gradient()", "TaoOptimizationSolver");
100 
101  PetscErrorCode ierr = 0;
102 
103  libmesh_assert(x);
104  libmesh_assert(g);
105  libmesh_assert(ctx);
106 
107  // ctx should be a pointer to the solver (it was passed in as void *)
108  TaoOptimizationSolver<Number> * solver =
109  static_cast<TaoOptimizationSolver<Number> *> (ctx);
110 
111  OptimizationSystem & sys = solver->system();
112 
113  // We'll use current_local_solution below, so let's ensure that it's consistent
114  // with the vector x that was passed in.
115  PetscVector<Number> & X_sys = *cast_ptr<PetscVector<Number> *>(sys.solution.get());
116  PetscVector<Number> X(x, sys.comm());
117 
118  // Perform a swap so that sys.solution points to the input vector
119  // "x", update sys.current_local_solution based on "x", then swap
120  // back.
121  X.swap(X_sys);
122  sys.update();
123  X.swap(X_sys);
124 
125  // We'll also pass the gradient in to the assembly routine
126  // so let's make a PETSc vector for that too.
127  PetscVector<Number> gradient(g, sys.comm());
128 
129  // Clear the gradient prior to assembly
130  gradient.zero();
131 
132  // Enforce constraints exactly on the current_local_solution.
133  sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
134 
135  if (solver->gradient_object != nullptr)
136  solver->gradient_object->gradient(*(sys.current_local_solution), gradient, sys);
137  else
138  libmesh_error_msg("Gradient function not defined in __libmesh_tao_gradient");
139 
140  gradient.close();
141 
142  return ierr;
143  }
PetscErrorCode ierr

◆ __libmesh_tao_hessian

template<typename T>
PetscErrorCode __libmesh_tao_hessian ( Tao  tao,
Vec  x,
Mat  h,
Mat  pc,
void *  ctx 
)
friend

Definition at line 148 of file tao_optimization_solver.C.

149  {
150  LOG_SCOPE("hessian()", "TaoOptimizationSolver");
151 
152  PetscErrorCode ierr = 0;
153 
154  libmesh_assert(x);
155  libmesh_assert(h);
156  libmesh_assert(pc);
157  libmesh_assert(ctx);
158 
159  // ctx should be a pointer to the solver (it was passed in as void *)
160  TaoOptimizationSolver<Number> * solver =
161  static_cast<TaoOptimizationSolver<Number> *> (ctx);
162 
163  OptimizationSystem & sys = solver->system();
164 
165  // We'll use current_local_solution below, so let's ensure that it's consistent
166  // with the vector x that was passed in.
167  PetscVector<Number> & X_sys = *cast_ptr<PetscVector<Number> *>(sys.solution.get());
168  PetscVector<Number> X(x, sys.comm());
169 
170  // Perform a swap so that sys.solution points to the input vector
171  // "x", update sys.current_local_solution based on "x", then swap
172  // back.
173  X.swap(X_sys);
174  sys.update();
175  X.swap(X_sys);
176 
177  // Let's also wrap pc and h in PetscMatrix objects for convenience
178  PetscMatrix<Number> PC(pc, sys.comm());
179  PetscMatrix<Number> hessian(h, sys.comm());
180  PC.attach_dof_map(sys.get_dof_map());
181  hessian.attach_dof_map(sys.get_dof_map());
182 
183  // Enforce constraints exactly on the current_local_solution.
184  sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
185 
186  if (solver->hessian_object != nullptr)
187  {
188  // Following PetscNonlinearSolver by passing in PC. It's not clear
189  // why we pass in PC and not hessian though?
190  solver->hessian_object->hessian(*(sys.current_local_solution), PC, sys);
191  }
192  else
193  libmesh_error_msg("Hessian function not defined in __libmesh_tao_hessian");
194 
195  PC.close();
196  hessian.close();
197 
198  return ierr;
199  }
PetscErrorCode ierr

◆ __libmesh_tao_inequality_constraints

template<typename T>
PetscErrorCode __libmesh_tao_inequality_constraints ( Tao  tao,
Vec  x,
Vec  cineq,
void *  ctx 
)
friend

Definition at line 306 of file tao_optimization_solver.C.

307  {
308  LOG_SCOPE("inequality_constraints()", "TaoOptimizationSolver");
309 
310  PetscErrorCode ierr = 0;
311 
312  libmesh_assert(x);
313  libmesh_assert(cineq);
314  libmesh_assert(ctx);
315 
316  // ctx should be a pointer to the solver (it was passed in as void *)
317  TaoOptimizationSolver<Number> * solver =
318  static_cast<TaoOptimizationSolver<Number> *> (ctx);
319 
320  OptimizationSystem & sys = solver->system();
321 
322  // We'll use current_local_solution below, so let's ensure that it's consistent
323  // with the vector x that was passed in.
324  PetscVector<Number> & X_sys = *cast_ptr<PetscVector<Number> *>(sys.solution.get());
325  PetscVector<Number> X(x, sys.comm());
326 
327  // Perform a swap so that sys.solution points to the input vector
328  // "x", update sys.current_local_solution based on "x", then swap
329  // back.
330  X.swap(X_sys);
331  sys.update();
332  X.swap(X_sys);
333 
334  // We'll also pass the constraints vector ce into the assembly routine
335  // so let's make a PETSc vector for that too.
336  PetscVector<Number> ineq_constraints(cineq, sys.comm());
337 
338  // Clear the gradient prior to assembly
339  ineq_constraints.zero();
340 
341  // Enforce constraints exactly on the current_local_solution.
342  sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
343 
344  if (solver->inequality_constraints_object != nullptr)
345  solver->inequality_constraints_object->inequality_constraints(*(sys.current_local_solution), ineq_constraints, sys);
346  else
347  libmesh_error_msg("Constraints function not defined in __libmesh_tao_inequality_constraints");
348 
349  ineq_constraints.close();
350 
351  return ierr;
352  }
PetscErrorCode ierr

◆ __libmesh_tao_inequality_constraints_jacobian

template<typename T>
PetscErrorCode __libmesh_tao_inequality_constraints_jacobian ( Tao  tao,
Vec  x,
Mat  J,
Mat  Jpre,
void *  ctx 
)
friend

Definition at line 358 of file tao_optimization_solver.C.

359  {
360  LOG_SCOPE("inequality_constraints_jacobian()", "TaoOptimizationSolver");
361 
362  PetscErrorCode ierr = 0;
363 
364  libmesh_assert(x);
365  libmesh_assert(J);
366  libmesh_assert(Jpre);
367 
368  // ctx should be a pointer to the solver (it was passed in as void *)
369  TaoOptimizationSolver<Number> * solver =
370  static_cast<TaoOptimizationSolver<Number> *> (ctx);
371 
372  OptimizationSystem & sys = solver->system();
373 
374  // We'll use current_local_solution below, so let's ensure that it's consistent
375  // with the vector x that was passed in.
376  PetscVector<Number> & X_sys = *cast_ptr<PetscVector<Number> *>(sys.solution.get());
377  PetscVector<Number> X(x, sys.comm());
378 
379  // Perform a swap so that sys.solution points to the input vector
380  // "x", update sys.current_local_solution based on "x", then swap
381  // back.
382  X.swap(X_sys);
383  sys.update();
384  X.swap(X_sys);
385 
386  // Let's also wrap J and Jpre in PetscMatrix objects for convenience
387  PetscMatrix<Number> J_petsc(J, sys.comm());
388  PetscMatrix<Number> Jpre_petsc(Jpre, sys.comm());
389 
390  // Enforce constraints exactly on the current_local_solution.
391  sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
392 
393  if (solver->inequality_constraints_jacobian_object != nullptr)
394  solver->inequality_constraints_jacobian_object->inequality_constraints_jacobian(*(sys.current_local_solution), J_petsc, sys);
395  else
396  libmesh_error_msg("Constraints function not defined in __libmesh_tao_inequality_constraints_jacobian");
397 
398  J_petsc.close();
399  Jpre_petsc.close();
400 
401  return ierr;
402  }
PetscErrorCode ierr

◆ __libmesh_tao_objective

template<typename T>
PetscErrorCode __libmesh_tao_objective ( Tao  tao,
Vec  x,
PetscReal *  objective,
void *  ctx 
)
friend

Definition at line 49 of file tao_optimization_solver.C.

50  {
51  LOG_SCOPE("objective()", "TaoOptimizationSolver");
52 
53  PetscErrorCode ierr = 0;
54 
55  libmesh_assert(x);
56  libmesh_assert(objective);
57  libmesh_assert(ctx);
58 
59  // ctx should be a pointer to the solver (it was passed in as void *)
60  TaoOptimizationSolver<Number> * solver =
61  static_cast<TaoOptimizationSolver<Number> *> (ctx);
62 
63  OptimizationSystem & sys = solver->system();
64 
65  // We'll use current_local_solution below, so let's ensure that it's consistent
66  // with the vector x that was passed in.
67  PetscVector<Number> & X_sys = *cast_ptr<PetscVector<Number> *>(sys.solution.get());
68  PetscVector<Number> X(x, sys.comm());
69 
70  // Perform a swap so that sys.solution points to the input vector
71  // "x", update sys.current_local_solution based on "x", then swap
72  // back.
73  X.swap(X_sys);
74  sys.update();
75  X.swap(X_sys);
76 
77  // Enforce constraints (if any) exactly on the
78  // current_local_solution. This is the solution vector that is
79  // actually used in the computation of the objective function
80  // below, and is not locked by debug-enabled PETSc the way that
81  // the solution vector is.
82  sys.get_dof_map().enforce_constraints_exactly(sys, sys.current_local_solution.get());
83 
84  if (solver->objective_object != nullptr)
85  (*objective) = solver->objective_object->objective(*(sys.current_local_solution), sys);
86  else
87  libmesh_error_msg("Objective function not defined in __libmesh_tao_objective");
88 
89  return ierr;
90  }
PetscErrorCode ierr

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::OptimizationSolver< T >::_is_initialized
protectedinherited

Flag indicating if the data structures have been initialized.

Definition at line 216 of file optimization_solver.h.

Referenced by libMesh::OptimizationSolver< T >::initialized().

◆ _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().

◆ _reason

template<typename T>
TaoConvergedReason libMesh::TaoOptimizationSolver< T >::_reason
protected

Store the reason for Tao convergence/divergence for use even after _tao has been cleared.

Note
print_converged_reason() will always try to get the current reason with TaoGetConvergedReason(), but if the Tao object has already been cleared, it will fall back on this stored value. This value is therefore necessarily not cleared by the clear() function.

Definition at line 136 of file tao_optimization_solver.h.

◆ _system

template<typename T >
sys_type& libMesh::OptimizationSolver< T >::_system
protectedinherited

A reference to the system we are solving.

Definition at line 211 of file optimization_solver.h.

Referenced by libMesh::OptimizationSolver< T >::system().

◆ _tao

template<typename T>
Tao libMesh::TaoOptimizationSolver< T >::_tao
protected

Optimization solver context

Definition at line 124 of file tao_optimization_solver.h.

Referenced by libMesh::TaoOptimizationSolver< T >::tao().

◆ equality_constraints_jacobian_object

template<typename T >
OptimizationSystem::ComputeEqualityConstraintsJacobian* libMesh::OptimizationSolver< T >::equality_constraints_jacobian_object
inherited

Object that computes the Jacobian of C_eq(X).

Definition at line 161 of file optimization_solver.h.

Referenced by libMesh::__libmesh_nlopt_equality_constraints(), and libMesh::__libmesh_tao_equality_constraints_jacobian().

◆ equality_constraints_object

template<typename T >
OptimizationSystem::ComputeEqualityConstraints* libMesh::OptimizationSolver< T >::equality_constraints_object
inherited

Object that computes the equality constraints vector C_eq(X). This will lead to the constraints C_eq(X) = 0 being imposed.

Definition at line 156 of file optimization_solver.h.

Referenced by libMesh::__libmesh_nlopt_equality_constraints(), and libMesh::__libmesh_tao_equality_constraints().

◆ gradient_object

template<typename T >
OptimizationSystem::ComputeGradient* libMesh::OptimizationSolver< T >::gradient_object
inherited

Object that computes the gradient grad_f(X) of the objective function at the input iterate X.

Definition at line 144 of file optimization_solver.h.

Referenced by libMesh::__libmesh_nlopt_objective(), and libMesh::__libmesh_tao_gradient().

◆ hessian_object

template<typename T >
OptimizationSystem::ComputeHessian* libMesh::OptimizationSolver< T >::hessian_object
inherited

Object that computes the Hessian H_f(X) of the objective function at the input iterate X.

Definition at line 150 of file optimization_solver.h.

Referenced by libMesh::__libmesh_tao_hessian().

◆ inequality_constraints_jacobian_object

template<typename T >
OptimizationSystem::ComputeInequalityConstraintsJacobian* libMesh::OptimizationSolver< T >::inequality_constraints_jacobian_object
inherited

Object that computes the Jacobian of C_ineq(X).

Definition at line 172 of file optimization_solver.h.

Referenced by libMesh::__libmesh_nlopt_inequality_constraints(), and libMesh::__libmesh_tao_inequality_constraints_jacobian().

◆ inequality_constraints_object

template<typename T >
OptimizationSystem::ComputeInequalityConstraints* libMesh::OptimizationSolver< T >::inequality_constraints_object
inherited

Object that computes the inequality constraints vector C_ineq(X). This will lead to the constraints C_ineq(X) >= 0 being imposed.

Definition at line 167 of file optimization_solver.h.

Referenced by libMesh::__libmesh_nlopt_inequality_constraints(), and libMesh::__libmesh_tao_inequality_constraints().

◆ lower_and_upper_bounds_object

template<typename T >
OptimizationSystem::ComputeLowerAndUpperBounds* libMesh::OptimizationSolver< T >::lower_and_upper_bounds_object
inherited

Object that computes the lower and upper bounds vectors.

Definition at line 177 of file optimization_solver.h.

◆ max_objective_function_evaluations

template<typename T >
unsigned int libMesh::OptimizationSolver< T >::max_objective_function_evaluations
inherited

Maximum number of objective function evaluations allowed.

Definition at line 194 of file optimization_solver.h.

◆ objective_function_relative_tolerance

template<typename T >
Real libMesh::OptimizationSolver< T >::objective_function_relative_tolerance
inherited

Required change in objective function which signals convergence.

Definition at line 199 of file optimization_solver.h.

◆ objective_object

template<typename T >
OptimizationSystem::ComputeObjective* libMesh::OptimizationSolver< T >::objective_object
inherited

Object that computes the objective function f(X) at the input iterate X.

Definition at line 138 of file optimization_solver.h.

Referenced by libMesh::__libmesh_nlopt_objective(), and libMesh::__libmesh_tao_objective().

◆ verbose

template<typename T >
bool libMesh::OptimizationSolver< T >::verbose
inherited

Control how much is output from the OptimizationSolver as it's running.

Definition at line 204 of file optimization_solver.h.

Referenced by libMesh::__libmesh_nlopt_objective().


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