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