37 #include "libmesh/numeric_vector.h" 38 #include "libmesh/equation_systems.h" 39 #include "libmesh/sparse_matrix.h" 40 #include "libmesh/dof_map.h" 41 #include "libmesh/nonlinear_solver.h" 42 #include "libmesh/petsc_linear_solver.h" 43 #include "libmesh/xdr_cxx.h" 44 #include "libmesh/mesh_tools.h" 45 #include "libmesh/utility.h" 46 #include "libmesh/libmesh_version.h" 47 #include "libmesh/generic_projector.h" 48 #include "libmesh/wrapped_functor.h" 49 #include "libmesh/fem_context.h" 50 #include "libmesh/parallel.h" 54 const std::string& name,
55 const unsigned int number):
56 libMesh::NonlinearImplicitSystem(es, name, number),
57 _initialize_B_matrix (false),
60 eigen_solver (nullptr),
61 _condensed_dofs_initialized (false),
62 _exchange_A_and_B (false),
63 _n_requested_eigenpairs (0),
64 _n_converged_eigenpairs (0),
66 _is_generalized_eigenproblem (false),
67 _eigen_problem_type (libMesh::NHEP),
100 libMesh::NonlinearImplicitSystem::clear();
118 libMesh::NonlinearImplicitSystem::init_data();
127 libMesh::DofMap& dof_map = this->get_dof_map();
130 matrix_A = libMesh::SparseMatrix<Real>::build(this->comm()).release();
137 matrix_B = libMesh::SparseMatrix<Real>::build(this->comm()).release();
144 if (libMesh::on_command_line(
"--solver_system_names")) {
147 std::string nm = this->name() +
"_";
148 EPSSetOptionsPrefix(eps, nm.c_str());
153 linear_solver.reset(
new libMesh::PetscLinearSolver<Real>(this->comm()));
154 if (libMesh::on_command_line(
"--solver_system_names")) {
156 std::string nm = this->name() +
"_";
166 libMesh::NonlinearImplicitSystem::reinit();
175 if (libMesh::on_command_line(
"--solver_system_names")) {
178 std::string nm = this->name() +
"_";
179 EPSSetOptionsPrefix(eps, nm.c_str());
184 linear_solver.reset(
new libMesh::PetscLinearSolver<Real>(this->comm()));
185 if (libMesh::on_command_line(
"--solver_system_names")) {
187 std::string nm = this->name() +
"_";
201 std::pair<unsigned int, Real>
204 this->set_solver_parameters();
205 return libMesh::NonlinearImplicitSystem::get_linear_solve_parameters();
218 libMesh::NonlinearImplicitSystem::ComputeResidualandJacobian
219 *old_ptr = this->nonlinear_solver->residual_and_jacobian_object;
221 this->nonlinear_solver->residual_and_jacobian_object = &assembly;
225 libMesh::NonlinearImplicitSystem::solve();
229 this->get_dof_map().enforce_constraints_exactly(*
this, this->solution.get());
231 this->nonlinear_solver->residual_and_jacobian_object = old_ptr;
247 START_LOG(
"eigensolve()",
"NonlinearSystem");
252 libMesh::EquationSystems& es = this->get_equation_systems();
264 tol = es.parameters.get<
Real>(
"linear solver tolerance");
267 maxits = es.parameters.get<
unsigned int>(
"linear solver maximum iterations"),
268 nev = es.parameters.get<
unsigned int>(
"eigenpairs"),
269 ncv = es.parameters.get<
unsigned int>(
"basis vectors");
271 std::pair<unsigned int, unsigned int> solve_data;
273 libMesh::SparseMatrix<Real>
323 std::unique_ptr<libMesh::SparseMatrix<Real> >
324 condensed_matrix_A(libMesh::SparseMatrix<Real>::build(this->comm()).release()),
325 condensed_matrix_B(libMesh::SparseMatrix<Real>::build(this->comm()).release());
328 matrix_A->create_submatrix(*condensed_matrix_A,
335 matrix_B->create_submatrix(*condensed_matrix_B,
347 eig_A = condensed_matrix_A.get();
348 eig_B = condensed_matrix_B.get();
351 eig_B = condensed_matrix_A.get();
352 eig_A = condensed_matrix_B.get();
367 solve_data =
eigen_solver->solve_standard (*condensed_matrix_A,
380 STOP_LOG(
"eigensolve()",
"NonlinearSystem");
390 std::pair<Real, Real>
398 Complex complex_val (val.first, val.second);
399 complex_val = 1./complex_val;
400 re = complex_val.real();
401 im = complex_val.imag();
411 libMesh::NumericVector<Real>& vec_re,
412 libMesh::NumericVector<Real>* vec_im) {
414 START_LOG(
"get_eigenpair()",
"NonlinearSystem");
415 std::pair<Real, Real>
421 libmesh_assert (vec_im ==
nullptr);
429 val = this->
eigen_solver->get_eigenpair (i, vec_re, vec_im);
436 Complex complex_val (val.first, val.second);
437 complex_val = 1./complex_val;
438 re = complex_val.real();
439 im = complex_val.imag();
447 std::unique_ptr< libMesh::NumericVector<Real> >
448 temp_re(libMesh::NumericVector<Real>::build(this->comm()).release()),
457 temp_re->init (n, n_local,
false, libMesh::PARALLEL);
462 temp_im.reset(libMesh::NumericVector<Real>::build(this->comm()).release());
463 temp_im->init (n, n_local,
false, libMesh::PARALLEL);
468 val = this->
eigen_solver->get_eigenpair (i, *temp_re, temp_im.get());
475 Complex complex_val (val.first, val.second);
476 complex_val = 1./complex_val;
477 re = complex_val.real();
478 im = complex_val.imag();
488 vec_re.set(index,(*temp_re)(temp_re->first_local_index()+j));
500 vec_im->set(index,(*temp_im)(temp_im->first_local_index()+j));
512 v = vec_re.dot(vec_re);
515 libmesh_assert_greater(v, 0.);
516 vec_re.scale(1./std::sqrt(v));
520 case libMesh::GHEP: {
522 std::unique_ptr<libMesh::NumericVector<Real> >
523 tmp(vec_re.zero_clone().release());
526 matrix_B->vector_mult(*tmp, vec_re);
529 v = tmp->dot(vec_re);
532 libmesh_assert_greater(v, 0.);
533 vec_re.scale(1./std::sqrt(v));
547 STOP_LOG(
"get_eigenpair()",
"NonlinearSystem");
559 std::vector<Real>& sens,
560 const std::vector<unsigned int>* indices) {
579 n_calc = indices?(
unsigned int)indices->size():nconv;
581 libmesh_assert_equal_to(n_calc, nconv);
583 std::vector<unsigned int> indices_to_calculate;
585 indices_to_calculate = *indices;
586 for (
unsigned int i=0; i<n_calc; i++) libmesh_assert_less(indices_to_calculate[i], nconv);
590 indices_to_calculate.resize(n_calc);
591 for (
unsigned int i=0; i<n_calc; i++) indices_to_calculate[i] = i;
596 sens.resize(n_calc, 0.);
598 std::vector<libMesh::NumericVector<Real>*>
602 std::unique_ptr<libMesh::NumericVector<Real> >
603 tmp (this->solution->zero_clone().release());
613 for (
unsigned int i=0; i<n_calc; i++) {
615 x_right[i] = (this->solution->zero_clone().release());
623 this->
get_eigenpair(indices_to_calculate[i], re, im, *x_right[i],
nullptr);
624 denom[i] = x_right[i]->dot(*x_right[i]);
629 case libMesh::GHEP: {
631 this->
get_eigenpair(indices_to_calculate[i], re, im, *x_right[i],
nullptr);
632 matrix_B->vector_mult(*tmp, *x_right[i]);
633 denom[i] = x_right[i]->dot(*tmp);
650 for (
unsigned int i=0; i<n_calc; i++) {
656 matrix_A->vector_mult(*tmp, *x_right[i]);
657 sens[i] = x_right[i]->dot(*tmp);
658 sens[i]-= eig[i] * x_right[i]->dot(*x_right[i]);
663 case libMesh::GHEP: {
665 matrix_A->vector_mult(*tmp, *x_right[i]);
666 sens[i] = x_right[i]->dot(*tmp);
667 matrix_B->vector_mult(*tmp, *x_right[i]);
668 sens[i]-= eig[i] * x_right[i]->dot(*tmp);
681 for (
unsigned int i=0; i<x_right.size(); i++)
695 const libMesh::DofMap & dof_map = this->get_dof_map();
697 std::set<unsigned int> global_dirichlet_dofs_set;
703 std::set<unsigned int> local_non_condensed_dofs_set;
704 for(
unsigned int i=this->get_dof_map().first_dof();
705 i<this->get_dof_map().end_dof();
707 if (!dof_map.is_constrained_dof(i))
708 local_non_condensed_dofs_set.insert(i);
711 std::set<unsigned int>::iterator
712 iter = global_dirichlet_dofs_set.begin(),
713 iter_end = global_dirichlet_dofs_set.end();
715 for ( ; iter != iter_end ; ++iter) {
717 unsigned int condensed_dof_index = *iter;
719 if ( (this->get_dof_map().first_dof() <= condensed_dof_index) &&
720 (condensed_dof_index < this->get_dof_map().end_dof()) )
721 local_non_condensed_dofs_set.erase(condensed_dof_index);
726 iter = local_non_condensed_dofs_set.begin();
727 iter_end = local_non_condensed_dofs_set.end();
731 for ( ; iter != iter_end; ++iter)
743 bool if_localize_sol,
747 bool if_assemble_jacobian) {
754 LOG_SCOPE(
"sensitivity_solve()",
"NonlinearSystem");
758 libMesh::NumericVector<Real>
759 &dsol = this->add_sensitivity_solution(),
760 &rhs = this->add_sensitivity_rhs();
762 if (if_assemble_jacobian)
771 std::pair<unsigned int, Real>
775 libMesh::SparseMatrix<Real> * pc = this->request_matrix(
"Preconditioner");
777 std::pair<unsigned int, Real> rval =
781 solver_params.second,
782 solver_params.first);
785 #ifdef LIBMESH_ENABLE_CONSTRAINTS 786 this->get_dof_map().enforce_constraints_exactly (*
this, &dsol,
true);
800 bool if_localize_sol,
804 bool if_assemble_jacobian) {
812 LOG_SCOPE(
"adjoint_solve()",
"NonlinearSystem");
814 libMesh::NumericVector<Real>
815 &dsol = this->add_adjoint_solution(),
816 &rhs = this->add_adjoint_rhs();
820 if (if_assemble_jacobian)
831 std::pair<unsigned int, Real>
834 const std::pair<unsigned int, Real> rval =
838 solver_params.second,
839 solver_params.first);
842 #ifdef LIBMESH_ENABLE_CONSTRAINTS 843 this->get_dof_map().enforce_adjoint_constraints_exactly(dsol, 0);
855 const std::string & directory_name,
856 const std::string & data_name,
857 const bool write_binary_vectors)
859 LOG_SCOPE(
"write_out_vector()",
"NonlinearSystem");
861 if (this->processor_id() == 0)
864 libMesh::Utility::mkdir(directory_name.c_str());
868 this->comm().barrier();
870 std::ostringstream file_name;
871 const std::string suffix = (write_binary_vectors ?
".xdr" :
".dat");
873 file_name << directory_name <<
"/" << data_name <<
"_data" << suffix;
874 libMesh::Xdr bf_data(file_name.str(),
875 write_binary_vectors ? libMesh::ENCODE : libMesh::WRITE);
877 std::string version(
"libMesh-" + libMesh::get_io_compatibility_version());
878 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 879 version +=
" with infinite elements";
881 bf_data.data(version ,
"# File Format Identifier");
883 this->write_header(bf_data, version,
false);
887 libMesh::MeshTools::Private::globally_renumber_nodes_and_elements(this->get_mesh());
892 std::vector<const libMesh::NumericVector<Real> *> bf_out(1);
894 this->write_serialized_vectors (bf_data, bf_out);
899 bf_data.set_version(LIBMESH_VERSION_ID(LIBMESH_MAJOR_VERSION,
900 LIBMESH_MINOR_VERSION,
901 LIBMESH_MICRO_VERSION));
905 this->get_mesh().fix_broken_node_and_element_numbering();
912 const std::string & directory_name,
913 const std::string & data_name,
914 const bool read_binary_vector) {
916 LOG_SCOPE(
"read_in_vector()",
"NonlinearSystem");
919 this->comm().barrier();
925 libMesh::MeshTools::Private::globally_renumber_nodes_and_elements(this->get_mesh());
927 std::ostringstream file_name;
928 const std::string suffix = (read_binary_vector ?
".xdr" :
".dat");
930 file_name << directory_name
932 <<
"_data" << suffix;
935 if (this->processor_id() == 0)
937 struct stat stat_info;
938 int stat_result = stat(file_name.str().c_str(), &stat_info);
940 if (stat_result != 0)
941 libmesh_error_msg(
"File does not exist: " + file_name.str());
944 if (!std::ifstream(file_name.str()))
945 libmesh_error_msg(
"File missing: " + file_name.str());
947 libMesh::Xdr vector_data(file_name.str(),
948 read_binary_vector ? libMesh::DECODE : libMesh::READ);
953 vector_data.data(version);
955 const std::string libMesh_label =
"libMesh-";
956 std::string::size_type lm_pos = version.find(libMesh_label);
957 if (lm_pos==std::string::npos)
959 libmesh_error_msg(
"version info missing in Xdr header");
962 std::istringstream iss(version.substr(lm_pos + libMesh_label.size()));
963 int ver_major = 0, ver_minor = 0, ver_patch = 0;
965 iss >> ver_major >> dot >> ver_minor >> dot >> ver_patch;
966 vector_data.set_version(LIBMESH_VERSION_ID(ver_major, ver_minor, ver_patch));
971 this->read_header(vector_data, version,
false,
false);
974 std::vector<libMesh::NumericVector<Real> *> bf_in(1);
976 this->read_serialized_vectors (vector_data, bf_in);
979 this->get_mesh().fix_broken_node_and_element_numbering();
987 libMesh::FunctionBase<Real>& f)
const {
989 LOG_SCOPE (
"project_vector_without_dirichlet()",
"NonlinearSystem");
991 libMesh::ConstElemRange active_local_range
992 (this->get_mesh().active_local_elements_begin(),
993 this->get_mesh().active_local_elements_end() );
995 libMesh::VectorSetAction<Real> setter(new_vector);
997 const unsigned int n_variables = this->n_vars();
999 std::vector<unsigned int> vars(n_variables);
1000 for (
unsigned int i=0; i != n_variables; ++i)
1005 libMesh::GenericProjector<libMesh::FEMFunctionWrapper<Real>, libMesh::FEMFunctionWrapper<libMesh::Gradient>,
1006 Real, libMesh::VectorSetAction<Real>> FEMProjector;
1008 libMesh::WrappedFunctor<Real> f_fem(f);
1009 libMesh::FEMFunctionWrapper<Real> fw(f_fem);
1011 #if (LIBMESH_MAJOR_VERSION == 1 && LIBMESH_MINOR_VERSION < 5) 1012 libMesh::Threads::parallel_for
1013 (active_local_range,
1014 FEMProjector(*
this, fw,
nullptr, setter, vars));
1016 FEMProjector projector(*
this, fw,
nullptr, setter, vars);
1017 projector.project(active_local_range);
1023 if (this->processor_id() == (this->n_processors()-1))
1026 libMesh::FEMContext context( *
this );
1028 const libMesh::DofMap & dof_map = this->get_dof_map();
1029 for (
unsigned int var=0; var<this->n_vars(); var++)
1030 if (this->variable(var).type().family == libMesh::SCALAR)
1035 libMesh::Elem * el =
const_cast<libMesh::Elem *
>(*(this->get_mesh().active_local_elements_begin()));
1036 context.pre_fe_reinit(*
this, el);
1038 std::vector<libMesh::dof_id_type> SCALAR_indices;
1039 dof_map.SCALAR_dof_indices (SCALAR_indices, var);
1040 const unsigned int n_SCALAR_dofs =
1041 libMesh::cast_int<unsigned int>(SCALAR_indices.size());
1043 for (
unsigned int i=0; i<n_SCALAR_dofs; i++)
1045 const libMesh::dof_id_type global_index = SCALAR_indices[i];
1046 const unsigned int component_index =
1047 this->variable_scalar_number(var,i);
1049 new_vector.set(global_index, f_fem.component(context, component_index, libMesh::Point(), this->time));
unsigned int _n_iterations
The number of iterations of the eigen solver algorithm.
virtual void eigenproblem_assemble(libMesh::SparseMatrix< Real > *A, libMesh::SparseMatrix< Real > *B)
assembles the matrices for eigenproblem depending on the analysis type
void initialize_condensed_dofs(MAST::PhysicsDisciplineBase &physics)
Loop over the dofs on each processor to initialize the list of non-condensed dofs.
void read_in_vector(libMesh::NumericVector< Real > &vec, const std::string &directory_name, const std::string &data_name, const bool read_binary_vectors)
reads the specified vector with the specified name in a directory.
This class implements a system for solution of nonlinear systems.
bool _is_generalized_eigenproblem
A boolean flag to indicate whether we are dealing with a generalized eigenvalue problem.
virtual void calculate_output_derivative(const libMesh::NumericVector< Real > &X, bool if_localize_sol, MAST::OutputAssemblyElemOperations &output, libMesh::NumericVector< Real > &dq_dX)
calculates
virtual void clear_elem_operation_object()
clears the association of this object with the assembly element operation object. ...
virtual void clear() libmesh_override
Clear all the data structures associated with the system.
void write_out_vector(libMesh::NumericVector< Real > &vec, const std::string &directory_name, const std::string &data_name, const bool write_binary_vectors)
writes the specified vector with the specified name in a directory.
unsigned int _n_requested_eigenpairs
The number of requested eigenpairs.
This provides the base class for definitin of element level contribution of output quantity in an ana...
std::unique_ptr< MAST::SlepcEigenSolver > eigen_solver
The EigenSolver, definig which interface, i.e solver package to use.
virtual void set_elem_operation_object(MAST::AssemblyElemOperations &elem_ops)
attaches a element operation to this object, and associated this with the element operation object...
bool _initialize_B_matrix
initialize the B matrix in addition to A, which might be needed for solution of complex system of equ...
This class inherits from libMesh::SlepcEigenSolver<Real> and implements a method for retriving the re...
void set_eigenproblem_type(libMesh::EigenProblemType ept)
Sets the type of the current eigen problem.
virtual void reinit() libmesh_override
Reinitializes the member data fields associated with the system, so that, e.g., assemble() may be use...
MAST::NonlinearSystem::Operation _operation
current operation of the system
virtual void eigenproblem_sensitivity_solve(MAST::AssemblyElemOperations &elem_ops, MAST::EigenproblemAssembly &assembly, const MAST::FunctionBase &f, std::vector< Real > &sens, const std::vector< unsigned int > *indices=nullptr)
Solves the sensitivity system, for the provided parameters.
virtual bool eigenproblem_sensitivity_assemble(const MAST::FunctionBase &f, libMesh::SparseMatrix< Real > *sensitivity_A, libMesh::SparseMatrix< Real > *sensitivity_B)
Assembly function.
void get_system_dirichlet_bc_dofs(libMesh::System &sys, std::set< unsigned int > &dof_ids) const
Prepares a list of the constrained dofs for system sys and returns in dof_ids.
Assembles the system of equations for an eigenproblem of type .
libMesh::SparseMatrix< Real > * matrix_A
The system matrix for standard eigenvalue problems.
void project_vector_without_dirichlet(libMesh::NumericVector< Real > &new_vector, libMesh::FunctionBase< Real > &f) const
bool _condensed_dofs_initialized
A private flag to indicate whether the condensed dofs have been initialized.
virtual void eigenproblem_solve(MAST::AssemblyElemOperations &elem_ops, MAST::EigenproblemAssembly &assembly)
Assembles & solves the eigen system.
virtual std::pair< unsigned int, Real > get_linear_solve_parameters()
calls NonlinearImplicitSystem::set_solver_parameters() before accessing the values.
virtual void init_data() libmesh_override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used...
libMesh::SparseMatrix< Real > * matrix_B
A second system matrix for generalized eigenvalue problems.
libMesh::EigenProblemType _eigen_problem_type
The type of the eigenvalue problem.
virtual void get_eigenvalue(unsigned int i, Real &re, Real &im)
gets the real and imaginary parts of the ith eigenvalue for the eigenproblem , and the associated eig...
virtual ~NonlinearSystem()
std::unique_ptr< libMesh::LinearSolver< Real > > linear_solver
The LinearSolver for solution of the linear equations.
unsigned int _n_converged_eigenpairs
The number of converged eigenpairs.
std::vector< libMesh::dof_id_type > _local_non_condensed_dofs_vector
Vector storing the local dof indices that will not be condensed.
virtual void get_eigenpair(unsigned int i, Real &re, Real &im, libMesh::NumericVector< Real > &vec_re, libMesh::NumericVector< Real > *vec_im=nullptr)
gets the real and imaginary parts of the ith eigenvalue for the eigenproblem , and the associated eig...
virtual void residual_and_jacobian(const libMesh::NumericVector< Real > &X, libMesh::NumericVector< Real > *R, libMesh::SparseMatrix< Real > *J, libMesh::NonlinearImplicitSystem &S)
function that assembles the matrices and vectors quantities for nonlinear solution ...
bool _exchange_A_and_B
flag to exchange the A and B matrices in the eigenproblem solution
virtual void sensitivity_solve(const libMesh::NumericVector< Real > &X, bool if_localize_sol, MAST::AssemblyElemOperations &elem_ops, MAST::AssemblyBase &assembly, const MAST::FunctionBase &p, bool if_assemble_jacobian=true)
Solves the sensitivity problem for the provided parameter.
NonlinearSystem(libMesh::EquationSystems &es, const std::string &name, const unsigned int number)
Default constructor.
virtual bool sensitivity_assemble(const libMesh::NumericVector< Real > &X, bool if_localize_sol, const MAST::FunctionBase &f, libMesh::NumericVector< Real > &sensitivity_rhs, bool close_vector=true)
Assembly function.
virtual void adjoint_solve(const libMesh::NumericVector< Real > &X, bool if_localize_sol, MAST::AssemblyElemOperations &elem_ops, MAST::OutputAssemblyElemOperations &output, MAST::AssemblyBase &assembly, bool if_assemble_jacobian=true)
solves the adjoint problem for the provided output function.
virtual void solve(MAST::AssemblyElemOperations &elem_ops, MAST::AssemblyBase &assembly)
solves the nonlinear problem with the specified assembly operation object