NProblemIPOPTc.cc 22.5 KB
Newer Older
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
//=============================================================================
//
//  CLASS IPOPTSolverLean - IMPLEMENTATION
//
//=============================================================================

//== INCLUDES =================================================================

//== COMPILE-TIME PACKAGE REQUIREMENTS ========================================
#include <CoMISo/Config/config.hh>
#if COMISO_IPOPT_AVAILABLE
//=============================================================================


#include "NProblemIPOPT.hh"
#include "NProblemGmmInterface.hh"
#include "NProblemInterface.hh"
#include "NConstraintInterface.hh"
#include "BoundConstraint.hh"
#include "CoMISo/Utils/CoMISoError.hh"

#include <Base/Debug/DebTime.hh>

24
#include <CoMISo/Utils/gmm.h>
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42

#include <IpTNLP.hpp>
#include <IpIpoptApplication.hpp>
#include <IpSolveStatistics.hpp>

//== NAMESPACES ===============================================================

namespace COMISO {



//== IMPLEMENTATION PROBLEM INSTANCE==========================================================


void
NProblemIPOPT::
split_constraints(const std::vector<NConstraintInterface*>& _constraints)
{
43
  DEB_enter_func;
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
  // split user-provided constraints into general-constraints and bound-constraints
  constraints_      .clear();       constraints_.reserve(_constraints.size());
  bound_constraints_.clear(); bound_constraints_.reserve(_constraints.size());

  for(unsigned int i=0; i<_constraints.size(); ++i)
  {
    BoundConstraint* bnd_ptr = dynamic_cast<BoundConstraint*>(_constraints[i]);

    if(bnd_ptr)
      bound_constraints_.push_back(bnd_ptr);
    else
      constraints_.push_back(_constraints[i]);
  }
}


//-----------------------------------------------------------------------------


void
NProblemIPOPT::
analyze_special_properties(const NProblemInterface* _problem, const std::vector<NConstraintInterface*>& _constraints)
{
  hessian_constant_ = true;
  jac_c_constant_   = true;
  jac_d_constant_   = true;

  if(!_problem->constant_hessian())
    hessian_constant_ = false;

  for(unsigned int i=0; i<_constraints.size(); ++i)
  {
    if(!_constraints[i]->constant_hessian())
      hessian_constant_ = false;

    if(!_constraints[i]->constant_gradient())
    {
      if(_constraints[i]->constraint_type() == NConstraintInterface::NC_EQUAL)
        jac_c_constant_ = false;
      else
        jac_d_constant_ = false;
    }

    // nothing else to check?
    if(!hessian_constant_ && !jac_c_constant_ && !jac_d_constant_)
      break;
  }

  //hessian of Lagrangian is only constant, if all hessians of the constraints are zero (due to lambda multipliers)
  if(!jac_c_constant_ || !jac_d_constant_)
    hessian_constant_ = false;
}


//-----------------------------------------------------------------------------


bool NProblemIPOPT::get_nlp_info(Index& n, Index& m, Index& nnz_jac_g,
                         Index& nnz_h_lag, IndexStyleEnum& index_style)
{
104
  DEB_enter_func;
105
  // number of variables
Andreas Fabri's avatar
Andreas Fabri committed
106
  n = static_cast<Index>(problem_->n_unknowns());
107
108

  // number of constraints
Andreas Fabri's avatar
Andreas Fabri committed
109
  m = static_cast<Index>(constraints_.size());
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

  // get non-zeros of hessian of lagrangian and jacobi of constraints
  nnz_jac_g = 0;
  nnz_h_lag = 0;

  // get nonzero structure
  std::vector<double> x(n);
  problem_->initial_x(P(x));

  // nonzeros in the jacobian of C_ and the hessian of the lagrangian
  SMatrixNP HP;
  SVectorNC g;
  SMatrixNC H;

  if (!hessian_approximation_)
  {
    problem_->eval_hessian(P(x), HP);

    // get nonzero structure of hessian of problem
    for(int i=0; i<HP.outerSize(); ++i)
      for (SMatrixNP::InnerIterator it(HP,i); it; ++it)
        if(it.row() >= it.col())
          ++nnz_h_lag;
  }

  // get nonzero structure of constraints
  for( int i=0; i<m; ++i)
  {
    constraints_[i]->eval_gradient(P(x),g);

140
    nnz_jac_g += static_cast<Index>(g.nonZeros());
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

    if(!hessian_approximation_)
    {
      // count lower triangular elements
      constraints_[i]->eval_hessian (P(x),H);

      SMatrixNC::iterator m_it = H.begin();
      for(; m_it != H.end(); ++m_it)
        if( m_it.row() >= m_it.col())
          ++nnz_h_lag;
    }
  }

  // We use the standard fortran index style for row/col entries
  index_style = C_STYLE;

  return true;
}


//-----------------------------------------------------------------------------


bool NProblemIPOPT::get_bounds_info(Index n, Number* x_l, Number* x_u,
                            Index m, Number* g_l, Number* g_u)
{
  DEB_enter_func;
  // check dimensions
  DEB_warning_if(( n != (Index)problem_->n_unknowns() ), 1,
    "IPOPT #unknowns != n " << n << problem_->n_unknowns() );
  DEB_warning_if(( m != (Index)constraints_.size() ), 1,
    "Warning: IPOPT #constraints != m " << m << constraints_.size() );


  // first clear all variable bounds
  for( int i=0; i<n; ++i)
  {
    // x_l[i] = Ipopt::nlp_lower_bound_inf;
    // x_u[i] = Ipopt::nlp_upper_bound_inf;

    x_l[i] = -1.0e19;
    x_u[i] =  1.0e19;
  }

  // iterate over bound constraints and set them
  for(unsigned int i=0; i<bound_constraints_.size(); ++i)
  {
    if((Index)(bound_constraints_[i]->idx()) < n)
    {
      switch(bound_constraints_[i]->constraint_type())
      {
      case NConstraintInterface::NC_LESS_EQUAL:
      {
        x_u[bound_constraints_[i]->idx()] = bound_constraints_[i]->bound();
      }break;

      case NConstraintInterface::NC_GREATER_EQUAL:
      {
        x_l[bound_constraints_[i]->idx()] = bound_constraints_[i]->bound();
      }break;

      case NConstraintInterface::NC_EQUAL:
      {
        x_l[bound_constraints_[i]->idx()] = bound_constraints_[i]->bound();
        x_u[bound_constraints_[i]->idx()] = bound_constraints_[i]->bound();
      }break;
      }
    }
    else
      DEB_warning(2, "invalid bound constraint in IPOPTSolverLean!!!")
  }

  // set bounds for constraints
  for( int i=0; i<m; ++i)
  {
    // enum ConstraintType {NC_EQUAL, NC_LESS_EQUAL, NC_GREATER_EQUAL};
    switch(constraints_[i]->constraint_type())
    {
      case NConstraintInterface::NC_EQUAL         : g_u[i] = 0.0   ; g_l[i] =  0.0   ; break;
      case NConstraintInterface::NC_LESS_EQUAL    : g_u[i] = 0.0   ; g_l[i] = -1.0e19; break;
      case NConstraintInterface::NC_GREATER_EQUAL : g_u[i] = 1.0e19; g_l[i] =  0.0   ; break;
      default                                     : g_u[i] = 1.0e19; g_l[i] = -1.0e19; break;
    }
  }

  return true;
}


//-----------------------------------------------------------------------------


bool NProblemIPOPT::get_starting_point(Index n, bool init_x, Number* x,
                               bool init_z, Number* z_L, Number* z_U,
                               Index m, bool init_lambda,
                               Number* lambda)
{
238
  DEB_enter_func;
239
240
241
242
243
244
245
246
247
248
249
250
  // get initial value of problem instance
  problem_->initial_x(x);

  return true;
}


//-----------------------------------------------------------------------------


bool NProblemIPOPT::eval_f(Index n, const Number* x, bool new_x, Number& obj_value)
{
251
  DEB_enter_func;
252
253
254
255
256
257
258
259
260
261
262
  // return the value of the objective function
  obj_value = problem_->eval_f(x);
  return true;
}


//-----------------------------------------------------------------------------


bool NProblemIPOPT::eval_grad_f(Index n, const Number* x, bool new_x, Number* grad_f)
{
263
  DEB_enter_func;
264
265
266
267
268
269
270
271
272
273
274
  problem_->eval_gradient(x, grad_f);

  return true;
}


//-----------------------------------------------------------------------------


bool NProblemIPOPT::eval_g(Index n, const Number* x, bool new_x, Index m, Number* g)
{
275
  DEB_enter_func;
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
  // evaluate all constraint functions
  for( int i=0; i<m; ++i)
    g[i] = constraints_[i]->eval_constraint(x);

  return true;
}


//-----------------------------------------------------------------------------


bool NProblemIPOPT::eval_jac_g(Index n, const Number* x, bool new_x,
                       Index m, Index nele_jac, Index* iRow, Index *jCol,
                       Number* values)
{
  DEB_enter_func;
  if (values == NULL)
  {
    // get x for evaluation (arbitrary position should be ok)
    std::vector<double> x_rnd(problem_->n_unknowns(), 0.0);

    int gi = 0;
    SVectorNC g;
    for( int i=0; i<m; ++i)
    {
      constraints_[i]->eval_gradient(&(x_rnd[0]), g);
      SVectorNC::InnerIterator v_it(g);
      for( ; v_it; ++v_it)
      {
        iRow[gi] = i;
        jCol[gi] = v_it.index();
        ++gi;
      }
    }
  }
  else
  {
    // return the values of the jacobian of the constraints

    // return the structure of the jacobian of the constraints
    // global index
    int gi = 0;
    SVectorNC g;

    for( int i=0; i<m; ++i)
    {
      constraints_[i]->eval_gradient(x, g);

      SVectorNC::InnerIterator v_it(g);

      for( ; v_it; ++v_it)
      {
        values[gi] = v_it.value();
        ++gi;
      }
    }

    DEB_warning_if((gi != nele_jac), 1,
      "number of non-zeros in Jacobian of C is incorrect: "
                << gi << " vs " << nele_jac)
  }

  return true;
}


//-----------------------------------------------------------------------------


bool NProblemIPOPT::eval_h(Index n, const Number* x, bool new_x,
                   Number obj_factor, Index m, const Number* lambda,
                   bool new_lambda, Index nele_hess, Index* iRow,
                   Index* jCol, Number* values)
{
  DEB_enter_func;
  if (values == NULL)
  {
    // return structure

    // get x for evaluation (arbitrary position should be ok)
    std::vector<double> x_rnd(problem_->n_unknowns(), 0.0);

     // global index
     int gi = 0;
     // get hessian of problem
     SMatrixNP HP;
     problem_->eval_hessian(&(x_rnd[0]), HP);

     for(int i=0; i<HP.outerSize(); ++i)
       for (SMatrixNP::InnerIterator it(HP,i); it; ++it)
       {
         // store lower triangular part only
         if(it.row() >= it.col())
         {
           //         it.value();
371
372
           iRow[gi] = static_cast<Index>(it.row());
           jCol[gi] = static_cast<Index>(it.col());
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
           ++gi;
         }
       }

    // Hessians of Constraints
    for(unsigned int j=0; j<constraints_.size(); ++j)
    {
      SMatrixNC H;
      constraints_[j]->eval_hessian(&(x_rnd[0]), H);

      SMatrixNC::iterator m_it  = H.begin();
      SMatrixNC::iterator m_end = H.end();

      for(; m_it != m_end; ++m_it)
      {
        // store lower triangular part only
        if( m_it.row() >= m_it.col())
        {
          iRow[gi] = m_it.row();
          jCol[gi] = m_it.col();
          ++gi;
        }
      }
    }

    // error check
    DEB_warning_if(( gi != nele_hess), 1,
      "number of non-zeros in Hessian of Lagrangian is incorrect while indexing: "
                << gi << " vs " << nele_hess )
  }
  else
  {
    // return values.

    // global index
    int gi = 0;
    // get hessian of problem
    SMatrixNP HP;
    problem_->eval_hessian(x, HP);

    for(int i=0; i<HP.outerSize(); ++i)
      for (SMatrixNP::InnerIterator it(HP,i); it; ++it)
      {
        // store lower triangular part only
        if(it.row() >= it.col())
        {
          values[gi] = obj_factor*it.value();
          ++gi;
        }
      }

    // Hessians of Constraints
    for(unsigned int j=0; j<constraints_.size(); ++j)
    {
      SMatrixNC H;
      constraints_[j]->eval_hessian(x, H);

      SMatrixNC::iterator m_it  = H.begin();
      SMatrixNC::iterator m_end = H.end();

      for(; m_it != m_end; ++m_it)
      {
        // store lower triangular part only
        if( m_it.row() >= m_it.col())
        {
          values[gi] = lambda[j]*(*m_it);
          ++gi;
        }
      }
    }

    // error check
    DEB_warning_if(( gi != nele_hess), 1,
      "number of non-zeros in Hessian of Lagrangian is incorrect2: "
                << gi << " vs " << nele_hess )
  }
  return true;
}


//-----------------------------------------------------------------------------


//inline double _QNT(double x)
//{
//	return double(float(x));
//}

//double _QNT(const double x)
//    {
//    // clear the 12 least significant mantissa bits to reduce noise
//    //const double fact = pow(2., 41);
//
//	const double fact = pow(2., 37);
//    int i;
//    double m = frexp(x, &i);
//    m *= fact;
//    int sgn_x = m < 0 ? -1 : 1;
//    m = sgn_x * floor(fabs(m));
//    m /= fact;
//    double xq = ldexp(m, i);
//    return xq;
//    }

double _QNT(const double x) { return x; }

void NProblemIPOPT::finalize_solution(SolverReturn status,
                              Index n, const Number* x, const Number* z_L, const Number* z_U,
                              Index m, const Number* g, const Number* lambda,
                              Number obj_value,
                              const IpoptData* ip_data,
                              IpoptCalculatedQuantities* ip_cq)
{
  DEB_enter_func;

  // problem knows what to do
489
  problem_->store_result(x);
490
491
492
493
494

  if(store_solution_)
  {
    x_.resize(n);
    for( Index i=0; i<n; ++i)
495
      x_[i] = x[i];
496
  }
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514


//  DEB_out(1, "Quantanizing the IPOPT solution\n");
//  std::vector<double> x_qnt(n);
//    for( Index i=0; i<n; ++i)
//      x_qnt[i] = _QNT(x[i]);
//
//
//
//  // problem knows what to do
//  problem_->store_result(&x_qnt[0]);
//
//  if(store_solution_)
//  {
//    x_.resize(n);
//    for( Index i=0; i<n; ++i)
//      x_[i] = x_qnt[i];
//  }
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
bool NProblemIPOPT::intermediate_callback(
  Ipopt::AlgorithmMode /*mode*/,
  Index /*iter*/, Number /*obj_value*/,
  Number /*inf_pr*/, Number /*inf_du*/,
  Number /*mu*/, Number /*d_norm*/,
  Number /*regularization_size*/,
  Number /*alpha_du*/, Number /*alpha_pr*/,
  Index /*ls_trials*/,
  const IpoptData* /*ip_data*/,
  IpoptCalculatedQuantities* /*ip_cq*/
) 
{
  PROGRESS_TICK;
  return true;
}


//-----------------------------------------------------------------------------


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
bool NProblemIPOPT::hessian_constant() const
{
  return hessian_constant_;
}


//-----------------------------------------------------------------------------


bool NProblemIPOPT::jac_c_constant() const
{
  return jac_c_constant_;
}


//-----------------------------------------------------------------------------


bool NProblemIPOPT::jac_d_constant() const
{
  return jac_d_constant_;
}


//== IMPLEMENTATION PROBLEM INSTANCE==========================================================


bool NProblemGmmIPOPT::get_nlp_info(Index& n, Index& m, Index& nnz_jac_g,
                         Index& nnz_h_lag, IndexStyleEnum& index_style)
{
571
  DEB_enter_func;
572
573
574
575
  // number of variables
  n = problem_->n_unknowns();

  // number of constraints
Andreas Fabri's avatar
Andreas Fabri committed
576
  m = static_cast<Index>(constraints_.size());
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

  // get nonzero structure
  std::vector<double> x(n);
  problem_->initial_x(&(x[0]));
  // ToDo: perturb x

  // nonzeros in the jacobian of C_ and the hessian of the lagrangian
  SMatrixNP HP;
  SVectorNC g;
  SMatrixNC H;
  problem_->eval_hessian(&(x[0]), HP);
  nnz_jac_g = 0;
  nnz_h_lag = 0;

  // clear old data
  jac_g_iRow_.clear();
  jac_g_jCol_.clear();
  h_lag_iRow_.clear();
  h_lag_jCol_.clear();

  // get non-zero structure of initial hessian
  // iterate over rows
  for( int i=0; i<n; ++i)
  {
    SVectorNP& ri = HP.row(i);

    SVectorNP_citer v_it  = gmm::vect_const_begin(ri);
    SVectorNP_citer v_end = gmm::vect_const_end  (ri);

    for(; v_it != v_end; ++v_it)
    {
      // store lower triangular part only
      if( i >= (int)v_it.index())
      {
        h_lag_iRow_.push_back(i);
Andreas Fabri's avatar
Andreas Fabri committed
612
        h_lag_jCol_.push_back(static_cast<int>(v_it.index()));
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
        ++nnz_h_lag;
      }
    }
  }


  // get nonzero structure of constraints
  for( int i=0; i<m; ++i)
  {
    constraints_[i]->eval_gradient(&(x[0]),g);
    constraints_[i]->eval_hessian (&(x[0]),H);

    // iterate over sparse vector
    SVectorNC::InnerIterator v_it(g);
    for(; v_it; ++v_it)
    {
      jac_g_iRow_.push_back(i);
      jac_g_jCol_.push_back(v_it.index());
      ++nnz_jac_g;
    }

    // iterate over superSparseMatrix
    SMatrixNC::iterator m_it  = H.begin();
    SMatrixNC::iterator m_end = H.end();
    for(; m_it != m_end; ++m_it)
      if( m_it.row() >= m_it.col())
      {
        h_lag_iRow_.push_back(m_it.row());
        h_lag_jCol_.push_back(m_it.col());
        ++nnz_h_lag;
      }
  }

  // store for error checking...
  nnz_jac_g_ = nnz_jac_g;
  nnz_h_lag_ = nnz_h_lag;

  // We use the standard fortran index style for row/col entries
  index_style = C_STYLE;

  return true;
}


//-----------------------------------------------------------------------------


bool NProblemGmmIPOPT::get_bounds_info(Index n, Number* x_l, Number* x_u,
                            Index m, Number* g_l, Number* g_u)
{
663
  DEB_enter_func;
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
  // first clear all variable bounds
  for( int i=0; i<n; ++i)
  {
    // x_l[i] = Ipopt::nlp_lower_bound_inf;
    // x_u[i] = Ipopt::nlp_upper_bound_inf;

    x_l[i] = -1.0e19;
    x_u[i] =  1.0e19;
  }

  // set bounds for constraints
  for( int i=0; i<m; ++i)
  {
    // enum ConstraintType {NC_EQUAL, NC_LESS_EQUAL, NC_GREATER_EQUAL};
    switch(constraints_[i]->constraint_type())
    {
      case NConstraintInterface::NC_EQUAL         : g_u[i] = 0.0   ; g_l[i] =  0.0   ; break;
      case NConstraintInterface::NC_LESS_EQUAL    : g_u[i] = 0.0   ; g_l[i] = -1.0e19; break;
      case NConstraintInterface::NC_GREATER_EQUAL : g_u[i] = 1.0e19; g_l[i] =  0.0   ; break;
      default                                     : g_u[i] = 1.0e19; g_l[i] = -1.0e19; break;
    }
  }

  return true;
}


//-----------------------------------------------------------------------------


bool NProblemGmmIPOPT::get_starting_point(Index n, bool init_x, Number* x,
                               bool init_z, Number* z_L, Number* z_U,
                               Index m, bool init_lambda,
                               Number* lambda)
{
699
  DEB_enter_func;
700
701
702
703
704
705
706
707
708
709
710
711
  // get initial value of problem instance
  problem_->initial_x(x);

  return true;
}


//-----------------------------------------------------------------------------


bool NProblemGmmIPOPT::eval_f(Index n, const Number* x, bool new_x, Number& obj_value)
{
712
  DEB_enter_func;
713
714
715
716
717
718
719
720
721
722
723
  // return the value of the objective function
  obj_value = problem_->eval_f(x);
  return true;
}


//-----------------------------------------------------------------------------


bool NProblemGmmIPOPT::eval_grad_f(Index n, const Number* x, bool new_x, Number* grad_f)
{
724
  DEB_enter_func;
725
726
727
728
729
730
731
732
733
734
735
  problem_->eval_gradient(x, grad_f);

  return true;
}


//-----------------------------------------------------------------------------


bool NProblemGmmIPOPT::eval_g(Index n, const Number* x, bool new_x, Index m, Number* g)
{
736
  DEB_enter_func;
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
  // evaluate all constraint functions
  for( int i=0; i<m; ++i)
    g[i] = constraints_[i]->eval_constraint(x);

  return true;
}


//-----------------------------------------------------------------------------


bool NProblemGmmIPOPT::eval_jac_g(Index n, const Number* x, bool new_x,
                       Index m, Index nele_jac, Index* iRow, Index *jCol,
                       Number* values)
{
  DEB_enter_func;
  if (values == NULL)
  {
    // return the (cached) structure of the jacobian of the constraints
    gmm::copy(jac_g_iRow_, VectorPTi(iRow, jac_g_iRow_.size()));
    gmm::copy(jac_g_jCol_, VectorPTi(jCol, jac_g_jCol_.size()));
  }
  else
  {
    // return the values of the jacobian of the constraints

    // return the structure of the jacobian of the constraints
    // global index
    int gi = 0;
    SVectorNC g;

    for( int i=0; i<m; ++i)
    {
      constraints_[i]->eval_gradient(x, g);

      SVectorNC::InnerIterator v_it(g);

      for( ; v_it; ++v_it)
      {
        if(gi < nele_jac)
          values[gi] = v_it.value();
        ++gi;
      }
    }

    DEB_warning_if(( gi != nele_jac), 1,
      "number of non-zeros in Jacobian of C is incorrect: "
       << gi << " vs " << nele_jac)
  }

  return true;
}


//-----------------------------------------------------------------------------


bool NProblemGmmIPOPT::eval_h(Index n, const Number* x, bool new_x,
                   Number obj_factor, Index m, const Number* lambda,
                   bool new_lambda, Index nele_hess, Index* iRow,
                   Index* jCol, Number* values)
{
  DEB_enter_func;
  if (values == NULL)
  {
    // return the (cached) structure of the hessian
    gmm::copy(h_lag_iRow_, VectorPTi(iRow, h_lag_iRow_.size()));
    gmm::copy(h_lag_jCol_, VectorPTi(jCol, h_lag_jCol_.size()));
  }
  else
  {
    // return values.

    // global index
    int gi = 0;

    // get hessian of problem
    problem_->eval_hessian(x, HP_);

    for( int i=0; i<n; ++i)
    {
      SVectorNP& ri = HP_.row(i);

      SVectorNP_citer v_it  = gmm::vect_const_begin(ri);
      SVectorNP_citer v_end = gmm::vect_const_end  (ri);

      for(; v_it != v_end; ++v_it)
      {
        // store lower triangular part only
        if( i >= (int)v_it.index())
        {
          if( gi < nele_hess)
            values[gi] = obj_factor*(*v_it);
          ++gi;
        }
      }
    }

    // Hessians of Constraints
    for(unsigned int j=0; j<constraints_.size(); ++j)
    {
      SMatrixNC H;
      constraints_[j]->eval_hessian(x, H);

      SMatrixNC::iterator m_it  = H.begin();
      SMatrixNC::iterator m_end = H.end();

      for(; m_it != m_end; ++m_it)
      {
        // store lower triangular part only
        if( m_it.row() >= m_it.col())
        {
          if( gi < nele_hess)
            values[gi] = lambda[j]*(*m_it);
          ++gi;
        }
      }
    }

    // error check
    DEB_warning_if(( gi != nele_hess), 1,
      "number of non-zeros in Hessian of Lagrangian is incorrect: "
        << gi << " vs " << nele_hess);
  }
  return true;
}


//-----------------------------------------------------------------------------


void NProblemGmmIPOPT::finalize_solution(SolverReturn status,
                              Index n, const Number* x, const Number* z_L, const Number* z_U,
                              Index m, const Number* g, const Number* lambda,
                              Number obj_value,
                              const IpoptData* ip_data,
                              IpoptCalculatedQuantities* ip_cq)
{
875
  DEB_enter_func;
876
877
878
879
  // problem knows what to do
  problem_->store_result(x);
}

880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
bool NProblemGmmIPOPT::intermediate_callback(
  Ipopt::AlgorithmMode /*mode*/,
  Index /*iter*/, Number /*obj_value*/,
  Number /*inf_pr*/, Number /*inf_du*/,
  Number /*mu*/, Number /*d_norm*/,
  Number /*regularization_size*/,
  Number /*alpha_du*/, Number /*alpha_pr*/,
  Index /*ls_trials*/,
  const IpoptData* /*ip_data*/,
  IpoptCalculatedQuantities* /*ip_cq*/
) 
{
  PROGRESS_TICK;
  return true;
}

896
897
898
899
900
901

//=============================================================================
} // namespace COMISO
//=============================================================================
#endif // COMISO_IPOPT_AVAILABLE
//=============================================================================