MoDeNa  1.0
Software framework facilitating sequential multi-scale modelling
module_solve_nonlinear.F90
1 
2 MODULE solve_nonlin
3 
4 ! Corrections to FUNCTION Enorm - 28 November 2003
5 
6 IMPLICIT NONE
7 INTEGER, PARAMETER :: dp = selected_real_kind(14, 60)
8 PRIVATE
9 PUBLIC :: hbrd, hybrd
10 
11 
12 CONTAINS
13 
14 
15 SUBROUTINE hbrd(fcn, n, x, fvec, epsfcn, tol, info, diag)
16 
17 ! Code converted using TO_F90 by Alan Miller
18 ! Date: 2003-07-15 Time: 13:27:42
19 
20 INTEGER, INTENT(IN) :: n
21 REAL (dp), INTENT(IN OUT) :: x(n)
22 REAL (dp), INTENT(IN OUT) :: fvec(n)
23 REAL (dp), INTENT(IN) :: epsfcn
24 REAL (dp), INTENT(IN) :: tol
25 INTEGER, INTENT(OUT) :: info
26 REAL (dp), INTENT(OUT) :: diag(n)
27 
28 ! EXTERNAL fcn
29 INTERFACE
30  SUBROUTINE fcn(N, X, FVEC, IFLAG)
31  IMPLICIT NONE
32  INTEGER, PARAMETER :: dp = selected_real_kind(14, 60)
33  INTEGER, INTENT(IN) :: n
34  REAL (dp), INTENT(IN) :: x(n)
35  REAL (dp), INTENT(OUT) :: fvec(n)
36  INTEGER, INTENT(IN OUT) :: iflag
37  END SUBROUTINE fcn
38 END INTERFACE
39 
100 INTEGER :: maxfev, ml, mode, mu, nfev, nprint
101 REAL (dp) :: xtol
102 REAL (dp), PARAMETER :: factor = 1.0_dp, zero = 0.0_dp
103 
104 info = 0
105 
106 ! CHECK THE INPUT PARAMETERS FOR ERRORS.
107 
108 
109 IF (n <= 0 .OR. epsfcn < zero .OR. tol < zero) GO TO 20
110 
111 ! CALL HYBRD.
112 
113 maxfev = 200*(n + 1)
114 xtol = tol
115 ml = n - 1
116 mu = n - 1
117 mode = 2
118 nprint = 0
119 CALL hybrd(fcn, n, x, fvec, xtol, maxfev, ml, mu, epsfcn, diag, mode, &
120  factor, nprint, info, nfev)
121 IF (info == 5) info = 4
122 20 RETURN
123 
124 ! LAST CARD OF SUBROUTINE HBRD.
125 
126 END SUBROUTINE hbrd
127 
128 
129 
130 SUBROUTINE hybrd(fcn, n, x, fvec, xtol, maxfev, ml, mu, epsfcn, diag, mode, &
131  factor, nprint, info, nfev)
132 
133 INTEGER, INTENT(IN) :: n
134 REAL (dp), INTENT(IN OUT) :: x(n)
135 REAL (dp), INTENT(IN OUT) :: fvec(n)
136 REAL (dp), INTENT(IN) :: xtol
137 INTEGER, INTENT(IN OUT) :: maxfev
138 INTEGER, INTENT(IN OUT) :: ml
139 INTEGER, INTENT(IN) :: mu
140 REAL (dp), INTENT(IN) :: epsfcn
141 REAL (dp), INTENT(OUT) :: diag(n)
142 INTEGER, INTENT(IN) :: mode
143 REAL (dp), INTENT(IN) :: factor
144 INTEGER, INTENT(IN OUT) :: nprint
145 INTEGER, INTENT(OUT) :: info
146 INTEGER, INTENT(OUT) :: nfev
147 
148 ! EXTERNAL fcn
149 INTERFACE
150  SUBROUTINE fcn(N, X, FVEC, IFLAG)
151  IMPLICIT NONE
152  INTEGER, PARAMETER :: dp = selected_real_kind(14, 60)
153  INTEGER, INTENT(IN) :: n
154  REAL (dp), INTENT(IN) :: x(n)
155  REAL (dp), INTENT(OUT) :: fvec(n)
156  INTEGER, INTENT(IN OUT) :: iflag
157  END SUBROUTINE fcn
158 END INTERFACE
159 
160 ! **********
161 
162 ! SUBROUTINE HYBRD
163 
164 ! THE PURPOSE OF HYBRD IS TO FIND A ZERO OF A SYSTEM OF N NONLINEAR
165 ! FUNCTIONS IN N VARIABLES BY A MODIFICATION OF THE POWELL HYBRID METHOD.
166 ! THE USER MUST PROVIDE A SUBROUTINE WHICH CALCULATES THE FUNCTIONS.
167 ! THE JACOBIAN IS THEN CALCULATED BY A FORWARD-DIFFERENCE APPROXIMATION.
168 
169 ! THE SUBROUTINE STATEMENT IS
170 
171 ! SUBROUTINE HYBRD(FCN, N, X, FVEC, XTOL, MAXFEV, ML, MU, EPSFCN,
172 ! DIAG, MODE, FACTOR, NPRINT, INFO, NFEV, FJAC,
173 ! LDFJAC, R, LR, QTF, WA1, WA2, WA3, WA4)
174 
175 ! WHERE
176 
177 ! FCN IS THE NAME OF THE USER-SUPPLIED SUBROUTINE WHICH CALCULATES
178 ! THE FUNCTIONS. FCN MUST BE DECLARED IN AN EXTERNAL STATEMENT IN
179 ! THE USER CALLING PROGRAM, AND SHOULD BE WRITTEN AS FOLLOWS.
180 
181 ! SUBROUTINE FCN(N, X, FVEC, IFLAG)
182 ! INTEGER N, IFLAG
183 ! REAL X(N), FVEC(N)
184 ! ----------
185 ! CALCULATE THE FUNCTIONS AT X AND
186 ! RETURN THIS VECTOR IN FVEC.
187 ! ---------
188 ! RETURN
189 ! END
190 
191 ! THE VALUE OF IFLAG SHOULD NOT BE CHANGED BY FCN UNLESS
192 ! THE USER WANTS TO TERMINATE EXECUTION OF HYBRD.
193 ! IN THIS CASE SET IFLAG TO A NEGATIVE INTEGER.
194 
195 ! N IS A POSITIVE INTEGER INPUT VARIABLE SET TO THE NUMBER
196 ! OF FUNCTIONS AND VARIABLES.
197 
198 ! X IS AN ARRAY OF LENGTH N. ON INPUT X MUST CONTAIN AN INITIAL
199 ! ESTIMATE OF THE SOLUTION VECTOR. ON OUTPUT X CONTAINS THE FINAL
200 ! ESTIMATE OF THE SOLUTION VECTOR.
201 
202 ! FVEC IS AN OUTPUT ARRAY OF LENGTH N WHICH CONTAINS
203 ! THE FUNCTIONS EVALUATED AT THE OUTPUT X.
204 
205 ! XTOL IS A NONNEGATIVE INPUT VARIABLE. TERMINATION OCCURS WHEN THE
206 ! RELATIVE ERROR BETWEEN TWO CONSECUTIVE ITERATES IS AT MOST XTOL.
207 
208 ! MAXFEV IS A POSITIVE INTEGER INPUT VARIABLE. TERMINATION OCCURS WHEN
209 ! THE NUMBER OF CALLS TO FCN IS AT LEAST MAXFEV BY THE END OF AN
210 ! ITERATION.
211 
212 ! ML IS A NONNEGATIVE INTEGER INPUT VARIABLE WHICH SPECIFIES THE
213 ! NUMBER OF SUBDIAGONALS WITHIN THE BAND OF THE JACOBIAN MATRIX.
214 ! IF THE JACOBIAN IS NOT BANDED, SET ML TO AT LEAST N - 1.
215 
216 ! MU IS A NONNEGATIVE INTEGER INPUT VARIABLE WHICH SPECIFIES THE NUMBER
217 ! OF SUPERDIAGONALS WITHIN THE BAND OF THE JACOBIAN MATRIX.
218 ! IF THE JACOBIAN IS NOT BANDED, SET MU TO AT LEAST N - 1.
219 
220 ! EPSFCN IS AN INPUT VARIABLE USED IN DETERMINING A SUITABLE STEP LENGTH
221 ! FOR THE FORWARD-DIFFERENCE APPROXIMATION. THIS APPROXIMATION
222 ! ASSUMES THAT THE RELATIVE ERRORS IN THE FUNCTIONS ARE OF THE ORDER
223 ! OF EPSFCN. IF EPSFCN IS LESS THAN THE MACHINE PRECISION,
224 ! IT IS ASSUMED THAT THE RELATIVE ERRORS IN THE FUNCTIONS ARE OF THE
225 ! ORDER OF THE MACHINE PRECISION.
226 
227 ! DIAG IS AN ARRAY OF LENGTH N. IF MODE = 1 (SEE BELOW),
228 ! DIAG IS INTERNALLY SET. IF MODE = 2, DIAG MUST CONTAIN POSITIVE
229 ! ENTRIES THAT SERVE AS MULTIPLICATIVE SCALE FACTORS FOR THE
230 ! VARIABLES.
231 
232 ! MODE IS AN INTEGER INPUT VARIABLE. IF MODE = 1, THE VARIABLES WILL BE
233 ! SCALED INTERNALLY. IF MODE = 2, THE SCALING IS SPECIFIED BY THE
234 ! INPUT DIAG. OTHER VALUES OF MODE ARE EQUIVALENT TO MODE = 1.
235 
236 ! FACTOR IS A POSITIVE INPUT VARIABLE USED IN DETERMINING THE
237 ! INITIAL STEP BOUND. THIS BOUND IS SET TO THE PRODUCT OF
238 ! FACTOR AND THE EUCLIDEAN NORM OF DIAG*X IF NONZERO, OR ELSE
239 ! TO FACTOR ITSELF. IN MOST CASES FACTOR SHOULD LIE IN THE
240 ! INTERVAL (.1,100.). 100. IS A GENERALLY RECOMMENDED VALUE.
241 
242 ! NPRINT IS AN INTEGER INPUT VARIABLE THAT ENABLES CONTROLLED
243 ! PRINTING OF ITERATES IF IT IS POSITIVE. IN THIS CASE,
244 ! FCN IS CALLED WITH IFLAG = 0 AT THE BEGINNING OF THE FIRST
245 ! ITERATION AND EVERY NPRINT ITERATIONS THEREAFTER AND
246 ! IMMEDIATELY PRIOR TO RETURN, WITH X AND FVEC AVAILABLE
247 ! FOR PRINTING. IF NPRINT IS NOT POSITIVE, NO SPECIAL CALLS
248 ! OF FCN WITH IFLAG = 0 ARE MADE.
249 
250 ! INFO IS AN INTEGER OUTPUT VARIABLE. IF THE USER HAS
251 ! TERMINATED EXECUTION, INFO IS SET TO THE (NEGATIVE)
252 ! VALUE OF IFLAG. SEE DESCRIPTION OF FCN. OTHERWISE,
253 ! INFO IS SET AS FOLLOWS.
254 
255 ! INFO = 0 IMPROPER INPUT PARAMETERS.
256 
257 ! INFO = 1 RELATIVE ERROR BETWEEN TWO CONSECUTIVE ITERATES
258 ! IS AT MOST XTOL.
259 
260 ! INFO = 2 NUMBER OF CALLS TO FCN HAS REACHED OR EXCEEDED MAXFEV.
261 
262 ! INFO = 3 XTOL IS TOO SMALL. NO FURTHER IMPROVEMENT IN
263 ! THE APPROXIMATE SOLUTION X IS POSSIBLE.
264 
265 ! INFO = 4 ITERATION IS NOT MAKING GOOD PROGRESS, AS
266 ! MEASURED BY THE IMPROVEMENT FROM THE LAST
267 ! FIVE JACOBIAN EVALUATIONS.
268 
269 ! INFO = 5 ITERATION IS NOT MAKING GOOD PROGRESS, AS MEASURED BY
270 ! THE IMPROVEMENT FROM THE LAST TEN ITERATIONS.
271 
272 ! NFEV IS AN INTEGER OUTPUT VARIABLE SET TO THE NUMBER OF CALLS TO FCN.
273 
274 ! FJAC IS AN OUTPUT N BY N ARRAY WHICH CONTAINS THE ORTHOGONAL MATRIX Q
275 ! PRODUCED BY THE QR FACTORIZATION OF THE FINAL APPROXIMATE JACOBIAN.
276 
277 ! LDFJAC IS A POSITIVE INTEGER INPUT VARIABLE NOT LESS THAN N
278 ! WHICH SPECIFIES THE LEADING DIMENSION OF THE ARRAY FJAC.
279 
280 ! R IS AN OUTPUT ARRAY OF LENGTH LR WHICH CONTAINS THE
281 ! UPPER TRIANGULAR MATRIX PRODUCED BY THE QR FACTORIZATION
282 ! OF THE FINAL APPROXIMATE JACOBIAN, STORED ROWWISE.
283 
284 ! LR IS A POSITIVE INTEGER INPUT VARIABLE NOT LESS THAN (N*(N+1))/2.
285 
286 ! QTF IS AN OUTPUT ARRAY OF LENGTH N WHICH CONTAINS
287 ! THE VECTOR (Q TRANSPOSE)*FVEC.
288 
289 ! WA1, WA2, WA3, AND WA4 ARE WORK ARRAYS OF LENGTH N.
290 
291 ! SUBPROGRAMS CALLED
292 
293 ! USER-SUPPLIED ...... FCN
294 
295 ! MINPACK-SUPPLIED ... DOGLEG,SPMPAR,ENORM,FDJAC1,
296 ! QFORM,QRFAC,R1MPYQ,R1UPDT
297 
298 ! FORTRAN-SUPPLIED ... ABS,MAX,MIN,MIN,MOD
299 
300 ! ARGONNE NATIONAL LABORATORY. MINPACK PROJECT. MARCH 1980.
301 ! BURTON S. GARBOW, KENNETH E. HILLSTROM, JORGE J. MORE
302 
303 ! **********
304 
305 INTEGER :: i, iflag, iter, j, jm1, l, lr, msum, ncfail, ncsuc, nslow1, &
306  nslow2
307 INTEGER :: iwa(1)
308 LOGICAL :: jeval, sing
309 REAL (dp) :: actred, delta, epsmch, fnorm, fnorm1, pnorm, prered, &
310  ratio, sum, temp, xnorm
311 REAL (dp), PARAMETER :: one = 1.0_dp, p1 = 0.1_dp, p5 = 0.5_dp, &
312  p001 = 0.001_dp, p0001 = 0.0001_dp, zero = 0.0_dp
313 
314 ! The following were workspace arguments
315 REAL (dp) :: fjac(n,n), r(n*(n+1)/2), qtf(n), wa1(n), wa2(n), &
316  wa3(n), wa4(n)
317 
318 ! EPSMCH IS THE MACHINE PRECISION.
319 
320 epsmch = epsilon(1.0_dp)
321 
322 info = 0
323 iflag = 0
324 nfev = 0
325 lr = n*(n+1)/2
326 
327 ! CHECK THE INPUT PARAMETERS FOR ERRORS.
328 
329 IF (n > 0 .AND. xtol >= zero .AND. maxfev > 0 .AND. ml >= 0 .AND. mu >= &
330  0 .AND. factor > zero ) THEN
331 IF (mode == 2) THEN
332  diag(1:n) = one
333 END IF
334 
335 ! EVALUATE THE FUNCTION AT THE STARTING POINT AND CALCULATE ITS NORM.
336 
337 iflag = 1
338 CALL fcn(n, x, fvec, iflag)
339 nfev = 1
340 IF (iflag >= 0) THEN
341  fnorm = enorm(n, fvec)
342 
343 ! DETERMINE THE NUMBER OF CALLS TO FCN NEEDED TO COMPUTE THE JACOBIAN MATRIX.
344 
345  msum = min(ml+mu+1,n)
346 
347 ! INITIALIZE ITERATION COUNTER AND MONITORS.
348 
349  iter = 1
350  ncsuc = 0
351  ncfail = 0
352  nslow1 = 0
353  nslow2 = 0
354 
355 ! BEGINNING OF THE OUTER LOOP.
356 
357  20 jeval = .true.
358 
359 ! CALCULATE THE JACOBIAN MATRIX.
360 
361  iflag = 2
362  CALL fdjac1(fcn, n, x, fvec, fjac, n, iflag, ml, mu, epsfcn, wa1, wa2)
363  nfev = nfev + msum
364  IF (iflag >= 0) THEN
365 
366 ! COMPUTE THE QR FACTORIZATION OF THE JACOBIAN.
367 
368  CALL qrfac(n, n, fjac, n, .false., iwa, 1, wa1, wa2, wa3)
369 
370 ! ON THE FIRST ITERATION AND IF MODE IS 1, SCALE ACCORDING
371 ! TO THE NORMS OF THE COLUMNS OF THE INITIAL JACOBIAN.
372 
373  IF (iter == 1) THEN
374  IF (mode /= 2) THEN
375  DO j = 1, n
376  diag(j) = wa2(j)
377  IF (wa2(j) == zero) diag(j) = one
378  END DO
379  END IF
380 
381 ! ON THE FIRST ITERATION, CALCULATE THE NORM OF THE SCALED X
382 ! AND INITIALIZE THE STEP BOUND DELTA.
383 
384  wa3(1:n) = diag(1:n) * x(1:n)
385  xnorm = enorm(n, wa3)
386  delta = factor * xnorm
387  IF (delta == zero) delta = factor
388  END IF
389 
390 ! FORM (Q TRANSPOSE)*FVEC AND STORE IN QTF.
391 
392  qtf(1:n) = fvec(1:n)
393  DO j = 1, n
394  IF (fjac(j,j) /= zero) THEN
395  sum = zero
396  DO i = j, n
397  sum = sum + fjac(i,j) * qtf(i)
398  END DO
399  temp = -sum / fjac(j,j)
400  DO i = j, n
401  qtf(i) = qtf(i) + fjac(i,j) * temp
402  END DO
403  END IF
404  END DO
405 
406 ! COPY THE TRIANGULAR FACTOR OF THE QR FACTORIZATION INTO R.
407 
408  sing = .false.
409  DO j = 1, n
410  l = j
411  jm1 = j - 1
412  IF (jm1 >= 1) THEN
413  DO i = 1, jm1
414  r(l) = fjac(i,j)
415  l = l + n - i
416  END DO
417  END IF
418  r(l) = wa1(j)
419  IF (wa1(j) == zero) sing = .true.
420  END DO
421 
422 ! ACCUMULATE THE ORTHOGONAL FACTOR IN FJAC.
423 
424  CALL qform(n, n, fjac, n, wa1)
425 
426 ! RESCALE IF NECESSARY.
427 
428  IF (mode /= 2) THEN
429  DO j = 1, n
430  diag(j) = max(diag(j), wa2(j))
431  END DO
432  END IF
433 
434 ! BEGINNING OF THE INNER LOOP.
435 
436 ! IF REQUESTED, CALL FCN TO ENABLE PRINTING OF ITERATES.
437 
438  120 IF (nprint > 0) THEN
439  iflag = 0
440  IF (mod(iter-1, nprint) == 0) CALL fcn(n, x, fvec, iflag)
441  IF (iflag < 0) GO TO 190
442  END IF
443 
444 ! DETERMINE THE DIRECTION P.
445 
446  CALL dogleg(n, r, lr, diag, qtf, delta, wa1, wa2, wa3)
447 
448 ! STORE THE DIRECTION P AND X + P. CALCULATE THE NORM OF P.
449 
450  DO j = 1, n
451  wa1(j) = -wa1(j)
452  wa2(j) = x(j) + wa1(j)
453  wa3(j) = diag(j) * wa1(j)
454  END DO
455  pnorm = enorm(n, wa3)
456 
457 ! ON THE FIRST ITERATION, ADJUST THE INITIAL STEP BOUND.
458 
459  IF (iter == 1) delta = min(delta, pnorm)
460 
461 ! EVALUATE THE FUNCTION AT X + P AND CALCULATE ITS NORM.
462 
463  iflag = 1
464  CALL fcn(n, wa2, wa4, iflag)
465  nfev = nfev + 1
466  IF (iflag >= 0) THEN
467  fnorm1 = enorm(n, wa4)
468 
469 ! COMPUTE THE SCALED ACTUAL REDUCTION.
470 
471  actred = -one
472  IF (fnorm1 < fnorm) actred = one - (fnorm1/fnorm) ** 2
473 
474 ! COMPUTE THE SCALED PREDICTED REDUCTION.
475 
476  l = 1
477  DO i = 1, n
478  sum = zero
479  DO j = i, n
480  sum = sum + r(l) * wa1(j)
481  l = l + 1
482  END DO
483  wa3(i) = qtf(i) + sum
484  END DO
485  temp = enorm(n, wa3)
486  prered = zero
487  IF (temp < fnorm) prered = one - (temp/fnorm) ** 2
488 
489 ! COMPUTE THE RATIO OF THE ACTUAL TO THE PREDICTED REDUCTION.
490 
491  ratio = zero
492  IF (prered > zero) ratio = actred / prered
493 
494 ! UPDATE THE STEP BOUND.
495 
496  IF (ratio < p1) THEN
497  ncsuc = 0
498  ncfail = ncfail + 1
499  delta = p5 * delta
500  ELSE
501  ncfail = 0
502  ncsuc = ncsuc + 1
503  IF (ratio >= p5 .OR. ncsuc > 1) delta = max(delta,pnorm/p5)
504  IF (abs(ratio-one) <= p1) delta = pnorm / p5
505  END IF
506 
507 ! TEST FOR SUCCESSFUL ITERATION.
508 
509  IF (ratio >= p0001) THEN
510 
511 ! SUCCESSFUL ITERATION. UPDATE X, FVEC, AND THEIR NORMS.
512 
513  DO j = 1, n
514  x(j) = wa2(j)
515  wa2(j) = diag(j) * x(j)
516  fvec(j) = wa4(j)
517  END DO
518  xnorm = enorm(n, wa2)
519  fnorm = fnorm1
520  iter = iter + 1
521  END IF
522 
523 ! DETERMINE THE PROGRESS OF THE ITERATION.
524 
525  nslow1 = nslow1 + 1
526  IF (actred >= p001) nslow1 = 0
527  IF (jeval) nslow2 = nslow2 + 1
528  IF (actred >= p1) nslow2 = 0
529 
530 ! TEST FOR CONVERGENCE.
531 
532  IF (delta <= xtol*xnorm .OR. fnorm == zero) info = 1
533  IF (info == 0) THEN
534 
535 ! TESTS FOR TERMINATION AND STRINGENT TOLERANCES.
536 
537  IF (nfev >= maxfev) info = 2
538  IF (p1*max(p1*delta, pnorm) <= epsmch*xnorm) info = 3
539  IF (nslow2 == 5) info = 4
540  IF (nslow1 == 10) info = 5
541  IF (info == 0) THEN
542 
543 ! CRITERION FOR RECALCULATING JACOBIAN APPROXIMATION
544 ! BY FORWARD DIFFERENCES.
545 
546  IF (ncfail /= 2) THEN
547 
548 ! CALCULATE THE RANK ONE MODIFICATION TO THE JACOBIAN
549 ! AND UPDATE QTF IF NECESSARY.
550 
551  DO j = 1, n
552  sum = zero
553  DO i = 1, n
554  sum = sum + fjac(i,j) * wa4(i)
555  END DO
556  wa2(j) = (sum-wa3(j)) / pnorm
557  wa1(j) = diag(j) * ((diag(j)*wa1(j))/pnorm)
558  IF (ratio >= p0001) qtf(j) = sum
559  END DO
560 
561 ! COMPUTE THE QR FACTORIZATION OF THE UPDATED JACOBIAN.
562 
563  CALL r1updt(n, n, r, lr, wa1, wa2, wa3, sing)
564  CALL r1mpyq(n, n, fjac, n, wa2, wa3)
565  CALL r1mpyq(1, n, qtf, 1, wa2, wa3)
566 
567 ! END OF THE INNER LOOP.
568 
569  jeval = .false.
570  GO TO 120
571  END IF
572 
573 ! END OF THE OUTER LOOP.
574 
575  GO TO 20
576  END IF
577  END IF
578  END IF
579  END IF
580 END IF
581 END IF
582 
583 ! TERMINATION, EITHER NORMAL OR USER IMPOSED.
584 
585 190 IF (iflag < 0) info = iflag
586 iflag = 0
587 IF (nprint > 0) CALL fcn(n, x, fvec, iflag)
588 RETURN
589 
590 ! LAST CARD OF SUBROUTINE HYBRD.
591 
592 END SUBROUTINE hybrd
593 
594 
595 
596 SUBROUTINE dogleg(n, r, lr, diag, qtb, delta, x, wa1, wa2)
597 
598 INTEGER, INTENT(IN) :: n
599 INTEGER, INTENT(IN) :: lr
600 REAL (dp), INTENT(IN) :: r(lr)
601 REAL (dp), INTENT(IN) :: diag(n)
602 REAL (dp), INTENT(IN) :: qtb(n)
603 REAL (dp), INTENT(IN) :: delta
604 REAL (dp), INTENT(IN OUT) :: x(n)
605 REAL (dp), INTENT(OUT) :: wa1(n)
606 REAL (dp), INTENT(OUT) :: wa2(n)
607 
608 
609 ! **********
610 
611 ! SUBROUTINE DOGLEG
612 
613 ! GIVEN AN M BY N MATRIX A, AN N BY N NONSINGULAR DIAGONAL
614 ! MATRIX D, AN M-VECTOR B, AND A POSITIVE NUMBER DELTA, THE
615 ! PROBLEM IS TO DETERMINE THE CONVEX COMBINATION X OF THE
616 ! GAUSS-NEWTON AND SCALED GRADIENT DIRECTIONS THAT MINIMIZES
617 ! (A*X - B) IN THE LEAST SQUARES SENSE, SUBJECT TO THE
618 ! RESTRICTION THAT THE EUCLIDEAN NORM OF D*X BE AT MOST DELTA.
619 
620 ! THIS SUBROUTINE COMPLETES THE SOLUTION OF THE PROBLEM
621 ! IF IT IS PROVIDED WITH THE NECESSARY INFORMATION FROM THE
622 ! QR FACTORIZATION OF A. THAT IS, IF A = Q*R, WHERE Q HAS
623 ! ORTHOGONAL COLUMNS AND R IS AN UPPER TRIANGULAR MATRIX,
624 ! THEN DOGLEG EXPECTS THE FULL UPPER TRIANGLE OF R AND
625 ! THE FIRST N COMPONENTS OF (Q TRANSPOSE)*B.
626 
627 ! THE SUBROUTINE STATEMENT IS
628 
629 ! SUBROUTINE DOGLEG(N,R,LR,DIAG,QTB,DELTA,X,WA1,WA2)
630 
631 ! WHERE
632 
633 ! N IS A POSITIVE INTEGER INPUT VARIABLE SET TO THE ORDER OF R.
634 
635 ! R IS AN INPUT ARRAY OF LENGTH LR WHICH MUST CONTAIN THE UPPER
636 ! TRIANGULAR MATRIX R STORED BY ROWS.
637 
638 ! LR IS A POSITIVE INTEGER INPUT VARIABLE NOT LESS THAN
639 ! (N*(N+1))/2.
640 
641 ! DIAG IS AN INPUT ARRAY OF LENGTH N WHICH MUST CONTAIN THE
642 ! DIAGONAL ELEMENTS OF THE MATRIX D.
643 
644 ! QTB IS AN INPUT ARRAY OF LENGTH N WHICH MUST CONTAIN THE FIRST
645 ! N ELEMENTS OF THE VECTOR (Q TRANSPOSE)*B.
646 
647 ! DELTA IS A POSITIVE INPUT VARIABLE WHICH SPECIFIES AN UPPER
648 ! BOUND ON THE EUCLIDEAN NORM OF D*X.
649 
650 ! X IS AN OUTPUT ARRAY OF LENGTH N WHICH CONTAINS THE DESIRED
651 ! CONVEX COMBINATION OF THE GAUSS-NEWTON DIRECTION AND THE
652 ! SCALED GRADIENT DIRECTION.
653 
654 ! WA1 AND WA2 ARE WORK ARRAYS OF LENGTH N.
655 
656 ! SUBPROGRAMS CALLED
657 
658 ! MINPACK-SUPPLIED ... SPMPAR,ENORM
659 
660 ! FORTRAN-SUPPLIED ... ABS,MAX,MIN,SQRT
661 
662 ! ARGONNE NATIONAL LABORATORY. MINPACK PROJECT. MARCH 1980.
663 ! BURTON S. GARBOW, KENNETH E. HILLSTROM, JORGE J. MORE
664 
665 ! **********
666 INTEGER :: i, j, jj, jp1, k, l
667 REAL (dp) :: alpha, bnorm, epsmch, gnorm, qnorm, sgnorm, sum, temp
668 
669 ! EPSMCH IS THE MACHINE PRECISION.
670 
671 epsmch = epsilon(1.0_dp)
672 
673 ! FIRST, CALCULATE THE GAUSS-NEWTON DIRECTION.
674 
675 jj = (n*(n+1)) / 2 + 1
676 DO k = 1, n
677  j = n - k + 1
678  jp1 = j + 1
679  jj = jj - k
680  l = jj + 1
681  sum = 0.0
682  IF (n >= jp1) THEN
683  DO i = jp1, n
684  sum = sum + r(l) * x(i)
685  l = l + 1
686  END DO
687  END IF
688  temp = r(jj)
689  IF (temp == 0.0) THEN
690  l = j
691  DO i = 1, j
692  temp = max(temp,abs(r(l)))
693  l = l + n - i
694  END DO
695  temp = epsmch * temp
696  IF (temp == 0.0) temp = epsmch
697  END IF
698  x(j) = (qtb(j)-sum) / temp
699 END DO
700 
701 ! TEST WHETHER THE GAUSS-NEWTON DIRECTION IS ACCEPTABLE.
702 
703 DO j = 1, n
704  wa1(j) = 0.0
705  wa2(j) = diag(j) * x(j)
706 END DO
707 qnorm = enorm(n, wa2)
708 IF (qnorm > delta) THEN
709 
710 ! THE GAUSS-NEWTON DIRECTION IS NOT ACCEPTABLE.
711 ! NEXT, CALCULATE THE SCALED GRADIENT DIRECTION.
712 
713  l = 1
714  DO j = 1, n
715  temp = qtb(j)
716  DO i = j, n
717  wa1(i) = wa1(i) + r(l) * temp
718  l = l + 1
719  END DO
720  wa1(j) = wa1(j) / diag(j)
721  END DO
722 
723 ! CALCULATE THE NORM OF THE SCALED GRADIENT AND TEST FOR
724 ! THE SPECIAL CASE IN WHICH THE SCALED GRADIENT IS ZERO.
725 
726  gnorm = enorm(n, wa1)
727  sgnorm = 0.0
728  alpha = delta / qnorm
729  IF (gnorm /= 0.0) THEN
730 
731 ! CALCULATE THE POINT ALONG THE SCALED GRADIENT
732 ! AT WHICH THE QUADRATIC IS MINIMIZED.
733 
734  DO j = 1, n
735  wa1(j) = (wa1(j)/gnorm) / diag(j)
736  END DO
737  l = 1
738  DO j = 1, n
739  sum = 0.0
740  DO i = j, n
741  sum = sum + r(l) * wa1(i)
742  l = l + 1
743  END DO
744  wa2(j) = sum
745  END DO
746  temp = enorm(n, wa2)
747  sgnorm = (gnorm/temp) / temp
748 
749 ! TEST WHETHER THE SCALED GRADIENT DIRECTION IS ACCEPTABLE.
750 
751  alpha = 0.0
752  IF (sgnorm < delta) THEN
753 
754 ! THE SCALED GRADIENT DIRECTION IS NOT ACCEPTABLE.
755 ! FINALLY, CALCULATE THE POINT ALONG THE DOGLEG
756 ! AT WHICH THE QUADRATIC IS MINIMIZED.
757 
758  bnorm = enorm(n, qtb)
759  temp = (bnorm/gnorm) * (bnorm/qnorm) * (sgnorm/delta)
760  temp = temp - (delta/qnorm) * (sgnorm/delta) ** 2 + sqrt(( &
761  temp-(delta/qnorm))**2+(1.0-(delta/qnorm)**2)*(1.0-( sgnorm/delta)**2))
762  alpha = ((delta/qnorm)*(1.0-(sgnorm/delta)**2)) / temp
763  END IF
764  END IF
765 
766 ! FORM APPROPRIATE CONVEX COMBINATION OF THE GAUSS-NEWTON
767 ! DIRECTION AND THE SCALED GRADIENT DIRECTION.
768 
769  temp = (1.0-alpha) * min(sgnorm,delta)
770  DO j = 1, n
771  x(j) = temp * wa1(j) + alpha * x(j)
772  END DO
773 END IF
774 RETURN
775 END SUBROUTINE dogleg
776 
777 
778 SUBROUTINE fdjac1(fcn, n, x, fvec, fjac, ldfjac, iflag, ml, mu, epsfcn, &
779  wa1, wa2)
780 
781 INTEGER, INTENT(IN) :: n
782 REAL (dp), INTENT(IN OUT) :: x(n)
783 REAL (dp), INTENT(IN) :: fvec(n)
784 INTEGER, INTENT(IN) :: ldfjac
785 REAL (dp), INTENT(OUT) :: fjac(ldfjac,n)
786 INTEGER, INTENT(IN OUT) :: iflag
787 INTEGER, INTENT(IN) :: ml
788 INTEGER, INTENT(IN) :: mu
789 REAL (dp), INTENT(IN) :: epsfcn
790 REAL (dp), INTENT(IN OUT) :: wa1(n)
791 REAL (dp), INTENT(OUT) :: wa2(n)
792 
793 ! EXTERNAL fcn
794 INTERFACE
795  SUBROUTINE fcn(N, X, FVEC, IFLAG)
796  IMPLICIT NONE
797  INTEGER, PARAMETER :: dp = selected_real_kind(14, 60)
798  INTEGER, INTENT(IN) :: n
799  REAL (dp), INTENT(IN) :: x(n)
800  REAL (dp), INTENT(OUT) :: fvec(n)
801  INTEGER, INTENT(IN OUT) :: iflag
802  END SUBROUTINE fcn
803 END INTERFACE
804 
805 ! **********
806 
807 ! SUBROUTINE FDJAC1
808 
809 ! THIS SUBROUTINE COMPUTES A FORWARD-DIFFERENCE APPROXIMATION TO THE N BY N
810 ! JACOBIAN MATRIX ASSOCIATED WITH A SPECIFIED PROBLEM OF N FUNCTIONS IN N
811 ! VARIABLES. IF THE JACOBIAN HAS A BANDED FORM, THEN FUNCTION EVALUATIONS
812 ! ARE SAVED BY ONLY APPROXIMATING THE NONZERO TERMS.
813 
814 ! THE SUBROUTINE STATEMENT IS
815 
816 ! SUBROUTINE FDJAC1(FCN,N,X,FVEC,FJAC,LDFJAC,IFLAG,ML,MU,EPSFCN,
817 ! WA1,WA2)
818 
819 ! WHERE
820 
821 ! FCN IS THE NAME OF THE USER-SUPPLIED SUBROUTINE WHICH CALCULATES
822 ! THE FUNCTIONS. FCN MUST BE DECLARED IN AN EXTERNAL STATEMENT IN
823 ! THE USER CALLING PROGRAM, AND SHOULD BE WRITTEN AS FOLLOWS.
824 
825 ! SUBROUTINE FCN(N,X,FVEC,IFLAG)
826 ! INTEGER N,IFLAG
827 ! REAL X(N),FVEC(N)
828 ! ----------
829 ! CALCULATE THE FUNCTIONS AT X AND
830 ! RETURN THIS VECTOR IN FVEC.
831 ! ----------
832 ! RETURN
833 ! END
834 
835 ! THE VALUE OF IFLAG SHOULD NOT BE CHANGED BY FCN UNLESS
836 ! THE USER WANTS TO TERMINATE EXECUTION OF FDJAC1.
837 ! IN THIS CASE SET IFLAG TO A NEGATIVE INTEGER.
838 
839 ! N IS A POSITIVE INTEGER INPUT VARIABLE SET TO THE NUMBER
840 ! OF FUNCTIONS AND VARIABLES.
841 
842 ! X IS AN INPUT ARRAY OF LENGTH N.
843 
844 ! FVEC IS AN INPUT ARRAY OF LENGTH N WHICH MUST CONTAIN THE
845 ! FUNCTIONS EVALUATED AT X.
846 
847 ! FJAC IS AN OUTPUT N BY N ARRAY WHICH CONTAINS THE
848 ! APPROXIMATION TO THE JACOBIAN MATRIX EVALUATED AT X.
849 
850 ! LDFJAC IS A POSITIVE INTEGER INPUT VARIABLE NOT LESS THAN N
851 ! WHICH SPECIFIES THE LEADING DIMENSION OF THE ARRAY FJAC.
852 
853 ! IFLAG IS AN INTEGER VARIABLE WHICH CAN BE USED TO TERMINATE
854 ! THE EXECUTION OF FDJAC1. SEE DESCRIPTION OF FCN.
855 
856 ! ML IS A NONNEGATIVE INTEGER INPUT VARIABLE WHICH SPECIFIES
857 ! THE NUMBER OF SUBDIAGONALS WITHIN THE BAND OF THE
858 ! JACOBIAN MATRIX. IF THE JACOBIAN IS NOT BANDED, SET
859 ! ML TO AT LEAST N - 1.
860 
861 ! EPSFCN IS AN INPUT VARIABLE USED IN DETERMINING A SUITABLE
862 ! STEP LENGTH FOR THE FORWARD-DIFFERENCE APPROXIMATION. THIS
863 ! APPROXIMATION ASSUMES THAT THE RELATIVE ERRORS IN THE
864 ! FUNCTIONS ARE OF THE ORDER OF EPSFCN. IF EPSFCN IS LESS
865 ! THAN THE MACHINE PRECISION, IT IS ASSUMED THAT THE RELATIVE
866 ! ERRORS IN THE FUNCTIONS ARE OF THE ORDER OF THE MACHINE PRECISION.
867 
868 ! MU IS A NONNEGATIVE INTEGER INPUT VARIABLE WHICH SPECIFIES
869 ! THE NUMBER OF SUPERDIAGONALS WITHIN THE BAND OF THE
870 ! JACOBIAN MATRIX. IF THE JACOBIAN IS NOT BANDED, SET
871 ! MU TO AT LEAST N - 1.
872 
873 ! WA1 AND WA2 ARE WORK ARRAYS OF LENGTH N. IF ML + MU + 1 IS AT
874 ! LEAST N, THEN THE JACOBIAN IS CONSIDERED DENSE, AND WA2 IS
875 ! NOT REFERENCED.
876 
877 ! SUBPROGRAMS CALLED
878 
879 ! MINPACK-SUPPLIED ... SPMPAR
880 
881 ! FORTRAN-SUPPLIED ... ABS,MAX,SQRT
882 
883 ! ARGONNE NATIONAL LABORATORY. MINPACK PROJECT. MARCH 1980.
884 ! BURTON S. GARBOW, KENNETH E. HILLSTROM, JORGE J. MORE
885 
886 ! **********
887 INTEGER :: i, j, k, msum
888 REAL (dp) :: eps, epsmch, h, temp
889 REAL (dp), PARAMETER :: zero = 0.0_dp
890 
891 ! EPSMCH IS THE MACHINE PRECISION.
892 
893 epsmch = epsilon(1.0_dp)
894 
895 eps = sqrt(max(epsfcn, epsmch))
896 msum = ml + mu + 1
897 IF (msum >= n) THEN
898 
899 ! COMPUTATION OF DENSE APPROXIMATE JACOBIAN.
900 
901  DO j = 1, n
902  temp = x(j)
903  h = eps * abs(temp)
904  IF (h == zero) h = eps
905  x(j) = temp + h
906  CALL fcn(n, x, wa1, iflag)
907  IF (iflag < 0) EXIT
908  x(j) = temp
909  DO i = 1, n
910  fjac(i,j) = (wa1(i)-fvec(i)) / h
911  END DO
912  END DO
913 ELSE
914 
915 ! COMPUTATION OF BANDED APPROXIMATE JACOBIAN.
916 
917  DO k = 1, msum
918  DO j = k, n, msum
919  wa2(j) = x(j)
920  h = eps * abs(wa2(j))
921  IF (h == zero) h = eps
922  x(j) = wa2(j) + h
923  END DO
924  CALL fcn(n, x, wa1, iflag)
925  IF (iflag < 0) EXIT
926  DO j = k, n, msum
927  x(j) = wa2(j)
928  h = eps * abs(wa2(j))
929  IF (h == zero) h = eps
930  DO i = 1, n
931  fjac(i,j) = zero
932  IF (i >= j-mu .AND. i <= j+ml) fjac(i,j) = (wa1(i)-fvec(i)) / h
933  END DO
934  END DO
935  END DO
936 END IF
937 RETURN
938 
939 ! LAST CARD OF SUBROUTINE FDJAC1.
940 
941 END SUBROUTINE fdjac1
942 
943 
944 
945 SUBROUTINE qform(m, n, q, ldq, wa)
946 
947 INTEGER, INTENT(IN) :: m
948 INTEGER, INTENT(IN) :: n
949 INTEGER, INTENT(IN) :: ldq
950 REAL (dp), INTENT(OUT) :: q(ldq,m)
951 REAL (dp), INTENT(OUT) :: wa(m)
952 
953 
954 ! **********
955 
956 ! SUBROUTINE QFORM
957 
958 ! THIS SUBROUTINE PROCEEDS FROM THE COMPUTED QR FACTORIZATION OF AN M BY N
959 ! MATRIX A TO ACCUMULATE THE M BY M ORTHOGONAL MATRIX Q FROM ITS FACTORED FORM.
960 
961 ! THE SUBROUTINE STATEMENT IS
962 
963 ! SUBROUTINE QFORM(M,N,Q,LDQ,WA)
964 
965 ! WHERE
966 
967 ! M IS A POSITIVE INTEGER INPUT VARIABLE SET TO THE NUMBER
968 ! OF ROWS OF A AND THE ORDER OF Q.
969 
970 ! N IS A POSITIVE INTEGER INPUT VARIABLE SET TO THE NUMBER OF COLUMNS OF A.
971 
972 ! Q IS AN M BY M ARRAY. ON INPUT THE FULL LOWER TRAPEZOID IN
973 ! THE FIRST MIN(M,N) COLUMNS OF Q CONTAINS THE FACTORED FORM.
974 ! ON OUTPUT Q HAS BEEN ACCUMULATED INTO A SQUARE MATRIX.
975 
976 ! LDQ IS A POSITIVE INTEGER INPUT VARIABLE NOT LESS THAN M
977 ! WHICH SPECIFIES THE LEADING DIMENSION OF THE ARRAY Q.
978 
979 ! WA IS A WORK ARRAY OF LENGTH M.
980 
981 ! SUBPROGRAMS CALLED
982 
983 ! FORTRAN-SUPPLIED ... MIN
984 
985 ! ARGONNE NATIONAL LABORATORY. MINPACK PROJECT. MARCH 1980.
986 ! BURTON S. GARBOW, KENNETH E. HILLSTROM, JORGE J. MORE
987 
988 ! **********
989 INTEGER :: i, j, jm1, k, l, minmn, np1
990 REAL (dp) :: sum, temp
991 REAL (dp), PARAMETER :: one = 1.0_dp, zero = 0.0_dp
992 
993 ! ZERO OUT UPPER TRIANGLE OF Q IN THE FIRST MIN(M,N) COLUMNS.
994 
995 minmn = min(m,n)
996 IF (minmn >= 2) THEN
997  DO j = 2, minmn
998  jm1 = j - 1
999  DO i = 1, jm1
1000  q(i,j) = zero
1001  END DO
1002  END DO
1003 END IF
1004 
1005 ! INITIALIZE REMAINING COLUMNS TO THOSE OF THE IDENTITY MATRIX.
1006 
1007 np1 = n + 1
1008 IF (m >= np1) THEN
1009  DO j = np1, m
1010  DO i = 1, m
1011  q(i,j) = zero
1012  END DO
1013  q(j,j) = one
1014  END DO
1015 END IF
1016 
1017 ! ACCUMULATE Q FROM ITS FACTORED FORM.
1018 
1019 DO l = 1, minmn
1020  k = minmn - l + 1
1021  DO i = k, m
1022  wa(i) = q(i,k)
1023  q(i,k) = zero
1024  END DO
1025  q(k,k) = one
1026  IF (wa(k) /= zero) THEN
1027  DO j = k, m
1028  sum = zero
1029  DO i = k, m
1030  sum = sum + q(i,j) * wa(i)
1031  END DO
1032  temp = sum / wa(k)
1033  DO i = k, m
1034  q(i,j) = q(i,j) - temp * wa(i)
1035  END DO
1036  END DO
1037  END IF
1038 END DO
1039 RETURN
1040 
1041 ! LAST CARD OF SUBROUTINE QFORM.
1042 
1043 END SUBROUTINE qform
1044 
1045 
1046 SUBROUTINE qrfac(m, n, a, lda, pivot, ipvt, lipvt, rdiag, acnorm, wa)
1047 
1048 INTEGER, INTENT(IN) :: m
1049 INTEGER, INTENT(IN) :: n
1050 INTEGER, INTENT(IN) :: lda
1051 REAL (dp), INTENT(IN OUT) :: a(lda,n)
1052 LOGICAL, INTENT(IN) :: pivot
1053 INTEGER, INTENT(IN) :: lipvt
1054 INTEGER, INTENT(OUT) :: ipvt(lipvt)
1055 REAL (dp), INTENT(OUT) :: rdiag(n)
1056 REAL (dp), INTENT(OUT) :: acnorm(n)
1057 REAL (dp), INTENT(OUT) :: wa(n)
1058 
1059 
1060 ! **********
1061 
1062 ! SUBROUTINE QRFAC
1063 
1064 ! THIS SUBROUTINE USES HOUSEHOLDER TRANSFORMATIONS WITH COLUMN PIVOTING
1065 ! (OPTIONAL) TO COMPUTE A QR FACTORIZATION OF THE M BY N MATRIX A.
1066 ! THAT IS, QRFAC DETERMINES AN ORTHOGONAL MATRIX Q, A PERMUTATION MATRIX P,
1067 ! AND AN UPPER TRAPEZOIDAL MATRIX R WITH DIAGONAL ELEMENTS OF NONINCREASING
1068 ! MAGNITUDE, SUCH THAT A*P = Q*R. THE HOUSEHOLDER TRANSFORMATION FOR
1069 ! COLUMN K, K = 1,2,...,MIN(M,N), IS OF THE FORM
1070 
1071 ! T
1072 ! I - (1/U(K))*U*U
1073 
1074 ! WHERE U HAS ZEROS IN THE FIRST K-1 POSITIONS. THE FORM OF THIS
1075 ! TRANSFORMATION AND THE METHOD OF PIVOTING FIRST APPEARED IN THE
1076 ! CORRESPONDING LINPACK SUBROUTINE.
1077 
1078 ! THE SUBROUTINE STATEMENT IS
1079 
1080 ! SUBROUTINE QRFAC(M,N,A,LDA,PIVOT,IPVT,LIPVT,RDIAG,ACNORM,WA)
1081 
1082 ! WHERE
1083 
1084 ! M IS A POSITIVE INTEGER INPUT VARIABLE SET TO THE NUMBER OF ROWS OF A.
1085 
1086 ! N IS A POSITIVE INTEGER INPUT VARIABLE SET TO THE NUMBER
1087 ! OF COLUMNS OF A.
1088 
1089 ! A IS AN M BY N ARRAY. ON INPUT A CONTAINS THE MATRIX FOR WHICH THE
1090 ! QR FACTORIZATION IS TO BE COMPUTED. ON OUTPUT THE STRICT UPPER
1091 ! TRAPEZOIDAL PART OF A CONTAINS THE STRICT UPPER TRAPEZOIDAL PART OF R,
1092 ! AND THE LOWER TRAPEZOIDAL PART OF A CONTAINS A FACTORED FORM OF Q
1093 ! (THE NON-TRIVIAL ELEMENTS OF THE U VECTORS DESCRIBED ABOVE).
1094 
1095 ! LDA IS A POSITIVE INTEGER INPUT VARIABLE NOT LESS THAN M
1096 ! WHICH SPECIFIES THE LEADING DIMENSION OF THE ARRAY A.
1097 
1098 ! PIVOT IS A LOGICAL INPUT VARIABLE. IF PIVOT IS SET TRUE,
1099 ! THEN COLUMN PIVOTING IS ENFORCED. IF PIVOT IS SET FALSE,
1100 ! THEN NO COLUMN PIVOTING IS DONE.
1101 
1102 ! IPVT IS AN INTEGER OUTPUT ARRAY OF LENGTH LIPVT. IPVT DEFINES THE
1103 ! PERMUTATION MATRIX P SUCH THAT A*P = Q*R.
1104 ! COLUMN J OF P IS COLUMN IPVT(J) OF THE IDENTITY MATRIX.
1105 ! IF PIVOT IS FALSE, IPVT IS NOT REFERENCED.
1106 
1107 ! LIPVT IS A POSITIVE INTEGER INPUT VARIABLE. IF PIVOT IS FALSE,
1108 ! THEN LIPVT MAY BE AS SMALL AS 1. IF PIVOT IS TRUE, THEN
1109 ! LIPVT MUST BE AT LEAST N.
1110 
1111 ! RDIAG IS AN OUTPUT ARRAY OF LENGTH N WHICH CONTAINS THE
1112 ! DIAGONAL ELEMENTS OF R.
1113 
1114 ! ACNORM IS AN OUTPUT ARRAY OF LENGTH N WHICH CONTAINS THE NORMS OF
1115 ! THE CORRESPONDING COLUMNS OF THE INPUT MATRIX A.
1116 ! IF THIS INFORMATION IS NOT NEEDED, THEN ACNORM CAN COINCIDE WITH RDIAG.
1117 
1118 ! WA IS A WORK ARRAY OF LENGTH N. IF PIVOT IS FALSE, THEN WA
1119 ! CAN COINCIDE WITH RDIAG.
1120 
1121 ! SUBPROGRAMS CALLED
1122 
1123 ! MINPACK-SUPPLIED ... SPMPAR,ENORM
1124 
1125 ! FORTRAN-SUPPLIED ... MAX,SQRT,MIN
1126 
1127 ! ARGONNE NATIONAL LABORATORY. MINPACK PROJECT. MARCH 1980.
1128 ! BURTON S. GARBOW, KENNETH E. HILLSTROM, JORGE J. MORE
1129 
1130 ! **********
1131 INTEGER :: i, j, jp1, k, kmax, minmn
1132 REAL (dp) :: ajnorm, epsmch, sum, temp
1133 REAL (dp), PARAMETER :: one = 1.0_dp, p05 = 0.05_dp, zero = 0.0_dp
1134 
1135 ! EPSMCH IS THE MACHINE PRECISION.
1136 
1137 epsmch = epsilon(1.0_dp)
1138 
1139 ! COMPUTE THE INITIAL COLUMN NORMS AND INITIALIZE SEVERAL ARRAYS.
1140 
1141 DO j = 1, n
1142  acnorm(j) = enorm(m, a(1:,j))
1143  rdiag(j) = acnorm(j)
1144  wa(j) = rdiag(j)
1145  IF (pivot) ipvt(j) = j
1146 END DO
1147 
1148 ! REDUCE A TO R WITH HOUSEHOLDER TRANSFORMATIONS.
1149 
1150 minmn = min(m,n)
1151 DO j = 1, minmn
1152  IF (pivot) THEN
1153 
1154 ! BRING THE COLUMN OF LARGEST NORM INTO THE PIVOT POSITION.
1155 
1156  kmax = j
1157  DO k = j, n
1158  IF (rdiag(k) > rdiag(kmax)) kmax = k
1159  END DO
1160  IF (kmax /= j) THEN
1161  DO i = 1, m
1162  temp = a(i,j)
1163  a(i,j) = a(i,kmax)
1164  a(i,kmax) = temp
1165  END DO
1166  rdiag(kmax) = rdiag(j)
1167  wa(kmax) = wa(j)
1168  k = ipvt(j)
1169  ipvt(j) = ipvt(kmax)
1170  ipvt(kmax) = k
1171  END IF
1172  END IF
1173 
1174 ! COMPUTE THE HOUSEHOLDER TRANSFORMATION TO REDUCE THE
1175 ! J-TH COLUMN OF A TO A MULTIPLE OF THE J-TH UNIT VECTOR.
1176 
1177  ajnorm = enorm(m-j+1, a(j:,j))
1178  IF (ajnorm /= zero) THEN
1179  IF (a(j,j) < zero) ajnorm = -ajnorm
1180  DO i = j, m
1181  a(i,j) = a(i,j) / ajnorm
1182  END DO
1183  a(j,j) = a(j,j) + one
1184 
1185 ! APPLY THE TRANSFORMATION TO THE REMAINING COLUMNS AND UPDATE THE NORMS.
1186 
1187  jp1 = j + 1
1188  IF (n >= jp1) THEN
1189  DO k = jp1, n
1190  sum = zero
1191  DO i = j, m
1192  sum = sum + a(i,j) * a(i,k)
1193  END DO
1194  temp = sum / a(j,j)
1195  DO i = j, m
1196  a(i,k) = a(i,k) - temp * a(i,j)
1197  END DO
1198  IF (.NOT.(.NOT.pivot.OR.rdiag(k) == zero)) THEN
1199  temp = a(j,k) / rdiag(k)
1200  rdiag(k) = rdiag(k) * sqrt(max(zero,one-temp**2))
1201  IF (p05*(rdiag(k)/wa(k))**2 <= epsmch) THEN
1202  rdiag(k) = enorm(m-j, a(jp1:,k))
1203  wa(k) = rdiag(k)
1204  END IF
1205  END IF
1206  END DO
1207  END IF
1208  END IF
1209  rdiag(j) = -ajnorm
1210 END DO
1211 RETURN
1212 
1213 ! LAST CARD OF SUBROUTINE QRFAC.
1214 
1215 END SUBROUTINE qrfac
1216 
1217 
1218 
1219 SUBROUTINE r1mpyq(m, n, a, lda, v, w)
1220 
1221 INTEGER, INTENT(IN) :: m
1222 INTEGER, INTENT(IN) :: n
1223 INTEGER, INTENT(IN) :: lda
1224 REAL (dp), INTENT(IN OUT) :: a(lda,n)
1225 REAL (dp), INTENT(IN) :: v(n)
1226 REAL (dp), INTENT(IN) :: w(n)
1227 
1228 
1229 ! **********
1230 
1231 ! SUBROUTINE R1MPYQ
1232 
1233 ! GIVEN AN M BY N MATRIX A, THIS SUBROUTINE COMPUTES A*Q WHERE
1234 ! Q IS THE PRODUCT OF 2*(N - 1) TRANSFORMATIONS
1235 
1236 ! GV(N-1)*...*GV(1)*GW(1)*...*GW(N-1)
1237 
1238 ! AND GV(I), GW(I) ARE GIVENS ROTATIONS IN THE (I,N) PLANE WHICH
1239 ! ELIMINATE ELEMENTS IN THE I-TH AND N-TH PLANES, RESPECTIVELY.
1240 ! Q ITSELF IS NOT GIVEN, RATHER THE INFORMATION TO RECOVER THE
1241 ! GV, GW ROTATIONS IS SUPPLIED.
1242 
1243 ! THE SUBROUTINE STATEMENT IS
1244 
1245 ! SUBROUTINE R1MPYQ(M, N, A, LDA, V, W)
1246 
1247 ! WHERE
1248 
1249 ! M IS A POSITIVE INTEGER INPUT VARIABLE SET TO THE NUMBER OF ROWS OF A.
1250 
1251 ! N IS A POSITIVE INTEGER INPUT VARIABLE SET TO THE NUMBER OF COLUMNS OF A.
1252 
1253 ! A IS AN M BY N ARRAY. ON INPUT A MUST CONTAIN THE MATRIX TO BE
1254 ! POSTMULTIPLIED BY THE ORTHOGONAL MATRIX Q DESCRIBED ABOVE.
1255 ! ON OUTPUT A*Q HAS REPLACED A.
1256 
1257 ! LDA IS A POSITIVE INTEGER INPUT VARIABLE NOT LESS THAN M
1258 ! WHICH SPECIFIES THE LEADING DIMENSION OF THE ARRAY A.
1259 
1260 ! V IS AN INPUT ARRAY OF LENGTH N. V(I) MUST CONTAIN THE INFORMATION
1261 ! NECESSARY TO RECOVER THE GIVENS ROTATION GV(I) DESCRIBED ABOVE.
1262 
1263 ! W IS AN INPUT ARRAY OF LENGTH N. W(I) MUST CONTAIN THE INFORMATION
1264 ! NECESSARY TO RECOVER THE GIVENS ROTATION GW(I) DESCRIBED ABOVE.
1265 
1266 ! SUBROUTINES CALLED
1267 
1268 ! FORTRAN-SUPPLIED ... ABS, SQRT
1269 
1270 ! ARGONNE NATIONAL LABORATORY. MINPACK PROJECT. MARCH 1980.
1271 ! BURTON S. GARBOW, KENNETH E. HILLSTROM, JORGE J. MORE
1272 
1273 ! **********
1274 INTEGER :: i, j, nmj, nm1
1275 REAL (dp) :: cos, sin, temp
1276 REAL (dp), PARAMETER :: one = 1.0_dp
1277 
1278 ! APPLY THE FIRST SET OF GIVENS ROTATIONS TO A.
1279 
1280 nm1 = n - 1
1281 IF (nm1 >= 1) THEN
1282  DO nmj = 1, nm1
1283  j = n - nmj
1284  IF (abs(v(j)) > one) cos = one / v(j)
1285  IF (abs(v(j)) > one) sin = sqrt(one-cos**2)
1286  IF (abs(v(j)) <= one) sin = v(j)
1287  IF (abs(v(j)) <= one) cos = sqrt(one-sin**2)
1288  DO i = 1, m
1289  temp = cos * a(i,j) - sin * a(i,n)
1290  a(i,n) = sin * a(i,j) + cos * a(i,n)
1291  a(i,j) = temp
1292  END DO
1293  END DO
1294 
1295 ! APPLY THE SECOND SET OF GIVENS ROTATIONS TO A.
1296 
1297  DO j = 1, nm1
1298  IF (abs(w(j)) > one) cos = one / w(j)
1299  IF (abs(w(j)) > one) sin = sqrt(one-cos**2)
1300  IF (abs(w(j)) <= one) sin = w(j)
1301  IF (abs(w(j)) <= one) cos = sqrt(one-sin**2)
1302  DO i = 1, m
1303  temp = cos * a(i,j) + sin * a(i,n)
1304  a(i,n) = -sin * a(i,j) + cos * a(i,n)
1305  a(i,j) = temp
1306  END DO
1307  END DO
1308 END IF
1309 RETURN
1310 
1311 ! LAST CARD OF SUBROUTINE R1MPYQ.
1312 
1313 END SUBROUTINE r1mpyq
1314 
1315 
1316 
1317 SUBROUTINE r1updt(m, n, s, ls, u, v, w, sing)
1318 
1319 INTEGER, INTENT(IN) :: m
1320 INTEGER, INTENT(IN) :: n
1321 INTEGER, INTENT(IN) :: ls
1322 REAL (dp), INTENT(IN OUT) :: s(ls)
1323 REAL (dp), INTENT(IN) :: u(m)
1324 REAL (dp), INTENT(IN OUT) :: v(n)
1325 REAL (dp), INTENT(OUT) :: w(m)
1326 LOGICAL, INTENT(OUT) :: sing
1327 
1328 
1329 ! **********
1330 
1331 ! SUBROUTINE R1UPDT
1332 
1333 ! GIVEN AN M BY N LOWER TRAPEZOIDAL MATRIX S, AN M-VECTOR U,
1334 ! AND AN N-VECTOR V, THE PROBLEM IS TO DETERMINE AN
1335 ! ORTHOGONAL MATRIX Q SUCH THAT
1336 
1337 ! T
1338 ! (S + U*V )*Q
1339 
1340 ! IS AGAIN LOWER TRAPEZOIDAL.
1341 
1342 ! THIS SUBROUTINE DETERMINES Q AS THE PRODUCT OF 2*(N - 1) TRANSFORMATIONS
1343 
1344 ! GV(N-1)*...*GV(1)*GW(1)*...*GW(N-1)
1345 
1346 ! WHERE GV(I), GW(I) ARE GIVENS ROTATIONS IN THE (I,N) PLANE
1347 ! WHICH ELIMINATE ELEMENTS IN THE I-TH AND N-TH PLANES, RESPECTIVELY.
1348 ! Q ITSELF IS NOT ACCUMULATED, RATHER THE INFORMATION TO RECOVER THE GV,
1349 ! GW ROTATIONS IS RETURNED.
1350 
1351 ! THE SUBROUTINE STATEMENT IS
1352 
1353 ! SUBROUTINE R1UPDT(M,N,S,LS,U,V,W,SING)
1354 
1355 ! WHERE
1356 
1357 ! M IS A POSITIVE INTEGER INPUT VARIABLE SET TO THE NUMBER OF ROWS OF S.
1358 
1359 ! N IS A POSITIVE INTEGER INPUT VARIABLE SET TO THE NUMBER
1360 ! OF COLUMNS OF S. N MUST NOT EXCEED M.
1361 
1362 ! S IS AN ARRAY OF LENGTH LS. ON INPUT S MUST CONTAIN THE LOWER
1363 ! TRAPEZOIDAL MATRIX S STORED BY COLUMNS. ON OUTPUT S CONTAINS
1364 ! THE LOWER TRAPEZOIDAL MATRIX PRODUCED AS DESCRIBED ABOVE.
1365 
1366 ! LS IS A POSITIVE INTEGER INPUT VARIABLE NOT LESS THAN
1367 ! (N*(2*M-N+1))/2.
1368 
1369 ! U IS AN INPUT ARRAY OF LENGTH M WHICH MUST CONTAIN THE VECTOR U.
1370 
1371 ! V IS AN ARRAY OF LENGTH N. ON INPUT V MUST CONTAIN THE VECTOR V.
1372 ! ON OUTPUT V(I) CONTAINS THE INFORMATION NECESSARY TO
1373 ! RECOVER THE GIVENS ROTATION GV(I) DESCRIBED ABOVE.
1374 
1375 ! W IS AN OUTPUT ARRAY OF LENGTH M. W(I) CONTAINS INFORMATION
1376 ! NECESSARY TO RECOVER THE GIVENS ROTATION GW(I) DESCRIBED ABOVE.
1377 
1378 ! SING IS A LOGICAL OUTPUT VARIABLE. SING IS SET TRUE IF ANY OF THE
1379 ! DIAGONAL ELEMENTS OF THE OUTPUT S ARE ZERO. OTHERWISE SING IS
1380 ! SET FALSE.
1381 
1382 ! SUBPROGRAMS CALLED
1383 
1384 ! MINPACK-SUPPLIED ... SPMPAR
1385 
1386 ! FORTRAN-SUPPLIED ... ABS,SQRT
1387 
1388 ! ARGONNE NATIONAL LABORATORY. MINPACK PROJECT. MARCH 1980.
1389 ! BURTON S. GARBOW, KENNETH E. HILLSTROM, JORGE J. MORE, JOHN L. NAZARETH
1390 
1391 ! **********
1392 INTEGER :: i, j, jj, l, nmj, nm1
1393 REAL (dp) :: cos, cotan, giant, sin, tan, tau, temp
1394 REAL (dp), PARAMETER :: one = 1.0_dp, p5 = 0.5_dp, p25 = 0.25_dp, zero = 0.0_dp
1395 
1396 ! GIANT IS THE LARGEST MAGNITUDE.
1397 
1398 giant = huge(1.0_dp)
1399 
1400 ! INITIALIZE THE DIAGONAL ELEMENT POINTER.
1401 
1402 jj = (n*(2*m-n+1)) / 2 - (m-n)
1403 
1404 ! MOVE THE NONTRIVIAL PART OF THE LAST COLUMN OF S INTO W.
1405 
1406 l = jj
1407 DO i = n, m
1408  w(i) = s(l)
1409  l = l + 1
1410 END DO
1411 
1412 ! ROTATE THE VECTOR V INTO A MULTIPLE OF THE N-TH UNIT VECTOR
1413 ! IN SUCH A WAY THAT A SPIKE IS INTRODUCED INTO W.
1414 
1415 nm1 = n - 1
1416 IF (nm1 >= 1) THEN
1417  DO nmj = 1, nm1
1418  j = n - nmj
1419  jj = jj - (m-j+1)
1420  w(j) = zero
1421  IF (v(j) /= zero) THEN
1422 
1423 ! DETERMINE A GIVENS ROTATION WHICH ELIMINATES THE J-TH ELEMENT OF V.
1424 
1425  IF (abs(v(n)) < abs(v(j))) THEN
1426  cotan = v(n) / v(j)
1427  sin = p5 / sqrt(p25+p25*cotan**2)
1428  cos = sin * cotan
1429  tau = one
1430  IF (abs(cos)*giant > one) tau = one / cos
1431  ELSE
1432  tan = v(j) / v(n)
1433  cos = p5 / sqrt(p25+p25*tan**2)
1434  sin = cos * tan
1435  tau = sin
1436  END IF
1437 
1438 ! APPLY THE TRANSFORMATION TO V AND STORE THE INFORMATION
1439 ! NECESSARY TO RECOVER THE GIVENS ROTATION.
1440 
1441  v(n) = sin * v(j) + cos * v(n)
1442  v(j) = tau
1443 
1444 ! APPLY THE TRANSFORMATION TO S AND EXTEND THE SPIKE IN W.
1445 
1446  l = jj
1447  DO i = j, m
1448  temp = cos * s(l) - sin * w(i)
1449  w(i) = sin * s(l) + cos * w(i)
1450  s(l) = temp
1451  l = l + 1
1452  END DO
1453  END IF
1454  END DO
1455 END IF
1456 
1457 ! ADD THE SPIKE FROM THE RANK 1 UPDATE TO W.
1458 
1459 DO i = 1, m
1460  w(i) = w(i) + v(n) * u(i)
1461 END DO
1462 
1463 ! ELIMINATE THE SPIKE.
1464 
1465 sing = .false.
1466 IF (nm1 >= 1) THEN
1467  DO j = 1, nm1
1468  IF (w(j) /= zero) THEN
1469 
1470 ! DETERMINE A GIVENS ROTATION WHICH ELIMINATES THE
1471 ! J-TH ELEMENT OF THE SPIKE.
1472 
1473  IF (abs(s(jj)) < abs(w(j))) THEN
1474  cotan = s(jj) / w(j)
1475  sin = p5 / sqrt(p25 + p25*cotan**2)
1476  cos = sin * cotan
1477  tau = one
1478  IF (abs(cos)*giant > one) tau = one / cos
1479  ELSE
1480  tan = w(j) / s(jj)
1481  cos = p5 / sqrt(p25+p25*tan**2)
1482  sin = cos * tan
1483  tau = sin
1484  END IF
1485 
1486 ! APPLY THE TRANSFORMATION TO S AND REDUCE THE SPIKE IN W.
1487 
1488  l = jj
1489  DO i = j, m
1490  temp = cos * s(l) + sin * w(i)
1491  w(i) = -sin * s(l) + cos * w(i)
1492  s(l) = temp
1493  l = l + 1
1494  END DO
1495 
1496 ! STORE THE INFORMATION NECESSARY TO RECOVER THE GIVENS ROTATION.
1497 
1498  w(j) = tau
1499  END IF
1500 
1501 ! TEST FOR ZERO DIAGONAL ELEMENTS IN THE OUTPUT S.
1502 
1503  IF (s(jj) == zero) sing = .true.
1504  jj = jj + (m-j+1)
1505  END DO
1506 END IF
1507 
1508 ! MOVE W BACK INTO THE LAST COLUMN OF THE OUTPUT S.
1509 
1510 l = jj
1511 DO i = n, m
1512  s(l) = w(i)
1513  l = l + 1
1514 END DO
1515 IF (s(jj) == zero) sing = .true.
1516 RETURN
1517 
1518 ! LAST CARD OF SUBROUTINE R1UPDT.
1519 
1520 END SUBROUTINE r1updt
1521 
1522 
1523 FUNCTION enorm(n, x) RESULT(fn_val)
1524 
1525 INTEGER, INTENT(IN) :: n
1526 REAL (dp), INTENT(IN) :: x(n)
1527 REAL (dp) :: fn_val
1528 
1529 ! **********
1530 
1531 ! FUNCTION ENORM
1532 
1533 ! GIVEN AN N-VECTOR X, THIS FUNCTION CALCULATES THE EUCLIDEAN NORM OF X.
1534 
1535 ! THE EUCLIDEAN NORM IS COMPUTED BY ACCUMULATING THE SUM OF SQUARES IN THREE
1536 ! DIFFERENT SUMS. THE SUMS OF SQUARES FOR THE SMALL AND LARGE COMPONENTS
1537 ! ARE SCALED SO THAT NO OVERFLOWS OCCUR. NON-DESTRUCTIVE UNDERFLOWS ARE
1538 ! PERMITTED. UNDERFLOWS AND OVERFLOWS DO NOT OCCUR IN THE COMPUTATION OF THE UNSCALED
1539 ! SUM OF SQUARES FOR THE INTERMEDIATE COMPONENTS.
1540 ! THE DEFINITIONS OF SMALL, INTERMEDIATE AND LARGE COMPONENTS DEPEND ON
1541 ! TWO CONSTANTS, RDWARF AND RGIANT. THE MAIN RESTRICTIONS ON THESE CONSTANTS
1542 ! ARE THAT RDWARF**2 NOT UNDERFLOW AND RGIANT**2 NOT OVERFLOW.
1543 ! THE CONSTANTS GIVEN HERE ARE SUITABLE FOR EVERY KNOWN COMPUTER.
1544 
1545 ! THE FUNCTION STATEMENT IS
1546 
1547 ! REAL FUNCTION ENORM(N, X)
1548 
1549 ! WHERE
1550 
1551 ! N IS A POSITIVE INTEGER INPUT VARIABLE.
1552 
1553 ! X IS AN INPUT ARRAY OF LENGTH N.
1554 
1555 ! SUBPROGRAMS CALLED
1556 
1557 ! FORTRAN-SUPPLIED ... ABS,SQRT
1558 
1559 ! ARGONNE NATIONAL LABORATORY. MINPACK PROJECT. MARCH 1980.
1560 ! BURTON S. GARBOW, KENNETH E. HILLSTROM, JORGE J. MORE
1561 
1562 ! **********
1563 INTEGER :: i
1564 REAL (dp) :: agiant, floatn, s1, s2, s3, xabs, x1max, x3max
1565 REAL (dp), PARAMETER :: rdwarf = 1.0d-100, rgiant = 1.0d+100
1566 
1567 s1 = 0.0_dp
1568 s2 = 0.0_dp
1569 s3 = 0.0_dp
1570 x1max = 0.0_dp
1571 x3max = 0.0_dp
1572 floatn = n
1573 agiant = rgiant / floatn
1574 DO i = 1, n
1575  xabs = abs(x(i))
1576  IF (xabs <= rdwarf .OR. xabs >= agiant) THEN
1577  IF (xabs > rdwarf) THEN
1578 
1579 ! SUM FOR LARGE COMPONENTS.
1580 
1581  IF (xabs > x1max) THEN
1582  s1 = 1.0_dp + s1 * (x1max/xabs) ** 2
1583  x1max = xabs
1584  ELSE
1585  s1 = s1 + (xabs/x1max) ** 2
1586  END IF
1587  ELSE
1588 
1589 ! SUM FOR SMALL COMPONENTS.
1590 
1591  IF (xabs > x3max) THEN
1592  s3 = 1.0_dp + s3 * (x3max/xabs) ** 2
1593  x3max = xabs
1594  ELSE
1595  IF (xabs /= 0.0_dp) s3 = s3 + (xabs/x3max) ** 2
1596  END IF
1597  END IF
1598  ELSE
1599 
1600 ! SUM FOR INTERMEDIATE COMPONENTS.
1601 
1602  s2 = s2 + xabs ** 2
1603  END IF
1604 END DO
1605 
1606 ! CALCULATION OF NORM.
1607 
1608 IF (s1 /= 0.0_dp) THEN
1609  fn_val = x1max * sqrt(s1 + (s2/x1max)/x1max)
1610 ELSE
1611  IF (s2 /= 0.0_dp) THEN
1612  IF (s2 >= x3max) fn_val = sqrt(s2*(1.0_dp + (x3max/s2)*(x3max*s3)))
1613  IF (s2 < x3max) fn_val = sqrt(x3max*((s2/x3max) + (x3max*s3)))
1614  ELSE
1615  fn_val = x3max * sqrt(s3)
1616  END IF
1617 END IF
1618 RETURN
1619 END FUNCTION enorm
1620 
1621 END MODULE solve_nonlin