-
Notifications
You must be signed in to change notification settings - Fork 0
/
krylov.hpp
1588 lines (1463 loc) · 67 KB
/
krylov.hpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/**
* \file krylov.hpp
* \brief templated Krylov-subspace methods
* \author Jason Hicken <[email protected]>
* \version 1.0
*/
#pragma once
#include <math.h>
#include <boost/numeric/ublas/vector.hpp>
#include <boost/numeric/ublas/vector_proxy.hpp>
#include <boost/numeric/ublas/matrix.hpp>
#include <boost/numeric/ublas/matrix_proxy.hpp>
#include <boost/scoped_ptr.hpp>
#include <boost/format.hpp>
#include <boost/property_tree/ptree.hpp>
#include <boost/tuple/tuple.hpp>
#include <ostream>
#include <iostream>
#include <iomanip>
#include <limits>
#include <string>
#include <vector>
namespace kona {
using std::cerr;
using std::cout;
using std::endl;
using std::ostream;
using std::ios;
using std::numeric_limits;
using std::string;
using std::vector;
namespace ublas = boost::numeric::ublas;
using ::boost::property_tree::ptree;
using ::boost::tuples::tuple;
using ::boost::tuples::tie;
using ::boost::tuples::make_tuple;
/// machine epsilon
const double kEpsilon = numeric_limits<double>::epsilon();
// ==============================================================================
/*!
* \class MatrixVectorProduct
* \brief abstract base class for defining matrix-vector products
* \tparam DomainVector - generic vector in the matrix domain space
* \tparam RangeVector - generic vector in the matrix range space
*
* The Krylov-subspace solvers require only matrix-vector products and
* not the actual matrix/Jacobian. We need some way to indicate which
* function will perform the product. However, sometimes the
* functions that define the product will require different numbers
* and types of inputs. For example, the forward-difference
* approximation to a Jacobian-vector product requires the vector that
* defines the Jacobian and a perturbation parameter. The
* MatrixVectorProduct class is used to derive child classes that can
* handle the different types of matrix-vector products and still be
* passed to a single implementation of the Krylov solvers.
*/
template <class DomainVector, class RangeVector = DomainVector>
class MatrixVectorProduct {
public:
virtual ~MatrixVectorProduct() {} ///< class destructor
virtual void MemoryRequired(ptree & num_required) {} ///< user vectors req.
virtual void Initialize(const ptree& prod_param) {} ///< post-constructor init
virtual void set_product_tol(const double & tol) {} ///< dynamic tolerance
virtual void operator()(const DomainVector & u,
RangeVector & v) = 0; ///< matrix-vector product
};
// ==============================================================================
/*!
* \class Preconditioner
* \brief abstract base class for defining preconditioning operation
* \tparam DomainVector - generic vector in the matrix domain space
* \tparam RangeVector - generic vector in the matrix range space
*
* See the remarks regarding the MatrixVectorProduct class. The same
* idea applies here to the preconditioning operation.
*/
template <class DomainVector, class RangeVector = DomainVector>
class Preconditioner {
public:
virtual ~Preconditioner() {} ///< class destructor
virtual void MemoryRequired(ptree & num_required) {} ///< user vectors req.
virtual void set_diagonal(const double & diag) {}; ///< add a diagonal matrix
virtual void operator()(DomainVector & u,
RangeVector & v) = 0; ///< preconditioning operation
};
// ==============================================================================
/*!
* \class IterativeSolver
* \brief abstract base class for defining iterative solvers
* \tparam DomainVector - generic vector in the matrix domain space
* \tparam RangeVector - generic vector in the matrix range space
*
* Hopefully, all Krylov solvers will transition to using this base class.
*/
template <class DomainVector, class RangeVector = DomainVector>
class IterativeSolver {
public:
virtual ~IterativeSolver() {} ///< class destructor
virtual void SubspaceSize(int m) = 0; ///< sets the solution subspace size
virtual void MemoryRequired(ptree& num_required) const = 0; ///< memory needed
virtual void Solve(const ptree & ptin, const RangeVector & b, DomainVector & x,
MatrixVectorProduct<DomainVector,RangeVector> & mat_vec,
Preconditioner<DomainVector,RangeVector> & precond,
ptree & ptout, ostream & fout = cout) = 0; ///< solve
virtual void ReSolve(const ptree & ptin, const RangeVector & b,
DomainVector & x, ptree & ptout,
ostream & fout) {} ///< resolve using subspace
};
// ==============================================================================
/*!
* \class FGMRESSolver
* \brief Flexible Generalized Minimal Residual (FGMRES) Krylov iterative method
* \tparam Vec - generic vector in the matrix domain and range
*/
template <class Vec>
class FGMRESSolver : public IterativeSolver<Vec,Vec> {
public:
/*!
* \brief default constructor
*/
FGMRESSolver();
/*!
* \brief destructor
*/
~FGMRESSolver() {}
/*!
* \brief sets the maximum number of subspace iterations, and sizes vectors
* \param[in] m - maximum number of Krylov subspace iterations
*/
void SubspaceSize(int m);
/*!
* \brief indicates how many user vectors are needed
* \param[out] num_required - stores the number of user vectors required
*
* num_required["num_vec"] stores number of Vec vectors required
*/
void MemoryRequired(ptree& num_required) const;
/*!
* \brief Solve a linear system using FGMRES
* \param[in] ptin - input parameters for solver
* \param[in] b - the rhs
* \param[in,out] x - the solution; on entry, the initial guess
* \param[in] mat_vec - object that defines matrix-vector product
* \param[in] precond - object that defines preconditioner
* \param[out] ptout - output information from solver
* \param[in] fout - object for writing the convergence history
*/
void Solve(const ptree & ptin, const Vec & b, Vec & x,
MatrixVectorProduct<Vec,Vec> & mat_vec,
Preconditioner<Vec,Vec> & precond,
ptree & ptout, ostream & fout = cout);
private:
int maxiter_; ///< maximum number of Krylov iterations
std::vector<Vec> W_; ///< vector space that spans the residual norm
std::vector<Vec> Z_; ///< vector space that the solution is built from
ublas::vector<double> g_; ///< reduced-space rhs
ublas::vector<double> y_; ///< reduced-space solution
ublas::vector<double> sn_; ///< sines for Givens rotations
ublas::vector<double> cs_; ///< cosines for Givens rotations
ublas::matrix<double> H_; ///< matrix from Arnoldi process
};
// ==============================================================================
/*!
* \class FFOMSolver
* \brief Flexible Full-Orthogonalization Method (FFOM) Krylov iterative method
* \tparam Vec - generic vector in the matrix domain and range
*/
template <class Vec>
class FFOMSolver : public IterativeSolver<Vec,Vec> {
public:
/*!
* \brief default constructor
*/
FFOMSolver();
/*!
* \brief destructor
*/
~FFOMSolver() {}
/*!
* \brief sets the maximum number of subspace iterations, and sizes vectors
* \param[in] m - maximum number of Krylov subspace iterations
*/
void SubspaceSize(int m);
/*!
* \brief indicates how many user vectors are needed
* \param[out] num_required - stores the number of user vectors required
*
* num_required["num_vec"] stores number of Vec vectors required
*/
void MemoryRequired(ptree& num_required) const;
/*!
* \brief Solve a linear system using FFOM
* \param[in] ptin - input parameters for solver
* \param[in] b - the rhs
* \param[in,out] x - the solution; on entry, the initial guess
* \param[in] mat_vec - object that defines matrix-vector product
* \param[in] precond - object that defines preconditioner
* \param[out] ptout - output information from solver
* \param[in] fout - object for writing the convergence history
*/
void Solve(const ptree & ptin, const Vec & b, Vec & x,
MatrixVectorProduct<Vec,Vec> & mat_vec,
Preconditioner<Vec,Vec> & precond,
ptree & ptout, ostream & fout = cout);
private:
int maxiter_; ///< maximum number of Krylov iterations
std::vector<Vec> W_; ///< vector space that spans the residual norm
std::vector<Vec> Z_; ///< vector space that the solution is built from
ublas::vector<double> g_; ///< reduced-space rhs
ublas::vector<double> y_; ///< reduced-space solution
ublas::vector<double> y_old_; ///< previous iteration reduced-space solution
ublas::vector<double> sn_; ///< sines for Givens rotations
ublas::vector<double> cs_; ///< cosines for Givens rotations
ublas::matrix<double> H_; ///< matrix from Arnoldi process
};
// ==============================================================================
/*!
* \class MINRESSolver
* \brief Minimum (MINRES) Krylov iterative method
* \tparam Vec - generic vector in the matrix domain and range
*
* This is an implementation of MINRES (see C. C. Paige and M. A. Saunders
* (1975), Solution of sparse indefinite systems of linear equations, SIAM
* J. Numer. Anal. 12(4), pp. 617-629). The code is based on the F90 version
* developed by the Stanford Systems Optimization Laboratory.
*/
template <class Vec>
class MINRESSolver : public IterativeSolver<Vec,Vec> {
public:
/*!
* \brief default constructor
*/
MINRESSolver();
/*!
* \brief destructor
*/
~MINRESSolver() {}
/*!
* \brief sets the maximum number of subspace iterations, and sizes vectors
* \param[in] m - maximum number of Krylov subspace iterations
*/
void SubspaceSize(int m);
/*!
* \brief indicates how many user vectors are needed
* \param[out] num_required - stores the number of user vectors required
*
* num_required["num_vec"] stores number of Vec vectors required
*/
void MemoryRequired(ptree& num_required) const;
/*!
* \brief Solve a linear system using MINRES
* \param[in] ptin - input parameters for solver
* \param[in] b - the rhs
* \param[in,out] x - the solution; on entry, the initial guess
* \param[in] mat_vec - object that defines matrix-vector product
* \param[in] precond - object that defines preconditioner
* \param[out] ptout - output information from solver
* \param[in] fout - object for writing the convergence history
*/
void Solve(const ptree & ptin, const Vec & b, Vec & x,
MatrixVectorProduct<Vec,Vec> & mat_vec,
Preconditioner<Vec,Vec> & precond,
ptree & ptout, ostream & fout = cout);
private:
int maxiter_; ///< maximum number of Krylov iterations
std::vector<Vec> work_; ///< "storage" for the user vectors needed
};
// ==============================================================================
/*!
* \class FITRSolver
* \brief Flexible Iterative Trust Region method
*/
template <class Vec>
class FITRSolver : public IterativeSolver<Vec,Vec> {
public:
/*!
* \brief default constructor
*/
FITRSolver();
/*!
* \brief class destructor
*/
~FITRSolver() {}
/*!
* \brief sets the maximum number of subspace iterations, and sizes vectors
* \param[in] m - maximum number of Krylov subspace iterations
*/
void SubspaceSize(int m);
/*!
* \brief indicates how many user vectors are needed
* \param[out] num_required - stores the number of user vectors required
*
* num_required["num_vec"] stores number of Vec vectors required
*/
void MemoryRequired(ptree& num_required) const;
/*!
* \brief solves the trust-region subproblem based on the inputs
* \param[in] inparam - input parameters for solver
* \param[in] b - the rhs (-gradient)
* \param[out] x - the solution (assumed that x = 0.0 initially)
* \param[in] mat_vec - object that defines matrix-vector product for Vec
* \param[in] precond - object that defines preconditioner for Vec
* \param[out] outparam - output information from solver
* \param[in] fout - object for writing the convergence history
*/
void Solve(const ptree & ptin, const Vec & b, Vec & x,
MatrixVectorProduct<Vec> & mat_vec, Preconditioner<Vec> & precond,
ptree & ptout, ostream & fout = cout);
/*!
* \brief recycles the subspace to resolve at a new radius
* \param[in] inparam - input parameters for solver
* \param[in] b - the rhs (-gradient)
* \param[out] x - the solution (assumed that x = 0.0 initially for this solver)
* \param[out] outparam - output information from solver
* \param[in] fout - object for writing the convergence history
*/
void ReSolve(const ptree & ptin, const Vec & b, Vec & x, ptree & ptout,
ostream & fout = cout);
private:
int maxiter_; ///< maximum number of Krylov iterations
vector<Vec> V_; ///< vector space that spans the residual norm
vector<Vec> Z_; ///< vector space that the solution is built from
boost::scoped_ptr<Vec> r_; ///< residual vector
ublas::vector<double> g_; ///< reduced-space rhs
ublas::vector<double> y_; ///< reduced-space solution
ublas::matrix<double> B_; ///< B = V^T Z
ublas::matrix<double> H_; ///< matrix from Arnoldi process
};
// ==============================================================================
/*!
* \class STCGSolver
* \brief Steihaug-Toint Conjugate Gradient (STCG) Krylov iterative method
* \tparam Vec - generic vector in the matrix domain and range
*/
template <class Vec>
class STCGSolver : public IterativeSolver<Vec,Vec> {
public:
/*!
* \brief default constructor
*/
STCGSolver();
/*!
* \brief destructor
*/
~STCGSolver() {}
/*!
* \brief sets the maximum number of subspace iterations, and sizes vectors
* \param[in] m - maximum number of Krylov subspace iterations
*/
void SubspaceSize(int m);
/*!
* \brief indicates how many user vectors are needed
* \param[out] num_required - stores the number of user vectors required
*
* num_required["num_vec"] stores number of Vec vectors required
*/
void MemoryRequired(ptree& num_required) const;
/*!
* \brief Solve a linear system using STCG
* \param[in] ptin - input parameters for solver
* \param[in] b - the rhs
* \param[in,out] x - the solution; on entry, the initial guess
* \param[in] mat_vec - object that defines matrix-vector product
* \param[in] precond - object that defines preconditioner
* \param[out] ptout - output information from solver
* \param[in] fout - object for writing the convergence history
*/
void Solve(const ptree & ptin, const Vec & b, Vec & x,
MatrixVectorProduct<Vec,Vec> & mat_vec,
Preconditioner<Vec,Vec> & precond,
ptree & ptout, ostream & fout = cout);
private:
int maxiter_; ///< maximum number of Krylov iterations
std::vector<Vec> work_; ///< "storage" for the user vectors needed
};
// ==============================================================================
/*!
* \class FFOMWithSMART
* \brief FFOM Krylov iterative method for KKT systems using SMART tests
* \tparam Vec - generic vector in the matrix domain and range
*/
template <class Vec, class PrimVec, class DualVec>
class FFOMWithSMART : public IterativeSolver<Vec,Vec> {
public:
/*!
* \brief default constructor
*/
FFOMWithSMART();
/*!
* \brief destructor
*/
~FFOMWithSMART() {}
/*!
* \brief sets the maximum number of subspace iterations, and sizes vectors
* \param[in] m - maximum number of Krylov subspace iterations
*/
void SubspaceSize(int m);
/*!
* \brief indicates how many user vectors are needed
* \param[out] num_required - stores number of vectors required
*
* num_required["num_vec"] stores number of Vec vectors required
* num_required["primal"] stores number of PrimVec vectors required
* num_required["dual"] stores number of DualVec vectors required
*/
void MemoryRequired(ptree& num_required) const;
/*!
* \brief Solve a linear system using FFOM
* \param[in] ptin - input parameters for solver
* \param[in] b - the rhs
* \param[in,out] x - the solution; on entry, x.primal() holds the obj grad!
* \param[in] mat_vec - object that defines matrix-vector product
* \param[in] precond - object that defines preconditioner
* \param[out] ptout - output information from solver
* \param[in] fout - object for writing the convergence history
*/
void Solve(const ptree & ptin, const Vec & b, Vec & x,
MatrixVectorProduct<Vec,Vec> & mat_vec,
Preconditioner<Vec,Vec> & precond,
ptree & ptout, ostream & fout = cout);
private:
int maxiter_; ///< maximum number of Krylov iterations
std::vector<Vec> V_; ///< vector space that spans the residual norm
std::vector<Vec> Z_; ///< vector space that the solution is built from
boost::scoped_ptr<Vec> res_; ///< residual, b - Ax
boost::scoped_ptr<DualVec> Ap_; ///< product of search dir. and cnstr. Jacobian
ublas::vector<double> g_; ///< reduced-space rhs
ublas::vector<double> y_; ///< reduced-space solution
ublas::vector<double> sn_; ///< sines for Givens rotations
ublas::vector<double> cs_; ///< cosines for Givens rotations
ublas::matrix<double> B_; ///< Hessenburg matrix from Arnoldi process
ublas::matrix<double> VtZ_; ///< V^T*Z
ublas::matrix<double> VtZprim_; ///< (V_primal)^T*(Z_primal)
ublas::matrix<double> VtZdual_; ///< (V_dual)^T(Z_dual)
/*!
* \brief lower and upper bounds for normal and tangential parts steps
* \param[in] p - the primal step
* \param[in] Ap - the product of the primal step with the Jacobian
* \param[in] A_norm - estimate of the norm of the Jacobian
* \param[out] upsilon - upper bound on the tangential part of the step
* \param[out] nu - lower bound on the normal part of the step
*/
void ComputeSMARTUpsilonAndNu(const PrimVec& p, const DualVec& Ap,
const double& A_norm, double& upsilon,
double& nu) const;
/*!
* \brief estimated bound for tangential part of step
* \param[in] i - current iteration/size of subspace
* \param[in] Z - subspace basis vectors
* \param[in] A - reduced KKT-matrix
* \param[in] gu - the reduced right-hand side for tangential part of step
* \param[out] work - work vector to store tangential part of step
* \param[out] upsilon - estimated bound on the tangential part of the step
*/
void ComputeSMARTUpsilon(int i, const std::vector<Vec>& Z,
const ublas::matrix<double>& A,
const ublas::vector<double>& gu,
PrimVec& work,
double& upsilon) const;
/*!
* \brief checks the Model Reduction Condition (3.12) from the SMART tests
* \param[in] grad_dot_search - (gradient)^T * (search dir)
* \param[in] pHp - the H-inner product, p^T * H * p, where H is Hessian
* \param[in] theta - 2*theta is an estimate of the smallest Rayleigh quotient
* \param[in] upsilon - upper bound on the tangential part of the step
* \param[in] sigma - merit function penalty parameter scaling
* \param[in] pi - merit function penalty parameter
* \param[in] dual_norm0 - initial norm of the constraints
* \param[in] dual_norm - norm of the linearized constraints using p
*/
bool CheckSMARTModelReduction(
const double& grad_dot_search, const double& pHp,
const double& theta, const double& upsilon, const double& sigma,
const double& pi, const double& dual_norm0, const double& dual_norm) const;
/*!
* \brief Calculates pi to satisfy model reduction when Test 2 is satisfied
* \param[in] grad_dot_search - (gradient)^T * (search dir)
* \param[in] pHp - the H-inner product, p^T * H * p, where H is Hessian
* \param[in] theta - 2*theta is an estimate of the smallest Rayleigh quotient
* \param[in] upsilon - upper bound on the tangential part of the step
* \param[in] tau - parameter that controls the growth of pi
* \param[in] dual_norm0 - initial norm of the constraints
* \param[in] dual_norm - norm of the linearized constraints using p
* \returns the new penalty parameter pi
*/
double UpdatePenaltyParameter(
const double& grad_dot_search, const double& pHp, const double& theta,
const double& upsilon, const double& tau, const double& dual_norm0,
const double& dual_norm) const;
};
// ==============================================================================
/*!
* \class FISQPSolver
* \brief Flexible Iterative Sequential Quadratic Programming method
*/
template <class Vec, class PrimVec, class DualVec>
class FISQPSolver : public IterativeSolver<Vec,Vec> {
public:
/*!
* \brief default constructor
*/
FISQPSolver();
/*!
* \brief class destructor
*/
~FISQPSolver() {}
/*!
* \brief sets the maximum number of subspace iterations, and sizes vectors
* \param[in] m - maximum number of Krylov subspace iterations
*/
void SubspaceSize(int m);
/*!
* \brief indicates how many user vectors are needed
* \param[out] num_required - stores the number of user vectors required
*
* num_required["num_vec"] stores number of Vec vectors required
* num_required["primal"] stores number of PrimVec vectors required
* num_required["dual"] stores number of DualVec vectors required
*/
void MemoryRequired(ptree& num_required) const;
/*!
* \brief solves the trust-region subproblem based on the inputs
* \param[in] ptin - input parameters for solver
* \param[in] b - the rhs (-gradient)
* \param[out] x - the solution (assumed that x = 0.0 initially)
* \param[in] mat_vec - object that defines matrix-vector product for Vec
* \param[in] precond - object that defines preconditioner for Vec
* \param[out] ptout - output information from solver
* \param[in] fout - object for writing the convergence history
*/
void Solve(const ptree & ptin, const Vec & b, Vec & x,
MatrixVectorProduct<Vec> & mat_vec, Preconditioner<Vec> & precond,
ptree & ptout, ostream & fout = cout);
/*!
* \brief resolves the trust-region problem using the existing subspace
* \param[in] ptin - input parameters for solver
* \param[in] b - the rhs (-gradient)
* \param[out] x - the solution (assumed that x = 0.0 initially)
* \param[out] ptout - output information from solver
* \param[in] fout - object for writing the convergence history
*/
void ReSolve(const ptree & ptin, const Vec & b, Vec & x, ptree & ptout,
ostream & fout);
private:
/*
* \brief solves the primal-dual and composite-step subspace problems
* \param[in] iter - present size of the subspace
* \param[in] radius - trust-region radius
* \param[in] H - upper Hessenberg matrix from Arnoldi
* \param[in] g - primal-dual subspace RHS
* \param[in] g_tang - tangential subspace RHS
* \param[in] ZtZ_prim - the primal part of the inner products Z^{T}*Z
* \param[in] VtZ - the inner products V^{T}*Z
* \param[in] VtZ_prim - the primal part of the inner products V^{T}*Z
* \param[in] VtZ_dual - the dual part of the inner products V^{T}*Z
* \param[in] VtV_dual - the dual part of the inner products V^{T}*V
* \param[out] y - subspace solution of the primal-dual problem
* \param[out] y_comp -subspace solution of the composite step problem
* \param[out] y_tang - tangential component of y_comp
* \param[out] y_norm - normal component of y_comp
* \returns (beta, gamma, gamma_comp, neg_curv, step_violate) (see below)
*
* The returned tuple holds the primal-dual residual (beta), the primal-dual
* constraint residual (gamma), the composite step constraint residual
* (gamma_comp), the null-vector quality (null_qual), a boolean for whether
* negative curvature was detected (neg_curv), and a boolean for whether the
* primal-dual step exceeded the trust radius (step_violate), the predicted
* reduction (pred), and the composite-step predicted reduction (pred_comp)
*/
boost::tuple<double, double, double, double, bool, bool, double, double>
SolveSubspaceProblems(
const int& iter, const double& radius, const ublas::matrix<double>& H,
const ublas::vector<double>& g, const ublas::vector<double>& g_tang,
const ublas::matrix<double>& ZtZ_prim, const ublas::matrix<double>& VtZ,
const ublas::matrix<double>& VtZ_prim,
const ublas::matrix<double>& VtZ_dual,
const ublas::matrix<double>& VtV_dual, ublas::vector<double>& y,
ublas::vector<double>& y_comp, ublas::vector<double>& y_tang,
ublas::vector<double>& y_norm);
/*
* \brief write header data for FISQP output
* \param[in,out] os - ostream class object for output
* \param[in] tol - target tolerance for primal-dual problem
* \param[in] res0 - the initial residual norm
* \param[in] feas0 - the initial constraint norm
*/
void WriteHeader(ostream& os, const double& tol, const double& res0,
const double& feas0);
/*
* \brief write data from a single FISQP iteration to output
* \param[in,out] os - ostream class object for output
* \param[in] iter - current iteration
* \param[in] res0 - the initial residual norm
* \param[in] feas0 - the initial constraint norm
* \param[in] res - relative residual norm of primal-dual solution
* \param[in] feas - relative constraint norm of primal-dual solution
* \param[in] feas_comp - relative constraint norm of the composite-step sol.
* \param[in] null_qual - measure of the quality of the null vector
*/
void WriteHistory(ostream& os, const int& iter, const double& res,
const double& feas, const double& feas_comp,
const double& null_qual);
int maxiter_; ///< maximum number of Krylov iterations
int iters_; ///< size of subspace used (needed for second-order correction)
std::vector<Vec> V_; ///< subspace basis for range of HZ
std::vector<Vec> Z_; ///< primal solution basis
boost::scoped_ptr<Vec> res_; ///< residual, b - Ax
ublas::vector<double> y_; /// primal solution in the reduced space
ublas::vector<double> y_old_; ///< previous iteration reduced-space solution
ublas::vector<double> sn_; ///< sines for Givens rotations
ublas::vector<double> cs_; ///< cosines for Givens rotations
ublas::vector<double> g_; /// initial residuals in the reduced space
ublas::vector<double> g_tang_; /// null-space projection problem rhs
ublas::matrix<double> H_; /// upper Hessenberg matrix
ublas::matrix<double> VtZ_; /// V^T * Z = VZprods_
ublas::matrix<double> ZtZ_prim_; // primal part of inner products Z^{T}*Z
ublas::matrix<double> VtV_dual_; // dual part of inner products V^{T}*V
ublas::matrix<double> VtZ_prim_; // primal part of inner products V^{T}*Z
ublas::matrix<double> VtZ_dual_; // dual part of inner products V^{T}*Z
};
// ==============================================================================
/*!
* \class FLECSSolver
* \brief Flexible Linearized-Equality-Constraint Subproblem Solver
*/
template <class Vec, class PrimVec, class DualVec>
class FLECSSolver : public IterativeSolver<Vec,Vec> {
public:
/*!
* \brief default constructor
*/
FLECSSolver();
/*!
* \brief class destructor
*/
~FLECSSolver() {}
/*!
* \brief sets the maximum number of subspace iterations, and sizes vectors
* \param[in] m - maximum number of Krylov subspace iterations
*/
void SubspaceSize(int m);
/*!
* \brief indicates how many user vectors are needed
* \param[out] num_required - stores the number of user vectors required
*
* num_required["num_vec"] stores number of Vec vectors required
* num_required["primal"] stores number of PrimVec vectors required
* num_required["dual"] stores number of DualVec vectors required
*/
void MemoryRequired(ptree& num_required) const;
/*!
* \brief solves the trust-region subproblem based on the inputs
* \param[in] ptin - input parameters for solver
* \param[in] b - the rhs (-gradient, -constraint)
* \param[out] x - the solution (assumed that x = 0.0 initially)
* \param[in] mat_vec - object that defines matrix-vector product for Vec
* \param[in] precond - object that defines preconditioner for Vec
* \param[out] ptout - output information from solver
* \param[in] fout - object for writing the convergence history
*/
void Solve(const ptree & ptin, const Vec & b, Vec & x,
MatrixVectorProduct<Vec> & mat_vec, Preconditioner<Vec> & precond,
ptree & ptout, ostream & fout = cout);
/*!
* \brief re-solves the trust-region problem using the existing subspace
* \param[in] ptin - input parameters for solver
* \param[in] b - the rhs (-gradient, -constraint)
* \param[out] x - the solution (assumed that x = 0.0 initially)
* \param[out] ptout - output information from solver
* \param[in] fout - object for writing the convergence history
*/
void ReSolve(const ptree & ptin, const Vec & b, Vec & x, ptree & ptout,
ostream & fout);
/*!
* \brief computes a second-order correction for the primal step
* \param[in] ptin - input parameters for solver
* \param[in] ceq - the (negative) constraint value
* \param[in] b - the rhs (-gradient, -constraint)
* \param[out] x - the solution (only x.primal() is updated)
* \param[out] ptout - output information from solver
* \param[in] fout - object for writing the convergence history
*/
void Correct2ndOrder(const ptree & ptin, const DualVec & ceq, const Vec & b,
Vec & x, ptree & ptout, ostream & fout);
private:
/*
* \brief solves the primal-dual and augmented-Lagrangian subspace problems
* \param[in] iter - present size of the subspace
* \param[in] radius - trust-region radius
* \param[in] H - upper Hessenberg matrix from Arnoldi
* \param[in] g - primal-dual subspace RHS
* \param[in] mu - augmented-Lagrangian penalty parameter
* \param[in] ZtZ_prim - the primal part of the inner products Z^{T}*Z
* \param[in] VtZ - the inner products V^{T}*Z
* \param[in] VtZ_prim - the primal part of the inner products V^{T}*Z
* \param[in] VtZ_dual - the dual part of the inner products V^{T}*Z
* \param[in] VtV_dual - the dual part of the inner products V^{T}*V
* \param[out] y - subspace solution of the primal-dual problem
* \param[out] y_aug -primal subspace solution of the composite step problem
* \param[out] y_mult - dual subspace solution of the composite step problem
* \returns (beta, gamma, beta_aug, gamma_aug, neg_curv, trust_active, pred, pred_aug) (see below)
*
* The returned tuple holds the primal-dual residual (beta), the primal-dual
* constraint residual (gamma), the primal-dual primal residual (omega), the
* aug-Lag step residual (beta_aug), the aug-Lag step constraint residual
* (gamma_aug), a boolean for whether negative curvature is likely (neg_curv),
* a boolean for whether the trust radius constraint is active (trust_active),
* the predicted reduction (pred), and augmented-Lagrangian predicted
* reduction (pred_aug)
*/
boost::tuple<double, double, double, double, double, bool, bool, double,
double>
SolveSubspaceProblems(
const int& iter, const double& radius, const ublas::matrix<double>& H,
const ublas::vector<double>& g, const double& mu,
const ublas::matrix<double>& ZtZ_prim, const ublas::matrix<double>& VtZ,
const ublas::matrix<double>& VtZ_prim,
const ublas::matrix<double>& VtZ_dual,
const ublas::matrix<double>& VtV_dual, ublas::vector<double>& y,
ublas::vector<double>& y_aug, ublas::vector<double>& y_mult);
/*
* \brief write header data for FLECS output
* \param[in,out] os - ostream class object for output
* \param[in] tol - target tolerance for primal-dual problem
* \param[in] res0 - the initial residual norm
* \param[in] grad0 - the initial gradient norm
* \param[in] feas0 - the initial constraint norm
*/
void WriteHeader(ostream& os, const double& tol, const double& res0,
const double& grad0, const double& feas0);
/*
* \brief write data from a single FLECS iteration to output
* \param[in,out] os - ostream class object for output
* \param[in] iter - current iteration
* \param[in] res - relative residual norm of primal-dual solution
* \param[in] grad - relative gradient norm of primal-dual solution
* \param[in] feas - relative constraint norm of primal-dual solution
* \param[in] feas_aug - relative constraint norm of the aug-Lag step
* \param[in] pred - predicted reduction in objective using FGMRES step
* \param[in] pred_aug - predicted reduction in objective using aug-Lag step
*/
void WriteHistory(ostream& os, const int& iter, const double& res,
const double& grad, const double& feas,
const double& feas_aug, const double& pred,
const double& pred_aug);
int maxiter_; ///< maximum number of Krylov iterations
int iters_; ///< size of subspace used (needed for second-order correction)
double mu_; ///< penalty parameter (needed for resolve)
std::vector<Vec> V_; ///< subspace basis for range of HZ
std::vector<Vec> Z_; ///< primal solution basis
boost::scoped_ptr<Vec> res_; ///< residual, b - Ax
ublas::vector<double> y_; /// primal solution in the reduced space
ublas::vector<double> g_; /// initial residuals in the reduced space
ublas::matrix<double> H_; /// upper Hessenberg matrix
ublas::matrix<double> VtZ_; /// V^T * Z = VZprods_
ublas::matrix<double> ZtZ_prim_; // primal part of inner products Z^{T}*Z
ublas::matrix<double> VtV_dual_; // dual part of inner products V^{T}*V
ublas::matrix<double> VtZ_prim_; // primal part of inner products V^{T}*Z
ublas::matrix<double> VtZ_dual_; // dual part of inner products V^{T}*Z
};
// ==============================================================================
/*!
* \class BPDSolver
* \brief Balanced Primal-Dual Solver
*/
template <class Vec, class PrimVec, class DualVec>
class BPDSolver : public IterativeSolver<Vec,Vec> {
public:
/*!
* \brief default constructor
*/
BPDSolver();
/*!
* \brief class destructor
*/
~BPDSolver() {}
/*!
* \brief sets the maximum number of subspace iterations, and sizes vectors
* \param[in] m - maximum number of Krylov subspace iterations
*/
void SubspaceSize(int m);
/*!
* \brief indicates how many user vectors are needed
* \param[out] num_required - stores the number of user vectors required
*
* num_required["num_vec"] stores number of Vec vectors required
* num_required["primal"] stores number of PrimVec vectors required
* num_required["dual"] stores number of DualVec vectors required
*/
void MemoryRequired(ptree& num_required) const;
/*!
* \brief solves a primal-dual system in a balanced way
* \param[in] ptin - input parameters for solver
* \param[in] b - the rhs (-gradient and -constraints)
* \param[out] x - the solution
* \param[in] mat_vec - object that defines primal-dual-matrix-vector product
* \param[in] precond - object that defines preconditioner
* \param[out] ptout - output information from solver
* \param[in] fout - object for writing the convergence history
*
* Given scalings provided in ptin, this solver insures that both the primal
* and the dual problems are inexactly solved (both have their residuals
* reduced equally).
*/
void Solve(const ptree & ptin, const Vec & b, Vec & x,
MatrixVectorProduct<Vec> & mat_vec, Preconditioner<Vec> & precond,
ptree & ptout, ostream & fout = cout);
private:
/*
* \brief write header data for BPD output
* \param[in,out] os - ostream class object for output
* \param[in] tol - target tolerance for primal-dual problem
* \param[in] res0 - the initial residual norm of problem
* \param[in] grad0 - the initial gradient norm
* \param[in] feas0 - the initial constraint norm
*/
void WriteHeader(ostream& os, const double& tol, const double& res0,
const double& grad0, const double& feas0);
/*
* \brief write data from a single BPD iteration to output
* \param[in,out] os - ostream class object for output
* \param[in] iter - current iteration
* \param[in] res - relative residual norm of primal-dual solution
* \param[in] grad - relative gradient norm
* \param[in] feas - relative constraint norm
*/
void WriteHistory(ostream& os, const int& iter, const double& res,
const double& grad, const double& feas);
int maxiter_; ///< maximum number of Krylov iterations
int iters_; ///< size of subspace used (needed for second-order correction)
std::vector<Vec> V_; ///< subspace basis for range of HZ
std::vector<Vec> Z_; ///< primal solution basis
ublas::vector<double> y_; /// primal solution in the reduced space
ublas::vector<double> g_; /// initial residuals in the reduced space
ublas::matrix<double> H_; /// upper Hessenberg matrix
ublas::matrix<double> VtV_prim_; // primal part of inner products V^{T}*V
ublas::matrix<double> VtV_dual_; // dual part of inner products V^{T}*V
};
// ==============================================================================
/*!
* \brief sign transfer function
* \param[in] x - value having sign prescribed
* \param[in] y - value that defined the sign
*/
double sign(const double & x, const double & y);
/*!
* \brief determines the perturbation parameter for Jacobian-free product
* \param[in] eval_at_norm - L2 norm of vector where Jacobian evaluated
* \param[in] mult_by_norm - L2 norm of vector that is being multiplied
* \returns the perturbation parameter for Jacobian-free product
*/
double CalcEpsilon(const double & eval_at_norm,
const double & mult_by_norm);
/*!
* \brief returns the eigenvalues of the symmetric part of a square matrix
* \param[in] n - size of the square matrix (actual size of A may be bigger)
* \param[in] A - matrix stored in dense format; not necessarily symmetric
* \param[out] eig - the eigenvalues in ascending order
*
* The matrix A is stored in dense format and is not assumed to be exactly
* symmetric. The eigenvalues are found by calling a LAPACK routine, which is
* given 0.5*(A^T + A) as the input matrix, not A itself.
*/
void eigenvalues(const int & n, const ublas::matrix<double> & A,
ublas::vector<double> & eig);
/*!
* \brief returns the eigenvalues and eigenvectors of a symmetric matrix
* \param[in] n - size of the square matrix (actual size of A may be bigger)
* \param[in] A - matrix stored in dense format; not necessarily symmetric
* \param[out] eig - the eigenvalues in ascending order
* \param[out] E - the eigenvectors corresponding to eig
*
* The matrix A is stored in dense format and is not assumed to be exactly
* symmetric. The eigenvalues are found by calling a LAPACK routine, which is
* given 0.5*(A^T + A) as the input matrix, not A itself.
*/
void eigenvaluesAndVectors(const int & n, const ublas::matrix<double> & A,
ublas::vector<double> & eig,
ublas::matrix<double> & E);
/*!
* \brief returns the Q-R factorization of the given matrix A
* \param[in] nrow - number of rows in matrix to consider: nrow <= A.size1()
* \param[in] ncol - number of columns in matrix to consider: ncol <= A.size2()
* \param[in] A - rectangular matrix stored in dense format
* \param[out] QR - QR decomposition (see storage details below)
*
* This routine is a front end for the LAPACK routine DGEQRF. As such, it
* follows the same storage scheme for QR. QR itself is a 1-d array. The first
* nrow*ncol part of QR stores Q and R in column-major (i.e. FORTRAN) ordering.
* R is stored in the upper triangular part of QR. Q is stored as elementary
* reflectors in the lower half of A and the factors for these reflectors are
* stored at the end of QR in QR(nrow*ncol:(nrow+1)*ncol-1). Users should not
* need to know these details, since operations with Q and R are provided by
* other routines.
*/
void factorQR(const int & nrow, const int & ncol,
const ublas::matrix<double> & A, ublas::vector<double> & QR);
/*!
* \brief solves the upper triangluar problem Rx=b, where R is from a QR-fac.
* \param[in] nrow - number of rows in the QR product
* \param[in] ncol - number of columns in the QR product (R is ncol X ncol)
* \param[in] QR - the QR factorization of some matrix
* \param[in] b - the rhs vector in the system
* \param[out] x - the solution vector
* \param[in] transpose - if true, solve R^{T}x = b
*/
void solveR(const int & nrow, const int & ncol, ublas::vector<double> & QR,
const ublas::vector<double> & b, ublas::vector<double> & x,
const bool & transpose = false);
/*!
* \brief applies the orthogonal matrix Q, from a QR factorization, to b
* \param[in] nrow - number of rows in the vector b
* \param[in] ncol - number of orthogonal vectors in Q
* \param[in] QR - the matrix Q stored as part of a QR factorization
* \param[in] b - the vector that Q is multiplying from the left
* \param[out] x - the result of the product: x = Q*b, or x = Q^T*b