44 #include "MoochoPack_ReducedHessianSerialization_Step.hpp"
45 #include "MoochoPack_moocho_algo_conversion.hpp"
46 #include "MoochoPack_Exceptions.hpp"
47 #include "IterationPack_print_algorithm_step.hpp"
48 #include "AbstractLinAlgPack_VectorSpace.hpp"
49 #include "AbstractLinAlgPack_VectorStdOps.hpp"
50 #include "AbstractLinAlgPack_VectorMutable.hpp"
51 #include "AbstractLinAlgPack_VectorOut.hpp"
52 #include "AbstractLinAlgPack_MatrixSymOpNonsing.hpp"
53 #include "AbstractLinAlgPack_MatrixSymInitDiag.hpp"
54 #include "AbstractLinAlgPack_MatrixOpOut.hpp"
55 #include "AbstractLinAlgPack_LinAlgOpPack.hpp"
56 #include "SerializationPack_Serializable.hpp"
57 #include "Teuchos_dyn_cast.hpp"
59 namespace MoochoPack {
62 const std::string &reduced_hessian_input_file_name
63 ,
const std::string &reduced_hessian_output_file_name
65 :reduced_hessian_input_file_name_(reduced_hessian_input_file_name)
66 ,reduced_hessian_output_file_name_(reduced_hessian_output_file_name)
79 NLPAlgo &algo = rsqp_algo(_algo);
82 EJournalOutputLevel olevel = algo.algo_cntr().journal_output_level();
83 EJournalOutputLevel ns_olevel = algo.algo_cntr().null_space_journal_output_level();
87 if( static_cast<int>(olevel) >= static_cast<int>(PRINT_ALGORITHM_STEPS) ) {
88 using IterationPack::print_algorithm_step;
89 print_algorithm_step( _algo, step_poss, type, assoc_step_poss, out );
92 IterQuantityAccess<MatrixSymOp> &rHL_iq = s.rHL();
94 if( !rHL_iq.updated_k(0) && reduced_hessian_input_file_name().length() ) {
95 int k_last_offset = rHL_iq.last_updated();
96 if( k_last_offset == IterQuantity::NONE_UPDATED ) {
97 if( (
int)olevel >= (
int)PRINT_ALGORITHM_STEPS ) {
99 <<
"\nNo previous matrix rHL was found!\n"
100 <<
"\nReading in the matrix rHL_k from the file \""<<reduced_hessian_input_file_name()<<
"\" ...\n";
102 MatrixSymOp &rHL_k = rHL_iq.set_k(0);
103 Serializable &rHL_serializable =
dyn_cast<Serializable>(rHL_k);
104 std::ifstream reduced_hessian_input_file(reduced_hessian_input_file_name().c_str());
106 !reduced_hessian_input_file, std::logic_error
107 ,
"ReducedHessianSerialization_Step::do_step(...): Error, the file \""<<reduced_hessian_input_file_name()<<
"\""
108 " could not be opened or contains no input!"
110 rHL_serializable.unserialize(reduced_hessian_input_file);
111 if( (
int)ns_olevel >= (
int)PRINT_ALGORITHM_STEPS ) {
112 out <<
"\nrHL_k.rows() = " << rHL_k.rows() << std::endl;
113 out <<
"\nrHL_k.cols() = " << rHL_k.cols() << std::endl;
114 if(algo.algo_cntr().calc_matrix_norms())
115 out <<
"\n||rHL_k||inf = " << rHL_k.calc_norm(MatrixOp::MAT_NORM_INF).value << std::endl;
116 if(algo.algo_cntr().calc_conditioning()) {
117 const MatrixSymOpNonsing
118 *rHL_ns_k =
dynamic_cast<const MatrixSymOpNonsing*
>(&rHL_k);
120 out <<
"\ncond_inf(rHL_k) = " << rHL_ns_k->calc_cond_num(MatrixOp::MAT_NORM_INF).value << std::endl;
123 if( (
int)ns_olevel >= (
int)PRINT_ITERATION_QUANTITIES )
124 out <<
"\nrHL_k = \n" << rHL_k;
126 const MatrixOp &Z_k = s.Z().get_k(0);
127 const VectorSpace &null_space = Z_k.space_rows();
129 !null_space.is_compatible(rHL_k.space_cols()) || !null_space.is_compatible(rHL_k.space_rows())
131 ,
"ReducedHessianSerialization_Step::do_step(...): Error, the read-in reduced Hessian of dimension "
132 << rHL_k.rows() <<
" x " << rHL_k.cols() <<
" is not compatible with the null space of dimension "
133 "Z_k.cols() = " << Z_k.cols() <<
"!"
150 const NLPAlgo &algo = rsqp_algo(_algo);
153 EJournalOutputLevel olevel = algo.algo_cntr().journal_output_level();
156 const IterQuantityAccess<MatrixSymOp> &rHL_iq = s.rHL();
158 int k_last_offset = rHL_iq.last_updated();
160 if( k_last_offset != IterQuantity::NONE_UPDATED && reduced_hessian_output_file_name().length() ) {
161 if( (
int)olevel >= (
int)PRINT_ALGORITHM_STEPS ) {
163 <<
"\nSerializing the matrix rHL_k("<<k_last_offset<<
") to the file "
164 <<
"\""<<reduced_hessian_output_file_name()<<
"\" ...\n";
166 const MatrixSymOp &rHL_k = rHL_iq.get_k(k_last_offset);
167 const Serializable &rHL_serializable =
dyn_cast<
const Serializable>(rHL_k);
168 std::ofstream reduced_hessian_output_file(reduced_hessian_output_file_name().c_str());
170 !reduced_hessian_output_file, std::logic_error
171 ,
"ReducedHessianSerialization_Step::finalize_step(...): Error, the file \""<<reduced_hessian_output_file_name()<<
"\""
172 " could not be opened!"
174 rHL_serializable.serialize(reduced_hessian_output_file);
180 ,std::ostream& out,
const std::string& L
184 << L <<
"*** Read in the reduced Hessian of the Lagrangian rHL from a file.\n"
185 << L <<
"if (rHL_k is not updated and reduced_hessian_input_file_name != \"\") then\n"
186 << L <<
" k_last_offset = last iteration rHL was updated for\n"
187 << L <<
" if k_last_offset==NONE_UPDATED then\n"
188 << L <<
" rHL_serializable = dyn_cast<Serializable>(rHL_k) *** Throws exception if fails!\n"
189 << L <<
" Unserialize into rHL_serializable from the file \""<<reduced_hessian_input_file_name()<<
"\"\n"
191 << L <<
" *** There is some reduced Hessian that exists for some past iteration so\n"
192 << L <<
" *** we will let some other step object initialize itQ\n"
195 << L <<
"*** Note: On finalization, this step object will serialize rHL_k to the file:\n"
196 << L <<
"*** \""<<reduced_hessian_output_file_name()<<
"\""
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
T_To & dyn_cast(T_From &from)
rSQP Algorithm control class.
void finalize_step(Algorithm &algo, poss_type step_poss, IterationPack::EDoStepType type, poss_type assoc_step_poss)
ReducedHessianSerialization_Step(const std::string &reduced_hessian_input_file_name="reduced_hessian.in", const std::string &reduced_hessian_output_file_name="reduced_hessian.out")
bool do_step(Algorithm &algo, poss_type step_poss, IterationPack::EDoStepType type, poss_type assoc_step_poss)
virtual std::ostream & journal_out() const
Reduced space SQP state encapsulation interface.
AlgorithmTracker & track()
void print_step(const Algorithm &algo, poss_type step_poss, IterationPack::EDoStepType type, poss_type assoc_step_poss, std::ostream &out, const std::string &leading_str) const
NLPAlgoState & rsqp_state()
<<std aggr>="">> members for algo_cntr