BlockIPPlatform
 All Classes Files Functions Variables
BlockIP.h
Go to the documentation of this file.
1 
6 //****************************************************************************
92 #ifndef BLOCKIP_H
93 #define BLOCKIP_H
94 
95 #include "MatrixBlockIP.h"
96 #include "StdForm.h"
97 
98 using namespace std;
99 
101 class BlockIP {
102 
103  // data representation
104 
105  public: // some types needed by users of the class
106  enum TYPE_PROBLEM {LINEAR, QUADRATIC, NONLINEAR};
107  enum TYPE_COMP_DY {CHOL_PWRS_PCG, FULL_CHOL, HYBRID_PCG, CHOL_THETA0_PCG};
108  // AUTOMATIC is NEWTON if NONLINEAR problem or PCG solver is used
109  enum TYPE_DIRECTION {NEWTON,SECOND_ORDER,AUTOMATIC};
110  enum TYPE_START_POINT {SIMPLE, QUAD_EQCONS_PROB};
111  // QUAD_REG: regularization term f(mu) * 1/2 x'Qx is added
112  // PROX_REG: proximal point regularization term f(mu) * 1/2* (x-xk)'Q(x-xk) is added
113  enum TYPE_REG {NO_REG, QUAD_REG, PROX_REG};
114  enum OUTPUT {NONE, SCREEN, FILE, BOTH};
115 
116  private:
117 
118  // constants
119  static constexpr double INF= 1.0e128;
120  static const int M_PW_PREC= 1;
121  static constexpr double SIGMA= 0.1;
122  static constexpr double RHO= 0.995;
123  static constexpr double OPTIM_GAP= 1.0e-6;
124  static constexpr double OPTIM_GAP_SAFEGUARD= 1.0e-5;
125  static constexpr double OPTIM_PFEAS= 1.0e-6;
126  static constexpr double OPTIM_DFEAS= 1.0e-6;
127  static const int OUTPUT_FREQ= 1;
128  static const int MAXITER= 200;
129  static constexpr double PCG_TOL_LIN= 1.0e-2;
130  static constexpr double PCG_TOL_QUAD= 1.0e-3;
131  static constexpr double MIN_PCGTOL= 1.0e-8;
132  static constexpr double RED_PCGTOL= 0.95;
133  static const TYPE_COMP_DY TYPE_COMP_DY_DEFAULT= CHOL_PWRS_PCG;
134  static const TYPE_DIRECTION TYPE_DIRECTION_DEFAULT= AUTOMATIC;
135  static const TYPE_START_POINT TYPE_START_POINT_DEFAULT= SIMPLE;
136  static const TYPE_REG TYPE_REG_DEFAULT= NO_REG;
137  static constexpr double FACTOR_REG_DEFAULT= 1.0e-6;
138  static const bool DEACTIVATELNK= true;
139 
140  // definition of problem, variables provided by the user and other not provided
141  // but required and computed
142  int input_problem; // -1 unchecked problem, 0 wrong problem, 1 ok problem
143 
144 //only to test
145 public:
146  int k_blocks;
147  int km;
148  int kn;
149  int m_cons;
150  int n_vars;
151  int l_link;
152  int *ini_n;
153  int *ini_m;
154  double *cost;
155  double *qcost;
156  double *lb;
157  double *ub;
158  double *lhs;
159  double *rhs;
160 
161  bool sameN;
162  bool sameL;
168 
169  void (*fobj) (int n, double x[], double &fx, double Gx[], double Hx[], void *params);
170  void *params;
171  double fx, *Gx, *Hx;
172  double dobj;
173 
174  bool inStdForm;
175 
176  string *blockNames;
177  string *varNames;
178  string *consNames;
179 
180  TYPE_PROBLEM type_objective;
182  double constantFObj;
183  int *origNVars;
184  int *origNCons;
185 
186 private:
187 //end only to test
188  bool initialized;
189 
191 
193 
195  WHO_PERMUTES whoperm;
196 
198  double inf;
199 
200  double optim_gap;
201  double optim_pfeas;
202  double optim_dfeas;
203  int maxiter;
204 
206  OUTPUT output;
207 
208  double epsmach;
209  double init_pcgtol;
210  double pcgtol;
211  double min_pcgtol;
212  double maxit_pcg;
213  double red_pcgtol;
214  int m_pw_prec;
215  int totcgit;
216  int cgit;
217  int it;
218  double est_specrad;
219  bool show_specrad;
221 
222  TYPE_COMP_DY type_comp_dy;
223  TYPE_DIRECTION type_direction;
224  TYPE_DIRECTION this_type_direction;
225  TYPE_START_POINT type_start_point;
226  TYPE_REG type_reg;
227  double factor_reg;
228 
229  double *x, *y, *z, *w, *s;
230  double *rb, *rc, *rxz, *rsw;
231  double *dxdz, *dsdw;
232  double normb, normc, normrb, normrc;
233  double alpha_p, alpha_d;
234  double rho;
235  double *dx, *dy, *dz, *dw;
236  double gap;
237  double mu, mu0;
238  double sigma;
239  double psi;
240  double *Theta;
241  double *Ax, *Aty, *r, *rhsdy, *rhspcg, *Cvaux, *Cvaux2;
242  double *r_cg, *z_cg, *p_cg;
243  double qreg;
244  double *valpha, *vbeta;
245 
250  bool *isActiveLnk;
253 
256  bool *isFreeVar;
259 
260  struct BackupLnk {
261  int lnk;
262  double rhs;
263  double **a;
264  };
265  int numBackupLnk;
266  BackupLnk *backupLnk;
267 
268  int numWithUb;
269  int numWithoutUb;
270  int *listWithoutUb;
271  int *listWithUb;
272 
273  /*// BlockIPminimize variables
274  double *b;
275  double *u;
276  double *c;
277  double *q;*/
278 
279  // member functions
280 
281  public: //interface routines
282 
283  // Constructor
284  BlockIP();
285 
286  // Create an instance with problem in standard form: 0 <= variables <= ub and constraints = rhs
287  BlockIP(TYPE_PROBLEM type_objective, double* &cost, double* &qcost, void (*fobj) (int n, double x[], double &fx, double Gx[], double Hx[], void *params), void *params, double* &ub, double* &rhs, int numBlocks, bool sameN, MatrixBlockIP N[], bool sameL, MatrixBlockIP L[]);
288 
289  // Create an instance with general problem: lb <= variables <= ub and lhs <= constraints <= rhs
290  BlockIP(TYPE_PROBLEM type_objective, double* &cost, double* &qcost, void (*fobj) (int n, double x[], double &fx, double Gx[], double Hx[], void *params), void *params, double* &lb, double* &ub, double* &lhs, double* &rhs, int numBlocks, bool sameN, MatrixBlockIP N[], bool sameL, MatrixBlockIP L[]);
291 
292  // Destructor
293  ~BlockIP();
294 
295  // Initializes class attributes
296  void initialize();
297  // Frees memory of class atrributes
298  void free_memory();
299 
300  // Create a problem in standard form: 0 <= variables <= ub and constraints = rhs
301  void create_problem(TYPE_PROBLEM type_objective, double* &cost, double* &qcost, void (*fobj) (int n, double x[], double &fx, double Gx[], double Hx[], void *params), void *params, double* &ub, double* &rhs, int numBlocks, bool sameN, MatrixBlockIP N[], bool sameL, MatrixBlockIP L[]);
302 
303  // Create a general problem: lb <= variables <= ub and lhs <= constraints <= rhs
304  void create_problem(TYPE_PROBLEM type_objective, double* &cost, double* &qcost, void (*fobj) (int n, double x[], double &fx, double Gx[], double Hx[], void *params), void *params, double* &lb, double* &ub, double* &lhs, double* &rhs, int numBlocks, bool sameN, MatrixBlockIP N[], bool sameL, MatrixBlockIP L[]);
305 
306 /* // Load a problem in standard form: 0 <= variables <= ub and constraints = rhs
307  void load_problem(double cost[], double qcost[], void (*fobj) (int n, double x[], double &fx, double Gx[], double Hx[], void *params), void *params, double ub[], double rhs[], int numBlocks, bool sameN, MatrixBlockIP N[], bool sameL, MatrixBlockIP L[], bool copy_vectors=true);
308 
309  // Load a general problem: lb <= variables <= ub and lhs <= constraints <= rhs
310  void load_problem(double cost[], double qcost[], void (*fobj) (int n, double x[], double &fx, double Gx[], double Hx[], void *params), void *params, double lb[], double ub[], double lhs[], double rhs[], int numBlocks, bool sameN, MatrixBlockIP N[], bool sameL, MatrixBlockIP L[], bool copy_vectors=true, bool convertToStd=true);*/
311 
312  // Convert a problem into standard form
313  void convert_to_standard();
314 
315  // Check the input problem and compute dimensions
316  void check_input_problem_and_compute_dimensions();
317 
318  // Check if the problems satisfies the conditions when using sameN
319  void check_sameN();
320 
321  // Problem minimization
322  int minimize();
323  // set to default values all the parameters of the minimization interior-point algorithm
324  void set_defaults_minimize();
325 
326  /* Retrieve results after optimization */
327 
328  // Return number of IP iterations
329  int get_it();
330 
331  // Return overall number of PCG iterations
332  int get_cgit();
333 
334  // Return objective function
335  double get_fobj();
336 
337  // Return dual objective function
338  double get_dobj();
339 
340  /* Change some parameters */
341  // Set infinity value to be considered
342  void set_inf(double inf= INF);
343 
344  // Get infinity value to be considered
345  double get_inf();
346 
347  // Set number of terms used as preconditioner of the power series expansion of (D-C'*B^{-1}*B)^{-1}
348  void set_m_pw_prec(int m_pw_prec=M_PW_PREC);
349 
350  // Get number of terms used as preconditioner of the power series expansion of (D-C'*B^{-1}*B)^{-1}
351  int get_m_pw_prec();
352 
353  // Set reduction of the centrality parameter at each IP iteration
354  void set_sigma(double sigma=SIGMA);
355 
356  // Get reduction of the centrality parameter at each IP iteration
357  double get_sigma();
358 
359  // Set reduction of the step-length for the primal and dual variables at each IP iteration
360  void set_rho(double rho=RHO);
361 
362  // Get reduction of the step-length for the primal and dual variables at each IP iteration
363  double get_rho();
364 
365  // Set optimality gap tolerance
366  void set_optim_gap(double optim_gap=OPTIM_GAP);
367 
368  // Get optimality gap tolerance
369  double get_optim_gap();
370 
371  // Set primal feasibility tolerance
372  void set_optim_pfeas(double optim_pfeas=OPTIM_PFEAS);
373 
374  // Get primal feasibility tolerance
375  double get_optim_pfeas();
376 
377  // Set dual feasibility tolerance
378  void set_optim_dfeas(double optim_dfeas=OPTIM_DFEAS);
379 
380  // Get dual feasibility tolerance
381  double get_optim_dfeas();
382 
383  // Set output information lines will be printed each output_freq IP iterations
384  void set_output_freq(int output_freq=OUTPUT_FREQ);
385 
386  // Get output information lines will be printed each output_freq IP iterations
387  int get_output_freq();
388 
389  // Set type of output
390  void set_output(OUTPUT output=SCREEN);
391 
392  // Get type of output
393  OUTPUT get_output();
394 
395  // Set maximum number of IP iterations
396  void set_maxiter(int maxiter=MAXITER);
397 
398  // Get maximum number of IP iterations
399  int get_maxiter();
400 
401  // Set initial tolerance for the conjugate gradient (by default PCG_TOL_LIN if linear problem, PCG_TOL_QUAD if quadratic or -convex- nonlinear)
402  void set_init_pcgtol(double init_pcgtol);
403 
404  // Get initial tolerance for the conjugate gradient (by default PCG_TOL_LIN if linear problem, PCG_TOL_QUAD if quadratic or -convex- nonlinear)
405  double get_init_pcgtol();
406 
407  // Set minimum tolerance for the conjugate gradient
408  void set_min_pcgtol(double min_pcgtol=MIN_PCGTOL);
409 
410  // Get minimum tolerance for the conjugate gradient
411  double get_min_pcgtol();
412 
413  // Set maximum number of pcg iterations (if 0, then the value will be computed by the code)
414  void set_maxit_pcg(double maxit_pcg=0);
415 
416  // Get maximum number of pcg iterations (if 0, then the value will be computed by the code)
417  double get_maxit_pcg();
418 
419  // Set reduction of pcg_tol at each IP iteration
420  void set_red_pcgtol(double red_pcgtol=RED_PCGTOL);
421 
422  // Get reduction of pcg_tol at each IP iteration
423  double get_red_pcgtol();
424 
425  // Set show_specrad at each IP iteration
426  void set_show_specrad(bool show_specrad=false);
427 
428  // Get show_specrad
429  bool get_show_specrad();
430 
431  // Set how dy direction is computed
432  void set_type_comp_dy(TYPE_COMP_DY type_comp_dy=TYPE_COMP_DY_DEFAULT);
433 
434  // Get how dy direction is computed
435  TYPE_COMP_DY get_type_comp_dy();
436 
437  // Set which direction is computed
438  void set_type_direction(TYPE_DIRECTION type_direction=TYPE_DIRECTION_DEFAULT);
439 
440  // Get which direction is computed
441  TYPE_DIRECTION get_type_direction();
442 
443  // Set how starting point is computed
444  void set_type_start_point(TYPE_START_POINT type_start_point=TYPE_START_POINT_DEFAULT);
445 
446  // Get how starting point is computed
447  TYPE_START_POINT get_type_start_point();
448 
449  // Set type of regularization
450  void set_type_reg(TYPE_REG type_reg=TYPE_REG_DEFAULT);
451 
452  // Get type of regularization
453  TYPE_REG get_type_reg();
454 
455  // Set factor of regularization
456  void set_factor_reg(double factor_reg=FACTOR_REG_DEFAULT);
457 
458  // Get factor of regularization
459  double get_factor_reg();
460 
461  // Set flag on deactivation of linking
462  void set_deactivateLnk(bool deactivateLnk=DEACTIVATELNK);
463 
464  // Get flag on deactivation of linking
465  bool get_deactivateLnk();
466 
467  // Get the number of blocks
468  int get_k_blocks();
469 
470  // Get the number of constraints without linking constraints
471  int get_km();
472 
473  // Get the number of variables without linking constraints
474  int get_kn();
475 
476  // Get the number of linking constraints
477  int get_l_link();
478 
479  // Get the number of variables of the problem to be optimized
480  int get_n_vars();
481 
482  // Get the number of variables of the problem to be optimized
483  int get_m_cons();
484 
485  // Get the number of variables of the original problem
486  int get_num_vars();
487 
488  // Get the number of constraints of the original problem
489  int get_num_cons();
490 
491  // Set who performs the permutation
492  void set_whoperm(WHO_PERMUTES whoperm_=CHOLESKY);
493 
494  // Get who performs the permutation
495  WHO_PERMUTES get_whoperm();
496 
497  // Get primal variables of problem optimized
498  const double* get_x();
499 
500  // Get dual variables of problem optimized (of equality constraints)
501  const double* get_y();
502 
503  // Get dual variable of problem optimized (of x>=0) i=1..n_vars
504  const double* get_z();
505 
506  // Get dual variable of problem optimized (of x<=u) i=1..n_vars
507  const double* get_w();
508 
509  // Get value of primal variable x(i) of problem optimized i=1..n_vars
510  double get_x(int i);
511 
512  // Get value of dual variable y(i) of problem optimized (of equality constraints) i=1..m_cons
513  double get_y(int i);
514 
515  // Get value of dual variable z(i) of problem optimized (of x>=0) i=1..n_vars
516  double get_z(int i);
517 
518  // Get value of dual variable w(i) of problem optimized (of x<=u) i=1..n_vars
519  double get_w(int i);
520 
521  // These two needed public for wrapper functions to call PCG
522  int min_PCG_Hv(int nn, double *v, double *Hv);
523  int min_PCG_Hz_eq_r(int nn, double *zz, double *rr);
524  int min_PCG_Theta0z_eq_r(int nn, double *zz, double *rr);
525 
526  // Set the constant to add to the objective function
527  void set_constant_fobj(double constant);
528 
529  // Get the constant to add to the objective function
530  double get_constant_fobj();
531 
532  // Set the names of the problem
533  void set_names(string *blockNames, string *varNames, string *consNames, bool copy_vectors=true);
534 
535  // Get the names of the problem
536  void get_names(const string* &blockNames, const string* &varNames, const string* &consNames);
537 
538  // Convert the problem to standard form (if it is not already converted) and write the problem in mps file
539  void convert_to_std_and_write_mps(const char *filename);
540 
541  // Write the problem in mps file
542  void write_mps(const char *filename);
543 
544  // Write the problem in mps file
545  void write_mps(const char *filename, TYPE_PROBLEM type_objective, double cost[], double qcost[], double lb[], double ub[], double lhs[], double rhs[], int numBlocks, bool sameN, MatrixBlockIP N[], bool sameL, MatrixBlockIP L[], string *blockNames=NULL, string *varNames=NULL, string *consNames=NULL, double constantFObj=0);
546 
547  // Create a problem from a mps file
548  void read_mps(const char *filename);
549 
550  // Write the problem in BlockIP format file
551  void write_BlockIP_format(const char *filename);
552 
553  // Create a problem from a BlockIP format file
554  void read_BlockIP_format(const char *filename);
555 
556  // Write all data related to the problem into a file
557  void write_problem(const char *filename);
558 
559  // Get the StdForm related to the problem
560  StdForm* get_std_form(bool keepStdFormAfterDelete=false);
561 
562  // Create names for blocks, variables and constraints if does not exist
563  void create_names();
564 
565  // Get the name of the variables including the block name
566  string* get_full_var_names();
567 
568  // Get the name of the constraints including the block name
569  string* get_full_cons_names();
570 
571  private:
572  void min_initializations();
573  void min_preprocess();
574  void min_normb_normc();
575  void min_Ax(double *Ax, const double *x);
576  void min_Aty(double *Aty, const double *y);
577  void min_compute_s();
578  void min_compute_mu();
579  void min_compute_sigma_psi();
580  void min_update_inactivelnk();
581  void min_KKT_residuals();
582  void min_starting_point();
583  void min_objective_functions_nonlin();
584  void min_objective_functions_linquad();
585  void eval_fobj(double x[], double &fx, double Gx[], double Hx[]);
586  void min_write_problem_information();
587  void min_write_current_point_info();
588  void min_steplengths(double &talpha_p, double &talpha_d, double tdx[], double tdz[], double tdw[]);
589  bool min_optimal();
590  void min_compute_Theta();
591  void min_compute_Theta0();
592  void min_compute_direction();
593  void min_Newton_direction();
594  void min_second_order_predictor_corrector_direction();
595  void min_update_newpoint();
596  void min_test_newpoint();
597  void min_backup_inactive(int lnk);
598  void min_restore_inactive(int lnk);
599  void min_compute_dy_cholpcg(double *dy, double *rhsdy, bool only_solve=false);
600  void min_compute_dy_fullchol(double *dy, double *rhsdy, bool only_solve=false);
601  void min_free_memory();
602  void min_check_Newton_direction_is_correct();
603  void min_check_predictor_direction_is_correct();
604  void min_check_predictor_corrector_direction_is_correct();
605  void min_solve_NThetaNt(double *v);
606  void min_Ctv(double *v, double *Ctv);
607  void min_Cv(double *v, double *Cv);
608  void min_debug_write_variables();
609  void min_debug_variables_of_inactivelnk();
610  void min_debug_variables_without_ub();
611  void min_update_qreg();
612 };
613 
614 inline void BlockIP::set_inf(double inf)
616 {
617  this->inf= inf;
618 }
619 
620 inline double BlockIP::get_inf()
622 {
623  return inf;
624 }
625 
626 inline void BlockIP::set_m_pw_prec(int m_pw_prec)
628 {
629  this->m_pw_prec= m_pw_prec;
630 }
631 
634 {
635  return m_pw_prec;
636 }
637 
638 inline void BlockIP::set_sigma(double sigma)
640 {
641  this->sigma= sigma;
642 }
643 
644 inline double BlockIP::get_sigma()
646 {
647  return sigma;
648 }
649 
650 inline void BlockIP::set_rho(double rho)
652 {
653  this->rho= rho;
654 }
655 
656 inline double BlockIP::get_rho()
658 {
659  return rho;
660 }
661 
662 inline void BlockIP::set_optim_gap(double optim_gap)
664 {
665  this->optim_gap= optim_gap;
666 }
667 
668 inline double BlockIP::get_optim_gap()
670 {
671  return optim_gap;
672 }
673 
674 inline void BlockIP::set_optim_pfeas(double optim_pfeas)
676 {
677  this->optim_pfeas= optim_pfeas;
678 }
679 
682 {
683  return optim_pfeas;
684 }
685 
686 inline void BlockIP::set_optim_dfeas(double optim_dfeas)
688 {
689  this->optim_dfeas= optim_dfeas;
690 }
691 
694 {
695  return optim_dfeas;
696 }
697 
698 inline void BlockIP::set_output_freq(int output_freq)
700 {
701  this->output_freq= output_freq;
702 }
703 
706 {
707  return output_freq;
708 }
709 
710 inline void BlockIP::set_output(OUTPUT output)
712 {
713  this->output= output;
714 }
715 
716 inline BlockIP::OUTPUT BlockIP::get_output()
718 {
719  return output;
720 }
721 
722 inline void BlockIP::set_maxiter(int maxiter)
724 {
725  this->maxiter= maxiter;
726 }
727 
730 {
731  return maxiter;
732 }
733 
734 inline void BlockIP::set_init_pcgtol(double init_pcgtol)
736 {
737  this->init_pcgtol= init_pcgtol;
738 }
739 
742 {
743  return init_pcgtol;
744 }
745 
746 inline void BlockIP::set_min_pcgtol(double min_pcgtol)
748 {
749  this->min_pcgtol= min_pcgtol;
750 }
751 
752 inline double BlockIP::get_min_pcgtol()
754 {
755  return min_pcgtol;
756 }
757 
758 inline void BlockIP::set_maxit_pcg(double maxit_pcg)
760 {
761  this->maxit_pcg= maxit_pcg;
762 }
763 
764 inline double BlockIP::get_maxit_pcg()
766 {
767  return maxit_pcg;
768 }
769 
770 inline void BlockIP::set_red_pcgtol(double red_pcgtol)
772 {
773  this->red_pcgtol= red_pcgtol;
774 }
775 
776 inline double BlockIP::get_red_pcgtol()
778 {
779  return red_pcgtol;
780 }
781 
782 inline void BlockIP::set_show_specrad(bool show_specrad)
784 {
785  this->show_specrad= show_specrad;
786 }
787 
790 {
791  return show_specrad;
792 }
793 
794 inline void BlockIP::set_type_comp_dy(TYPE_COMP_DY type_comp_dy)
796 {
797  this->type_comp_dy= type_comp_dy;
798 }
799 
800 inline BlockIP::TYPE_COMP_DY BlockIP::get_type_comp_dy()
802 {
803  return type_comp_dy;
804 }
805 
806 inline void BlockIP::set_type_direction(TYPE_DIRECTION type_direction)
808 {
809  this->type_direction= type_direction;
810 }
811 
812 inline BlockIP::TYPE_DIRECTION BlockIP::get_type_direction()
814 {
815  return type_direction;
816 }
817 
818 inline void BlockIP::set_type_start_point(TYPE_START_POINT type_start_point)
820 {
821  this->type_start_point= type_start_point;
822 }
823 
824 inline BlockIP::TYPE_START_POINT BlockIP::get_type_start_point()
826 {
827  return type_start_point;
828 }
829 
830 inline void BlockIP::set_type_reg(TYPE_REG type_reg)
832 {
833  this->type_reg= type_reg;
834 }
835 
836 inline BlockIP::TYPE_REG BlockIP::get_type_reg()
838 {
839  return type_reg;
840 }
841 
842 inline void BlockIP::set_factor_reg(double factor_reg)
844 {
845  this->factor_reg= factor_reg;
846 }
847 
848 inline double BlockIP::get_factor_reg()
850 {
851  return factor_reg;
852 }
853 
854 inline void BlockIP::set_deactivateLnk(bool deactivateLnk)
856 {
857  this->deactivateLnk= deactivateLnk;
858 }
859 
862 {
863  return deactivateLnk;
864 }
865 
868 
871 {
872  return k_blocks;
873 }
874 
875 inline int BlockIP::get_km()
877 
880 {
881  return km;
882 }
883 
884 inline int BlockIP::get_kn()
886 
889 {
890  return kn;
891 }
892 
895 
898 {
899  return l_link;
900 }
901 
904 
907 {
908  return n_vars;
909 }
910 
913 
916 {
917  return m_cons;
918 }
919 
922 
925 {
926  if (stdForm)
927  return stdForm->numOrigVars;
928  else
929  return n_vars;
930 }
931 
934 
937 {
938  if (stdForm)
939  return m_cons + stdForm->deletedRows->size();
940  else
941  return m_cons;
942 }
943 
944 inline void BlockIP::set_whoperm(WHO_PERMUTES whoperm_)
946 {
947  whoperm= whoperm_;
948 }
949 
950 inline WHO_PERMUTES BlockIP::get_whoperm()
952 {
953  return whoperm;
954 }
955 
956 inline const double* BlockIP::get_x()
958 {
959  if (stdForm)
960  return stdForm->xOrig;
961  else
962  return x;
963 }
964 
965 inline const double* BlockIP::get_y()
967 {
968  if (stdForm)
969  return stdForm->yOrig;
970  else
971  return y;
972 }
973 
974 inline const double* BlockIP::get_z()
976 {
977  if (stdForm)
978  return stdForm->zOrig;
979  else
980  return z;
981 }
982 
983 inline const double* BlockIP::get_w()
985 {
986  if (stdForm)
987  return stdForm->wOrig;
988  else
989  return w;
990 }
991 
992 inline double BlockIP::get_x(int i)
994 {
995  if (stdForm)
996  return stdForm->xOrig[i];
997  else
998  return x[i];
999 }
1000 
1001 inline double BlockIP::get_y(int i)
1003 {
1004  if (stdForm)
1005  return stdForm->yOrig[i];
1006  else
1007  return y[i];
1008 }
1009 
1010 inline double BlockIP::get_z(int i)
1012 {
1013  if (stdForm)
1014  return stdForm->zOrig[i];
1015  else
1016  return z[i];
1017 }
1018 
1019 inline double BlockIP::get_w(int i)
1021 {
1022  if (stdForm)
1023  return stdForm->wOrig[i];
1024  else
1025  return w[i];
1026 }
1027 
1028 
1029 
1030 inline int BlockIP::get_it()
1032 {
1033  return(it);
1034 }
1035 
1036 inline int BlockIP::get_cgit()
1038 {
1039  return(totcgit);
1040 }
1041 
1042 inline double BlockIP::get_fobj()
1044 {
1045  if (stdForm)
1046  return fx + stdForm->constantFObj;
1047  else
1048  return(fx);
1049 }
1050 
1051 inline double BlockIP::get_dobj()
1053 {
1054  if (stdForm)
1055  return dobj + stdForm->constantFObj;
1056  else
1057  return(dobj);
1058 }
1059 
1060 inline StdForm* BlockIP::get_std_form(bool keepStdFormAfterDelete)
1062 
1068 {
1069  this->keepStdFormAfterDelete = keepStdFormAfterDelete;
1070  return stdForm;
1071 }
1072 
1073 inline void BlockIP::set_constant_fobj(double constant)
1075 {
1076  constantFObj = constant;
1077 }
1078 
1081 {
1082  return constantFObj;
1083 }
1084 
1085 inline void BlockIP::eval_fobj(double x[], double &fx, double Gx[], double Hx[])
1086 {
1087  if (stdForm)
1088  stdForm->fobj_stdform(n_vars, x, fx, Gx, Hx, params);
1089  else
1090  fobj(n_vars, x, fx, Gx, Hx, params);
1091 }
1092 
1093 
1094 #endif //BLOCKIP_H
1095