Purpose
To compute a set of parameters for approximating a Wiener system in a least-squares sense, using a neural network approach and a MINPACK-like Levenberg-Marquardt algorithm. The Wiener system consists of a linear part and a static nonlinearity, and it is represented as x(t+1) = A*x(t) + B*u(t) z(t) = C*x(t) + D*u(t), y(t) = f(z(t),wb(1:L)), where t = 1, 2, ..., NSMP, and f is a nonlinear function, evaluated by the SLICOT Library routine NF01AY. The parameter vector X is partitioned as X = ( wb(1), ..., wb(L), theta ), where theta corresponds to the linear part, and wb(i), i = 1 : L, correspond to the nonlinear part. See SLICOT Library routine NF01AD for further details. The sum of squares of the error functions, defined by e(t) = y(t) - Y(t), t = 1, 2, ..., NSMP, is minimized, where Y(t) is the measured output vector. The functions and their Jacobian matrices are evaluated by SLICOT Library routine NF01BF (the FCN routine in the call of MD03BD).Specification
SUBROUTINE IB03BD( INIT, NOBR, M, L, NSMP, N, NN, ITMAX1, ITMAX2, $ NPRINT, U, LDU, Y, LDY, X, LX, TOL1, TOL2, $ IWORK, DWORK, LDWORK, IWARN, INFO ) C .. Scalar Arguments .. CHARACTER INIT INTEGER INFO, ITMAX1, ITMAX2, IWARN, L, LDU, LDWORK, $ LDY, LX, M, N, NN, NOBR, NPRINT, NSMP DOUBLE PRECISION TOL1, TOL2 C .. Array Arguments .. DOUBLE PRECISION DWORK(*), U(LDU, *), X(*), Y(LDY, *) INTEGER IWORK(*)Arguments
Mode Parameters
INIT CHARACTER*1 Specifies which parts have to be initialized, as follows: = 'L' : initialize the linear part only, X already contains an initial approximation of the nonlinearity; = 'S' : initialize the static nonlinearity only, X already contains an initial approximation of the linear part; = 'B' : initialize both linear and nonlinear parts; = 'N' : do not initialize anything, X already contains an initial approximation. If INIT = 'S' or 'B', the error functions for the nonlinear part, and their Jacobian matrices, are evaluated by SLICOT Library routine NF01BE (used as a second FCN routine in the MD03BD call for the initialization step, see METHOD).Input/Output Parameters
NOBR (input) INTEGER If INIT = 'L' or 'B', NOBR is the number of block rows, s, in the input and output block Hankel matrices to be processed for estimating the linear part. NOBR > 0. (In the MOESP theory, NOBR should be larger than n, the estimated dimension of state vector.) This parameter is ignored if INIT is 'S' or 'N'. M (input) INTEGER The number of system inputs. M >= 0. L (input) INTEGER The number of system outputs. L >= 0, and L > 0, if INIT = 'L' or 'B'. NSMP (input) INTEGER The number of input and output samples, t. NSMP >= 0, and NSMP >= 2*(M+L+1)*NOBR - 1, if INIT = 'L' or 'B'. N (input/output) INTEGER The order of the linear part. If INIT = 'L' or 'B', and N < 0 on entry, the order is assumed unknown and it will be found by the routine. Otherwise, the input value will be used. If INIT = 'S' or 'N', N must be non-negative. The values N >= NOBR, or N = 0, are not acceptable if INIT = 'L' or 'B'. NN (input) INTEGER The number of neurons which shall be used to approximate the nonlinear part. NN >= 0. ITMAX1 (input) INTEGER The maximum number of iterations for the initialization of the static nonlinearity. This parameter is ignored if INIT is 'N' or 'L'. Otherwise, ITMAX1 >= 0. ITMAX2 (input) INTEGER The maximum number of iterations. ITMAX2 >= 0. NPRINT (input) INTEGER This parameter enables controlled printing of iterates if it is positive. In this case, FCN is called with IFLAG = 0 at the beginning of the first iteration and every NPRINT iterations thereafter and immediately prior to return, and the current error norm is printed. Other intermediate results could be printed by modifying the corresponding FCN routine (NF01BE and/or NF01BF). If NPRINT <= 0, no special calls of FCN with IFLAG = 0 are made. U (input) DOUBLE PRECISION array, dimension (LDU, M) The leading NSMP-by-M part of this array must contain the set of input samples, U = ( U(1,1),...,U(1,M); ...; U(NSMP,1),...,U(NSMP,M) ). LDU INTEGER The leading dimension of array U. LDU >= MAX(1,NSMP). Y (input) DOUBLE PRECISION array, dimension (LDY, L) The leading NSMP-by-L part of this array must contain the set of output samples, Y = ( Y(1,1),...,Y(1,L); ...; Y(NSMP,1),...,Y(NSMP,L) ). LDY INTEGER The leading dimension of array Y. LDY >= MAX(1,NSMP). X (input/output) DOUBLE PRECISION array dimension (LX) On entry, if INIT = 'L', the leading (NN*(L+2) + 1)*L part of this array must contain the initial parameters for the nonlinear part of the system. On entry, if INIT = 'S', the elements lin1 : lin2 of this array must contain the initial parameters for the linear part of the system, corresponding to the output normal form, computed by SLICOT Library routine TB01VD, where lin1 = (NN*(L+2) + 1)*L + 1; lin2 = (NN*(L+2) + 1)*L + N*(L+M+1) + L*M. On entry, if INIT = 'N', the elements 1 : lin2 of this array must contain the initial parameters for the nonlinear part followed by the initial parameters for the linear part of the system, as specified above. This array need not be set on entry if INIT = 'B'. On exit, the elements 1 : lin2 of this array contain the optimal parameters for the nonlinear part followed by the optimal parameters for the linear part of the system, as specified above. LX (input/output) INTEGER On entry, this parameter must contain the intended length of X. If N >= 0, then LX >= NX := lin2 (see parameter X). If N is unknown (N < 0 on entry), a large enough estimate of N should be used in the formula of lin2. On exit, if N < 0 on entry, but LX is not large enough, then this parameter contains the actual length of X, corresponding to the computed N. Otherwise, its value is unchanged.Tolerances
TOL1 DOUBLE PRECISION If INIT = 'S' or 'B' and TOL1 >= 0, TOL1 is the tolerance which measures the relative error desired in the sum of squares, as well as the relative error desired in the approximate solution, for the initialization step of nonlinear part. Termination occurs when either both the actual and predicted relative reductions in the sum of squares, or the relative error between two consecutive iterates are at most TOL1. If the user sets TOL1 < 0, then SQRT(EPS) is used instead TOL1, where EPS is the machine precision (see LAPACK Library routine DLAMCH). This parameter is ignored if INIT is 'N' or 'L'. TOL2 DOUBLE PRECISION If TOL2 >= 0, TOL2 is the tolerance which measures the relative error desired in the sum of squares, as well as the relative error desired in the approximate solution, for the whole optimization process. Termination occurs when either both the actual and predicted relative reductions in the sum of squares, or the relative error between two consecutive iterates are at most TOL2. If the user sets TOL2 < 0, then SQRT(EPS) is used instead TOL2. This default value could require many iterations, especially if TOL1 is larger. If INIT = 'S' or 'B', it is advisable that TOL2 be larger than TOL1, and spend more time with cheaper iterations.Workspace
IWORK INTEGER array, dimension (MAX( LIW1, LIW2, LIW3 )), where LIW1 = LIW2 = 0, if INIT = 'S' or 'N'; otherwise, LIW1 = M+L; LIW2 = MAX(M*NOBR+N,M*(N+L)); LIW3 = 3+MAX(NN*(L+2)+2,NX+L), if INIT = 'S' or 'B'; LIW3 = 3+NX+L, if INIT = 'L' or 'N'. On output, if INFO = 0, IWORK(1) and IWORK(2) return the (total) number of function and Jacobian evaluations, respectively (including the initialization step, if it was performed), and if INIT = 'L' or INIT = 'B', IWORK(3) specifies how many locations of DWORK contain reciprocal condition number estimates (see below); otherwise, IWORK(3) = 0. If INFO = 0, the entries 4 to 3+NX of IWORK define a permutation matrix P such that J*P = Q*R, where J is the final calculated Jacobian, Q is an orthogonal matrix (not stored), and R is upper triangular with diagonal elements of nonincreasing magnitude (possibly for each block column of J). Column j of P is column IWORK(3+j) of the identity matrix. Moreover, the entries 4+NX:3+NX+L of this array contain the ranks of the final submatrices S_k (see description of LMPARM in MD03BD). DWORK DOUBLE PRECISION array, dimension (LDWORK) On entry, if desired, and if INIT = 'S' or 'B', the entries DWORK(1:4) are set to initialize the random numbers generator for the nonlinear part parameters (see the description of the argument XINIT of SLICOT Library routine MD03BD); this enables to obtain reproducible results. The same seed is used for all outputs. On exit, if INFO = 0, DWORK(1) returns the optimal value of LDWORK, DWORK(2) returns the residual error norm (the sum of squares), DWORK(3) returns the number of iterations performed, and DWORK(4) returns the final Levenberg factor, for optimizing the parameters of both the linear part and the static nonlinearity part. If INIT = 'S' or INIT = 'B' and INFO = 0, then the elements DWORK(5) to DWORK(8) contain the corresponding four values for the initialization step (see METHOD). (If L > 1, DWORK(8) contains the maximum of the Levenberg factors for all outputs.) If INIT = 'L' or INIT = 'B', and INFO = 0, DWORK(9) to DWORK(8+IWORK(3)) contain reciprocal condition number estimates set by SLICOT Library routines IB01AD, IB01BD, and IB01CD. On exit, if INFO = -21, DWORK(1) returns the minimum value of LDWORK. LDWORK INTEGER The length of the array DWORK. In the formulas below, N should be taken not larger than NOBR - 1, if N < 0 on entry. LDWORK = MAX( LW1, LW2, LW3, LW4 ), where LW1 = 0, if INIT = 'S' or 'N'; otherwise, LW1 = MAX( 2*(M+L)*NOBR*(2*(M+L)*(NOBR+1)+3) + L*NOBR, 4*(M+L)*NOBR*(M+L)*NOBR + (N+L)*(N+M) + MAX( LDW1, LDW2 ), (N+L)*(N+M) + N + N*N + 2 + N*(N+M+L) + MAX( 5*N, 2, MIN( LDW3, LDW4 ), LDW5, LDW6 ), where, LDW1 >= MAX( 2*(L*NOBR-L)*N+2*N, (L*NOBR-L)*N+N*N+7*N, L*NOBR*N + MAX( (L*NOBR-L)*N+2*N + (2*M+L)*NOBR+L, 2*(L*NOBR-L)*N+N*N+8*N, N+4*(M*NOBR+N)+1, M*NOBR+3*N+L ) ) LDW2 >= 0, if M = 0; LDW2 >= L*NOBR*N + M*NOBR*(N+L)*(M*(N+L)+1) + MAX( (N+L)**2, 4*M*(N+L)+1 ), if M > 0; LDW3 = NSMP*L*(N+1) + 2*N + MAX( 2*N*N, 4*N ), LDW4 = N*(N+1) + 2*N + MAX( N*L*(N+1) + 2*N*N + L*N, 4*N ); LDW5 = NSMP*L + (N+L)*(N+M) + 3*N+M+L; LDW6 = NSMP*L + (N+L)*(N+M) + N + MAX(1, N*N*L + N*L + N, N*N + MAX(N*N + N*MAX(N,L) + 6*N + MIN(N,L), N*M)); LW2 = LW3 = 0, if INIT = 'L' or 'N'; otherwise, LW2 = NSMP*L + BSN + MAX( 4, NSMP + MAX( NSMP*BSN + MAX( 2*NN, 5*BSN + 1 ), BSN**2 + BSN + MAX( NSMP + 2*NN, 5*BSN ) ) ); LW3 = MAX( LDW7, NSMP*L + (N+L)*(2*N+M) + 2*N ); LDW7 = NSMP*L + (N+L)*(N+M) + 3*N+M+L, if M > 0; LDW7 = NSMP*L + (N+L)*N + 2*N+L, if M = 0; LW4 = NSMP*L + NX + MAX( 4, NSMP*L + MAX( NSMP*L*( BSN + LTHS ) + MAX( NSMP*L + L1, L2 + NX ), NX*( BSN + LTHS ) + NX + MAX( NSMP*L + L1, NX + L3 ) ) ), L0 = MAX( N*(N+L), N+M+L ), if M > 0; L0 = MAX( N*(N+L), L ), if M = 0; L1 = NSMP*L + MAX( 2*NN, (N+L)*(N+M) + 2*N + L0); L2 = 4*NX + 1, if L <= 1 or BSN = 0; otherwise, L2 = BSN + MAX(3*BSN+1,LTHS); L2 = MAX(L2,4*LTHS+1), if NSMP > BSN; L2 = MAX(L2,(NSMP-BSN)*(L-1)), if BSN < NSMP < 2*BSN; L3 = 4*NX, if L <= 1 or BSN = 0; L3 = LTHS*BSN + 2*NX + 2*MAX(BSN,LTHS), if L > 1 and BSN > 0, with BSN = NN*( L + 2 ) + 1, LTHS = N*( L + M + 1 ) + L*M. For optimum performance LDWORK should be larger.Warning Indicator
IWARN INTEGER < 0: the user set IFLAG = IWARN in (one of) the subroutine(s) FCN, i.e., NF01BE, if INIT = 'S' or 'B', and/or NF01BF; this value cannot be returned without changing the FCN routine(s); otherwise, IWARN has the value k*100 + j*10 + i, where k is defined below, i refers to the whole optimization process, and j refers to the initialization step (j = 0, if INIT = 'L' or 'N'), and the possible values for i and j have the following meaning (where TOL* denotes TOL1 or TOL2, and similarly for ITMAX*): = 1: both actual and predicted relative reductions in the sum of squares are at most TOL*; = 2: relative error between two consecutive iterates is at most TOL*; = 3: conditions for i or j = 1 and i or j = 2 both hold; = 4: the cosine of the angle between the vector of error function values and any column of the Jacobian is at most EPS in absolute value; = 5: the number of iterations has reached ITMAX* without satisfying any convergence condition; = 6: TOL* is too small: no further reduction in the sum of squares is possible; = 7: TOL* is too small: no further improvement in the approximate solution X is possible; = 8: the vector of function values e is orthogonal to the columns of the Jacobian to machine precision. The digit k is normally 0, but if INIT = 'L' or 'B', it can have a value in the range 1 to 6 (see IB01AD, IB01BD and IB01CD). In all these cases, the entries DWORK(1:4), DWORK(5:8) (if INIT = 'S' or 'B'), and DWORK(9:8+IWORK(3)) (if INIT = 'L' or 'B'), are set as described above.Error Indicator
INFO INTEGER = 0: successful exit; < 0: if INFO = -i, the i-th argument had an illegal value; otherwise, INFO has the value k*100 + j*10 + i, where k is defined below, i refers to the whole optimization process, and j refers to the initialization step (j = 0, if INIT = 'L' or 'N'), and the possible values for i and j have the following meaning: = 1: the routine FCN returned with INFO <> 0 for IFLAG = 1; = 2: the routine FCN returned with INFO <> 0 for IFLAG = 2; = 3: the routine QRFACT returned with INFO <> 0; = 4: the routine LMPARM returned with INFO <> 0. In addition, if INIT = 'L' or 'B', i could also be = 5: if a Lyapunov equation could not be solved; = 6: if the identified linear system is unstable; = 7: if the QR algorithm failed on the state matrix of the identified linear system. QRFACT and LMPARM are generic names for SLICOT Library routines NF01BS and NF01BP, respectively, for the whole optimization process, and MD03BA and MD03BB, respectively, for the initialization step (if INIT = 'S' or 'B'). The digit k is normally 0, but if INIT = 'L' or 'B', it can have a value in the range 1 to 10 (see IB01AD/IB01BD).Method
If INIT = 'L' or 'B', the linear part of the system is approximated using the combined MOESP and N4SID algorithm. If necessary, this algorithm can also choose the order, but it is advantageous if the order is already known. If INIT = 'S' or 'B', the output of the approximated linear part is computed and used to calculate an approximation of the static nonlinearity using the Levenberg-Marquardt algorithm [1,3]. This step is referred to as the (nonlinear) initialization step. As last step, the Levenberg-Marquardt algorithm is used again to optimize the parameters of the linear part and the static nonlinearity as a whole. Therefore, it is necessary to parametrise the matrices of the linear part. The output normal form [2] parameterisation is used. The Jacobian is computed analytically, for the nonlinear part, and numerically, for the linear part.References
[1] More, J.J., Garbow, B.S, and Hillstrom, K.E. User's Guide for MINPACK-1. Applied Math. Division, Argonne National Laboratory, Argonne, Illinois, Report ANL-80-74, 1980. [2] Peeters, R.L.M., Hanzon, B., and Olivi, M. Balanced realizations of discrete-time stable all-pass systems and the tangential Schur algorithm. Proceedings of the European Control Conference, 31 August - 3 September 1999, Karlsruhe, Germany. Session CP-6, Discrete-time Systems, 1999. [3] More, J.J. The Levenberg-Marquardt algorithm: implementation and theory. In Watson, G.A. (Ed.), Numerical Analysis, Lecture Notes in Mathematics, vol. 630, Springer-Verlag, Berlin, Heidelberg and New York, pp. 105-116, 1978.Numerical Aspects
The Levenberg-Marquardt algorithm described in [3] is scaling invariant and globally convergent to (maybe local) minima. The convergence rate near a local minimum is quadratic, if the Jacobian is computed analytically, and linear, if the Jacobian is computed numerically.Further Comments
NoneExample
Program Text
* IB03BD EXAMPLE PROGRAM TEXT * Copyright (c) 2002-2017 NICONET e.V. * * .. Parameters .. INTEGER NIN, NOUT PARAMETER ( NIN = 5, NOUT = 6 ) INTEGER BSNM, LDU, LDY, LIWORK, LMAX, LTHS, LXM, MMAX, $ NMAX, NNMAX, NOBRMX, NSMPMX PARAMETER ( LMAX = 2, MMAX = 3, NOBRMX = 10, NNMAX = 12, $ NMAX = 4, NSMPMX = 1024, $ BSNM = NNMAX*( LMAX + 2 ) + 1, $ LTHS = NMAX*( LMAX + MMAX + 1 ) + LMAX*MMAX, $ LDU = NSMPMX, LDY = NSMPMX, $ LXM = BSNM*LMAX + LTHS, $ LIWORK = MAX( MMAX + LMAX, MMAX*NOBRMX + NMAX, $ MMAX*( NMAX + LMAX ), 3 + $ MAX( BSNM + 1, LXM + LMAX ) ) ) INTEGER L0, L1M, L2M, L3M, LDW1, LDW2, LDW3, LDW4, LDW5, $ LDW6, LDW7, LDWORK, LW1, LW2, LW3, LW4 PARAMETER ( L0 = MAX( NMAX*( NMAX + LMAX ), $ NMAX + MMAX + LMAX ), $ L1M = NSMPMX*LMAX + $ MAX( 2*NNMAX, $ ( NMAX + LMAX )*( NMAX + MMAX ) + $ 2*NMAX + L0 ), $ L2M = MAX( 4*LXM + 1, BSNM + $ MAX( 3*BSNM + 1, LTHS ), $ NSMPMX*( LMAX - 1 ) ), $ L3M = MAX( 4*LXM, LTHS*BSNM + 2*LXM + $ 2*MAX( BSNM, LTHS ) ), $ LDW1 = MAX( 2*( LMAX*NOBRMX - LMAX )*NMAX + $ 2*NMAX, $ ( LMAX*NOBRMX - LMAX )*NMAX + $ NMAX*NMAX + 7*NMAX, $ LMAX*NOBRMX*NMAX + $ MAX( ( LMAX*NOBRMX - LMAX )*NMAX + $ 2*NMAX + LMAX + $ ( 2*MMAX + LMAX )*NOBRMX, $ 2*( LMAX*NOBRMX - LMAX )*NMAX $ + NMAX*NMAX + 8*NMAX, $ NMAX + 4*( MMAX*NOBRMX + $ NMAX ) + 1, $ MMAX*NOBRMX + 3*NMAX + LMAX ) $ ), $ LDW2 = LMAX*NOBRMX*NMAX + $ MMAX*NOBRMX*( NMAX + LMAX )* $ ( MMAX*( NMAX + LMAX ) + 1 ) + $ MAX( ( NMAX + LMAX )**2, $ 4*MMAX*( NMAX + LMAX ) + 1 ), $ LDW3 = NSMPMX*LMAX*( NMAX + 1 ) + 2*NMAX + $ MAX( 2*NMAX*NMAX, 4*NMAX ), $ LDW4 = NMAX*( NMAX + 1 ) + 2*NMAX + $ MAX( NMAX*LMAX*( NMAX + 1 ) + $ 2*NMAX*NMAX + LMAX*NMAX, 4*NMAX ), $ LDW5 = NSMPMX*LMAX + ( NMAX + LMAX )* $ ( NMAX + MMAX ) + 3*NMAX + MMAX + LMAX, $ LDW6 = NSMPMX*LMAX + ( NMAX + LMAX )* $ ( NMAX + MMAX ) + NMAX + $ MAX( 1, NMAX*NMAX*LMAX + NMAX*LMAX + $ NMAX, NMAX*NMAX + $ MAX( NMAX*NMAX + $ NMAX*MAX( NMAX, LMAX ) + $ 6*NMAX + MIN( NMAX, LMAX ), $ NMAX*MMAX ) ), $ LDW7 = NSMPMX*LMAX + ( NMAX + LMAX )* $ ( NMAX + MMAX ) + 3*NMAX + MMAX + LMAX, $ LW1 = MAX( 2*( MMAX + LMAX )*NOBRMX* $ ( 2*( MMAX + LMAX )*( NOBRMX + 1 ) $ + 3 ) + LMAX*NOBRMX, $ 4*( MMAX + LMAX )*NOBRMX* $ ( MMAX + LMAX )*NOBRMX + $ ( NMAX + LMAX )*( NMAX + MMAX ) + $ MAX( LDW1, LDW2 ), $ ( NMAX + LMAX )*( NMAX + MMAX ) + $ NMAX + NMAX*NMAX + 2 + $ NMAX*( NMAX + MMAX + LMAX ) + $ MAX( 5*NMAX, 2, MIN( LDW3, LDW4 ), $ LDW5, LDW6 ) ), $ LW2 = NSMPMX*LMAX + BSNM + $ MAX( 4, NSMPMX + $ MAX( NSMPMX*BSNM + $ MAX( 2*NNMAX, 5*BSNM + 1 ), $ BSNM**2 + BSNM + $ MAX( NSMPMX + 2*NNMAX, $ 5*BSNM ) ) ), $ LW3 = MAX( LDW7, NSMPMX*LMAX + $ ( NMAX + LMAX )*( 2*NMAX + MMAX )+ $ 2*NMAX ), $ LW4 = NSMPMX*LMAX + LXM + $ MAX( 4, NSMPMX*LMAX + $ MAX( NSMPMX*LMAX*( BSNM + LTHS ) + $ MAX( NSMPMX*LMAX + L1M, $ L2M + LXM ), $ LXM*( BSNM + LTHS ) + $ LXM + $ MAX( NSMPMX*LMAX + L1M, $ LXM + L3M ) ) ), $ LDWORK = MAX( LW1, LW2, LW3, LW4 ) ) * .. Local Scalars .. LOGICAL INIT1, INITB, INITL, INITN, INITS CHARACTER*1 INIT INTEGER BSN, I, INFO, INI, ITER, ITMAX1, ITMAX2, IWARN, $ J, L, L1, L2, LPAR, LX, M, N, NN, NOBR, NPRINT, $ NS, NSMP DOUBLE PRECISION TOL1, TOL2 * .. Array Arguments .. INTEGER IWORK(LIWORK) DOUBLE PRECISION DWORK(LDWORK), U(LDU,MMAX), X(LXM), Y(LDY,LMAX) * .. External Functions .. LOGICAL LSAME EXTERNAL LSAME * .. External Subroutines .. EXTERNAL IB03BD * .. Intrinsic Functions .. INTRINSIC MAX, MIN * .. Executable Statements .. * WRITE ( NOUT, FMT = 99999 ) * Skip the heading in the data file and read the data. READ ( NIN, FMT = '()' ) READ ( NIN, FMT = * ) NOBR, M, L, NSMP, N, NN, ITMAX1, ITMAX2, $ NPRINT, TOL1, TOL2, INIT INITL = LSAME( INIT, 'L' ) INITS = LSAME( INIT, 'S' ) INITB = LSAME( INIT, 'B' ) INITN = LSAME( INIT, 'N' ) INIT1 = INITL .OR. INITB IF( M.LE.0 .OR. M.GT.MMAX ) THEN WRITE ( NOUT, FMT = 99993 ) M ELSE IF( L.LE.0 .OR. L.GT.LMAX ) THEN WRITE ( NOUT, FMT = 99992 ) L ELSE NS = N IF( INIT1 ) THEN IF( NOBR.LE.0 .OR. NOBR.GT.NOBRMX ) THEN WRITE ( NOUT, FMT = 99991 ) NOBR STOP ELSEIF( NSMP.LT.2*( M + L + 1 )*NOBR - 1 ) THEN WRITE ( NOUT, FMT = 99990 ) NSMP STOP ELSEIF( N.EQ.0 .OR. N.GE.NOBR ) THEN WRITE ( NOUT, FMT = 99989 ) N STOP END IF IF ( N.LT.0 ) $ N = NOBR - 1 ELSE IF( NSMP.LT.0 ) THEN WRITE ( NOUT, FMT = 99990 ) NSMP STOP ELSEIF( N.LT.0 .OR. N.GT.NMAX ) THEN WRITE ( NOUT, FMT = 99989 ) N STOP END IF END IF IF( NN.LT.0 .OR. NN.GT.NNMAX ) THEN WRITE ( NOUT, FMT = 99988 ) NN ELSE BSN = NN*( L + 2 ) + 1 L1 = BSN*L L2 = N*( L + M + 1 ) + L*M LX = L1 + L2 INI = 1 IF ( INITL ) THEN LPAR = L1 ELSEIF ( INITS ) THEN INI = L1 + 1 LPAR = L2 ELSEIF ( INITN ) THEN LPAR = LX END IF IF( INIT1 ) $ N = NS * Read the input-output data, initial parameters, and seed. READ ( NIN, FMT = * ) ( ( U(I,J), J = 1,M ), I = 1,NSMP ) READ ( NIN, FMT = * ) ( ( Y(I,J), J = 1,L ), I = 1,NSMP ) IF ( .NOT.INITB ) $ READ ( NIN, FMT = * ) ( X(I), I = INI,INI+LPAR-1 ) IF ( INITS .OR. INITB ) $ READ ( NIN, FMT = * ) ( DWORK(I), I = 1,4 ) * Solve a Wiener system identification problem. CALL IB03BD( INIT, NOBR, M, L, NSMP, N, NN, ITMAX1, $ ITMAX2, NPRINT, U, LDU, Y, LDY, X, LX, TOL1, $ TOL2, IWORK, DWORK, LDWORK, IWARN, INFO ) * IF ( INFO.NE.0 ) THEN WRITE ( NOUT, FMT = 99998 ) INFO ELSE IF( IWARN.NE.0 ) WRITE ( NOUT, FMT = 99987 ) IWARN ITER = DWORK(3) WRITE ( NOUT, FMT = 99997 ) DWORK(2) WRITE ( NOUT, FMT = 99996 ) ITER, IWORK(1), IWORK(2) * Recompute LX is necessary. IF ( INIT1 .AND. NS.LT.0 ) $ LX = L1 + N*( L + M + 1 ) + L*M WRITE ( NOUT, FMT = 99994 ) WRITE ( NOUT, FMT = 99995 ) ( X(I), I = 1, LX ) END IF END IF END IF END IF STOP * 99999 FORMAT (' IB03BD EXAMPLE PROGRAM RESULTS',/1X) 99998 FORMAT (' INFO on exit from IB03BD = ',I4) 99997 FORMAT (/' Final 2-norm of the residuals = ',D15.7) 99996 FORMAT (/' Number of iterations = ', I7, $ /' Number of function evaluations = ', I7, $ /' Number of Jacobian evaluations = ', I7) 99995 FORMAT (10(1X,F9.4)) 99994 FORMAT (/' Final approximate solution is ' ) 99993 FORMAT (/' M is out of range.',/' M = ',I5) 99992 FORMAT (/' L is out of range.',/' L = ',I5) 99991 FORMAT (/' NOBR is out of range.',/' NOBR = ',I5) 99990 FORMAT (/' NSMP is out of range.',/' NSMP = ',I5) 99989 FORMAT (/' N is out of range.',/' N = ',I5) 99988 FORMAT (/' NN is out of range.',/' NN = ',I5) 99987 FORMAT (' IWARN on exit from IB03BD = ',I4) ENDProgram Data
IB03BD EXAMPLE PROGRAM DATA 10 1 1 1024 4 12 500 1000 0 .00001 .00001 B 2.2183165e-01 3.9027807e-02 -5.0295887e-02 8.5386224e-03 7.2431159e-02 -1.7082198e-03 -1.7176287e-01 -2.6198104e-01 -1.7194108e-01 1.8566868e-02 1.5625362e-01 1.7463811e-01 1.1564450e-01 2.8779248e-02 -8.4265993e-02 -2.0978501e-01 -2.6591828e-01 -1.7268680e-01 2.1525013e-02 1.4363602e-01 7.3101431e-02 -1.0259212e-01 -1.6380473e-01 -1.0021167e-02 2.0263451e-01 2.1983417e-01 -2.1636523e-02 -3.0986057e-01 -3.8521982e-01 -2.1785179e-01 -1.4761096e-02 3.7005180e-02 -2.8119028e-02 -4.2167901e-02 5.2117694e-02 1.2023747e-01 1.8863385e-02 -1.9506434e-01 -3.0192175e-01 -1.7000747e-01 8.0740471e-02 2.0188076e-01 8.5108288e-02 -1.3270970e-01 -2.3646822e-01 -1.6505385e-01 -4.7448014e-02 -2.7886815e-02 -1.0152026e-01 -1.4155374e-01 -6.1650823e-02 8.3519614e-02 1.5926650e-01 8.6142760e-02 -9.4385381e-02 -2.6609066e-01 -3.2883874e-01 -2.5908050e-01 -1.1648940e-01 -3.0653766e-03 1.0326675e-02 -5.3445909e-02 -9.2412724e-02 -3.0279541e-02 8.4846832e-02 1.1133075e-01 -3.2135250e-02 -2.5308181e-01 -3.5670882e-01 -2.4458860e-01 -2.5254261e-02 9.3714332e-02 1.8643667e-02 -1.4592119e-01 -2.2730880e-01 -1.7140060e-01 -7.4131665e-02 -3.9669515e-02 -5.1266129e-02 -1.1752833e-02 1.0785565e-01 2.0665525e-01 1.6117322e-01 -2.6938653e-02 -2.1941152e-01 -2.7753567e-01 -1.8805912e-01 -4.6845025e-02 5.8585698e-02 1.2218407e-01 1.7838638e-01 2.2169815e-01 1.9825589e-01 8.0215288e-02 -7.2135308e-02 -1.4381520e-01 -6.8724371e-02 1.0191205e-01 2.3766633e-01 2.3876101e-01 1.1678077e-01 -2.0428168e-02 -5.8973233e-02 3.1326900e-02 1.7391495e-01 2.4558570e-01 1.7650262e-01 1.2444292e-02 -1.1538234e-01 -9.5917970e-02 6.4762165e-02 2.4258524e-01 3.0102251e-01 2.1222960e-01 7.8706189e-02 3.1500466e-02 1.0297577e-01 1.9875173e-01 1.9434906e-01 5.8146667e-02 -1.1941921e-01 -2.1038478e-01 -1.5594967e-01 1.8552198e-03 1.6878529e-01 2.5937416e-01 2.2516346e-01 6.6144472e-02 -1.5623019e-01 -3.3161105e-01 -3.6695732e-01 -2.6565333e-01 -1.3254832e-01 -8.0101064e-02 -1.2531889e-01 -1.8843171e-01 -1.9038956e-01 -1.3230055e-01 -7.0889306e-02 -3.9679280e-02 -2.6286077e-02 -2.3630770e-02 -6.0652834e-02 -1.4929250e-01 -2.2155095e-01 -1.7331044e-01 5.2693564e-03 1.7683919e-01 1.8244690e-01 2.5118458e-02 -1.1051051e-01 -5.1764984e-02 1.6342054e-01 3.1563281e-01 2.3808751e-01 -4.4871135e-03 -1.8778679e-01 -1.6017584e-01 2.3481991e-02 1.9209185e-01 2.4281065e-01 2.1224192e-01 1.8825017e-01 1.9811718e-01 2.0202486e-01 1.6812825e-01 1.1444796e-01 7.2452475e-02 4.0090973e-02 -6.7139529e-03 -6.8721730e-02 -1.1460099e-01 -1.1914168e-01 -8.9852521e-02 -4.5942222e-02 1.0932686e-02 8.1900393e-02 1.3092374e-01 9.0790221e-02 -6.3538148e-02 -2.5119963e-01 -3.2585173e-01 -2.0850925e-01 1.7922009e-02 1.6783753e-01 1.2518317e-01 -4.3517162e-02 -1.5783138e-01 -1.0686847e-01 4.4782565e-02 1.3893172e-01 9.8691579e-02 2.6311282e-03 -1.6073049e-02 7.8512306e-02 1.9453537e-01 2.2504627e-01 1.6121235e-01 7.8124056e-02 2.9774586e-02 -5.3899280e-03 -6.5745322e-02 -1.2329059e-01 -9.5096521e-02 5.5471394e-02 2.5017082e-01 3.4773286e-01 2.6656242e-01 5.3705965e-02 -1.6135006e-01 -2.7310977e-01 -2.6814818e-01 -2.1074926e-01 -1.7743213e-01 -1.9796482e-01 -2.4059041e-01 -2.4663820e-01 -1.8780129e-01 -9.8317382e-02 -4.7848155e-02 -7.3425069e-02 -1.3529842e-01 -1.4739094e-01 -6.2482366e-02 6.8729554e-02 1.3251322e-01 6.1482940e-02 -8.5065014e-02 -1.6074078e-01 -6.7974104e-02 1.3976672e-01 2.9838081e-01 2.8233998e-01 1.1391411e-01 -7.1966946e-02 -1.5876983e-01 -1.3805556e-01 -8.2998592e-02 -5.7864811e-02 -6.5300733e-02 -7.0590592e-02 -5.5847027e-02 -4.1219301e-02 -6.1578267e-02 -1.3176243e-01 -2.2968907e-01 -3.0193311e-01 -2.8770451e-01 -1.5729276e-01 5.4414593e-02 2.5362617e-01 3.4482230e-01 3.0119122e-01 1.8534835e-01 9.6712488e-02 9.3385279e-02 1.6057572e-01 2.4424680e-01 3.0164891e-01 3.1693510e-01 2.8441517e-01 1.9948758e-01 7.3600888e-02 -5.4291337e-02 -1.3721320e-01 -1.5626045e-01 -1.3464149e-01 -1.1510541e-01 -1.2587072e-01 -1.6605420e-01 -2.1242088e-01 -2.3059410e-01 -1.8785957e-01 -7.8188380e-02 5.0484398e-02 1.0697957e-01 2.7421051e-02 -1.4419852e-01 -2.5888039e-01 -1.8018121e-01 7.8519535e-02 3.4009981e-01 4.0793257e-01 2.3842529e-01 -2.7029751e-02 -1.9919385e-01 -2.0420528e-01 -1.1389043e-01 -3.5602606e-02 5.7385906e-04 3.8759790e-02 1.0691941e-01 1.6303496e-01 1.4314046e-01 4.7786789e-02 -4.1030659e-02 -3.5960232e-02 7.0498851e-02 2.0120383e-01 2.6638170e-01 2.3249669e-01 1.2937468e-01 1.3309043e-02 -6.2770099e-02 -5.8936178e-02 3.4143049e-02 1.6425689e-01 2.2228910e-01 1.2062705e-01 -1.0832755e-01 -3.0711352e-01 -3.2002334e-01 -1.4072879e-01 7.6263091e-02 1.6385270e-01 1.0093887e-01 1.7269577e-02 4.3458474e-02 1.6769625e-01 2.4967945e-01 1.7314220e-01 -2.7519776e-02 -1.9806822e-01 -2.1140982e-01 -7.2758850e-02 1.1057470e-01 2.3440218e-01 2.5956640e-01 1.9629970e-01 7.2200120e-02 -6.6390448e-02 -1.4805958e-01 -1.1487691e-01 1.3561014e-02 1.3146288e-01 1.3205007e-01 1.5159726e-02 -9.9141126e-02 -7.9831031e-02 8.4487631e-02 2.6348526e-01 2.9617209e-01 1.3322758e-01 -1.1642178e-01 -2.7289866e-01 -2.2996687e-01 -3.5143323e-02 1.5983180e-01 2.3035457e-01 1.7179773e-01 7.3333592e-02 1.1653452e-02 -1.8499701e-02 -6.7962911e-02 -1.4361094e-01 -1.7665147e-01 -9.1259528e-02 9.8323111e-02 2.6912800e-01 2.8047779e-01 9.9377687e-02 -1.5436535e-01 -2.9569363e-01 -2.3017874e-01 -4.1007324e-02 8.2484352e-02 2.1760384e-02 -1.5212456e-01 -2.4257965e-01 -1.2641528e-01 1.0676585e-01 2.2865135e-01 1.0211687e-01 -1.6408728e-01 -3.0761461e-01 -1.7309336e-01 1.2302931e-01 3.0157576e-01 1.9992664e-01 -6.5766948e-02 -2.2490680e-01 -1.3209725e-01 9.1452627e-02 1.9707770e-01 7.0972862e-02 -1.6016460e-01 -2.7859962e-01 -2.0288880e-01 -4.9817844e-02 1.3587087e-02 -5.2447125e-02 -1.4164147e-01 -1.3776729e-01 -3.9470574e-02 5.4688171e-02 5.9780155e-02 -2.0666265e-02 -1.2306679e-01 -1.9150051e-01 -1.9953793e-01 -1.3072099e-01 1.7129752e-02 1.9139299e-01 2.8015628e-01 1.9737258e-01 -1.0273734e-02 -1.6921879e-01 -1.2914132e-01 8.3866166e-02 2.8290870e-01 3.0288568e-01 1.5939055e-01 1.4121758e-02 -8.0309556e-03 5.7046152e-02 7.8808779e-02 -4.0300321e-04 -9.3021531e-02 -6.6955916e-02 1.0073094e-01 2.8905786e-01 3.4946321e-01 2.4220689e-01 5.3331283e-02 -1.0609621e-01 -1.9358889e-01 -2.2728166e-01 -2.1680862e-01 -1.4144032e-01 -5.2173696e-03 1.1701944e-01 1.2668247e-01 4.8375112e-03 -1.4889224e-01 -1.9905951e-01 -9.9563224e-02 6.4580042e-02 1.5505008e-01 9.7617503e-02 -6.4905019e-02 -2.1769152e-01 -2.6787937e-01 -2.0919394e-01 -1.1033568e-01 -4.3266567e-02 -1.8066266e-02 1.3641281e-02 9.0806946e-02 1.8645977e-01 2.3150216e-01 1.9334856e-01 1.1238648e-01 4.9498545e-02 1.3155560e-02 -3.5876844e-02 -1.0537074e-01 -1.2612890e-01 -1.8934023e-02 1.8850628e-01 3.4290627e-01 3.0108912e-01 9.0554124e-02 -9.4812468e-02 -8.8842381e-02 6.3160674e-02 1.4646977e-01 1.7441277e-02 -2.2104173e-01 -3.1862778e-01 -1.5530235e-01 1.1291463e-01 2.1663682e-01 7.1521680e-02 -1.2722266e-01 -1.3147084e-01 6.8036453e-02 2.2914846e-01 1.4875917e-01 -8.5725554e-02 -1.9280127e-01 -3.7053987e-02 1.9484616e-01 2.0627194e-01 -5.0290692e-02 -2.9703694e-01 -2.4262627e-01 7.3980280e-02 3.1209111e-01 2.0500085e-01 -1.4678863e-01 -3.9620361e-01 -3.3299784e-01 -8.5315346e-02 7.0026906e-02 3.1783466e-02 -5.6224174e-02 -3.8238612e-02 4.1162402e-02 1.4020902e-02 -1.6267337e-01 -3.2229719e-01 -2.8405914e-01 -8.0208074e-02 7.7279407e-02 5.2461001e-02 -5.6931255e-02 -5.7081867e-02 8.4722273e-02 1.8989091e-01 9.1251490e-02 -1.4913841e-01 -3.0047660e-01 -2.2924644e-01 -4.5027749e-02 4.5847665e-02 -1.0582268e-02 -7.0165157e-02 8.8253349e-03 1.7968871e-01 2.6336655e-01 1.6274839e-01 -3.4038513e-02 -1.6866975e-01 -1.7822821e-01 -1.1212378e-01 -2.2511191e-02 9.2633595e-02 2.2273027e-01 2.8312792e-01 1.8855450e-01 -1.3339719e-02 -1.4451328e-01 -7.9411873e-02 9.5243626e-02 1.5825934e-01 8.6924573e-03 -1.9762612e-01 -2.0963986e-01 3.0881541e-02 3.1088543e-01 3.7605990e-01 2.0371110e-01 3.1659734e-03 -4.2255731e-02 2.7937777e-02 4.3768827e-02 -5.0975761e-02 -1.2013869e-01 -1.9514056e-02 1.9409077e-01 3.0061057e-01 1.6772761e-01 -8.4377993e-02 -2.0596833e-01 -8.8137439e-02 1.3053768e-01 2.3231724e-01 1.5592782e-01 3.3546556e-02 1.2609146e-02 8.8143918e-02 1.3076425e-01 5.2445727e-02 -9.1540218e-02 -1.6532665e-01 -8.9700956e-02 9.2256458e-02 2.6287064e-01 3.2206114e-01 2.4782579e-01 1.0180547e-01 -1.2653507e-02 -2.4053903e-02 4.5165362e-02 9.2697417e-02 3.9645255e-02 -7.0244568e-02 -9.7812594e-02 4.0489353e-02 2.5706426e-01 3.5970764e-01 2.4838839e-01 2.8758245e-02 -9.2051146e-02 -1.8531616e-02 1.4540527e-01 2.2483594e-01 1.6366159e-01 6.0613849e-02 2.6700790e-02 4.8805007e-02 2.4088984e-02 -8.7776563e-02 -1.9182802e-01 -1.5875230e-01 2.1332672e-02 2.1574747e-01 2.8121193e-01 1.9605244e-01 5.2140821e-02 -6.0594054e-02 -1.3111027e-01 -1.9003660e-01 -2.3031943e-01 -1.9896872e-01 -7.1576527e-02 8.7126470e-02 1.5966083e-01 8.0700885e-02 -9.6050487e-02 -2.3768453e-01 -2.4174619e-01 -1.1781079e-01 2.4058534e-02 6.3114157e-02 -3.4924911e-02 -1.8708629e-01 -2.5777811e-01 -1.7457598e-01 2.3256558e-03 1.2615984e-01 9.1298660e-02 -7.2869748e-02 -2.3064584e-01 -2.6487668e-01 -1.7896622e-01 -8.1019614e-02 -7.2160218e-02 -1.5109102e-01 -2.2270453e-01 -1.9311631e-01 -5.5949947e-02 1.0558527e-01 1.9015867e-01 1.5010510e-01 9.3491571e-03 -1.6206410e-01 -2.7872156e-01 -2.6789883e-01 -1.0908763e-01 1.3219241e-01 3.2581004e-01 3.6597785e-01 2.5860903e-01 1.1593033e-01 5.3232658e-02 8.9253999e-02 1.5038178e-01 1.6325136e-01 1.2516262e-01 8.1000365e-02 5.6249003e-02 4.1260796e-02 3.6021307e-02 7.0909773e-02 1.5431016e-01 2.1909293e-01 1.6946538e-01 1.3913978e-03 -1.5472276e-01 -1.5445369e-01 -6.5114694e-03 1.1511921e-01 5.3537688e-02 -1.4926948e-01 -2.8563000e-01 -2.0489020e-01 2.2256191e-02 1.8089745e-01 1.3686717e-01 -4.3194077e-02 -1.9185844e-01 -2.2260927e-01 -1.8688905e-01 -1.7299493e-01 -1.9552456e-01 -2.0311384e-01 -1.6521655e-01 -1.1035364e-01 -7.5596967e-02 -5.2167223e-02 -5.0648414e-03 6.7754101e-02 1.2412118e-01 1.2838133e-01 9.0308482e-02 4.0708671e-02 -1.2463102e-02 -7.6325303e-02 -1.2432208e-01 -9.0380523e-02 5.7426602e-02 2.4318485e-01 3.1839858e-01 2.0029814e-01 -2.6893656e-02 -1.7351791e-01 -1.2458940e-01 4.6580380e-02 1.5624992e-01 9.9382689e-02 -5.1882624e-02 -1.4100610e-01 -1.0040874e-01 -1.2845131e-02 -3.6737447e-03 -9.7637188e-02 -2.0172142e-01 -2.1938378e-01 -1.5223806e-01 -7.5818447e-02 -3.6932476e-02 -8.3361793e-03 4.9321106e-02 1.0828653e-01 8.6261922e-02 -5.6487106e-02 -2.4839500e-01 -3.5078033e-01 -2.7598256e-01 -6.2963150e-02 1.5901166e-01 2.7685307e-01 2.7164897e-01 2.1079033e-01 1.7714997e-01 2.0086813e-01 2.4438441e-01 2.4570310e-01 1.8078261e-01 9.0365447e-02 4.4844498e-02 7.6311118e-02 1.4103984e-01 1.5313326e-01 6.6678933e-02 -6.7720328e-02 -1.3565971e-01 -6.6316159e-02 8.3832277e-02 1.6588475e-01 7.6147385e-02 -1.3444251e-01 -2.9759248e-01 -2.8274479e-01 -1.1318459e-01 7.1421886e-02 1.5414324e-01 1.3182338e-01 8.0829372e-02 6.0814130e-02 6.6565578e-02 6.1490382e-02 3.4525574e-02 1.4709018e-02 3.9340413e-02 1.1733787e-01 2.1846966e-01 2.8684125e-01 2.6688313e-01 1.3632576e-01 -6.7370697e-02 -2.5502586e-01 -3.3949317e-01 -3.0013913e-01 -1.9871892e-01 -1.2610649e-01 -1.2941580e-01 -1.8923457e-01 -2.5813995e-01 -3.0533743e-01 -3.1970649e-01 -2.8788006e-01 -1.9500297e-01 -5.4155345e-02 8.1116905e-02 1.5269009e-01 1.4976106e-01 1.1681611e-01 1.0728712e-01 1.3670700e-01 1.8344060e-01 2.2041268e-01 2.2972773e-01 1.9334746e-01 9.8734288e-02 -2.6231283e-02 -9.9070456e-02 -4.1644202e-02 1.2360480e-01 2.5212308e-01 1.9060093e-01 -6.5066267e-02 -3.3581971e-01 -4.0871250e-01 -2.3222990e-01 4.0796545e-02 2.0553146e-01 1.9047036e-01 8.7982654e-02 2.1078714e-02 1.1947834e-02 -7.4158796e-03 -8.0649898e-02 -1.5932177e-01 -1.5963498e-01 -6.7654645e-02 3.3754864e-02 4.5488264e-02 -5.1656648e-02 -1.8439778e-01 -2.5821552e-01 -2.3168258e-01 -1.3075945e-01 -1.4319768e-02 6.0276859e-02 5.2808278e-02 -4.2009846e-02 -1.6857834e-01 -2.1862301e-01 -1.0815610e-01 1.2758494e-01 3.3007803e-01 3.4236071e-01 1.5606744e-01 -7.3906241e-02 -1.7487103e-01 -1.1779263e-01 -2.8797157e-02 -4.2649366e-02 -1.5603253e-01 -2.3465677e-01 -1.6213440e-01 3.1155521e-02 1.9455902e-01 2.0308035e-01 6.4105637e-02 -1.1373221e-01 -2.2912186e-01 -2.4930244e-01 -1.8794162e-01 -6.9023299e-02 6.6894859e-02 1.4860950e-01 1.1319286e-01 -2.1622177e-02 -1.4430675e-01 -1.4139382e-01 -1.4679189e-02 1.0606471e-01 8.3987908e-02 -8.6549724e-02 -2.6473902e-01 -2.8787546e-01 -1.1665499e-01 1.3032718e-01 2.7649250e-01 2.2886289e-01 4.1972959e-02 -1.4166947e-01 -2.1351821e-01 -1.7294568e-01 -9.5242426e-02 -3.9988034e-02 6.0215518e-04 6.4278100e-02 1.4411085e-01 1.7008073e-01 7.6346726e-02 -1.1397897e-01 -2.7942868e-01 -2.8837790e-01 -1.1356283e-01 1.2995490e-01 2.6791352e-01 2.1050936e-01 3.2758432e-02 -8.8492035e-02 -3.6187051e-02 1.3102808e-01 2.2789768e-01 1.2664599e-01 -9.9240525e-02 -2.3008477e-01 -1.1958430e-01 1.3943384e-01 2.8863442e-01 1.6130336e-01 -1.3747854e-01 -3.2522857e-01 -2.2524885e-01 5.3864511e-02 2.3305883e-01 1.5177574e-01 -7.4373920e-02 -1.8870441e-01 -6.7093573e-02 1.6495747e-01 2.8369836e-01 2.0511206e-01 5.1011236e-02 -6.5929875e-03 6.8964562e-02 1.6340844e-01 1.5740112e-01 5.4023734e-02 -4.3471011e-02 -5.1346211e-02 2.3145779e-02 1.1745308e-01 1.8212689e-01 1.9584070e-01 1.4022670e-01 5.9022790e-03 -1.6079919e-01 -2.4935419e-01 -1.7100378e-01 3.1256057e-02 1.8605482e-01 1.4297623e-01 -7.3243962e-02 -2.7593402e-01 -2.9797544e-01 -1.5307840e-01 -4.0914832e-03 2.1269662e-02 -4.1497170e-02 -5.9046655e-02 2.7976789e-02 1.2846949e-01 1.0303296e-01 -7.5938937e-02 -2.8392411e-01 -3.6123552e-01 -2.5664252e-01 -5.3262494e-02 1.2879625e-01 2.3255706e-01 2.6842403e-01 2.5122050e-01 1.7087253e-01 3.4014290e-02 -9.3227815e-02 -1.2001867e-01 -2.1139059e-02 1.2023890e-01 1.7758447e-01 9.6606085e-02 -5.2792108e-02 -1.3892628e-01 -8.4350032e-02 7.1620365e-02 2.1524576e-01 2.5910116e-01 2.0627091e-01 1.2532985e-01 7.1727643e-02 3.8319163e-02 -1.9240088e-02 -1.1662856e-01 -2.1107703e-01 -2.4258539e-01 -1.9809090e-01 -1.2271124e-01 -6.5266079e-02 -2.6001544e-02 2.6587042e-02 8.9979857e-02 1.0112134e-01 -1.6495775e-03 -1.8712095e-01 -3.2285436e-01 -2.8769737e-01 -1.0373843e-01 6.3283390e-02 6.4192144e-02 -6.9141383e-02 -1.4546154e-01 -2.2743165e-02 2.1671482e-01 3.3495240e-01 1.9730942e-01 -6.4245098e-02 -1.8430371e-01 -5.9313975e-02 1.3285821e-01 1.3988590e-01 -6.3313853e-02 -2.3781208e-01 -1.6565753e-01 7.8634007e-02 2.0643470e-01 6.3051903e-02 -1.7337120e-01 -1.9553447e-01 5.8877424e-02 3.1320739e-01 2.6455767e-01 -5.6738794e-02 -3.0614673e-01 -2.0738949e-01 1.4261991e-01 3.9321755e-01 3.3131011e-01 8.6485026e-02 -6.3943179e-02 -2.3354764e-02 5.9552949e-02 3.1845636e-02 -5.2189216e-02 -1.8514555e-02 1.7050716e-01 3.3649462e-01 2.9310084e-01 7.8582244e-02 -8.5200138e-02 -5.9242022e-02 5.3629257e-02 5.3919799e-02 -9.1290610e-02 -1.9983794e-01 -1.0236954e-01 1.3831631e-01 2.9035137e-01 -1.7703630e-01 -1.1470789e-01 -1.7257803e-02 7.3360924e-02 1.2806267e-01 1.3650217e-01 1.0539571e-01 5.4901306e-02 1.0347593e-02 -1.4210364e-02 -2.9316079e-02 -5.9818410e-02 -1.1287079e-01 -1.5651256e-01 -1.3759239e-01 -3.1325918e-02 1.2118952e-01 2.2925439e-01 2.1688928e-01 8.3280850e-02 -9.0968958e-02 -1.9863421e-01 -1.7919413e-01 -5.4874063e-02 9.1323774e-02 1.7241745e-01 1.4973591e-01 5.1202694e-02 -5.0722214e-02 -8.6474562e-02 -3.6675604e-02 5.0794719e-02 9.2852996e-02 3.5475423e-02 -9.8019853e-02 -2.1560266e-01 -2.2054921e-01 -8.4207430e-02 1.2773783e-01 2.9411889e-01 3.1432928e-01 1.7183620e-01 -5.3673166e-02 -2.3087548e-01 -2.5206313e-01 -9.9556443e-02 1.3579254e-01 3.0302360e-01 2.8345210e-01 6.9698019e-02 -2.2311064e-01 -4.2606792e-01 -4.1979542e-01 -2.0235411e-01 1.1680679e-01 3.8269042e-01 4.7499251e-01 3.6130151e-01 1.0698485e-01 -1.5666457e-01 -2.9684785e-01 -2.5130444e-01 -6.7456399e-02 1.2329504e-01 1.8968350e-01 8.9456729e-02 -1.0185072e-01 -2.4339863e-01 -2.2562726e-01 -4.5215735e-02 1.9190737e-01 3.3930982e-01 3.0360010e-01 1.0486525e-01 -1.3364785e-01 -2.6276635e-01 -2.0355127e-01 -1.0514338e-03 2.0109829e-01 2.5410141e-01 1.0538640e-01 -1.6182684e-01 -3.7724711e-01 -3.8906986e-01 -1.6075631e-01 2.0065197e-01 5.0030087e-01 5.6260189e-01 3.3306758e-01 -8.1981699e-02 -4.6637054e-01 -6.1157444e-01 -4.3578631e-01 -3.4787751e-02 3.6943357e-01 5.5331393e-01 4.1651911e-01 3.8203811e-02 -3.6624642e-01 -5.6531588e-01 -4.4111547e-01 -5.7977077e-02 3.6800859e-01 5.8749279e-01 4.6334166e-01 5.9154789e-02 -3.8817476e-01 -6.0585734e-01 -4.5438072e-01 -2.1770889e-02 4.2269933e-01 5.9388393e-01 3.7277877e-01 -1.1367643e-01 -5.6785416e-01 -7.0538273e-01 -4.3261293e-01 9.5667577e-02 5.7311674e-01 7.2849359e-01 4.8697304e-01 9.0040534e-03 -4.1643634e-01 -5.5375692e-01 -3.6053568e-01 1.0675442e-03 2.8391467e-01 3.2050851e-01 1.2014875e-01 -1.5499683e-01 -3.0636590e-01 -2.2845450e-01 3.0168597e-02 3.0447079e-01 4.1814633e-01 2.9408146e-01 3.3795396e-03 -2.8043536e-01 -3.9163122e-01 -2.7524621e-01 -1.6330862e-02 2.2338646e-01 3.1163298e-01 2.1884631e-01 2.0034460e-02 -1.6244160e-01 -2.3122765e-01 -1.5928083e-01 4.5460308e-03 1.6378113e-01 2.2566835e-01 1.5187573e-01 -1.8633628e-02 -1.8835877e-01 -2.5597784e-01 -1.7568160e-01 1.6144538e-02 2.1796548e-01 3.1334397e-01 2.3350541e-01 9.9054075e-04 -2.7139443e-01 -4.3349329e-01 -3.8409180e-01 -1.3941008e-01 1.6850242e-01 3.6865127e-01 3.5669633e-01 1.5962938e-01 -8.6421861e-02 -2.2603591e-01 -1.7879992e-01 1.5608870e-02 2.2316774e-01 2.9540664e-01 1.5777130e-01 -1.3932674e-01 -4.3707134e-01 -5.5308393e-01 -3.9056636e-01 -6.9866596e-03 4.0342788e-01 6.1470960e-01 5.0478901e-01 1.3556472e-01 -2.7661265e-01 -4.8754120e-01 -3.7410263e-01 -1.0933935e-02 3.7332700e-01 5.3265415e-01 3.5296792e-01 -7.5112937e-02 -5.0630963e-01 -6.8543131e-01 -5.0254861e-01 -6.3204556e-02 3.7616490e-01 5.6861420e-01 4.2839911e-01 7.7256895e-02 -2.4286013e-01 -3.2974149e-01 -1.4621212e-01 1.6396591e-01 3.7227253e-01 3.1398669e-01 -1.5203951e-03 -3.8826155e-01 -5.9422715e-01 -4.6290884e-01 -4.4082503e-02 4.2614489e-01 6.6944646e-01 5.4057059e-01 1.1914310e-01 -3.4186097e-01 -5.7361170e-01 -4.5144665e-01 -6.3037624e-02 3.5015696e-01 5.3940241e-01 3.9354970e-01 6.6063109e-05 -4.0735798e-01 -5.8396114e-01 -4.1610263e-01 1.0313382e-02 4.5449701e-01 6.5638620e-01 4.8903578e-01 3.8482894e-02 -4.3952337e-01 -6.6436421e-01 -4.9492372e-01 -1.7915270e-02 4.9445240e-01 7.3828446e-01 5.5772875e-01 4.3827397e-02 -5.1216643e-01 -7.8827423e-01 -6.2373284e-01 -1.1577453e-01 4.4053448e-01 7.3121649e-01 6.0691719e-01 1.6037942e-01 -3.4101558e-01 -6.1837622e-01 -5.3898039e-01 -1.7955555e-01 2.3296574e-01 4.6098842e-01 3.9204767e-01 9.4586522e-02 -2.3425494e-01 -3.9383077e-01 -2.9901136e-01 -2.1727093e-02 2.6290754e-01 3.8667642e-01 2.8641038e-01 3.4299620e-02 -2.1199530e-01 -3.0703990e-01 -2.0539827e-01 1.3733625e-02 1.9989717e-01 2.2856610e-01 8.0442398e-02 -1.4924794e-01 -3.1635143e-01 -3.2043874e-01 -1.6226330e-01 6.7449386e-02 2.5253008e-01 3.1855044e-01 2.6051993e-01 1.2699840e-01 -1.6342455e-02 -1.1750854e-01 -1.5094063e-01 -1.1699324e-01 -3.6407066e-02 5.7070826e-02 1.2470744e-01 1.3295525e-01 6.7237676e-02 -5.6199791e-02 -1.8928499e-01 -2.6860491e-01 -2.4751370e-01 -1.2546869e-01 4.7269068e-02 1.9379936e-01 2.5012057e-01 1.9757699e-01 6.9603172e-02 -6.6884197e-02 -1.4260360e-01 -1.1800895e-01 -4.5690911e-03 1.3505757e-01 2.1176910e-01 1.5667518e-01 -2.9715225e-02 -2.6058872e-01 -4.0072162e-01 -3.4636170e-01 -1.0002597e-01 2.1522385e-01 4.2116592e-01 3.9178740e-01 1.3552073e-01 -2.0194672e-01 -4.2193015e-01 -3.9351670e-01 -1.3365470e-01 2.0423921e-01 4.2544835e-01 4.1162219e-01 1.8730580e-01 -1.0283670e-01 -2.8986993e-01 -2.8756628e-01 -1.3866788e-01 2.8290398e-02 9.5513335e-02 3.5118646e-02 -8.2724881e-02 -1.5147446e-01 -1.0799938e-01 2.6949604e-02 1.6959254e-01 2.3358015e-01 1.8482066e-01 5.6424609e-02 -7.8806247e-02 -1.5583364e-01 -1.5299245e-01 -9.3729273e-02 -1.9708548e-02 3.8600307e-02 7.1469845e-02 7.8472613e-02 5.5625386e-02 -1.0621857e-03 -8.0782039e-02 -1.5057837e-01 -1.6705428e-01 -1.0304932e-01 2.9389143e-02 1.7801990e-01 2.7318425e-01 2.6234323e-01 1.3834554e-01 -5.4215912e-02 -2.3593270e-01 -3.2392000e-01 -2.6898405e-01 -8.5844039e-02 1.4215609e-01 2.9652172e-01 2.8801270e-01 1.1683545e-01 -1.1688760e-01 -2.6947626e-01 -2.4573958e-01 -6.4329645e-02 1.5353975e-01 2.6653313e-01 2.0755588e-01 2.4602079e-02 -1.5772495e-01 -2.2567844e-01 -1.4875573e-01 9.9414396e-03 1.4397851e-01 1.7486115e-01 9.6314112e-02 -3.2169687e-02 -1.2887854e-01 -1.3861783e-01 -5.9693947e-02 6.1826068e-02 1.6117670e-01 1.8758542e-01 1.2643056e-01 4.7038639e-03 -1.2089033e-01 -1.8936563e-01 -1.6676448e-01 -6.8240952e-02 4.6702545e-02 1.0911959e-01 8.7135042e-02 1.1538006e-02 -4.4789930e-02 -2.4262269e-02 6.5437901e-02 1.5116338e-01 1.4886934e-01 3.3820535e-02 -1.3097789e-01 -2.3522600e-01 -2.0099760e-01 -4.2018915e-02 1.4060900e-01 2.2430878e-01 1.4698003e-01 -4.9334401e-02 -2.4015379e-01 -2.9449301e-01 -1.5978257e-01 9.9469238e-02 3.3553927e-01 4.0432846e-01 2.5275189e-01 -4.8157255e-02 -3.4363559e-01 -4.8101858e-01 -3.9093124e-01 -1.2065446e-01 1.9561509e-01 4.0816957e-01 4.2449571e-01 2.4947873e-01 -2.2290220e-02 -2.5535821e-01 -3.3965313e-01 -2.4442241e-01 -3.2717407e-02 1.7386538e-01 2.6131002e-01 1.8344736e-01 -1.4617105e-02 -2.2004617e-01 -3.0989410e-01 -2.1648361e-01 2.9614296e-02 3.0600899e-01 4.6010027e-01 3.9585763e-01 1.3407054e-01 -1.9445050e-01 -4.2254041e-01 -4.4190341e-01 -2.6148822e-01 2.4561144e-03 1.9639531e-01 2.2058130e-01 8.8618067e-02 -8.2771773e-02 -1.5145974e-01 -4.8116921e-02 1.7081593e-01 3.5448643e-01 3.5655964e-01 1.3834184e-01 -1.9528570e-01 -4.5613811e-01 -4.9089820e-01 -2.7873232e-01 5.5837539e-02 3.2156811e-01 3.7683870e-01 2.1007687e-01 -6.1195486e-02 -2.6670692e-01 -2.8529736e-01 -1.1252984e-01 1.4069959e-01 3.1548805e-01 3.0070613e-01 1.0177110e-01 -1.6096596e-01 -3.2711612e-01 -2.9842835e-01 -9.9492033e-02 1.4305421e-01 2.8418081e-01 2.4879424e-01 7.0440776e-02 -1.3708347e-01 -2.5105923e-01 -2.1001593e-01 -4.5285982e-02 1.4155737e-01 2.4209754e-01 2.0725941e-01 7.3959838e-02 -6.6466455e-02 -1.3533231e-01 -1.1722667e-01 -5.6247689e-02 -8.2151160e-03 4.6646596e-03 -5.3013327e-05 6.4836935e-03 3.4885521e-02 7.2093769e-02 9.6085499e-02 9.0621414e-02 5.0063443e-02 -1.9216694e-02 -9.5194586e-02 -1.4177512e-01 -1.2554939e-01 -4.1561203e-02 7.4612994e-02 1.6458119e-01 1.8370169e-01 1.2694288e-01 2.5574339e-02 -7.6209464e-02 -1.4292208e-01 -1.5717793e-01 -1.2150507e-01 -5.7465582e-02 3.0433319e-03 3.8135050e-02 5.3444515e-02 7.4126764e-02 1.1232692e-01 1.4266966e-01 1.1713381e-01 1.2919877e-02 -1.3094351e-01 -2.2903887e-01 -2.1083457e-01 -7.7741149e-02 9.2251468e-02 1.9732652e-01 1.8027267e-01 6.1530912e-02 -8.1015797e-02 -1.6435623e-01 -1.4922825e-01 -5.8874212e-02 3.9408110e-02 7.8379546e-02 3.6886774e-02 -4.2241134e-02 -8.1505612e-02 -2.9557008e-02 9.2798034e-02 2.0055247e-01 2.0414883e-01 7.6944227e-02 -1.2029199e-01 -2.7519345e-01 -2.9408814e-01 -1.6081545e-01 5.1070794e-02 2.1840144e-01 2.3874816e-01 9.4335060e-02 -1.2904879e-01 -2.8774773e-01 -2.6899028e-01 -6.6408095e-02 2.1071698e-01 4.0356249e-01 3.9994180e-01 1.9633323e-01 -1.0730235e-01 -3.6601054e-01 -4.6248715e-01 -3.5922221e-01 -1.1354600e-01 1.4870456e-01 2.9521055e-01 2.5966678e-01 8.3040302e-02 -1.0914113e-01 -1.8742442e-01 -1.0478464e-01 7.3317409e-02 2.1546569e-01 2.1382067e-01 5.6531581e-02 -1.6427012e-01 -3.1183656e-01 -2.9186150e-01 -1.1383004e-01 1.1231696e-01 2.4506533e-01 2.0292544e-01 1.9811075e-02 -1.7391062e-01 -2.3677906e-01 -1.1242105e-01 1.2953875e-01 3.3467916e-01 3.5946938e-01 1.6169418e-01 -1.6880410e-01 -4.5538345e-01 -5.3000472e-01 -3.2991559e-01 5.7588162e-02 4.3386984e-01 5.9508457e-01 4.4813661e-01 6.8860243e-02 -3.3635714e-01 -5.4527976e-01 -4.4370745e-01 -8.9647493e-02 3.1753702e-01 5.4673805e-01 4.6318145e-01 1.0733728e-01 -3.1949400e-01 -5.6446899e-01 -4.7269412e-01 -8.8269356e-02 3.6150197e-01 5.9965309e-01 4.7275161e-01 5.2712510e-02 -4.0097128e-01 -6.0010920e-01 -4.1032807e-01 6.1089052e-02 5.2877389e-01 7.0388838e-01 4.7272792e-01 -3.2841140e-02 -5.1806125e-01 -7.0615746e-01 -5.0443062e-01 -5.3964611e-02 3.6781621e-01 5.2531916e-01 3.6514315e-01 3.1895267e-02 -2.4276338e-01 -2.9561167e-01 -1.2568333e-01 1.2380832e-01 2.6979551e-01 2.0920891e-01 -2.0179145e-02 -2.6980104e-01 -3.7620139e-01 -2.6519009e-01 -1.4966321e-04 2.5905182e-01 3.5875119e-01 2.4783584e-01 5.4317821e-03 -2.1770753e-01 -2.9814845e-01 -2.0810260e-01 -1.7395596e-02 1.5890290e-01 2.2758901e-01 1.6085463e-01 3.3576307e-03 -1.5297196e-01 -2.1737064e-01 -1.5023570e-01 1.2479222e-02 1.7606639e-01 2.4089523e-01 1.6216345e-01 -2.3230254e-02 -2.1504218e-01 -3.0098784e-01 -2.1779026e-01 8.8067567e-03 2.6812984e-01 4.1695437e-01 3.6159556e-01 1.2203070e-01 -1.7147580e-01 -3.5437470e-01 -3.3058973e-01 -1.3341351e-01 9.9954914e-02 2.1969740e-01 1.5589313e-01 -4.1996520e-02 -2.3771826e-01 -2.9083527e-01 -1.4002506e-01 1.5548285e-01 4.3862419e-01 5.3769302e-01 3.6811228e-01 -6.9569482e-03 -3.9769165e-01 -5.8956799e-01 -4.7193386e-01 -1.1138894e-01 2.8025332e-01 4.6943948e-01 3.4372376e-01 -1.6555081e-02 -3.8429530e-01 -5.2185674e-01 -3.2705351e-01 1.0055685e-01 5.1629500e-01 6.7570174e-01 4.8204840e-01 4.6679399e-02 -3.7892485e-01 -5.5799051e-01 -4.1189337e-01 -6.3130989e-02 2.4927425e-01 3.2624429e-01 1.3391859e-01 -1.7899014e-01 -3.7999275e-01 -3.0718591e-01 1.9919795e-02 4.0587411e-01 5.9872071e-01 4.5200311e-01 2.6827172e-02 -4.3774484e-01 -6.7014857e-01 -5.3423365e-01 -1.1312830e-01 3.4367827e-01 5.7281717e-01 4.5156693e-01 6.5481027e-02 -3.4683106e-01 -5.3783781e-01 -3.9562633e-01 -5.2304328e-03 4.0256826e-01 5.8408144e-01 4.2300297e-01 -1.8218267e-04 -4.4833216e-01 -6.5943295e-01 -5.0033881e-01 -5.1578103e-02 4.3192551e-01 6.6545648e-01 5.0237264e-01 2.6477477e-02 -4.8897549e-01 -7.3697545e-01 -5.5960739e-01 -4.7597748e-02 5.0867228e-01 7.8911527e-01 6.3269313e-01 1.3197226e-01 -4.2464681e-01 -7.2603682e-01 -6.1784801e-01 -1.8264666e-01 3.2014735e-01 6.1135123e-01 5.4895999e-01 1.9768580e-01 -2.2062099e-01 -4.6220719e-01 -4.0211731e-01 -9.9950534e-02 2.4465654e-01 4.1872319e-01 3.2500596e-01 3.2810917e-02 -2.7440750e-01 -4.1536442e-01 -3.1832701e-01 -5.5989066e-02 2.0726049e-01 3.1798239e-01 2.2484797e-01 5.1703651e-03 -1.8889751e-01 -2.2927380e-01 -9.1914974e-02 1.3314428e-01 3.0513495e-01 3.2224987e-01 1.7778028e-01 -4.7100451e-02 -2.4007922e-01 -3.2145867e-01 -2.7615883e-01 -1.4545755e-01 4.2822900e-03 1.1399372e-01 1.5138712e-01 1.1530153e-01 3.0234280e-02 -6.4234624e-02 -1.2615802e-01 -1.2407054e-01 -4.9317670e-02 7.5619816e-02 2.0015044e-01 2.6472178e-01 2.3118708e-01 1.0699863e-01 -5.5412012e-02 -1.8550876e-01 -2.3096135e-01 -1.8218227e-01 -7.2615500e-02 4.0881922e-02 1.0372451e-01 8.6362391e-02 -1.1351454e-03 -1.0889033e-01 -1.6548976e-01 -1.1405709e-01 4.6560657e-02 2.4386985e-01 3.6111476e-01 3.0662373e-01 8.1468123e-02 -2.0497551e-01 -3.9165036e-01 -3.6309524e-01 -1.2535574e-01 1.8954273e-01 3.9793935e-01 3.7486538e-01 1.3124068e-01 -1.9174474e-01 -4.0848802e-01 -4.0149539e-01 -1.8960477e-01 9.0301438e-02 2.7507284e-01 2.7972729e-01 1.4341274e-01 -1.2566755e-02 -7.8032703e-02 -2.7425697e-02 7.5351759e-02 1.3487633e-01 9.5488652e-02 -2.4590018e-02 -1.5233210e-01 -2.1189289e-01 -1.7248897e-01 -6.2455423e-02 5.4933614e-02 1.2398028e-01 1.2778044e-01 8.7386392e-02 3.4966577e-02 -1.0850501e-02 -4.6716543e-02 -6.9020828e-02 -6.3681635e-02 -1.6203206e-02 6.7394491e-02 1.5127737e-01 1.8399090e-01 1.2920707e-01 -7.0434827e-03 -1.7216342e-01 -2.8937677e-01 -2.9509198e-01 -1.7314710e-01 3.2745183e-02 2.3542177e-01 3.4097958e-01 2.9247721e-01 1.0411948e-01 -1.3495077e-01 -2.9868629e-01 -2.9240849e-01 -1.1517683e-01 1.2871323e-01 2.8803761e-01 2.6146766e-01 6.7234759e-02 -1.6729947e-01 -2.9180077e-01 -2.3297675e-01 -3.8493954e-02 1.6188055e-01 2.4607750e-01 1.7580193e-01 1.0770499e-02 -1.3917580e-01 -1.8630712e-01 -1.1496682e-01 1.8120146e-02 1.2605380e-01 1.4532251e-01 6.9056099e-02 -5.5814690e-02 -1.6001831e-01 -1.8912751e-01 -1.2778372e-01 -4.4698128e-03 1.2208903e-01 1.8963074e-01 1.6384408e-01 6.0799128e-02 -5.7339158e-02 -1.1860919e-01 -9.0086196e-02 -4.5798607e-03 6.0280807e-02 4.1676388e-02 -5.5180320e-02 -1.5518201e-01 -1.6828578e-01 -6.2049884e-02 1.0561621e-01 2.2337555e-01 2.0643187e-01 5.9839911e-02 -1.2043322e-01 -2.1083864e-01 -1.4415945e-01 4.3538937e-02 2.3203364e-01 2.9044234e-01 1.6171416e-01 -9.5674666e-02 -3.3749265e-01 -4.1795872e-01 -2.7746809e-01 2.0648626e-02 3.2603206e-01 4.8410918e-01 4.1672303e-01 1.5905611e-01 -1.6318595e-01 -3.9931562e-01 -4.4568803e-01 -2.9169291e-01 -2.0960934e-02 2.3175866e-01 3.4693819e-01 2.7877641e-01 7.7125945e-02 -1.4069530e-01 -2.5367798e-01 -2.0150506e-01 -1.6778161e-02 1.9116819e-01 2.9409556e-01 2.1593628e-01 -1.9610708e-02 -2.9401135e-01 -4.5512990e-01 -4.0311941e-01 -1.5075705e-01 1.7921653e-01 4.2153577e-01 4.6143206e-01 2.9688389e-01 3.5275834e-02 -1.7206796e-01 -2.2040717e-01 -1.1280250e-01 4.6014479e-02 1.2005000e-01 3.5297082e-02 -1.6459920e-01 -3.4121448e-01 -3.5130088e-01 -1.4787707e-01 1.7615712e-01 4.3972643e-01 4.8949447e-01 2.9899548e-01 -1.6059656e-02 -2.7414987e-01 -3.4124596e-01 -2.0476598e-01 3.1287353e-02 2.1535118e-01 2.3693813e-01 8.7039128e-02 -1.3914592e-01 -2.9731202e-01 -2.8057123e-01 -8.9244625e-02 1.6445576e-01 3.2621002e-01 2.9949560e-01 1.0678193e-01 -1.3016725e-01 -2.7225661e-01 -2.4687907e-01 -8.3173776e-02 1.1381888e-01 2.2819642e-01 1.9830143e-01 4.8505476e-02 -1.2763594e-01 -2.2560309e-01 -1.9560311e-01 -7.1212054e-02 6.0380807e-02 1.2445307e-01 1.0835168e-01 5.5609724e-02 1.7269294e-02 9.3997346e-03 1.1223045e-02 -4.3543819e-03 -4.2668837e-02 -8.5657964e-02 -1.0909342e-01 -9.7154374e-02 -4.6781850e-02 3.1101930e-02 1.0973840e-01 1.5122945e-01 1.2531404e-01 3.3620966e-02 -8.3194568e-02 -1.6716420e-01 1998. 1999. 2000. 2001.Program Results
IB03BD EXAMPLE PROGRAM RESULTS IWARN on exit from IB03BD = 12 Final 2-norm of the residuals = 0.2995840D+00 Number of iterations = 42 Number of function evaluations = 898 Number of Jacobian evaluations = 295 Final approximate solution is 14.1294 1.1232 6.4322 -11.2418 7.6380 -33.4730 -64.7203 747.1515 -0.4623 -92.6092 6.1682 -0.7672 0.1194 0.3558 0.9091 0.2948 1.3465 0.0093 0.0560 -0.0035 -0.4179 -0.0455 -2.0871 -0.9196 1.0777 0.9213 0.5373 1.0412 -0.3978 7.6832 -6.8614 -31.6119 -0.1092 -9.8984 0.1257 0.4056 0.0472 7.5819 -13.3969 2.4869 -66.0727 -0.8411 -0.7040 1.9641 1.3059 -0.2046 -0.9326 0.0040 0.4032 0.1479