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