Title: | Nonparametric Instrumental Variables Estimation and Inference |
---|---|
Description: | Implements methods introduced in Chen, Christensen, and Kankanala (2024) <doi:10.1093/restud/rdae025> for estimating and constructing uniform confidence bands for nonparametric structural functions using instrumental variables, including data-driven choice of tuning parameters. All methods in this package apply to nonparametric regression as a special case. |
Authors: | Jeffrey S. Racine [aut], Timothy Christensen [aut, cre], Patrick Alken [ctb], Rhys Ulerich [ctb], Simon N. Wood [ctb] |
Maintainer: | Timothy Christensen <[email protected]> |
License: | GPL (>= 3) |
Version: | 0.1.3 |
Built: | 2025-01-08 17:14:16 UTC |
Source: | CRAN |
This package implements the nonparametric instrumental variables estimation and inference methods described in Chen, Christensen, and Kankanala (2024) and Chen and Christensen (2018). The function npiv
estimates the nonparametric structural function h0
using B-splines and constructs uniform confidence bands for h0
. The function npiv_choose_J
performs data-driven choice of sieve dimension. All methods in this package apply to estimation and inference for nonparametric regression as a special case.
This package provides a function npiv(...)
with a simple interface for performing nonparametric instrumental variable estimation and inference.
Given a dependent variable vector Y
, matrix of endogenous regressors X
, and matrix of instruments W
, npiv
nonparametrically estimates the structural function h0
and its derivative using B-splines. npiv
can also be used for estimting the conditional mean h0
of Y
given X
, as as well as the derivative of the conditional mean function, by nonparametric regression.
The function npiv
also constructs uniform confidence bands for h0
and its derivative.
Sieve dimensions are determined in a data-dependent way if not provided by the user via the function npiv_choose_J
, which implements the methods described in Chen, Christensen, and Kankanala (2024). This data-driven choice of sieve dimension ensures estimators of h0
and its derivatives converge at the optimal sup-norm rate. The resulting uniform confidence bands for h0
and its derivative contract within a logarithmic factor of the optimal rate. In this way, npiv
facilitates fully data-driven estimation and uniform inference on h0
and its derivative.
If sieve dimensions are provided by the user, npiv
implements the bootstrap-based procedure of Chen and Christensen (2018) to construct uniform confidence bands for h0
and its derivative.
Jeffrey S. Racine <[email protected]>, Timothy Christensen <[email protected]>
Maintainer: Timothy Christensen <[email protected]>
Chen, X. and T. Christensen (2018). “Optimal Sup-norm Rates and Uniform Inference on Nonlinear Functionals of Nonparametric IV Regression.” Quantitative Economics, 9(1), 39-85. doi:10.3982/QE722
Chen, X., T. Christensen and S. Kankanala (2024). “Adaptive Estimation and Uniform Confidence Bands for Nonparametric Structural Functions and Elasticities.” Review of Economic Studies, forthcoming. doi:10.1093/restud/rdae025
This dataset is based on a sample taken from the British Family Expenditure Survey for 1995. It includes households consisting of married or cohabiting couples with an employed head of household, aged between 25 and 55 years, and with at most two children. There are 1655 household-level observations in total.
data("Engel95")
data("Engel95")
A data frame with 10 columns, and 1655 rows.
expenditure share on food, of type numeric
expenditure share on catering, of type numeric
expenditure share on alcohol, of type numeric
expenditure share on fuel, of type numeric
expenditure share on motor, of type numeric
expenditure share on fares, of type numeric
expenditure share on leisure, of type numeric
logarithm of total expenditure, of type numeric
logarithm of total earnings, of type numeric
'0' indicates no children, '1' indicates 1-2 children, of type numeric
Richard Blundell and Dennis Kristensen
Blundell, R., X. Chen and D. Kristensen (2007). “Semi-Nonparametric IV Estimation of Shape-Invariant Engel Curves.” Econometrica, 75(6), 1613-1669. doi:10.1111/j.1468-0262.2007.00808.x
Chen, X. and T. Christensen (2018). “Optimal Sup-norm Rates and Uniform Inference on Nonlinear Functionals of Nonparametric IV Regression.” Quantitative Economics, 9(1), 39-85. doi:10.3982/QE722
Chen, X., T. Christensen and S. Kankanala (2024). “Adaptive Estimation and Uniform Confidence Bands for Nonparametric Structural Functions and Elasticities.” Review of Economic Studies, forthcoming. doi:10.1093/restud/rdae025
## Load data data("Engel95", package = "npiv") ## Sort on logexp (the regressor) for plotting purposes Engel95 <- Engel95[order(Engel95$logexp),] attach(Engel95) logexp.eval <- seq(4.5,6.5,length=100) ## Estimate the Engel curve for food using logwages as an instrument food_engel <- npiv(food, logexp, logwages, X.eval = logexp.eval) ## Plot the estimated function and uniform confidence bands plot(food_engel, showdata = TRUE)
## Load data data("Engel95", package = "npiv") ## Sort on logexp (the regressor) for plotting purposes Engel95 <- Engel95[order(Engel95$logexp),] attach(Engel95) logexp.eval <- seq(4.5,6.5,length=100) ## Estimate the Engel curve for food using logwages as an instrument food_engel <- npiv(food, logexp, logwages, X.eval = logexp.eval) ## Plot the estimated function and uniform confidence bands plot(food_engel, showdata = TRUE)
gsl.bs
generates the B-spline basis matrix for a
polynomial spline and (optionally) the B-spline basis matrix
derivative of a specified order with respect to each predictor
gsl.bs(...) ## Default S3 method: gsl.bs(x, degree = 3, nbreak = 2, deriv = 0, x.min = NULL, x.max = NULL, intercept = FALSE, knots = NULL, ...)
gsl.bs(...) ## Default S3 method: gsl.bs(x, degree = 3, nbreak = 2, deriv = 0, x.min = NULL, x.max = NULL, intercept = FALSE, knots = NULL, ...)
x |
the predictor variable. Missing values are not allowed |
degree |
degree of the piecewise polynomial - default is ‘3’ (cubic spline) |
nbreak |
number of breaks in each interval - default is ‘2’ |
deriv |
the order of the derivative to be computed-default if
|
x.min |
the lower bound on which to construct the spline -
defaults to |
x.max |
the upper bound on which to construct the spline -
defaults to |
intercept |
if ‘TRUE’, an intercept is included in the basis; default is ‘FALSE’ |
knots |
a vector (default |
... |
optional arguments |
Typical usages are (see below for a list of options and also the examples at the end of this help file)
B <- gsl.bs(x,degree=10) B.predict <- predict(gsl.bs(x,degree=10),newx=xeval)
gsl.bs
returns a gsl.bs
object. A matrix of dimension
‘c(length(x), degree+nbreak-1)’. The generic function
predict
extracts (or generates) predictions from the
returned object.
A primary use is in modelling formulas to directly specify a piecewise polynomial term in a model. See https://www.gnu.org/software/gsl/ for further details.
Jeffrey S. Racine [email protected]
Chen, X., T. Christensen and S. Kankanala (2024). “Adaptive Estimation and Uniform Confidence Bands for Nonparametric Structural Functions and Elasticities.” Review of Economic Studies, forthcoming. doi:10.1093/restud/rdae025
## Plot the spline bases and their first order derivatives x <- seq(0,1,length=100) matplot(x,gsl.bs(x,degree=5),type="l") matplot(x,gsl.bs(x,degree=5,deriv=1),type="l") ## Regression example n <- 1000 x <- sort(runif(n)) y <- cos(2*pi*x) + rnorm(n,sd=.25) B <- gsl.bs(x,degree=5,intercept=FALSE) plot(x,y,cex=.5,col="grey") lines(x,fitted(lm(y~B)))
## Plot the spline bases and their first order derivatives x <- seq(0,1,length=100) matplot(x,gsl.bs(x,degree=5),type="l") matplot(x,gsl.bs(x,degree=5,deriv=1),type="l") ## Regression example n <- 1000 x <- sort(runif(n)) y <- cos(2*pi*x) + rnorm(n,sd=.25) B <- gsl.bs(x,degree=5,intercept=FALSE) plot(x,y,cex=.5,col="grey") lines(x,fitted(lm(y~B)))
npiv
performs nonparametric a structural function h0
and its derivatives using a B-spline sieve. It also constructs uniform confidence bands for h0
and its derivative.
Sieve dimensions are determined in a data-dependent way if not provided by the user, via the methods described in Chen, Christensen, and Kankanala (2024). This data-driven choice of sieve dimension ensures estimators of h0
and its derivatives converge at the optimal sup-norm rate. The resulting uniform confidence bands for h0
and its derivatives also converge at the minimax rate up to log factors; see Chen, Christensen, and Kankanala (2024).
If sieve dimensions are provided by the user, npiv
implements the bootstrap-based procedure of Chen and Christensen (2018) to construct uniform confidence bands based on undersmoothing for h0
and its derivatives.
The methods in npiv
apply to estimation and inference on a nonparametric regression function as a special case.
npiv(...) ## S3 method for class 'formula' npiv(formula, data=NULL, newdata=NULL, subset=NULL, na.action="na.omit", call, ...) ## Default S3 method: npiv(Y, X, W, X.eval=NULL, X.grid=NULL, alpha=0.05, basis=c("tensor","additive","glp"), boot.num=99, check.is.fullrank=FALSE, deriv.index=1, deriv.order=1, grid.num=50, J.x.degree=3, J.x.segments=NULL, K.w.degree=4, K.w.segments=NULL, K.w.smooth=2, knots=c("uniform","quantiles"), progress=TRUE, ucb.h=TRUE, ucb.deriv=TRUE, W.max=NULL, W.min=NULL, X.min=NULL, X.max=NULL, ...)
npiv(...) ## S3 method for class 'formula' npiv(formula, data=NULL, newdata=NULL, subset=NULL, na.action="na.omit", call, ...) ## Default S3 method: npiv(Y, X, W, X.eval=NULL, X.grid=NULL, alpha=0.05, basis=c("tensor","additive","glp"), boot.num=99, check.is.fullrank=FALSE, deriv.index=1, deriv.order=1, grid.num=50, J.x.degree=3, J.x.segments=NULL, K.w.degree=4, K.w.segments=NULL, K.w.smooth=2, knots=c("uniform","quantiles"), progress=TRUE, ucb.h=TRUE, ucb.deriv=TRUE, W.max=NULL, W.min=NULL, X.min=NULL, X.max=NULL, ...)
formula |
a symbolic description of the model to be fit. |
data |
an optional data frame containing the variables in the model. |
newdata |
an optional data frame in which to look for variables with which to predict (i.e., predictors in |
subset |
an optional vector specifying a subset of observations to be used in the fitting process (see additional details about how this argument interacts with data-dependent bases in the ‘Details’ section of the |
na.action |
a function which indicates what should happen when the data contain NAs. The default is set by the |
call |
the original function call (this is passed internally by |
Y |
dependent variable vector. |
X |
matrix of endogenous regressors. |
W |
matrix of instrumental variables. Set |
X.eval |
optional matrix of evaluation data for the endogenous regressors. |
X.grid |
optional vector of grid points for |
alpha |
nominal size of the uniform confidence bands. Default is |
basis |
basis type (if
|
boot.num |
number of bootstrap replications. |
check.is.fullrank |
check that |
deriv.index |
integer indicating the column of |
deriv.order |
integer indicating the order of derivative to be computed. |
grid.num |
number of grid points for each |
J.x.degree |
B-spline degree (integer or vector of integers of length |
J.x.segments |
B-spline number of segments (integer or vector of integers of length |
K.w.degree |
B-spline degree (integer or vector of integers of lenth |
K.w.segments |
B-spline number of segments (integer or vector of integers of length |
K.w.smooth |
non-negative integer. Basis for the nonparametric first-stage uses |
knots |
knots type, a character string. Options are:
|
progress |
whether to display progress bar or not. Default is |
ucb.h |
whether to compute a uniform confidence band for the structural function. Default is |
ucb.deriv |
whether to compute a uniform confidence band for the derivative of the structural function. Default is |
W.min |
lower bound on the support of each |
W.max |
upper bound on the support of each |
X.min |
lower bound on the support of each |
X.max |
upper bound on the support of each |
... |
optional arguments |
npiv
estimates and constructs uniform confidence bands for a nonparametric structural function and its derivatives in the model
Estimation is performed using nonparametric two-stage least-squares with a B-spline sieve. The key tuning parameter is the dimension
of the sieve used to approximate
. The dimension is tuned via modifying the number and placement of interior knots in the B-spline basis (equivalently, the number of segments of the basis). Sieve dimensions can be user-provided or data-determined using the procedure of Chen, Christensen, and Kankanala (2024).
Typical usages mirror ivreg
(see above and below for a list of options and the example at the bottom of this document)
foo <- npiv(y~x|w) foo <- npiv(y~x1+x2|w1+w2) foo <- npiv(Y=y,X=x,W=w)
npiv
can be used in two ways:
1. Data-driven sieve dimension is invoked if either K.w.segments
or J.x.segments
are unspecified or NULL
(the default). Sieve dimensions are chosen automatically using npiv_choose_J
. Uniform confidence bands for and its derivatives are constructed using the data-driven method of Chen, Christensen, and Kankanala (2024).
2. The user may specify the sieve dimensions of both bases by specifying values for K.w.segments
and J.x.segments
. Uniform confidence bands for and its derivatives are constructed using the method of Chen and Christensen (2018).
npiv
can also be used for estimation and inference on a nonparametric regression function by setting W=X
.
npiv
returns a npiv
object. The generic function fitted
extracts the estimated values for the sample (or evaluation data, if provided), while the generic function residuals
extracts the sample residuals. The generic function summary
provides a simple model summary. The generic function plot
also plots the estimated function and derivative, together with uniform confidence bands.
The function npiv
returns a list with the following components:
h |
estimated structural function evaluated at the sample data (or evaluation data, if provided). |
residuals |
residuals for the sample data. |
deriv |
estimated derivative of the structural function evaluated at the sample data (or evaluation data, if provided). |
asy.se |
pre-asymptotic standard errors for the estimator of the structural function evaluated at the sample data (or evaluation data, if provided) |
deriv.asy.se |
pre-asymptotic standard errors for the estimator of the derivative of the structural function evaluated at the sample data (or evaluation data, if provided). |
deriv.index |
index for the estimated derivative. |
deriv.order |
order of the estimated derivative. |
K.w.degree |
value of |
K.w.segments |
value of |
J.x.degree |
value of |
J.x.segments |
value of |
beta |
vector of estimated spline coefficients. |
Jeffrey S. Racine <[email protected]>, Timothy Christensen <[email protected]>
Chen, X. and T. Christensen (2018). “Optimal Sup-norm Rates and Uniform Inference on Nonlinear Functionals of Nonparametric IV Regression.” Quantitative Economics, 9(1), 39-85. doi:10.3982/QE722
Chen, X., T. Christensen and S. Kankanala (2024). “Adaptive Estimation and Uniform Confidence Bands for Nonparametric Structural Functions and Elasticities.” Review of Economic Studies, forthcoming. doi:10.1093/restud/rdae025
## load data data("Engel95", package = "npiv") ## sort on logexp (the regressor) for plotting purposes Engel95 <- Engel95[order(Engel95$logexp),] attach(Engel95) ## Estimate the Engel curve for food using logwages as an instrument fm1 <- npiv(food ~ logexp | logwages) ## Plot the estimated Engel curve and data-driven uniform confidence bands plot(logexp,food, ylab="Food Budget Share", xlab="log(Total Household Expenditure)", xlim=c(4.75, 6.25), ylim=c(0, 0.4), main="", type="p", cex=.5, col="lightgrey") lines(logexp,fm1$h,col="blue",lwd=2,lty=1) lines(logexp,fm1$h.upper,col="blue",lwd=2,lty=2) lines(logexp,fm1$h.lower,col="blue",lwd=2,lty=2) ## Estimate the Engel curve using pre-specified sieve dimension ## (dimension 5 for logexp, dimension 9 for logwages) fm2 <- npiv(food ~ logexp | logwages, J.x.segments = 2, K.w.segments = 5) ## Plot uniform confidence bands based on undersmoothing lines(logexp,fm2$h.upper,col="red",lwd=2,lty=2) lines(logexp,fm2$h.lower,col="red",lwd=2,lty=2) ## Plot pointwise confidence bands based on pre-asymptotic standard errors lines(logexp,fm2$h+1.96*fm2$asy.se,col="red",lwd=2,lty=3) lines(logexp,fm2$h-1.96*fm2$asy.se,col="red",lwd=2,lty=3) legend("topright", legend=c("Data-driven Estimate", "Data-driven UCBs", "Undersmoothed UCBs", "Pointwise CBs"), col=c("blue","blue","red","red"), lty=c(1,2,2,3), lwd=c(2,2,2,2)) ## Plot the data-driven estimate of the derivative of the Engel curve plot(logexp,fm1$deriv,col="blue",lwd=2,lty=1,type="l", ylab="Derivative of Food Budget Share", xlab="log(Total Household Expenditure)", xlim=c(4.75, 6.25), ylim=c(-1,1)) ## Plot data-driven uniform confidence bands for the derivative lines(logexp,fm1$h.upper.deriv,col="blue",lwd=2,lty=2) lines(logexp,fm1$h.lower.deriv,col="blue",lwd=2,lty=2) ## Plot uniform confidence bands based on undersmoothing lines(logexp,fm2$h.upper.deriv,col="red",lwd=2,lty=2) lines(logexp,fm2$h.lower.deriv,col="red",lwd=2,lty=2) ## Plot pointwise confidence bands based on pre-asymptotic standard errors lines(logexp,fm2$deriv+1.96*fm2$deriv.asy.se,col="red",lwd=2,lty=3) lines(logexp,fm2$deriv-1.96*fm2$deriv.asy.se,col="red",lwd=2,lty=3) legend("topright", legend=c("Data-driven Estimate", "Data-driven UCBs", "Undersmoothed UCBs", "Pointwise CBs"), col=c("blue","blue","red","red"), lty=c(1,2,2,3), lwd=c(2,2,2,2))
## load data data("Engel95", package = "npiv") ## sort on logexp (the regressor) for plotting purposes Engel95 <- Engel95[order(Engel95$logexp),] attach(Engel95) ## Estimate the Engel curve for food using logwages as an instrument fm1 <- npiv(food ~ logexp | logwages) ## Plot the estimated Engel curve and data-driven uniform confidence bands plot(logexp,food, ylab="Food Budget Share", xlab="log(Total Household Expenditure)", xlim=c(4.75, 6.25), ylim=c(0, 0.4), main="", type="p", cex=.5, col="lightgrey") lines(logexp,fm1$h,col="blue",lwd=2,lty=1) lines(logexp,fm1$h.upper,col="blue",lwd=2,lty=2) lines(logexp,fm1$h.lower,col="blue",lwd=2,lty=2) ## Estimate the Engel curve using pre-specified sieve dimension ## (dimension 5 for logexp, dimension 9 for logwages) fm2 <- npiv(food ~ logexp | logwages, J.x.segments = 2, K.w.segments = 5) ## Plot uniform confidence bands based on undersmoothing lines(logexp,fm2$h.upper,col="red",lwd=2,lty=2) lines(logexp,fm2$h.lower,col="red",lwd=2,lty=2) ## Plot pointwise confidence bands based on pre-asymptotic standard errors lines(logexp,fm2$h+1.96*fm2$asy.se,col="red",lwd=2,lty=3) lines(logexp,fm2$h-1.96*fm2$asy.se,col="red",lwd=2,lty=3) legend("topright", legend=c("Data-driven Estimate", "Data-driven UCBs", "Undersmoothed UCBs", "Pointwise CBs"), col=c("blue","blue","red","red"), lty=c(1,2,2,3), lwd=c(2,2,2,2)) ## Plot the data-driven estimate of the derivative of the Engel curve plot(logexp,fm1$deriv,col="blue",lwd=2,lty=1,type="l", ylab="Derivative of Food Budget Share", xlab="log(Total Household Expenditure)", xlim=c(4.75, 6.25), ylim=c(-1,1)) ## Plot data-driven uniform confidence bands for the derivative lines(logexp,fm1$h.upper.deriv,col="blue",lwd=2,lty=2) lines(logexp,fm1$h.lower.deriv,col="blue",lwd=2,lty=2) ## Plot uniform confidence bands based on undersmoothing lines(logexp,fm2$h.upper.deriv,col="red",lwd=2,lty=2) lines(logexp,fm2$h.lower.deriv,col="red",lwd=2,lty=2) ## Plot pointwise confidence bands based on pre-asymptotic standard errors lines(logexp,fm2$deriv+1.96*fm2$deriv.asy.se,col="red",lwd=2,lty=3) lines(logexp,fm2$deriv-1.96*fm2$deriv.asy.se,col="red",lwd=2,lty=3) legend("topright", legend=c("Data-driven Estimate", "Data-driven UCBs", "Undersmoothed UCBs", "Pointwise CBs"), col=c("blue","blue","red","red"), lty=c(1,2,2,3), lwd=c(2,2,2,2))
npiv_choose_J
implements the data-driven choice of sieve dimension developed in Chen, Christensen, and Kankanala (2024) for nonparametric instrumental variables estimation using a B-spline sieve. It applies to nonparametric regression as a special case.
npiv_choose_J(Y, X, W, X.grid = NULL, J.x.degree = 3, K.w.degree = 4, K.w.smooth = 2, knots = c("uniform", "quantiles"), basis = c("tensor", "additive", "glp"), X.min = NULL, X.max = NULL, W.min = NULL, W.max = NULL, grid.num = 50, boot.num = 99, check.is.fullrank = FALSE, progress = TRUE)
npiv_choose_J(Y, X, W, X.grid = NULL, J.x.degree = 3, K.w.degree = 4, K.w.smooth = 2, knots = c("uniform", "quantiles"), basis = c("tensor", "additive", "glp"), X.min = NULL, X.max = NULL, W.min = NULL, W.max = NULL, grid.num = 50, boot.num = 99, check.is.fullrank = FALSE, progress = TRUE)
Y |
dependent variable vector. |
X |
matrix of endogenous regressors. |
W |
matrix of instrumental variables. Set |
X.grid |
vector of grid point(s). Default uses 50 equally spaced points over the support of each |
J.x.degree |
B-spline degree (integer or vector of integers of length |
K.w.degree |
B-spline degree (integer or vector of integers of lenth |
K.w.smooth |
non-negative integer. Basis for the nonparametric first-stage uses |
knots |
knots type, a character string. Options are:
|
basis |
basis type (if
|
X.min |
lower bound on the support of each |
X.max |
upper bound on the support of each |
W.min |
lower bound on the support of each |
W.max |
upper bound on the support of each |
grid.num |
number of grid points for each |
boot.num |
number of bootstrap replications. |
check.is.fullrank |
check that |
progress |
whether to display progress bar or not. Default is |
J.hat.max |
largest element of candidate set of sieve dimensions searched over. |
J.hat.n |
second largest element of candidate set of sieve dimensions searched over. |
J.hat |
bootstrap-based Lepski choice of sieve dimension. |
J.tilde |
data-driven choice of sieve dimension using the method of Chen, Christensen, and Kankanala (2024). Minimum of |
J.x.seg |
data-driven number of segments for |
K.w.seg |
data-driven number of segments for |
theta.star |
Lepski critical value used in determination of |
Jeffrey S. Racine <[email protected]>, Timothy Christensen <[email protected]>
Chen, X., T. Christensen and S. Kankanala (2024). “Adaptive Estimation and Uniform Confidence Bands for Nonparametric Structural Functions and Elasticities.” Review of Economic Studies, forthcoming. doi:10.1093/restud/rdae025
library(MASS) ## Simulate the data n <- 10000 cov.ux <- 0.5 var.u <- 0.1 mu <- c(1,1,0) Sigma <- matrix(c(1.0,0.85,cov.ux, 0.85,1.0,0.0, cov.ux,0.0,1.0), 3,3, byrow=TRUE) foo <- mvrnorm(n = n, mu, Sigma) X <- 2*pnorm(foo[,1],mean=mu[1],sd=sqrt(Sigma[1,1])) -1 W <- 2*pnorm(foo[,2],mean=mu[2],sd=sqrt(Sigma[2,2])) -1 U <- foo[,3] ## Cosine structural function h0 <- sin(pi*X) Y <- h0 + sqrt(var.u)*U npiv_choose_J(Y,X,W)
library(MASS) ## Simulate the data n <- 10000 cov.ux <- 0.5 var.u <- 0.1 mu <- c(1,1,0) Sigma <- matrix(c(1.0,0.85,cov.ux, 0.85,1.0,0.0, cov.ux,0.0,1.0), 3,3, byrow=TRUE) foo <- mvrnorm(n = n, mu, Sigma) X <- 2*pnorm(foo[,1],mean=mu[1],sd=sqrt(Sigma[1,1])) -1 W <- 2*pnorm(foo[,2],mean=mu[2],sd=sqrt(Sigma[2,2])) -1 U <- foo[,3] ## Cosine structural function h0 <- sin(pi*X) Y <- h0 + sqrt(var.u)*U npiv_choose_J(Y,X,W)