class
BicopA class for bivariate copula models.
The copula model is fully characterized by the family, rotation, and parameters.
Constructors, destructors, conversion operators
-
Bicop(const BicopFamily family = BicopFamily::
indep, const int rotation = 0, const Eigen::MatrixXd& parameters = Eigen::MatrixXd(), const std::vector<std::string>& var_types = { "c", "c" }) - Instantiates a specific bivariate copula model.
- Bicop(const Eigen::MatrixXd& data, const FitControlsBicop& controls = FitControlsBicop(), const std::vector<std::string>& var_types = { "c", "c" }) explicit
- Instantiates from data.
- Bicop(const Bicop& other)
- Copy constructor (deep copy)
- Bicop(const std::string& filename) explicit
- Instantiates from a JSON file.
- Bicop(const nlohmann::json& input) explicit
- Instantiates from a nlohmann::json object.
Public functions
- auto operator=(Bicop other) -> Bicop&
- Copy assignment operator (deep copy)
- auto to_json() const -> nlohmann::json
- Convert the copula into a nlohmann::json object.
- void to_file(const std::string& filename) const
- Write the copula object into a JSON file.
- auto get_npars() const -> double
- Returns the actual number of parameters for parameteric families.
- auto pdf(const Eigen::MatrixXd& u) const -> Eigen::VectorXd
- Evaluates the copula density.
- auto cdf(const Eigen::MatrixXd& u) const -> Eigen::VectorXd
- Evaluates the copula distribution.
- auto hfunc1(const Eigen::MatrixXd& u) const -> Eigen::VectorXd
- Evaluates the first h-function.
- auto hfunc2(const Eigen::MatrixXd& u) const -> Eigen::VectorXd
- Evaluates the second h-function.
- auto hinv1(const Eigen::MatrixXd& u) const -> Eigen::VectorXd
- Evaluates the inverse of the first h-function.
- auto hinv2(const Eigen::MatrixXd& u) const -> Eigen::VectorXd
- Evaluates the inverse of the second h-function.
- auto simulate(const size_t& n, const bool qrng = false, const std::vector<int>& seeds = std::vector<int>()) const -> Eigen::MatrixXd
- Simulates from a bivariate copula.
- void fit(const Eigen::MatrixXd& data, const FitControlsBicop& controls = FitControlsBicop())
- void select(const Eigen::MatrixXd& data, FitControlsBicop controls = FitControlsBicop())
- Selects the best fitting model.
- auto loglik(const Eigen::MatrixXd& u = Eigen::MatrixXd()) const -> double
- Evaluates the log-likelihood.
- auto aic(const Eigen::MatrixXd& u = Eigen::MatrixXd()) const -> double
- Evaluates the Akaike information criterion (AIC).
- auto bic(const Eigen::MatrixXd& u = Eigen::MatrixXd()) const -> double
- Evaluates the Bayesian information criterion (BIC).
- auto mbic(const Eigen::MatrixXd& u = Eigen::MatrixXd(), const double psi0 = 0.9) const -> double
- Evaluates the modified Bayesian information criterion (mBIC).
- auto parameters_to_tau(const Eigen::MatrixXd& parameters) const -> double
- Converts the copula parameters to Kendall's .
- auto tau_to_parameters(const double& tau) const -> Eigen::MatrixXd
- Converts a Kendall's into copula parameters.
Getters and setters
- auto get_family() const -> BicopFamily
- Gets the copula family.
- auto get_family_name() const -> std::string
- Gets the copula family as a string.
- auto get_rotation() const -> int
- Gets the rotation.
- auto get_parameters() const -> Eigen::MatrixXd
- Gets the parameters.
- auto get_tau() const -> double
- Gets the Kendall's tau.
- auto get_loglik() const -> double
- Gets the log-likelihood (only for fitted objects).
- auto get_nobs() const -> size_t
- Gets the number of observations (only for fitted objects).
- auto get_aic() const -> double
- Gets the aic (only for fitted objects).
- auto get_bic() const -> double
- Gets the bic (only for fitted objects).
- auto get_mbic(const double psi0 = 0.9) const -> double
- Gets the modified bic (only for fitted objects).
- void set_rotation(const int rotation)
- Sets the rotation.
- void set_var_types(const std::vector<std::string>& var_types = { "c", "c" })
- Sets variable types.
- auto get_var_types() const -> std::vector<std::string>
- Gets variable types.
Utilities
useful functions for bivariate copulas
- auto str() const -> std::string
- Summarizes the model into a string (can be used for printing).
- void flip()
- Adjusts the copula model to a change in the variable order.
- auto get_parameters_lower_bounds() const -> Eigen::MatrixXd
- Gets lower bounds for copula parameters.
- auto get_parameters_upper_bounds() const -> Eigen::MatrixXd
- Gets upper bounds for copula parameters.
Function documentation
vinecopulib:: Bicop:: Bicop(const BicopFamily family = BicopFamily:: indep,
const int rotation = 0,
const Eigen::MatrixXd& parameters = Eigen::MatrixXd(),
const std::vector<std::string>& var_types = { "c", "c" })
Instantiates a specific bivariate copula model.
Parameters | |
---|---|
family | The copula family. |
rotation | The rotation of the copula; one of 0, 90, 180, or 270 (for Independence, Gaussian, Student, Frank, and nonparametric families, only 0 is allowed). |
parameters | The copula parameters. |
var_types | Two strings specifying the types of the variables, e.g., ("c", "d") means first variable continuous, second discrete. |
vinecopulib:: Bicop:: Bicop(const Eigen::MatrixXd& data,
const FitControlsBicop& controls = FitControlsBicop(),
const std::vector<std::string>& var_types = { "c", "c" }) explicit
Instantiates from data.
Parameters | |
---|---|
data | See select() . |
controls | See select() . |
var_types | Two strings specifying the types of the variables, e.g., ("c", "d") means first variable continuous, second discrete. |
Equivalent to creating a default Bicop()
and then selecting the model using select()
.
vinecopulib:: Bicop:: Bicop(const std::string& filename) explicit
Instantiates from a JSON file.
Parameters | |
---|---|
filename | The name of the JSON file to read. |
The input file contains four attributes: "fam"
, "rot"
, "par"
, "vt"
respectively a string for the family name, an integer for the rotation, and a numeric matrix for the parameters, and a list of two strings for the variable types.
nlohmann::json vinecopulib:: Bicop:: to_json() const
Convert the copula into a nlohmann::json object.
Returns | the nlohmann::json object containing the copula. |
---|
The nlohmann::json is contains of three values named "fam"
, "rot"
, "par"
, "vt"
, respectively a string for the family name, an integer for the rotation, a numeric matrix for the parameters and a list of two strings for the variables types.
void vinecopulib:: Bicop:: to_file(const std::string& filename) const
Write the copula object into a JSON file.
Parameters | |
---|---|
filename | The name of the file to write. |
The written file contains four attributes: "fam"
, "rot"
, "par"
, "vt"
, "nobs"
, "ll"
, "npars"
respectively a string for the family name, an integer for the rotation, and a numeric matrix for the parameters, a list of two strings for the variable types, an integer for the number of observations (if fitted), a double for the log-likelihood (if fitted), and a double for the number of parameters (can be non-integer in nonparametric models).
double vinecopulib:: Bicop:: get_npars() const
Returns the actual number of parameters for parameteric families.
For nonparametric families, there is a conceptually similar definition in the sense that it can be used in the calculation of fit statistics.
Eigen::VectorXd vinecopulib:: Bicop:: pdf(const Eigen::MatrixXd& u) const
Evaluates the copula density.
Parameters | |
---|---|
u | An matrix of observations contained in , where is the number of discrete variables. |
Returns | The copula density evaluated at u . |
The copula density is defined as joint density divided by marginal densities, irrespective of variable types.
Eigen::VectorXd vinecopulib:: Bicop:: cdf(const Eigen::MatrixXd& u) const
Evaluates the copula distribution.
Parameters | |
---|---|
u | An matrix of observations contained in , where is the number of discrete variables. |
Returns | The copula distribution evaluated at u . |
Eigen::VectorXd vinecopulib:: Bicop:: hfunc1(const Eigen::MatrixXd& u) const
Evaluates the first h-function.
Parameters | |
---|---|
u | An matrix of observations contained in , where is the number of discrete variables. |
The first h-function is .
Eigen::VectorXd vinecopulib:: Bicop:: hfunc2(const Eigen::MatrixXd& u) const
Evaluates the second h-function.
Parameters | |
---|---|
u | An matrix of observations contained in , where is the number of discrete variables. |
The second h-function is .
Eigen::VectorXd vinecopulib:: Bicop:: hinv1(const Eigen::MatrixXd& u) const
Evaluates the inverse of the first h-function.
Parameters | |
---|---|
u | An matrix of observations contained in , where is the number of discrete variables. |
The first h-function is . The inverse is calulated w.r.t. the second argument.
Eigen::VectorXd vinecopulib:: Bicop:: hinv2(const Eigen::MatrixXd& u) const
Evaluates the inverse of the second h-function.
Parameters | |
---|---|
u | An matrix of observations contained in , where is the number of discrete variables. |
The second h-function is . The inverse is calculated w.r.t. the first argument.
Eigen::MatrixXd vinecopulib:: Bicop:: simulate(const size_t& n,
const bool qrng = false,
const std::vector<int>& seeds = std::vector<int>()) const
Simulates from a bivariate copula.
Parameters | |
---|---|
n | Number of observations. |
qrng | Set to true for quasi-random numbers. |
seeds | Seeds of the (quasi-)random number generator; if empty (default), the (quasi-)random number generator is seeded randomly. |
Returns | An matrix of samples from the copula model. |
If qrng = TRUE
, generalized Halton sequences are used. For more information on Generalized Halton sequences, see Faure, H., Lemieux, C. (2009). Generalized Halton Sequences in 2008: A Comparative Study. ACM-TOMACS 19(4), Article 15.
void vinecopulib:: Bicop:: fit(const Eigen::MatrixXd& data,
const FitControlsBicop& controls = FitControlsBicop())
Fits a bivariate copula (with fixed family) to data.
For parametric models, two different methods are available. "mle"
fits the parameters by maximum-likelihood. "itau"
uses inversion of Kendall's , but is only available for one-parameter families and the Student t copula. For the latter, there is a one-to-one transformation for the first parameter, the second is found by profile likelihood optimization (with accuracy of at least 0.5). Nonparametric families have specialized methods, no specification is required.
When at least one variable is discrete, two types of "observations" are required: the first block contains realizations of . Let denote the number of discrete variables (either one or two). Then the second block contains realizations of . The minus indicates a left-sided limit of the cdf. For continuous variables the left limit and the cdf itself coincide. For, e.g., an integer-valued variable, it holds .
Incomplete observations (i.e., ones with a NaN value) are discarded. @param data An \f$ n \times (2 + k) \f$ matrix of observations contained in \f$(0, 1) \f$, where \f$ k \f$ is the number of discrete variables. @param controls The controls (see FitControlsBicop).
void vinecopulib:: Bicop:: select(const Eigen::MatrixXd& data,
FitControlsBicop controls = FitControlsBicop())
Selects the best fitting model.
The function calls fit()
for all families in family_set
and selecting the best fitting model by either BIC or AIC, see bic()
and aic()
.
When at least one variable is discrete, two types of "observations" are required: the first block contains realizations of . Let denote the number of discrete variables (either one or two). Then the second block contains realizations of . The minus indicates a left-sided limit of the cdf. For continuous variables the left limit and the cdf itself coincide. For, e.g., an integer-valued variable, it holds .
Incomplete observations (i.e., ones with a NaN value) are discarded. @param data An \f$ n \times (2 + k) \f$ matrix of observations contained in \f$(0, 1) \f$, where \f$ k \f$ is the number of discrete variables. @param controls The controls (see FitControlsBicop).
double vinecopulib:: Bicop:: loglik(const Eigen::MatrixXd& u = Eigen::MatrixXd()) const
Evaluates the log-likelihood.
Parameters | |
---|---|
u | An matrix of observations contained in , where is the number of discrete variables. |
The log-likelihood is defined as
where is the copula density, see Bicop::
.
double vinecopulib:: Bicop:: aic(const Eigen::MatrixXd& u = Eigen::MatrixXd()) const
Evaluates the Akaike information criterion (AIC).
Parameters | |
---|---|
u | An matrix of observations contained in , where is the number of discrete variables. |
The AIC is defined as
where is the log-liklihood (see loglik()
) and is the (effective) number of parameters of the model. The AIC is a consistent model selection criterion even for nonparametric models.
double vinecopulib:: Bicop:: bic(const Eigen::MatrixXd& u = Eigen::MatrixXd()) const
Evaluates the Bayesian information criterion (BIC).
Parameters | |
---|---|
u | An matrix of observations contained in , where is the number of discrete variables. |
The BIC is defined as
where is the log-liklihood (see loglik()
) and is the (effective) number of parameters of the model. The BIC is a consistent model selection criterion for parametric models.
double vinecopulib:: Bicop:: mbic(const Eigen::MatrixXd& u = Eigen::MatrixXd(),
const double psi0 = 0.9) const
Evaluates the modified Bayesian information criterion (mBIC).
Parameters | |
---|---|
u | An matrix of observations contained in , where is the number of discrete variables. |
psi0 | Prior probability of a non-independence copula. |
The mBIC is defined as
where is the \log-liklihood (see loglik()
), is the (effective) number of parameters of the model, and is the prior probability of having a non-independence copula and is an indicator for the family being non-independence.
double vinecopulib:: Bicop:: parameters_to_tau(const Eigen::MatrixXd& parameters) const
Converts the copula parameters to Kendall's .
Parameters | |
---|---|
parameters | The parameters (must be a valid parametrization of the current family). |
Eigen::MatrixXd vinecopulib:: Bicop:: tau_to_parameters(const double& tau) const
Converts a Kendall's into copula parameters.
Parameters | |
---|---|
tau | A value in . |
It only works for one-parameter families.
void vinecopulib:: Bicop:: set_var_types(const std::vector<std::string>& var_types = { "c", "c" })
Sets variable types.
Parameters | |
---|---|
var_types | A vector of size two specifying the types of the variables, e.g., {"c", "d"} means first variable continuous, second discrete. |