MLSE

   more...
Computes least‑squares estimates of a linear regression model for a multichannel time series with a specified base channel.
Required Arguments
XNOBSX by NCHANX matrix containing the time series. (Input)
Each row of X corresponds to an observation of a multivariate time series, and each column of X corresponds to a univariate time series. The base time series or output channel is contained in the first column.
NDIFF — Vector of length NCHANX containing the order of differencing for each channel of X. (Input)
The elements of NDIFF must be greater than or equal to zero.
NDPREG — Vector of length NCHANX containing the number of regression parameters in the differenced form of the model for each channel of X. (Input) The elements of NDPREG must be greater than or equal to zero.
LAG — Vector of length NCHANX containing the amount of time that each channel of X is to lag the base series. (Input)
The elements of LAG must be greater than or equal to zero.
CNST — Estimate of the overall constant. (Output)
NPREG — Number of regression parameters in the undifferenced model. (Output)
NPREG = IADD + (NDPREG(1) + NDIFF(1)) +  + (NDPREG(NCHANX) +  DIFF(NCHANX)
where
IADD = NDIFF(1), if NDPREG(1) = 0
IADD = max (0, min(LAG(1)  1, NDIFF(1))), if NDPREG(1) > 0.
PREG — Vector of length NPREG containing the regression parameters in the undifferenced model. (Output)
The parameter estimates are concatenated as follows.
Channel 1:                  REG(i), i = 1, 2, …, IADD + NDPREG(1) + NDIFF(1)
Channel j:      PREG(i), i = I(j) + 1, I(j) + 2, …, I(j) + NDPREG(j) + NDIFF(j)
where
I(j) = IADD + NDPREG(1) + NDIFF(1) +  + NDPREG(j  1) + NDIFF(j  1)
for j = 2, 3, …, NCHANX.
Optional Arguments
NOBSX — Number of observations in each channel of the time series X. (Input)
NOBSX must be less than or equal LDX and greater than max(NDPREG(i) + LAG(i)) for
i = 1, 2, …, NCHANX.
Default: NOBSX = size (X,1).
NCHANX — Number of channels in the time series X. (Input)
NCHANX must be greater than or equal to one.
Default: NCHANX = size (X,2).
LDX — Leading dimension of X exactly as specified in the dimension statement of the calling program. (Input)
Default: LDX = size (X,1).
IMEAN — Option for computation of the means of the channels of X. (Input)
Default: IMEAN = 1.
IMEAN
Action
0
XMEAN is user specified.
1
XMEAN is set to the vector of arithmetic means of the channels of X.
XMEAN — Vector of length NCHANX containing the means of the channels of X. (Input, if IMEAN = 0; output, if IMEAN = 1)
FORTRAN 90 Interface
Generic: CALL MLSE (X, NDIFF, NDPREG, LAG, CNST, NPREG, PREG [])
Specific: The specific interface names are S_MLSE and D_MLSE.
FORTRAN 77 Interface
Single: CALL MLSE (NOBSX, NCHANX, X, LDX, IMEAN, XMEAN, NDIFF, NDPREG, LAG, CNST, NPREG, PREG)
Double: The double precision name is DMLSE.
Description
Routine MLSE performs least‑squares estimation of a linear regression model for a multichannel time series with a specified base channel.
Define the multichannel time series X by
X = (X1, X2, …, Xm)
where
Xj = (X1j, X2j, …, Xnj)T      j = 1, 2, …, m
with n = NOBSX and m = NCHANX. The columns of X correspond to individual channels of a multichannel time series and may be examined from a univariate perspective. The rows of X correspond to observations of an m‑variate time series and may be examined from a multivariate perspective. Note that an alternative characterization of the multivariate time series X considers the columns of X to be observations of an m‑variate time series with the rows of X containing univariate time series. For example, see Priestley (1981, page 692) and Fuller (1976, page 14).
The model is formed by regressing the base series X1 on its previous values and on the remaining channels X2Xm. The differenced form of the model is given by
where θ0 = CNST is the overall constant, dk = NDIFF(k) is the order of differencing Xk, lk = LAG(k) is the amount Xk lags X1,
and pk = NDPREG(k) for k = 1, 2, …, m.
The undifferenced form of the model is given by
where the undifferenced parameters are defined by
for k = 1, 2, …, m. Note that if l1 d1 0, the base series terms Xtj,1 at lags j = 1, …, (l1  1) are omitted from the right‑hand side of the above model when d1 1. In the actual computations, these terms are included.
Comments
1. Workspace may be explicitly provided, if desired, by use of M2SE/DM2SE. The reference is:
CALL M2SE (NOBSX, NCHANX, X, LDX, IMEAN, XMEAN, NDIFF, NDPREG, LAG, CNST, NPREG, PREG, XWK, IWK)
The additional arguments are as follows:
XWK — Work vector of length NOBSX * NCHANX+ 2 * NSUM2+ max(IADD, NCHANX + NSUM), where NSUM = NDPREG(1) +  + NDPREG(NCHANX).
IWK — Work vector of length NSUM.
2. Prior to parameter estimation, the channels of X are centered and/or differenced according to XMEAN and NDIFF, respectively.
3. The undifferenced predictive form of the model is
X(t, 1) = CNST + PREG(1) * X(t  1, 1) +  + PREG(IADD) * X(t  IADD, 1) +
PREG(IADD + 1) * X(t  LAG(1), 1) +  + PREG(IADD + NDPREG(1) + NDIFF(1)) *
X(t  LAG(1) + 1  NDPREG(1)  NDIFF(1), 1) +  + PREG(I(j) + 1) * X(t  LAG(j), j)
+  + PREG(I(j)+NDPREG(j)+NDIFF(j)) * X(t  LAG(j) + 1  NDPREG(j NDIFF(j),j)
 + 
where
I(j) = IADD + NDPREG(1) + NDIFF(1) +  + NDPREG(j  1)
+ NDIFF(j  1)
for j = 2, 3, …, NCHANX.
Examples
Example 1
Consider the Wölfer Sunspot Data (Box and Jenkins 1976, page 530) along with data on northern light activity and earthquake activity (Robinson 1967, page 204) to be a three‑channel time series. Routine MLSE is applied to these data to examine the regressive relationship between the channels.
 
USE GDATA_INT
USE MLSE_INT
USE UMACH_INT
 
IMPLICIT NONE
INTEGER LDX, NCHANX, NOBSX
PARAMETER (NCHANX=3, NOBSX=100, LDX=NOBSX)
!
INTEGER I, IMEAN, LAG(NCHANX), NCOL, NDIFF(NCHANX), &
NDPREG(NCHANX), NOUT, NPREG, NROW
REAL CNST, PREG(20), RDATA(100,4), X(LDX,NCHANX), &
XMEAN(NCHANX)
!
EQUIVALENCE (X(1,1), RDATA(1,2)), (X(1,2), RDATA(1,3)), &
(X(1,3), RDATA(1,4))
!
DATA NDIFF(1), NDIFF(2), NDIFF(3)/1, 1, 0/
DATA LAG(1), LAG(2), LAG(3)/1, 2, 1/
DATA NDPREG(1), NDPREG(2), NDPREG(3)/2, 1, 3/
!
CALL GDATA (8, RDATA, NROW, NCOL)
! USE Default Option to estimate channel means
! Compute regression parameters
CALL MLSE (X, NDIFF, NDPREG, LAG, CNST, NPREG, PREG, XMEAN=XMEAN)
!
! Print results
CALL UMACH (2, NOUT)
WRITE (NOUT,99993)
99993 FORMAT (//, 1X, ' Results of MLSE/M2SE')
WRITE (NOUT,99994)
99994 FORMAT (1X, ' I NDIFF(I) LAG(I) NDPREG(I) XMEAN(I)')
DO 10 I=1, NCHANX
WRITE (NOUT,99995) I, NDIFF(I), LAG(I), NDPREG(I), XMEAN(I)
99995 FORMAT (1X, 4(I3,6X), F12.4)
10 CONTINUE
WRITE (NOUT,99996) CNST
99996 FORMAT (1X, 'Overall constant, CNST = ', F12.4)
WRITE (NOUT,99997) NPREG
99997 FORMAT (//, 1X, 'Total number of parameters, NPREG = ', I2)
WRITE (NOUT,99998)
99998 FORMAT (//, 1X, ' I PREG(I)')
DO 20 I=1, NPREG
WRITE (NOUT,99999) I, PREG(I)
99999 FORMAT (1X, I2, 5X, F12.4)
20 CONTINUE
!
END
Output
 
Results of MLSE/M2SE
I NDIFF(I) LAG(I) NDPREG(I) XMEAN(I)
1 1 1 2 46.9400
2 1 2 1 63.4300
3 0 1 3 97.9700
Overall constant, CNST = -7.2698
 
Total number of parameters, NPREG = 8
 
I PREG(I)
1 -0.1481
2 -1.3444
3 0.4925
4 -0.0302
5 0.0302
6 -0.0210
7 0.0187
8 0.0765
 
Example 2
Consider the Gas Furnace Data (Box and Jenkins 1976, pages 532–533) where X1 is the percent CO2 in the outlet gas and X2 is the input gas rate in cubic feet/minute. Application of routine MLSE to these data produces the following results:
 
USE GDATA_INT
USE SCOPY_INT
USE MLSE_INT
USE UMACH_INT
 
IMPLICIT NONE
INTEGER LDX, NCHANX, NOBSX
PARAMETER (NCHANX=2, NOBSX=296, LDX=NOBSX)
!
INTEGER I, IMEAN, LAG(NCHANX), NCOL, NDIFF(NCHANX), &
NDPREG(NCHANX), NOUT, NPREG, NROW
REAL CNST, PREG(20), RDATA(296,2), X(LDX,NCHANX), &
XMEAN(NCHANX)
!
DATA NDIFF(1), NDIFF(2)/0, 0/
DATA LAG(1), LAG(2)/1, 3/
DATA NDPREG(1), NDPREG(2)/2, 3/
! Gas Furnace Data
CALL GDATA (7, RDATA, NROW, NCOL)
! Multichannel X consists of
! Column 1: Output percent CO2
! Column 2: Input gas rate
CALL SCOPY (NOBSX, RDATA(1:,2), 1, X(1:,1), 1)
CALL SCOPY (NOBSX, RDATA(1:,1), 1, X(1:,2), 1)
! Option to estimate channel means
IMEAN = 1
! Compute regression parameters
CALL MLSE (X, NDIFF, NDPREG, LAG, CNST, NPREG, PREG, XMEAN=XMEAN)
!
! Print results
CALL UMACH (2, NOUT)
WRITE (NOUT,99993)
99993 FORMAT (1X, 'Results of MLSE/M2SE on Gas Furnace Data')
WRITE (NOUT,99994)
99994 FORMAT (1X, ' I NDIFF(I) LAG(I) NDPREG(I) XMEAN(I)')
DO 10 I=1, NCHANX
WRITE (NOUT,99995) I, NDIFF(I), LAG(I), NDPREG(I), XMEAN(I)
99995 FORMAT (1X, 4(I3,6X), F12.4)
10 CONTINUE
WRITE (NOUT,99996) CNST
99996 FORMAT (1X, 'Overall constant, CNST = ', F12.4)
WRITE (NOUT,99997) NPREG
99997 FORMAT (1X, 'Total number of parameters, NPREG = ', I2)
WRITE (NOUT,99998)
99998 FORMAT (1X, ' I PREG(I)')
DO 20 I=1, NPREG
WRITE (NOUT,99999) I, PREG(I)
99999 FORMAT (1X, I2, 5X, F12.4)
20 CONTINUE
!
END
Output
 
Results of MLSE/M2SE on Gas Furnace Data
I NDIFF(I) LAG(I) NDPREG(I) XMEAN(I)
1 0 1 2 53.5091
2 0 3 3 -0.0568
Overall constant, CNST = 2.6562
Total number of parameters, NPREG = 5
I PREG(I)
1 1.6063
2 -0.6561
3 -0.4837
4 -0.1653
5 0.5052
 
Published date: 03/19/2020
Last modified date: 03/19/2020