Function documentation
Problem construction and manipulation
NPDemand.define_problem
— Functiondefine_problem(df::DataFrame; exchange = [], index_vars = ["prices"], FE = [], constraints = [], bO = 2, tol = 1e-5)
Constructs a problem
::NPDProblem using the provided problem characteristics. Inputs:
exchange
: A vector of groups of products which are exchangeable. E.g., with 4 goods, if the first
and second are exchangeable and so are the third and fourth, set exchange
= [[1 2], [3 4]].
index_vars
: String array listing column names indf
which represent variables that enter the inverted index.- "prices" must be the first element of
index_vars
- "prices" must be the first element of
FE
: String array listing column names indf
which should be included as fixed effects.- Note: All fixed effects are estimated as parameters by the minimizer, so be careful adding fixed effects for variables that take many values.
bO
: Order of the univariate Bernstein polynomials in market shares. Default is 2.constraint_tol
: Tolerance specifying tightness of constraintschunk_size
: Controls chunk size in ForwardDiff.Gradient autodiff calculation for nonlinear constraints. Only used if :subs_in_group specified.constraints
: A list of symbols of accepted constraints. Currently supported constraints are:- :monotone
- :all_substitutes
- :diagonal_dominance_group
- :diagonal_dominance_all
- :subs_in_group (Note: this constraint is the only available nonlinear constraint and will slow down estimation considerably)
verbose
: iffalse
, will not print updates as problem is generated
NPDemand.update_constraints!
— Functionupdate_constraints!(problem::NPDProblem, new_constraints::Vector{Symbol})
Re-calculates constraint matrices under new_constraints
, ignoring previous constraints used to define the problem. The exchangeability
constraint is an exception. To change anything about the structure of exchangeability for the problem changes, define a new problem.
NPDemand.list_constraints
— Functionlist_constraints()
Returns a dictionary of constraints supported by the define_problem
function in the NPDemand package. Each constraint is represented by a symbol and its corresponding description.
Estimation
NPDemand.estimate!
— Functionestimate!(problem::NPDProblem;
verbose = true,
linear_solver = "Ipopt",
quasi_bayes = false,
sampler = [],
n_samples::Int = 50_000,
burn_in::Real = 0.25,
skip::Int = 5,
n_attempts = 0,
penalty = 0,
step::Union{Real, Symbol} = 0.01)
Estimates the problem using the specified parameters.
Arguments
problem::NPDProblem
: The problem to be estimated.verbose::Bool
: Whether to print verbose output. Default istrue
.linear_solver::String
: The linear solver to use. Must be either "Ipopt" or "OSQP". Default is "Ipopt".quasi_bayes::Bool
: Whether to use quasi-bayes sampling. Default isfalse
.sampler
: The sampler to use for quasi-bayes sampling. Default is an empty array.n_samples::Int
: The number of samples to draw. Default is 50,000.burn_in::Real
: The fraction of samples to drop as burn-in. Must be less than 1. Default is 0.25.skip::Int
: The number of samples to skip between saved samples. Default is 5.n_attempts
: The number of attempts to find a valid starting point for the sampler. Default is 0.penalty
: The penalty value for the objective function. Default is 0.step::Union{Real, Symbol}
: The step size for the sampler. Can be a real number or the symbol:auto
to automatically calculate the step size. Default is 0.01.
NPDemand.smc!
— Functionsmc!(problem::NPDemand.NPDProblem;
grid_points::Int = 50,
max_penalty::Real = 5,
ess_threshold::Real = 100,
step::Real = 0.1,
skip::Int = 5,
burn_in::Real = 0.25,
mh_steps = max(5, floor(size(problem.results.filtered_chain, 2))/10),
seed = 4132,
smc_method = :adaptive,
max_iter = 1000,
adaptive_tolerance = false,
max_violations = 0.01)
Run sequentially constrained Monte Carlo (SMC) on the problem.
Arguments
problem::NPDemand.NPDProblem
: The problem object on which SMC will be run.
Optional Arguments
grid_points::Int
: The number of grid points for the SMC grid. Default is 50.max_penalty::Real
: The maximum penalty value for the SMC algorithm. Default is 5.ess_threshold::Real
: The effective sample size threshold for the SMC algorithm. Default is 100.step::Real
: The step size for the SMC algorithm. Default is 0.1.skip::Int
: The thinning factor for the SMC chain. Default is 5.burn_in::Real
: The fraction of samples to be discarded as burn-in. Default is 0.25.mh_steps
: The number of Metropolis-Hastings steps per iteration. Default is calculated based on the size of the filtered chain.seed
: The random seed for the SMC algorithm. Default is 4132.smc_method
: The method for choosing the SMC grid. Default is :adaptive. Other options are [:linear_grid, :geometric_grid, and :logit_grid], which specify grids of each form between zero and the maximum penalty.max_iter
: The maximum number of iterations for the SMC algorithm. Default is 1000.adaptive_tolerance
: Whether to use adaptive tolerance for the SMC algorithm. Default is false.max_violations
: The maximum allowed fraction of markets with violations. Default is 0.01.
The function will overwrite the results in the problem object with the resulting chain.
For harder or slower problems, it may be necessary to increase the number of Metropolis-Hastings steps per iteration (mh_steps
), the number of iterations (max_iter
), or the maximum allowed fraction markets with violations (max_violations
).
burn_in
and skip
control the number of samples to drop and the thinning of the chain, respectively.
Post-estimation tools
NPDemand.price_elasticities!
— Functionprice_elasticities!(problem;
CI::Union{Vector{Any}, Real} = [],
n_draws::Union{Vector{Any}, Int} = [])
Takes the solved problem
as first argument, a DataFrame
as the second argument, and evaluates all price elasticities in-sample. Currently does not calculate out-of-sample price elasticities. For this, use the function compute_demand_function!
.
Results of this function are stored as a DataFrame
in problem.allelasticities. Results can be summarized by hand or using the `summarizeelasticities` function. We also store the Jacobian of the demand function with respect to prices, which can be used to calculate other quantities of interest.
NPDemand.elasticity_quantiles
— Functionelasticity_quantiles(problem::NPDProblem, ind1::Int, ind2::Int;
quantiles = collect(0.01:0.01:0.99),
n_draws::Int = 100)
Convenience function for calculating quantiles of price elasticities. problem
should be a solved NPDProblem
after running price_elasticities!
. ind1
and ind2
are the indices of the products for which you want to calculate the quantiles. E.g., ind1=1, ind2=1 returns quantiles of the own-price elasticity for product 1. The user can control the set of quantiles to return and the number of draws (n_draws
) to use in the posterior integration.
NPDemand.report_constraint_violations
— Functionreport_constraint_violations(problem;
verbose = true,
params = [],
output = "dict",
n_draws::Int = 0)
This function reports the constraint violations for a given problem.
Arguments
problem
: The problem for which constraint violations need to be reported.verbose
: A boolean indicating whether to print verbose output. Default istrue
.params
: An optional parameter vector. If not provided, the function uses the GMM result if there's no chain, or the filtered chain if available. Default is[]
.output
: The output format for the violations. Default is"dict"
.n_draws
: The number of draws to use from the filtered chain. Default is0
.
Returns
violations
: A dictionary or a single violation value, depending on theoutput
parameter.
Example output:
Dict{Symbol, Float64} with 8 entries:
:monotone => 0.0
:all_substitutes => 0.0
:diagonal_dominance_all => 0.0
NPDemand.compute_demand_function!
— Functioncompute_demand_function!(problem, df; max_iter = 1000, show_trace = false)
compute_demand_function!
estimates the demand function/curve using NPD estimates calculated via estimate!.
The function takes in an estimated problem::NPDProblem and a dataframe with counterfactual values of the covariates in the utility function. One must specify all fields that were used in estimation (including shares). The function will change the values of df[!,r"shares"] to take on the value of the estimated demand function.
Options:
max_iter
: controls the number of iterations for the nonlinear solver calculating market shares. Default is 1000 but well-estimated problems should converge faster.show_trace
: iftrue
, Optim will print the trace for each iteration of the nonlinear solver.compute_own_elasticities
: NOT yet implemented– iftrue
, will also generate columns calledown_elast
which will include estimated own-price elasticities at all counterfactual points.
NPDemand.summarize_elasticities
— Functionsummarize_elasticities(problem::NPDProblem, which_elasticities::String, stat::String;
q = [],
integrate = false,
n_draws::Int = 100,
CI::Real = 0.95)
Convenience function for summarizing price elasticities. problem
should be a solved NPDProblem
after running price_elasticities!
. stat
should be in ["mean", "quantile"], and if stat
=="quantile", the option q
should include the quantile of interest (e.g., 0.75 for the 75th percentile price elasticities).
When a problem has been estimated via quasi-Bayesian methods, the function can integrate over the posterior distribution of the parameters to provide a posterior distribution of the summarized value. n_draws
controls the number of draws from the posterior to use in the integration, and CI
controls the with of the credible interval to use in the integration (0.95 for a 95% credible interval).
Output is a NamedTuple: (;Statistic, PosteriorMean, PosteriorMedian, Posterior_CI)
Back-end functions
NPDemand.NPD_parameters
— TypeNPD_parameters
Custom struct to store estimated parameters specifically. This can be used to replace the candidate parameters in an NPDProblem struct. The two key fields
are `minimizer` and `filtered_chain`. The `minimizer` field stores the estimated parameters, while the `filtered_chain` field stores the Markov chain for quasi-Bayes
methods, after filtering out burn-in and thinning but before reformatting into the full parameter sieve.
NPDemand.simulate_logit
— Functionsimulate_logit(J,T, beta, v)
Simulates logit demand for J
products in T
markets, with price preference parameter beta
and market shocks with standard deviation v
.
NPDemand.toDataFrame
— FunctiontoDataFrame(s::Matrix, p::Matrix, z::Matrix, x::Matrix = zeros(size(s)))
Converts the results of simulate_logit
to a DataFrame with column names that can be processed by NPDemand.jl
NPDemand.bern
— Functionbern(t, order)
Returns a univariate Bernstein polynomial of order order
constructed from array/matrix t
NPDemand.dbern
— Functiondbern(t, order)
Returns the derivative of a univariate Bernstein polynomial of order order
constructed from array/matrix t