forked from royshil/OpenHPE
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathtnc.c
1928 lines (1692 loc) · 47.4 KB
/
tnc.c
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
/* tnc : truncated newton bound constrained minimization
using gradient information, in C */
/*
* Copyright (c) 2002-2005, Jean-Sebastien Roy ([email protected])
*
* Permission is hereby granted, free of charge, to any person obtaining a
* copy of this software and associated documentation files (the
* "Software"), to deal in the Software without restriction, including
* without limitation the rights to use, copy, modify, merge, publish,
* distribute, sublicense, and/or sell copies of the Software, and to
* permit persons to whom the Software is furnished to do so, subject to
* the following conditions:
*
* The above copyright notice and this permission notice shall be included
* in all copies or substantial portions of the Software.
*
* THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS
* OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF
* MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT.
* IN NO EVENT SHALL THE AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY
* CLAIM, DAMAGES OR OTHER LIABILITY, WHETHER IN AN ACTION OF CONTRACT,
* TORT OR OTHERWISE, ARISING FROM, OUT OF OR IN CONNECTION WITH THE
* SOFTWARE OR THE USE OR OTHER DEALINGS IN THE SOFTWARE.
*/
/*
* This software is a C implementation of TNBC, a truncated newton minimization
* package originally developed by Stephen G. Nash in Fortran.
*
* The original source code can be found at :
* http://iris.gmu.edu/~snash/nash/software/software.html
*
* Copyright for the original TNBC fortran routines:
*
* TRUNCATED-NEWTON METHOD: SUBROUTINES
* WRITTEN BY: STEPHEN G. NASH
* SCHOOL OF INFORMATION TECHNOLOGY & ENGINEERING
* GEORGE MASON UNIVERSITY
* FAIRFAX, VA 22030
*/
/*
* Conversion into C by Elisabeth Nguyen & Jean-Sebastien Roy
* Modifications by Jean-Sebastien Roy, 2001-2002
*/
static char const rcsid[] =
"@(#) $Jeannot: tnc.c,v 1.205 2005/01/28 18:27:31 js Exp $";
static char const copyright[] =
"(c) 2002-2003, Jean-Sebastien Roy ([email protected])";
#include <math.h>
#include <stdio.h>
#include <stdlib.h>
#include "tnc.h"
typedef enum
{
TNC_FALSE = 0,
TNC_TRUE
} logical;
/*
* Return code strings
*/
char *tnc_rc_string[11] =
{
"Memory allocation failed",
"Invalid parameters (n<0)",
"Infeasible (low bound > up bound)",
"Local minima reach (|pg| ~= 0)",
"Converged (|f_n-f_(n-1)| ~= 0)",
"Converged (|x_n-x_(n-1)| ~= 0)",
"Maximum number of function evaluations reached",
"Linear search failed",
"All lower bounds are equal to the upper bounds",
"Unable to progress",
"User requested end of minimization"
};
/*
* getptc return codes
*/
typedef enum
{
GETPTC_OK = 0, /* Suitable point found */
GETPTC_EVAL = 1, /* Function evaluation required */
GETPTC_EINVAL = 2, /* Bad input values */
GETPTC_FAIL = 3 /* No suitable point found */
} getptc_rc;
/*
* linearSearch return codes
*/
typedef enum
{
LS_OK = 0, /* Suitable point found */
LS_MAXFUN = 1, /* Max. number of function evaluations reach */
LS_FAIL = 2, /* No suitable point found */
LS_USERABORT = 3, /* User requested end of minimization */
LS_ENOMEM = 4 /* Memory allocation failed */
} ls_rc;
/*
* Prototypes
*/
static tnc_rc tnc_minimize(int n, double x[], double *f, double g[],
tnc_function *function, void *state,
double xscale[], double xoffset[], double *fscale,
double low[], double up[], tnc_message messages,
int maxCGit, int maxnfeval, int *nfeval,
double eta, double stepmx, double accuracy,
double fmin, double ftol, double xtol, double pgtol, double rescale);
static getptc_rc getptcInit(double *reltol, double *abstol, double tnytol,
double eta, double rmu, double xbnd,
double *u, double *fu, double *gu, double *xmin,
double *fmin, double *gmin, double *xw, double *fw,
double *gw, double *a, double *b, double *oldf,
double *b1, double *scxbnd, double *e, double *step,
double *factor, logical *braktd, double *gtest1,
double *gtest2, double *tol);
static getptc_rc getptcIter(double big, double
rtsmll, double *reltol, double *abstol, double tnytol,
double fpresn, double xbnd,
double *u, double *fu, double *gu, double *xmin,
double *fmin, double *gmin, double *xw, double *fw,
double *gw, double *a, double *b, double *oldf,
double *b1, double *scxbnd, double *e, double *step,
double *factor, logical *braktd, double *gtest1,
double *gtest2, double *tol);
static void printCurrentIteration(int n, double f, double g[], int niter,
int nfeval, int pivot[]);
static double initialStep(double fnew, double fmin, double gtp, double smax);
static ls_rc linearSearch(int n, tnc_function *function, void *state,
double low[], double up[],
double xscale[], double xoffset[], double fscale, int pivot[],
double eta, double ftol, double xbnd,
double p[], double x[], double *f,
double *alpha, double gfull[], int maxnfeval, int *nfeval);
static int tnc_direction(double *zsol, double *diagb,
double *x, double *g, int n,
int maxCGit, int maxnfeval, int *nfeval,
logical upd1, double yksk, double yrsr,
double *sk, double *yk, double *sr, double *yr,
logical lreset, tnc_function *function, void *state,
double xscale[], double xoffset[], double fscale,
int *pivot, double accuracy,
double gnorm, double xnorm, double *low, double *up);
static double stepMax(double step, int n, double x[], double p[], int pivot[],
double low[], double up[], double xscale[], double xoffset[]);
/* Active set of constraints */
static void setConstraints(int n, double x[], int pivot[], double xscale[],
double xoffset[], double low[], double up[]);
static logical addConstraint(int n, double x[], double p[], int pivot[],
double low[], double up[], double xscale[], double xoffset[]);
static logical removeConstraint(double gtpnew, double gnorm, double pgtolfs,
double f, double fLastConstraint, double g[], int pivot[], int n);
static void project(int n, double x[], int pivot[]);
static int hessianTimesVector(double v[], double gv[], int n,
double x[], double g[], tnc_function *function, void *state,
double xscale[], double xoffset[], double fscale,
double accuracy, double xnorm, double low[], double up[]);
static int msolve(double g[], double *y, int n,
double sk[], double yk[], double diagb[], double sr[],
double yr[], logical upd1, double yksk, double yrsr,
logical lreset);
static void diagonalScaling(int n, double e[], double v[], double gv[],
double r[]);
static void ssbfgs(int n, double gamma, double sj[], double *hjv,
double hjyj[], double yjsj,
double yjhyj, double vsj, double vhyj, double hjp1v[]);
static int initPreconditioner(double diagb[], double emat[], int n,
logical lreset, double yksk, double yrsr,
double sk[], double yk[], double sr[], double yr[],
logical upd1);
/* Scaling */
static void coercex(int n, double x[], double low[], double up[]);
static void unscalex(int n, double x[], double xscale[], double xoffset[]);
static void scaleg(int n, double g[], double xscale[], double fscale);
static void scalex(int n, double x[], double xscale[], double xoffset[]);
static void projectConstants(int n, double x[], double xscale[]);
/* Machine precision */
static double mchpr1(void);
/* Special blas for incx=incy=1 */
static double ddot1(int n, double dx[], double dy[]);
static void dxpy1(int n, double dx[], double dy[]);
static void daxpy1(int n, double da, double dx[], double dy[]);
static void dcopy1(int n, double dx[], double dy[]);
static double dnrm21(int n, double dx[]);
/* additionnal blas-like functions */
static void dneg1(int n, double v[]);
/*
* This routine solves the optimization problem
*
* minimize f(x)
* x
* subject to low <= x <= up
*
* where x is a vector of n real variables. The method used is
* a truncated-newton algorithm (see "newton-type minimization via
* the lanczos algorithm" by s.g. nash (technical report 378, math.
* the lanczos method" by s.g. nash (siam j. numer. anal. 21 (1984),
* pp. 770-778). this algorithm finds a local minimum of f(x). It does
* not assume that the function f is convex (and so cannot guarantee a
* global solution), but does assume that the function is bounded below.
* it can solve problems having any number of variables, but it is
* especially useful when the number of variables (n) is large.
*
*/
extern int tnc(int n, double x[], double *f, double g[], tnc_function *function,
void *state, double low[], double up[], double scale[], double offset[],
int messages, int maxCGit, int maxnfeval, double eta, double stepmx,
double accuracy, double fmin, double ftol, double xtol, double pgtol,
double rescale, int *nfeval)
{
int rc, frc, i, nc, nfeval_local,
free_low = TNC_FALSE, free_up = TNC_FALSE,
free_g = TNC_FALSE;
double *xscale = NULL, fscale, epsmch, rteps, *xoffset = NULL;
if(nfeval==NULL)
{
/* Ignore nfeval */
nfeval = &nfeval_local;
}
*nfeval = 0;
/* Version info */
if (messages & TNC_MSG_VERS)
{
fprintf(stderr, "tnc: Version %s, %s\n",TNC_VERSION,copyright);
fprintf(stderr, "tnc: RCS ID: %s\n",rcsid);
}
/* Check for errors in the input parameters */
if (n == 0)
{
rc = TNC_CONSTANT;
goto cleanup;
}
if (n < 0)
{
rc = TNC_EINVAL;
goto cleanup;
}
/* Check bounds arrays */
if (low == NULL)
{
low = malloc(n*sizeof(*low));
if (low == NULL)
{
rc = TNC_ENOMEM;
goto cleanup;
}
free_low = TNC_TRUE;
for (i = 0 ; i < n ; i++) low[i] = -HUGE_VAL;
}
if (up == NULL)
{
up = malloc(n*sizeof(*up));
if (up == NULL)
{
rc = TNC_ENOMEM;
goto cleanup;
}
free_up = TNC_TRUE;
for (i = 0 ; i < n ; i++) up[i] = HUGE_VAL;
}
/* Coherency check */
for (i = 0 ; i < n ; i++)
{
if (low[i] > up [i])
{
rc = TNC_INFEASIBLE;
goto cleanup;
}
}
/* Coerce x into bounds */
coercex(n, x, low, up);
if (maxnfeval < 1)
{
rc = TNC_MAXFUN;
goto cleanup;
}
/* Allocate g if necessary */
if(g == NULL)
{
g = malloc(n*sizeof(*g));
if (g == NULL)
{
rc = TNC_ENOMEM;
goto cleanup;
}
free_g = TNC_TRUE;
}
/* Initial function evaluation */
frc = function(x, f, g, state);
(*nfeval) ++;
if (frc)
{
rc = TNC_USERABORT;
goto cleanup;
}
/* Constant problem ? */
for (nc = 0, i = 0 ; i < n ; i++)
if ((low[i] == up[i]) || (scale != NULL && scale[i] == 0.0))
nc ++;
if (nc == n)
{
rc = TNC_CONSTANT;
goto cleanup;
}
/* Scaling parameters */
xscale = malloc(sizeof(*xscale)*n);
if (xscale == NULL)
{
rc = TNC_ENOMEM;
goto cleanup;
}
xoffset = malloc(sizeof(*xoffset)*n);
if (xoffset == NULL)
{
rc = TNC_ENOMEM;
goto cleanup;
}
fscale = 1.0;
for (i = 0 ; i < n ; i++)
{
if (scale != NULL)
{
xscale[i] = fabs(scale[i]);
if (xscale[i] == 0.0)
xoffset[i] = low[i] = up[i] = x[i];
}
else if (low[i] != -HUGE_VAL && up[i] != HUGE_VAL)
{
xscale[i] = up[i] - low[i];
xoffset[i] = (up[i]+low[i])*0.5;
}
else
{
xscale[i] = 1.0+fabs(x[i]);
xoffset[i] = x[i];
}
if (offset != NULL)
xoffset[i] = offset[i];
}
/* Default values for parameters */
epsmch = mchpr1();
rteps = sqrt(epsmch);
if (stepmx < rteps * 10.0) stepmx = 1.0e1;
if (eta < 0.0 || eta >= 1.0) eta = 0.25;
if (rescale < 0) rescale = 1.3;
if (maxCGit < 0) /* maxCGit == 0 is valid */
{
maxCGit = n / 2;
if (maxCGit < 1) maxCGit = 1;
else if (maxCGit > 50) maxCGit = 50;
}
if (maxCGit > n) maxCGit = n;
if (accuracy <= epsmch) accuracy = rteps;
if (ftol < 0.0) ftol = accuracy;
if (pgtol < 0.0) pgtol = 1e-2 * sqrt(accuracy);
if (xtol < 0.0) xtol = rteps;
/* Optimisation */
rc = tnc_minimize(n, x, f, g, function, state,
xscale, xoffset, &fscale, low, up, messages,
maxCGit, maxnfeval, nfeval, eta, stepmx, accuracy, fmin, ftol, xtol, pgtol,
rescale);
cleanup:
if (messages & TNC_MSG_EXIT)
fprintf(stderr, "tnc: %s\n", tnc_rc_string[rc - TNC_MINRC]);
if (xscale) free(xscale);
if (free_low) free(low);
if (free_up) free(up);
if (free_g) free(g);
if (xoffset) free(xoffset);
return rc;
}
/* Coerce x into bounds */
static void coercex(int n, double x[], double low[], double up[])
{
int i;
for (i = 0 ; i < n ; i++)
{
if (x[i]<low[i]) x[i] = low[i];
else if (x[i]>up[i]) x[i] = up[i];
}
}
/* Unscale x */
static void unscalex(int n, double x[], double xscale[], double xoffset[])
{
int i;
for (i = 0 ; i < n ; i++)
x[i] = x[i]*xscale[i]+xoffset[i];
}
/* Scale x */
static void scalex(int n, double x[], double xscale[], double xoffset[])
{
int i;
for (i = 0 ; i < n ; i++)
if (xscale[i]>0.0)
x[i] = (x[i]-xoffset[i])/xscale[i];
}
/* Scale g */
static void scaleg(int n, double g[], double xscale[], double fscale)
{
int i;
for (i = 0 ; i < n ; i++)
g[i] *= xscale[i]*fscale;
}
/* Caculate the pivot vector */
static void setConstraints(int n, double x[], int pivot[], double xscale[],
double xoffset[], double low[], double up[])
{
int i;
double epsmch;
epsmch = mchpr1();
for (i = 0; i < n; i++)
{
/* tolerances should be better ajusted */
if (xscale[i] == 0.0)
{
pivot[i] = 2;
}
else
{
if (low[i] != - HUGE_VAL &&
(x[i]*xscale[i]+xoffset[i] - low[i] <= epsmch * 10.0 * (fabs(low[i]) + 1.0)))
pivot[i] = -1;
else
{
if (up[i] != HUGE_VAL &&
(x[i]*xscale[i]+xoffset[i] - up[i] >= epsmch * 10.0 * (fabs(up[i]) + 1.0)))
pivot[i] = 1;
else
pivot[i] = 0;
}
}
}
}
/*
* This routine is a bounds-constrained truncated-newton method.
* the truncated-newton method is preconditioned by a limited-memory
* quasi-newton method (this preconditioning strategy is developed
* in this routine) with a further diagonal scaling
* (see routine diagonalscaling).
*/
static tnc_rc tnc_minimize(int n, double x[],
double *f, double gfull[], tnc_function *function, void *state,
double xscale[], double xoffset[], double *fscale,
double low[], double up[], tnc_message messages,
int maxCGit, int maxnfeval, int *nfeval, double eta, double stepmx,
double accuracy, double fmin, double ftol, double xtol, double pgtol,
double rescale)
{
double fLastReset, difnew, epsmch, epsred, oldgtp,
difold, oldf, xnorm, newscale,
gnorm, ustpmax, fLastConstraint, spe, yrsr, yksk,
*temp = NULL, *sk = NULL, *yk = NULL, *diagb = NULL, *sr = NULL,
*yr = NULL, *oldg = NULL, *pk = NULL, *g = NULL;
double alpha = 0.0; /* Default unused value */
int i, icycle, niter = 0, oldnfeval, *pivot = NULL, frc;
logical lreset, newcon, upd1, remcon;
tnc_rc rc = TNC_ENOMEM; /* Default error */
/* Allocate temporary vectors */
oldg = malloc(sizeof(*oldg)*n);
if (oldg == NULL) goto cleanup;
g = malloc(sizeof(*g)*n);
if (g == NULL) goto cleanup;
temp = malloc(sizeof(*temp)*n);
if (temp == NULL) goto cleanup;
diagb = malloc(sizeof(*diagb)*n);
if (diagb == NULL) goto cleanup;
pk = malloc(sizeof(*pk)*n);
if (pk == NULL) goto cleanup;
sk = malloc(sizeof(*sk)*n);
if (sk == NULL) goto cleanup;
yk = malloc(sizeof(*yk)*n);
if (yk == NULL) goto cleanup;
sr = malloc(sizeof(*sr)*n);
if (sr == NULL) goto cleanup;
yr = malloc(sizeof(*yr)*n);
if (yr == NULL) goto cleanup;
pivot = malloc(sizeof(*pivot)*n);
if (pivot == NULL) goto cleanup;
/* Initialize variables */
epsmch = mchpr1();
difnew = 0.0;
epsred = 0.05;
upd1 = TNC_TRUE;
icycle = n - 1;
newcon = TNC_TRUE;
/* Uneeded initialisations */
lreset = TNC_FALSE;
yrsr = 0.0;
yksk = 0.0;
/* Initial scaling */
scalex(n, x, xscale, xoffset);
(*f) *= *fscale;
/* initial pivot calculation */
setConstraints(n, x, pivot, xscale, xoffset, low, up);
dcopy1(n, gfull, g);
scaleg(n, g, xscale, *fscale);
/* Test the lagrange multipliers to see if they are non-negative. */
for (i = 0; i < n; i++)
if (-pivot[i] * g[i] < 0.0)
pivot[i] = 0;
project(n, g, pivot);
/* Set initial values to other parameters */
gnorm = dnrm21(n, g);
fLastConstraint = *f; /* Value at last constraint */
fLastReset = *f; /* Value at last reset */
if (messages & TNC_MSG_ITER) fprintf(stderr,
" NIT NF F GTG\n");
if (messages & TNC_MSG_ITER) printCurrentIteration(n, *f / *fscale, gfull,
niter, *nfeval, pivot);
/* Set the diagonal of the approximate hessian to unity. */
for (i = 0; i < n; i++) diagb[i] = 1.0;
/* Start of main iterative loop */
while(TNC_TRUE)
{
/* Local minimum test */
if (dnrm21(n, g) <= pgtol * (*fscale))
{
/* |PG| == 0.0 => local minimum */
dcopy1(n, gfull, g);
project(n, g, pivot);
if (messages & TNC_MSG_INFO) fprintf(stderr,
"tnc: |pg| = %g -> local minimum\n", dnrm21(n, g) / (*fscale));
rc = TNC_LOCALMINIMUM;
break;
}
/* Terminate if more than maxnfeval evaluations have been made */
if (*nfeval >= maxnfeval)
{
rc = TNC_MAXFUN;
break;
}
/* Rescale function if necessary */
newscale = dnrm21(n, g);
if ((newscale > epsmch) && (fabs(log10(newscale)) > rescale))
{
newscale = 1.0/newscale;
*f *= newscale;
*fscale *= newscale;
gnorm *= newscale;
fLastConstraint *= newscale;
fLastReset *= newscale;
difnew *= newscale;
for (i = 0; i < n; i++) g[i] *= newscale;
for (i = 0; i < n; i++) diagb[i] = 1.0;
upd1 = TNC_TRUE;
icycle = n - 1;
newcon = TNC_TRUE;
if (messages & TNC_MSG_INFO) fprintf(stderr,
"tnc: fscale = %g\n", *fscale);
}
dcopy1(n, x, temp);
project(n, temp, pivot);
xnorm = dnrm21(n, temp);
oldnfeval = *nfeval;
/* Compute the new search direction */
frc = tnc_direction(pk, diagb, x, g, n, maxCGit, maxnfeval, nfeval,
upd1, yksk, yrsr, sk, yk, sr, yr,
lreset, function, state, xscale, xoffset, *fscale,
pivot, accuracy, gnorm, xnorm, low, up);
if (frc == -1)
{
rc = TNC_ENOMEM;
break;
}
if (frc)
{
rc = TNC_USERABORT;
break;
}
if (!newcon)
{
if (!lreset)
{
/* Compute the accumulated step and its corresponding gradient
difference. */
dxpy1(n, sk, sr);
dxpy1(n, yk, yr);
icycle++;
}
else
{
/* Initialize the sum of all the changes */
dcopy1(n, sk, sr);
dcopy1(n, yk, yr);
fLastReset = *f;
icycle = 1;
}
}
dcopy1(n, g, oldg);
oldf = *f;
oldgtp = ddot1(n, pk, g);
/* Maximum unconstrained step length */
ustpmax = stepmx / (dnrm21(n, pk) + epsmch);
/* Maximum constrained step length */
spe = stepMax(ustpmax, n, x, pk, pivot, low, up, xscale, xoffset);
if (spe > 0.0)
{
ls_rc lsrc;
/* Set the initial step length */
alpha = initialStep(*f, fmin / (*fscale), oldgtp, spe);
/* Perform the linear search */
lsrc = linearSearch(n, function, state, low, up,
xscale, xoffset, *fscale, pivot,
eta, ftol, spe, pk, x, f, &alpha, gfull, maxnfeval, nfeval);
if (lsrc == LS_ENOMEM)
{
rc = TNC_ENOMEM;
break;
}
if (lsrc == LS_USERABORT)
{
rc = TNC_USERABORT;
break;
}
if (lsrc == LS_FAIL)
{
rc = TNC_LSFAIL;
break;
}
/* If we went up to the maximum unconstrained step, increase it */
if (alpha >= 0.9 * ustpmax)
{
stepmx *= 1e2;
if (messages & TNC_MSG_INFO) fprintf(stderr,
"tnc: stepmx = %g\n", stepmx);
}
/* If we went up to the maximum constrained step,
a new constraint was encountered */
if (alpha - spe >= -epsmch * 10.0)
{
newcon = TNC_TRUE;
}
else
{
/* Break if the linear search has failed to find a lower point */
if (lsrc != LS_OK)
{
if (lsrc == LS_MAXFUN) rc = TNC_MAXFUN;
else rc = TNC_LSFAIL;
break;
}
newcon = TNC_FALSE;
}
}
else
{
/* Maximum constrained step == 0.0 => new constraint */
newcon = TNC_TRUE;
}
if (newcon)
{
if(!addConstraint(n, x, pk, pivot, low, up, xscale, xoffset))
{
if(*nfeval == oldnfeval)
{
rc = TNC_NOPROGRESS;
break;
}
}
fLastConstraint = *f;
}
niter++;
/* Set up parameters used in convergence and resetting tests */
difold = difnew;
difnew = oldf - *f;
/* If this is the first iteration of a new cycle, compute the
percentage reduction factor for the resetting test */
if (icycle == 1)
{
if (difnew > difold * 2.0) epsred += epsred;
if (difnew < difold * 0.5) epsred *= 0.5;
}
dcopy1(n, gfull, g);
scaleg(n, g, xscale, *fscale);
dcopy1(n, g, temp);
project(n, temp, pivot);
gnorm = dnrm21(n, temp);
/* Reset pivot */
remcon = removeConstraint(oldgtp, gnorm, pgtol * (*fscale), *f,
fLastConstraint, g, pivot, n);
/* If a constraint is removed */
if (remcon)
{
/* Recalculate gnorm and reset fLastConstraint */
dcopy1(n, g, temp);
project(n, temp, pivot);
gnorm = dnrm21(n, temp);
fLastConstraint = *f;
}
if (!remcon && !newcon)
{
/* No constraint removed & no new constraint : tests for convergence */
if (fabs(difnew) <= ftol * (*fscale))
{
if (messages & TNC_MSG_INFO) fprintf(stderr,
"tnc: |fn-fn-1] = %g -> convergence\n", fabs(difnew) / (*fscale));
rc = TNC_FCONVERGED;
break;
}
if (alpha * dnrm21(n, pk) <= xtol)
{
if (messages & TNC_MSG_INFO) fprintf(stderr,
"tnc: |xn-xn-1] = %g -> convergence\n", alpha * dnrm21(n, pk));
rc = TNC_XCONVERGED;
break;
}
}
project(n, g, pivot);
if (messages & TNC_MSG_ITER) printCurrentIteration(n, *f / *fscale, gfull,
niter, *nfeval, pivot);
/* Compute the change in the iterates and the corresponding change in the
gradients */
if (!newcon)
{
for (i = 0; i < n; i++)
{
yk[i] = g[i] - oldg[i];
sk[i] = alpha * pk[i];
}
/* Set up parameters used in updating the preconditioning strategy */
yksk = ddot1(n, yk, sk);
if (icycle == (n - 1) || difnew < epsred * (fLastReset - *f))
lreset = TNC_TRUE;
else
{
yrsr = ddot1(n, yr, sr);
if (yrsr <= 0.0) lreset = TNC_TRUE;
else lreset = TNC_FALSE;
}
upd1 = TNC_FALSE;
}
}
if (messages & TNC_MSG_ITER) printCurrentIteration(n, *f / *fscale, gfull,
niter, *nfeval, pivot);
/* Unscaling */
unscalex(n, x, xscale, xoffset);
coercex(n, x, low, up);
(*f) /= *fscale;
cleanup:
if (oldg) free(oldg);
if (g) free(g);
if (temp) free(temp);
if (diagb) free(diagb);
if (pk) free(pk);
if (sk) free(sk);
if (yk) free(yk);
if (sr) free(sr);
if (yr) free(yr);
if (pivot) free(pivot);
return rc;
}
/* Print the results of the current iteration */
static void printCurrentIteration(int n, double f, double g[], int niter,
int nfeval, int pivot[])
{
int i;
double gtg;
gtg = 0.0;
for (i = 0; i < n; i++)
if (pivot[i] == 0)
gtg += g[i] * g[i];
fprintf(stderr, " %4d %4d %22.15E %15.8E\n", niter, nfeval, f, gtg);
}
/*
* Set x[i] = 0.0 if direction i is currently constrained
*/
static void project(int n, double x[], int pivot[])
{
int i;
for (i = 0; i < n; i++)
if (pivot[i] != 0)
x[i] = 0.0;
}
/*
* Set x[i] = 0.0 if direction i is constant
*/
static void projectConstants(int n, double x[], double xscale[])
{
int i;
for (i = 0; i < n; i++)
if (xscale[i] == 0.0)
x[i] = 0.0;
}
/*
* Compute the maximum allowable step length
*/
static double stepMax(double step, int n, double x[], double dir[],
int pivot[], double low[], double up[], double xscale[], double xoffset[])
{
int i;
double t;
/* Constrained maximum step */
for (i = 0; i < n; i++)
{
if ((pivot[i] == 0) && (dir[i] != 0.0))
{
if (dir[i] < 0.0)
{
t = (low[i]-xoffset[i])/xscale[i] - x[i];
if (t > step * dir[i]) step = t / dir[i];
}
else
{
t = (up[i]-xoffset[i])/xscale[i] - x[i];
if (t < step * dir[i]) step = t / dir[i];
}
}
}
return step;
}
/*
* Update the constraint vector pivot if a new constraint is encountered
*/
static logical addConstraint(int n, double x[], double p[], int pivot[],
double low[], double up[], double xscale[], double xoffset[])
{
int i, newcon = TNC_FALSE;
double tol, epsmch;
epsmch = mchpr1();
for (i = 0; i < n; i++)
{
if ((pivot[i] == 0) && (p[i] != 0.0))
{
if (p[i] < 0.0 && low[i] != - HUGE_VAL)
{
tol = epsmch * 10.0 * (fabs(low[i]) + 1.0);
if (x[i]*xscale[i]+xoffset[i] - low[i] <= tol)
{
pivot[i] = -1;
x[i] = (low[i]-xoffset[i])/xscale[i];
newcon = TNC_TRUE;
}
}
else if (up[i] != HUGE_VAL)
{
tol = epsmch * 10.0 * (fabs(up[i]) + 1.0);
if (up[i] - (x[i]*xscale[i]+xoffset[i]) <= tol)
{
pivot[i] = 1;
x[i] = (up[i]-xoffset[i])/xscale[i];
newcon = TNC_TRUE;
}
}
}
}
return newcon;
}
/*
* Check if a constraint is no more active
*/
static logical removeConstraint(double gtpnew, double gnorm, double pgtolfs,
double f, double fLastConstraint, double g[], int pivot[], int n)
{
double cmax, t;
int imax, i;
if (((fLastConstraint - f) <= (gtpnew * -0.5)) && (gnorm > pgtolfs))
return TNC_FALSE;
imax = -1;
cmax = 0.0;
for (i = 0; i < n; i++)
{
if (pivot[i] == 2)
continue;
t = -pivot[i] * g[i];
if (t < cmax)
{
cmax = t;
imax = i;