trstlp.f 14 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438
  1. SUBROUTINE TRSTLP (N,M,A,B,RHO,DX,IFULL,IACT,Z,ZDOTA,VMULTC,
  2. 1 SDIRN,DXNEW,VMULTD)
  3. IMPLICIT DOUBLE PRECISION (A-H,O-Z)
  4. DIMENSION A(N,*),B(*),DX(*),IACT(*),Z(N,*),ZDOTA(*),
  5. 1 VMULTC(*),SDIRN(*),DXNEW(*),VMULTD(*)
  6. C
  7. C This subroutine calculates an N-component vector DX by applying the
  8. C following two stages. In the first stage, DX is set to the shortest
  9. C vector that minimizes the greatest violation of the constraints
  10. C A(1,K)*DX(1)+A(2,K)*DX(2)+...+A(N,K)*DX(N) .GE. B(K), K=2,3,...,M,
  11. C subject to the Euclidean length of DX being at most RHO. If its length is
  12. C strictly less than RHO, then we use the resultant freedom in DX to
  13. C minimize the objective function
  14. C -A(1,M+1)*DX(1)-A(2,M+1)*DX(2)-...-A(N,M+1)*DX(N)
  15. C subject to no increase in any greatest constraint violation. This
  16. C notation allows the gradient of the objective function to be regarded as
  17. C the gradient of a constraint. Therefore the two stages are distinguished
  18. C by MCON .EQ. M and MCON .GT. M respectively. It is possible that a
  19. C degeneracy may prevent DX from attaining the target length RHO. Then the
  20. C value IFULL=0 would be set, but usually IFULL=1 on return.
  21. C
  22. C In general NACT is the number of constraints in the active set and
  23. C IACT(1),...,IACT(NACT) are their indices, while the remainder of IACT
  24. C contains a permutation of the remaining constraint indices. Further, Z is
  25. C an orthogonal matrix whose first NACT columns can be regarded as the
  26. C result of Gram-Schmidt applied to the active constraint gradients. For
  27. C J=1,2,...,NACT, the number ZDOTA(J) is the scalar product of the J-th
  28. C column of Z with the gradient of the J-th active constraint. DX is the
  29. C current vector of variables and here the residuals of the active
  30. C constraints should be zero. Further, the active constraints have
  31. C nonnegative Lagrange multipliers that are held at the beginning of
  32. C VMULTC. The remainder of this vector holds the residuals of the inactive
  33. C constraints at DX, the ordering of the components of VMULTC being in
  34. C agreement with the permutation of the indices of the constraints that is
  35. C in IACT. All these residuals are nonnegative, which is achieved by the
  36. C shift RESMAX that makes the least residual zero.
  37. C
  38. C Initialize Z and some other variables. The value of RESMAX will be
  39. C appropriate to DX=0, while ICON will be the index of a most violated
  40. C constraint if RESMAX is positive. Usually during the first stage the
  41. C vector SDIRN gives a search direction that reduces all the active
  42. C constraint violations by one simultaneously.
  43. C
  44. IFULL=1
  45. MCON=M
  46. NACT=0
  47. RESMAX=0.0d0
  48. DO 20 I=1,N
  49. DO 10 J=1,N
  50. 10 Z(I,J)=0.0d0
  51. Z(I,I)=1.0d0
  52. 20 DX(I)=0.0d0
  53. IF (M .GE. 1) THEN
  54. DO 30 K=1,M
  55. IF (B(K) .GT. RESMAX) THEN
  56. RESMAX=B(K)
  57. ICON=K
  58. END IF
  59. 30 CONTINUE
  60. DO 40 K=1,M
  61. IACT(K)=K
  62. 40 VMULTC(K)=RESMAX-B(K)
  63. END IF
  64. IF (RESMAX .EQ. 0.0d0) GOTO 480
  65. DO 50 I=1,N
  66. 50 SDIRN(I)=0.0d0
  67. C
  68. C End the current stage of the calculation if 3 consecutive iterations
  69. C have either failed to reduce the best calculated value of the objective
  70. C function or to increase the number of active constraints since the best
  71. C value was calculated. This strategy prevents cycling, but there is a
  72. C remote possibility that it will cause premature termination.
  73. C
  74. 60 OPTOLD=0.0d0
  75. ICOUNT=0
  76. 70 IF (MCON .EQ. M) THEN
  77. OPTNEW=RESMAX
  78. ELSE
  79. OPTNEW=0.0d0
  80. DO 80 I=1,N
  81. 80 OPTNEW=OPTNEW-DX(I)*A(I,MCON)
  82. END IF
  83. IF (ICOUNT .EQ. 0 .OR. OPTNEW .LT. OPTOLD) THEN
  84. OPTOLD=OPTNEW
  85. NACTX=NACT
  86. ICOUNT=3
  87. ELSE IF (NACT .GT. NACTX) THEN
  88. NACTX=NACT
  89. ICOUNT=3
  90. ELSE
  91. ICOUNT=ICOUNT-1
  92. IF (ICOUNT .EQ. 0) GOTO 490
  93. END IF
  94. C
  95. C If ICON exceeds NACT, then we add the constraint with index IACT(ICON) to
  96. C the active set. Apply Givens rotations so that the last N-NACT-1 columns
  97. C of Z are orthogonal to the gradient of the new constraint, a scalar
  98. C product being set to zero if its nonzero value could be due to computer
  99. C rounding errors. The array DXNEW is used for working space.
  100. C
  101. IF (ICON .LE. NACT) GOTO 260
  102. KK=IACT(ICON)
  103. DO 90 I=1,N
  104. 90 DXNEW(I)=A(I,KK)
  105. TOT=0.0d0
  106. K=N
  107. 100 IF (K .GT. NACT) THEN
  108. SP=0.0d0
  109. SPABS=0.0d0
  110. DO 110 I=1,N
  111. TEMP=Z(I,K)*DXNEW(I)
  112. SP=SP+TEMP
  113. 110 SPABS=SPABS+ABS(TEMP)
  114. ACCA=SPABS+0.1d0*ABS(SP)
  115. ACCB=SPABS+0.2d0*ABS(SP)
  116. IF (SPABS .GE. ACCA .OR. ACCA .GE. ACCB) SP=0.0d0
  117. IF (TOT .EQ. 0.0d0) THEN
  118. TOT=SP
  119. ELSE
  120. KP=K+1
  121. TEMP=SQRT(SP*SP+TOT*TOT)
  122. ALPHA=SP/TEMP
  123. BETA=TOT/TEMP
  124. TOT=TEMP
  125. DO 120 I=1,N
  126. TEMP=ALPHA*Z(I,K)+BETA*Z(I,KP)
  127. Z(I,KP)=ALPHA*Z(I,KP)-BETA*Z(I,K)
  128. 120 Z(I,K)=TEMP
  129. END IF
  130. K=K-1
  131. GOTO 100
  132. END IF
  133. C
  134. C Add the new constraint if this can be done without a deletion from the
  135. C active set.
  136. C
  137. IF (TOT .NE. 0.0d0) THEN
  138. NACT=NACT+1
  139. ZDOTA(NACT)=TOT
  140. VMULTC(ICON)=VMULTC(NACT)
  141. VMULTC(NACT)=0.0d0
  142. GOTO 210
  143. END IF
  144. C
  145. C The next instruction is reached if a deletion has to be made from the
  146. C active set in order to make room for the new active constraint, because
  147. C the new constraint gradient is a linear combination of the gradients of
  148. C the old active constraints. Set the elements of VMULTD to the multipliers
  149. C of the linear combination. Further, set IOUT to the index of the
  150. C constraint to be deleted, but branch if no suitable index can be found.
  151. C
  152. RATIO=-1.0d0
  153. K=NACT
  154. 130 ZDOTV=0.0d0
  155. ZDVABS=0.0d0
  156. DO 140 I=1,N
  157. TEMP=Z(I,K)*DXNEW(I)
  158. ZDOTV=ZDOTV+TEMP
  159. 140 ZDVABS=ZDVABS+ABS(TEMP)
  160. ACCA=ZDVABS+0.1d0*ABS(ZDOTV)
  161. ACCB=ZDVABS+0.2d0*ABS(ZDOTV)
  162. IF (ZDVABS .LT. ACCA .AND. ACCA .LT. ACCB) THEN
  163. TEMP=ZDOTV/ZDOTA(K)
  164. IF (TEMP .GT. 0.0d0 .AND. IACT(K) .LE. M) THEN
  165. TEMPA=VMULTC(K)/TEMP
  166. IF (RATIO .LT. 0.0d0 .OR. TEMPA .LT. RATIO) THEN
  167. RATIO=TEMPA
  168. IOUT=K
  169. END IF
  170. END IF
  171. IF (K .GE. 2) THEN
  172. KW=IACT(K)
  173. DO 150 I=1,N
  174. 150 DXNEW(I)=DXNEW(I)-TEMP*A(I,KW)
  175. END IF
  176. VMULTD(K)=TEMP
  177. ELSE
  178. VMULTD(K)=0.0d0
  179. END IF
  180. K=K-1
  181. IF (K .GT. 0) GOTO 130
  182. IF (RATIO .LT. 0.0d0) GOTO 490
  183. C
  184. C Revise the Lagrange multipliers and reorder the active constraints so
  185. C that the one to be replaced is at the end of the list. Also calculate the
  186. C new value of ZDOTA(NACT) and branch if it is not acceptable.
  187. C
  188. DO 160 K=1,NACT
  189. 160 VMULTC(K)=DMAX1(0.0d0,VMULTC(K)-RATIO*VMULTD(K))
  190. IF (ICON .LT. NACT) THEN
  191. ISAVE=IACT(ICON)
  192. VSAVE=VMULTC(ICON)
  193. K=ICON
  194. 170 KP=K+1
  195. KW=IACT(KP)
  196. SP=0.0d0
  197. DO 180 I=1,N
  198. 180 SP=SP+Z(I,K)*A(I,KW)
  199. TEMP=SQRT(SP*SP+ZDOTA(KP)**2)
  200. ALPHA=ZDOTA(KP)/TEMP
  201. BETA=SP/TEMP
  202. ZDOTA(KP)=ALPHA*ZDOTA(K)
  203. ZDOTA(K)=TEMP
  204. DO 190 I=1,N
  205. TEMP=ALPHA*Z(I,KP)+BETA*Z(I,K)
  206. Z(I,KP)=ALPHA*Z(I,K)-BETA*Z(I,KP)
  207. 190 Z(I,K)=TEMP
  208. IACT(K)=KW
  209. VMULTC(K)=VMULTC(KP)
  210. K=KP
  211. IF (K .LT. NACT) GOTO 170
  212. IACT(K)=ISAVE
  213. VMULTC(K)=VSAVE
  214. END IF
  215. TEMP=0.0d0
  216. DO 200 I=1,N
  217. 200 TEMP=TEMP+Z(I,NACT)*A(I,KK)
  218. IF (TEMP .EQ. 0.0d0) GOTO 490
  219. ZDOTA(NACT)=TEMP
  220. VMULTC(ICON)=0.0d0
  221. VMULTC(NACT)=RATIO
  222. C
  223. C Update IACT and ensure that the objective function continues to be
  224. C treated as the last active constraint when MCON>M.
  225. C
  226. 210 IACT(ICON)=IACT(NACT)
  227. IACT(NACT)=KK
  228. IF (MCON .GT. M .AND. KK .NE. MCON) THEN
  229. K=NACT-1
  230. SP=0.0d0
  231. DO 220 I=1,N
  232. 220 SP=SP+Z(I,K)*A(I,KK)
  233. TEMP=SQRT(SP*SP+ZDOTA(NACT)**2)
  234. ALPHA=ZDOTA(NACT)/TEMP
  235. BETA=SP/TEMP
  236. ZDOTA(NACT)=ALPHA*ZDOTA(K)
  237. ZDOTA(K)=TEMP
  238. DO 230 I=1,N
  239. TEMP=ALPHA*Z(I,NACT)+BETA*Z(I,K)
  240. Z(I,NACT)=ALPHA*Z(I,K)-BETA*Z(I,NACT)
  241. 230 Z(I,K)=TEMP
  242. IACT(NACT)=IACT(K)
  243. IACT(K)=KK
  244. TEMP=VMULTC(K)
  245. VMULTC(K)=VMULTC(NACT)
  246. VMULTC(NACT)=TEMP
  247. END IF
  248. C
  249. C If stage one is in progress, then set SDIRN to the direction of the next
  250. C change to the current vector of variables.
  251. C
  252. IF (MCON .GT. M) GOTO 320
  253. KK=IACT(NACT)
  254. TEMP=0.0d0
  255. DO 240 I=1,N
  256. 240 TEMP=TEMP+SDIRN(I)*A(I,KK)
  257. TEMP=TEMP-1.0d0
  258. TEMP=TEMP/ZDOTA(NACT)
  259. DO 250 I=1,N
  260. 250 SDIRN(I)=SDIRN(I)-TEMP*Z(I,NACT)
  261. GOTO 340
  262. C
  263. C Delete the constraint that has the index IACT(ICON) from the active set.
  264. C
  265. 260 IF (ICON .LT. NACT) THEN
  266. ISAVE=IACT(ICON)
  267. VSAVE=VMULTC(ICON)
  268. K=ICON
  269. 270 KP=K+1
  270. KK=IACT(KP)
  271. SP=0.0d0
  272. DO 280 I=1,N
  273. 280 SP=SP+Z(I,K)*A(I,KK)
  274. TEMP=SQRT(SP*SP+ZDOTA(KP)**2)
  275. ALPHA=ZDOTA(KP)/TEMP
  276. BETA=SP/TEMP
  277. ZDOTA(KP)=ALPHA*ZDOTA(K)
  278. ZDOTA(K)=TEMP
  279. DO 290 I=1,N
  280. TEMP=ALPHA*Z(I,KP)+BETA*Z(I,K)
  281. Z(I,KP)=ALPHA*Z(I,K)-BETA*Z(I,KP)
  282. 290 Z(I,K)=TEMP
  283. IACT(K)=KK
  284. VMULTC(K)=VMULTC(KP)
  285. K=KP
  286. IF (K .LT. NACT) GOTO 270
  287. IACT(K)=ISAVE
  288. VMULTC(K)=VSAVE
  289. END IF
  290. NACT=NACT-1
  291. C
  292. C If stage one is in progress, then set SDIRN to the direction of the next
  293. C change to the current vector of variables.
  294. C
  295. IF (MCON .GT. M) GOTO 320
  296. TEMP=0.0d0
  297. DO 300 I=1,N
  298. 300 TEMP=TEMP+SDIRN(I)*Z(I,NACT+1)
  299. DO 310 I=1,N
  300. 310 SDIRN(I)=SDIRN(I)-TEMP*Z(I,NACT+1)
  301. GO TO 340
  302. C
  303. C Pick the next search direction of stage two.
  304. C
  305. 320 TEMP=1.0d0/ZDOTA(NACT)
  306. DO 330 I=1,N
  307. 330 SDIRN(I)=TEMP*Z(I,NACT)
  308. C
  309. C Calculate the step to the boundary of the trust region or take the step
  310. C that reduces RESMAX to zero. The two statements below that include the
  311. C factor 1.0E-6 prevent some harmless underflows that occurred in a test
  312. C calculation. Further, we skip the step if it could be zero within a
  313. C reasonable tolerance for computer rounding errors.
  314. C
  315. 340 DD=RHO*RHO
  316. SD=0.0d0
  317. SS=0.0d0
  318. DO 350 I=1,N
  319. IF (ABS(DX(I)) .GE. 1.0E-6*RHO) DD=DD-DX(I)**2
  320. SD=SD+DX(I)*SDIRN(I)
  321. 350 SS=SS+SDIRN(I)**2
  322. IF (DD .LE. 0.0d0) GOTO 490
  323. TEMP=SQRT(SS*DD)
  324. IF (ABS(SD) .GE. 1.0E-6*TEMP) TEMP=SQRT(SS*DD+SD*SD)
  325. STPFUL=DD/(TEMP+SD)
  326. STEP=STPFUL
  327. IF (MCON .EQ. M) THEN
  328. ACCA=STEP+0.1d0*RESMAX
  329. ACCB=STEP+0.2d0*RESMAX
  330. IF (STEP .GE. ACCA .OR. ACCA .GE. ACCB) GOTO 480
  331. STEP=DMIN1(STEP,RESMAX)
  332. END IF
  333. C
  334. C Set DXNEW to the new variables if STEP is the steplength, and reduce
  335. C RESMAX to the corresponding maximum residual if stage one is being done.
  336. C Because DXNEW will be changed during the calculation of some Lagrange
  337. C multipliers, it will be restored to the following value later.
  338. C
  339. DO 360 I=1,N
  340. 360 DXNEW(I)=DX(I)+STEP*SDIRN(I)
  341. IF (MCON .EQ. M) THEN
  342. RESOLD=RESMAX
  343. RESMAX=0.0d0
  344. DO 380 K=1,NACT
  345. KK=IACT(K)
  346. TEMP=B(KK)
  347. DO 370 I=1,N
  348. 370 TEMP=TEMP-A(I,KK)*DXNEW(I)
  349. RESMAX=DMAX1(RESMAX,TEMP)
  350. 380 CONTINUE
  351. END IF
  352. C
  353. C Set VMULTD to the VMULTC vector that would occur if DX became DXNEW. A
  354. C device is included to force VMULTD(K)=0.0 if deviations from this value
  355. C can be attributed to computer rounding errors. First calculate the new
  356. C Lagrange multipliers.
  357. C
  358. K=NACT
  359. 390 ZDOTW=0.0d0
  360. ZDWABS=0.0d0
  361. DO 400 I=1,N
  362. TEMP=Z(I,K)*DXNEW(I)
  363. ZDOTW=ZDOTW+TEMP
  364. 400 ZDWABS=ZDWABS+ABS(TEMP)
  365. ACCA=ZDWABS+0.1d0*ABS(ZDOTW)
  366. ACCB=ZDWABS+0.2d0*ABS(ZDOTW)
  367. IF (ZDWABS .GE. ACCA .OR. ACCA .GE. ACCB) ZDOTW=0.0d0
  368. VMULTD(K)=ZDOTW/ZDOTA(K)
  369. IF (K .GE. 2) THEN
  370. KK=IACT(K)
  371. DO 410 I=1,N
  372. 410 DXNEW(I)=DXNEW(I)-VMULTD(K)*A(I,KK)
  373. K=K-1
  374. GOTO 390
  375. END IF
  376. IF (MCON .GT. M) VMULTD(NACT)=DMAX1(0.0d0,VMULTD(NACT))
  377. C
  378. C Complete VMULTC by finding the new constraint residuals.
  379. C
  380. DO 420 I=1,N
  381. 420 DXNEW(I)=DX(I)+STEP*SDIRN(I)
  382. IF (MCON .GT. NACT) THEN
  383. KL=NACT+1
  384. DO 440 K=KL,MCON
  385. KK=IACT(K)
  386. SUM=RESMAX-B(KK)
  387. SUMABS=RESMAX+ABS(B(KK))
  388. DO 430 I=1,N
  389. TEMP=A(I,KK)*DXNEW(I)
  390. SUM=SUM+TEMP
  391. 430 SUMABS=SUMABS+ABS(TEMP)
  392. ACCA=SUMABS+0.1*ABS(SUM)
  393. ACCB=SUMABS+0.2*ABS(SUM)
  394. IF (SUMABS .GE. ACCA .OR. ACCA .GE. ACCB) SUM=0.0
  395. 440 VMULTD(K)=SUM
  396. END IF
  397. C
  398. C Calculate the fraction of the step from DX to DXNEW that will be taken.
  399. C
  400. RATIO=1.0d0
  401. ICON=0
  402. DO 450 K=1,MCON
  403. IF (VMULTD(K) .LT. 0.0d0) THEN
  404. TEMP=VMULTC(K)/(VMULTC(K)-VMULTD(K))
  405. IF (TEMP .LT. RATIO) THEN
  406. RATIO=TEMP
  407. ICON=K
  408. END IF
  409. END IF
  410. 450 CONTINUE
  411. C
  412. C Update DX, VMULTC and RESMAX.
  413. C
  414. TEMP=1.0d0-RATIO
  415. DO 460 I=1,N
  416. 460 DX(I)=TEMP*DX(I)+RATIO*DXNEW(I)
  417. DO 470 K=1,MCON
  418. 470 VMULTC(K)=DMAX1(0.0d0,TEMP*VMULTC(K)+RATIO*VMULTD(K))
  419. IF (MCON .EQ. M) RESMAX=RESOLD+RATIO*(RESMAX-RESOLD)
  420. C
  421. C If the full step is not acceptable then begin another iteration.
  422. C Otherwise switch to stage two or end the calculation.
  423. C
  424. IF (ICON .GT. 0) GOTO 70
  425. IF (STEP .EQ. STPFUL) GOTO 500
  426. 480 MCON=M+1
  427. ICON=MCON
  428. IACT(MCON)=MCON
  429. VMULTC(MCON)=0.0d0
  430. GOTO 60
  431. C
  432. C We employ any freedom that may be available to reduce the objective
  433. C function before returning a DX whose length is less than RHO.
  434. C
  435. 490 IF (MCON .EQ. M) GOTO 480
  436. IFULL=0
  437. 500 RETURN
  438. END