//' Construct a copula using uniform sampling from the unit simplex
//'
//' Given two families of parallel hyperplanes intersecting the canonical simplex, this function uniformly samples from the canonical simplex and construct an approximation of the bivariate probability distribution, called copula.
//'
//' @param h1 A \eqn{d}-dimensional vector that describes the direction of the first family of parallel hyperplanes.
//' @param h2 A \eqn{d}-dimensional vector that describes the direction of the second family of parallel hyperplanes.
//' @param numSlices The number of the slices for the copula. Default value is 100.
//' @param N The number of points to sample. Default value is \eqn{4\cdot 10^6}.
//'
//' @references \cite{L. Cales, A. Chalkis, I.Z. Emiris, V. Fisikopoulos,
//' \dQuote{Practical volume computation of structured convex bodies, and an application to modeling portfolio dependencies and financial crises,} \emph{Proc. of Symposium on Computational Geometry, Budapest, Hungary,} 2018.}
//'
//' @return A \eqn{numSlices\times numSlices} numerical matrix that corresponds to a copula.
//' @examples
//' # compute a copula for two random families of parallel hyperplanes
//' h1 = runif(n = 10, min = 1, max = 1000)
//' h1 = h1 / 1000
//' h2=runif(n = 10, min = 1, max = 1000)
//' h2 = h2 / 1000
//' cop = copula1(h1=h1, h2=h2, numSlices = 10, N = 100000)
//' @export
// [[Rcpp::export]]
Rcpp::NumericMatrix copula1 (Rcpp::NumericVector h1, Rcpp::NumericVector h2,
                             Rcpp::Nullable<unsigned int> numSlices, Rcpp::Nullable<unsigned int> N){

    typedef double NT;
    typedef Cartesian<NT>    Kernel;
    typedef typename Kernel::Point    Point;
    typedef boost::mt19937    RNGType;
    unsigned int num_slices = 100, numpoints = 4000000;

    if (numSlices.isNotNull()) {
        num_slices = Rcpp::as<unsigned int>(numSlices);
    }

    if (N.isNotNull()) {
        numpoints = Rcpp::as<unsigned int>(N);
    }

    Rcpp::NumericMatrix copula(num_slices, num_slices);
    std::vector<std::vector<NT> > StdCopula;
    unsigned int dim = h1.size(), i, j;

    std::vector<NT> hyp1 = Rcpp::as<std::vector<NT> >(h1);
    std::vector<NT> hyp2 = Rcpp::as<std::vector<NT> >(h2);

    StdCopula = twoParHypFam<Point, RNGType >(dim, numpoints, num_slices, hyp1, hyp2);

    for(i=0; i<num_slices; i++) {
        for(j=0; j<num_slices; j++){
            copula(i,j) = StdCopula[i][j];
        }
    }

    return copula;
}
Exemple #2
0
	void initAlias(const distancetype* posWeights,
                 const distancetype* negWeights,
                Rcpp::Nullable<Rcpp::NumericVector> seed) {
		negAlias.initialize(negWeights);
		posAlias.initialize(posWeights);

		if (seed.isNotNull()) {
#ifdef _OPENMP
			storedThreads = omp_get_max_threads();
			omp_set_num_threads(1);
#endif
			vertexidxtype innerSeed = Rcpp::NumericVector(seed)[0];
			innerSeed = negAlias.initRandom(innerSeed);
			posAlias.initRandom(innerSeed);
		} else {
			negAlias.initRandom();
			posAlias.initRandom();
		}
	}
Exemple #3
0
// [[Rcpp::export]]
arma::mat sgd(arma::mat& coords,
              arma::ivec& targets_i, // vary randomly
              arma::ivec& sources_j, // ordered
              arma::ivec& ps, // N+1 length vector of indices to start of each row j in vector is
              arma::vec& weights, // w{ij}
              const double& gamma,
              const double& rho,
              const arma::uword& n_samples,
              const int& M,
              const double& alpha,
              const Rcpp::Nullable<Rcpp::NumericVector> momentum,
              const bool& useDegree,
              const Rcpp::Nullable<Rcpp::NumericVector> seed,
              const Rcpp::Nullable<Rcpp::NumericVector> threads,
              const bool verbose) {
#ifdef _OPENMP
	checkCRAN(threads);
#endif
	const dimidxtype D = coords.n_rows;
	const vertexidxtype N = coords.n_cols;
	const edgeidxtype E = targets_i.n_elem;

	Visualizer* v;
	if (momentum.isNull()) v = new Visualizer(
			sources_j.memptr(), targets_i.memptr(), coords.memptr(),
     	D, N, E,
     	rho, n_samples,
     	M, alpha, gamma);
	else {
		float moment = NumericVector(momentum)[0];
		if (moment < 0) throw Rcpp::exception("Momentum cannot be negative.");
		if (moment > 0.95) throw Rcpp::exception("Bad things happen when momentum is > 0.95.");
		v = new MomentumVisualizer(
			 sources_j.memptr(), targets_i.memptr(), coords.memptr(),
	     D, N, E,
	     rho, n_samples, moment,
	     M, alpha, gamma);
	}

	distancetype* negweights = new distancetype[N];
	std::fill(negweights, negweights + N, 0);
	if (useDegree) {
		std::for_each(targets_i.begin(), targets_i.end(), [&negweights](const sword& e) {negweights[e]++;});
	} else {
		for (vertexidxtype p = 0; p < N; ++p) {
			for (edgeidxtype e = ps[p]; e != ps[p + 1]; ++e) {
				negweights[p] += weights[e];
			}
		}
	}
	std::for_each(negweights, negweights + N, [](distancetype& weight) {weight = pow(weight, 0.75);});
	v -> initAlias(weights.memptr(), negweights, seed);
	delete[] negweights;

	const uword batchSize = BATCHSIZE;
#ifdef _OPENMP
	const unsigned int ts = omp_get_max_threads();
#else
	const unsigned int ts = 2;
#endif
	Progress progress(max((uword) ts, n_samples / BATCHSIZE), verbose);
#ifdef _OPENMP
#pragma omp parallel for
#endif
	for (unsigned int t = 0; t < ts; ++t) {
		v->thread(progress, batchSize);
	}
	delete v;
	return coords;
}
Exemple #4
0
// [[Rcpp::export]]
Rcpp::List europeanOptionEngine(std::string type,
                                double underlying,
                                double strike,
                                double dividendYield,
                                double riskFreeRate,
                                double maturity,
                                double volatility,
                                Rcpp::Nullable<Rcpp::NumericVector> discreteDividends,
                                Rcpp::Nullable<Rcpp::NumericVector> discreteDividendsTimeUntil) {

#ifdef QL_HIGH_RESOLUTION_DATE    
    // in minutes
    boost::posix_time::time_duration length = boost::posix_time::minutes(boost::uint64_t(maturity * 360 * 24 * 60)); 
#else
    int length           = int(maturity*360 + 0.5); // FIXME: this could be better
#endif
    
    QuantLib::Option::Type optionType = getOptionType(type);
    QuantLib::Date today = QuantLib::Date::todaysDate();
    QuantLib::Settings::instance().evaluationDate() = today;

    // new framework as per QuantLib 0.3.5
    QuantLib::DayCounter dc = QuantLib::Actual360();
    boost::shared_ptr<QuantLib::SimpleQuote> spot(new QuantLib::SimpleQuote( underlying ));
    boost::shared_ptr<QuantLib::SimpleQuote> vol(new QuantLib::SimpleQuote( volatility ));
    boost::shared_ptr<QuantLib::BlackVolTermStructure> volTS = flatVol(today, vol, dc);
    boost::shared_ptr<QuantLib::SimpleQuote> qRate(new QuantLib::SimpleQuote( dividendYield ));
    boost::shared_ptr<QuantLib::YieldTermStructure> qTS = flatRate(today, qRate, dc);
    boost::shared_ptr<QuantLib::SimpleQuote> rRate(new QuantLib::SimpleQuote( riskFreeRate ));
    boost::shared_ptr<QuantLib::YieldTermStructure> rTS = flatRate(today, rRate, dc);

    bool withDividends = discreteDividends.isNotNull() && discreteDividendsTimeUntil.isNotNull();
    
#ifdef QL_HIGH_RESOLUTION_DATE
    QuantLib::Date exDate(today.dateTime() + length);
#else 
    QuantLib::Date exDate = today + length;
#endif

    boost::shared_ptr<QuantLib::Exercise> exercise(new QuantLib::EuropeanExercise(exDate));
    boost::shared_ptr<QuantLib::StrikedTypePayoff> payoff(new QuantLib::PlainVanillaPayoff(optionType, strike));
    
    if (withDividends) {
        Rcpp::NumericVector divvalues(discreteDividends), divtimes(discreteDividendsTimeUntil);
        int n = divvalues.size();
        std::vector<QuantLib::Date> discDivDates(n);
        std::vector<double> discDividends(n);
        for (int i = 0; i < n; i++) {
#ifdef QL_HIGH_RESOLUTION_DATE
            boost::posix_time::time_duration discreteDividendLength = boost::posix_time::minutes(boost::uint64_t(divtimes[i] * 360 * 24 * 60));
            discDivDates[i] = QuantLib::Date(today.dateTime() + discreteDividendLength);
#else
            discDivDates[i] = today + int(divtimes[i] * 360 + 0.5); 
#endif    
            discDividends[i] = divvalues[i];
        }
        
        boost::shared_ptr<QuantLib::BlackScholesMertonProcess> 
            stochProcess(new QuantLib::BlackScholesMertonProcess(QuantLib::Handle<QuantLib::Quote>(spot),
                                                                 QuantLib::Handle<QuantLib::YieldTermStructure>(qTS),
                                                                 QuantLib::Handle<QuantLib::YieldTermStructure>(rTS),
                                                                 QuantLib::Handle<QuantLib::BlackVolTermStructure>(volTS)));
        
        boost::shared_ptr<QuantLib::PricingEngine> engine(new QuantLib::AnalyticDividendEuropeanEngine(stochProcess));
        
        QuantLib::DividendVanillaOption option(payoff, exercise, discDivDates, discDividends);
        option.setPricingEngine(engine);
        
        return Rcpp::List::create(Rcpp::Named("value") = option.NPV(),
                                  Rcpp::Named("delta") = option.delta(),
                                  Rcpp::Named("gamma") = option.gamma(),
                                  Rcpp::Named("vega") = option.vega(),
                                  Rcpp::Named("theta") = option.theta(),
                                  Rcpp::Named("rho") = option.rho(),
                                  Rcpp::Named("divRho") = R_NaReal);
    }
    else {
        
        boost::shared_ptr<QuantLib::VanillaOption> option = makeOption(payoff, exercise, spot, qTS, rTS, volTS);
        
        return Rcpp::List::create(Rcpp::Named("value") = option->NPV(),
                                  Rcpp::Named("delta") = option->delta(),
                                  Rcpp::Named("gamma") = option->gamma(),
                                  Rcpp::Named("vega") = option->vega(),
                                  Rcpp::Named("theta") = option->theta(),
                                  Rcpp::Named("rho") = option->rho(),
                                  Rcpp::Named("divRho") = option->dividendRho());
    }
    
}
Exemple #5
0
// [[Rcpp::export]]
Rcpp::List americanOptionEngine(std::string type,
                                double underlying,
                                double strike,
                                double dividendYield,
                                double riskFreeRate,
                                double maturity,
                                double volatility,
                                int timeSteps,
                                int gridPoints,
                                std::string engine,
                                Rcpp::Nullable<Rcpp::NumericVector> discreteDividends,
                                Rcpp::Nullable<Rcpp::NumericVector> discreteDividendsTimeUntil) {

#ifdef QL_HIGH_RESOLUTION_DATE    
    // in minutes
    boost::posix_time::time_duration length = boost::posix_time::minutes(boost::uint64_t(maturity * 360 * 24 * 60)); 
#else
    int length = int(maturity * 360 + 0.5); // FIXME: this could be better
    
#endif
    QuantLib::Option::Type optionType = getOptionType(type);

    // new framework as per QuantLib 0.3.5, updated for 0.3.7
    // updated again for 0.9.0, see eg test-suite/americanoption.cpp
    QuantLib::Date today = QuantLib::Date::todaysDate();
    QuantLib::Settings::instance().evaluationDate() = today;
    QuantLib::DayCounter dc = QuantLib::Actual360();
    boost::shared_ptr<QuantLib::SimpleQuote> spot(new QuantLib::SimpleQuote(underlying));
    boost::shared_ptr<QuantLib::SimpleQuote> qRate(new QuantLib::SimpleQuote(dividendYield));
    boost::shared_ptr<QuantLib::YieldTermStructure> qTS = flatRate(today,qRate,dc);
    boost::shared_ptr<QuantLib::SimpleQuote> rRate(new QuantLib::SimpleQuote(riskFreeRate));
    boost::shared_ptr<QuantLib::YieldTermStructure> rTS = flatRate(today,rRate,dc);
    boost::shared_ptr<QuantLib::SimpleQuote> vol(new QuantLib::SimpleQuote(volatility));
    boost::shared_ptr<QuantLib::BlackVolTermStructure> volTS = flatVol(today, vol, dc);
    
    bool withDividends = discreteDividends.isNotNull() && discreteDividendsTimeUntil.isNotNull();
    
#ifdef QL_HIGH_RESOLUTION_DATE
    QuantLib::Date exDate(today.dateTime() + length);
#else 
    QuantLib::Date exDate = today + length;
#endif

    boost::shared_ptr<QuantLib::StrikedTypePayoff> payoff(new QuantLib::PlainVanillaPayoff(optionType, strike));
    boost::shared_ptr<QuantLib::Exercise> exercise(new QuantLib::AmericanExercise(today, exDate));

    boost::shared_ptr<QuantLib::BlackScholesMertonProcess> 
        stochProcess(new QuantLib::BlackScholesMertonProcess(QuantLib::Handle<QuantLib::Quote>(spot),
                                                             QuantLib::Handle<QuantLib::YieldTermStructure>(qTS),
                                                             QuantLib::Handle<QuantLib::YieldTermStructure>(rTS),
                                                             QuantLib::Handle<QuantLib::BlackVolTermStructure>(volTS)));
    
    if (withDividends) {
        Rcpp::NumericVector divvalues(discreteDividends), divtimes(discreteDividendsTimeUntil);
        int n = divvalues.size();
        std::vector<QuantLib::Date> discDivDates(n);
        std::vector<double> discDividends(n);
        for (int i = 0; i < n; i++) {
#ifdef QL_HIGH_RESOLUTION_DATE
            boost::posix_time::time_duration discreteDividendLength = boost::posix_time::minutes(boost::uint64_t(divtimes[i] * 360 * 24 * 60));
            discDivDates[i] = QuantLib::Date(today.dateTime() + discreteDividendLength);
#else
            discDivDates[i] = today + int(divtimes[i] * 360 + 0.5); 
#endif    
            discDividends[i] = divvalues[i];
        }

        QuantLib::DividendVanillaOption option(payoff, exercise, discDivDates, discDividends);
        if (engine=="BaroneAdesiWhaley") { 
            Rcpp::warning("Discrete dividends, engine switched to CrankNicolson");
            engine = "CrankNicolson";
        }
       
        if (engine=="CrankNicolson") { // FDDividendAmericanEngine only works with CrankNicolson
            // suggestion by Bryan Lewis: use CrankNicolson for greeks
            boost::shared_ptr<QuantLib::PricingEngine> 
            fdcnengine(new QuantLib::FDDividendAmericanEngine<QuantLib::CrankNicolson>(stochProcess, timeSteps, gridPoints));
            option.setPricingEngine(fdcnengine);
            return Rcpp::List::create(Rcpp::Named("value") = option.NPV(),
                                      Rcpp::Named("delta") = option.delta(),
                                      Rcpp::Named("gamma") = option.gamma(),
                                      Rcpp::Named("vega") = R_NaReal,
                                      Rcpp::Named("theta") = R_NaReal,
                                      Rcpp::Named("rho") = R_NaReal,
                                      Rcpp::Named("divRho") = R_NaReal);
        } else {
            throw std::range_error("Unknown engine " + engine);
        }
        
    } else {
        QuantLib::VanillaOption option(payoff, exercise);
        
        if (engine=="BaroneAdesiWhaley") {
            // new from 0.3.7 BaroneAdesiWhaley
            
            boost::shared_ptr<QuantLib::PricingEngine> engine(new QuantLib::BaroneAdesiWhaleyApproximationEngine(stochProcess));
            option.setPricingEngine(engine);
            return Rcpp::List::create(Rcpp::Named("value") = option.NPV(),
                                      Rcpp::Named("delta") = R_NaReal,
                                      Rcpp::Named("gamma") = R_NaReal,
                                      Rcpp::Named("vega") = R_NaReal,
                                      Rcpp::Named("theta") = R_NaReal,
                                      Rcpp::Named("rho") = R_NaReal,
                                      Rcpp::Named("divRho") = R_NaReal);
        } else if (engine=="CrankNicolson") {
            // suggestion by Bryan Lewis: use CrankNicolson for greeks
            boost::shared_ptr<QuantLib::PricingEngine> 
            fdcnengine(new QuantLib::FDAmericanEngine<QuantLib::CrankNicolson>(stochProcess, timeSteps, gridPoints));
            option.setPricingEngine(fdcnengine);
            return Rcpp::List::create(Rcpp::Named("value") = option.NPV(),
                                      Rcpp::Named("delta") = option.delta(),
                                      Rcpp::Named("gamma") = option.gamma(),
                                      Rcpp::Named("vega") = R_NaReal,
                                      Rcpp::Named("theta") = R_NaReal,
                                      Rcpp::Named("rho") = R_NaReal,
                                      Rcpp::Named("divRho") = R_NaReal);
        } else {
            throw std::range_error("Unknown engine " + engine);
        }
    }
    
    
}
//' Sample points from a convex Polytope (H-polytope, V-polytope or a zonotope) or use direct methods for uniform sampling from the unit or the canonical or an arbitrary \eqn{d}-dimensional simplex and the boundary or the interior of a \eqn{d}-dimensional hypersphere
//'
//' Sample N points with uniform or multidimensional spherical gaussian -centered in an internal point- target distribution.
//' The \eqn{d}-dimensional unit simplex is the set of points \eqn{\vec{x}\in \R^d}, s.t.: \eqn{\sum_i x_i\leq 1}, \eqn{x_i\geq 0}. The \eqn{d}-dimensional canonical simplex is the set of points \eqn{\vec{x}\in \R^d}, s.t.: \eqn{\sum_i x_i = 1}, \eqn{x_i\geq 0}.
//'
//' @param P A convex polytope. It is an object from class (a) Hpolytope or (b) Vpolytope or (c) Zonotope.
//' @param N The number of points that the function is going to sample from the convex polytope. The default value is \eqn{100}.
//' @param distribution Optional. A string that declares the target distribution: a) \code{'uniform'} for the uniform distribution or b) \code{'gaussian'} for the multidimensional spherical distribution. The default target distribution is uniform.
//' @param WalkType Optional. A string that declares the random walk method: a) \code{'CDHR'} for Coordinate Directions Hit-and-Run, b) \code{'RDHR'} for Random Directions Hit-and-Run or c) \code{'BW'} for Ball Walk. The default walk is \code{'CDHR'}.
//' @param walk_step Optional. The number of the steps for the random walk. The default value is \eqn{\lfloor 10 + d/10\rfloor}, where \eqn{d} implies the dimension of the polytope.
//' @param exact A boolean parameter. It should be used for the uniform sampling from the boundary or the interior of a hypersphere centered at the origin or from the unit or the canonical or an arbitrary simplex. The arbitrary simplex has to be given as a V-polytope. For the rest well known convex bodies the dimension has to be declared and the type of body as well as the radius of the hypersphere.
//' @param body A string that declares the type of the body for the exact sampling: a) \code{'unit simplex'} for the unit simplex, b) \code{'canonical simplex'} for the canonical simplex, c) \code{'hypersphere'} for the boundary of a hypersphere centered at the origin, d) \code{'ball'} for the interior of a hypersphere centered at the origin.
//' @param Parameters A list for the parameters of the methods:
//' \itemize{
//' \item{\code{variance} }{ The variance of the multidimensional spherical gaussian. The default value is 1.}
//' \item{\code{dimension} }{ An integer that declares the dimension when exact sampling is enabled for a simplex or a hypersphere.}
//' \item{\code{radius} }{ The radius of the \eqn{d}-dimensional hypersphere. The default value is \eqn{1}.}
//' \item{\code{BW_rad} }{ The radius for the ball walk.}
//' }
//' @param InnerPoint A \eqn{d}-dimensional numerical vector that defines a point in the interior of polytope P.
//'
//' @references \cite{R.Y. Rubinstein and B. Melamed,
//' \dQuote{Modern simulation and modeling} \emph{ Wiley Series in Probability and Statistics,} 1998.}
//' @references \cite{A Smith, Noah and W Tromble, Roy,
//' \dQuote{Sampling Uniformly from the Unit Simplex,} \emph{ Center for Language and Speech Processing Johns Hopkins University,} 2004.}
//' @references \cite{Art B. Owen,
//' \dQuote{Monte Carlo theory, methods and examples,} \emph{ Copyright Art Owen,} 2009-2013.}
//'
//' @return A \eqn{d\times N} matrix that contains, column-wise, the sampled points from the convex polytope P.
//' @examples
//' # uniform distribution from the 3d unit cube in V-representation using ball walk
//' P = GenCube(3, 'V')
//' points = sample_points(P, WalkType = "BW", walk_step = 5)
//'
//' # gaussian distribution from the 2d unit simplex in H-representation with variance = 2
//' A = matrix(c(-1,0,0,-1,1,1), ncol=2, nrow=3, byrow=TRUE)
//' b = c(0,0,1)
//' P = Hpolytope$new(A,b)
//' points = sample_points(P, distribution = "gaussian", Parameters = list("variance" = 2))
//'
//' # uniform points from the boundary of a 10-dimensional hypersphere
//' points = sample_points(exact = TRUE, body = "hypersphere", Parameters = list("dimension" = 10))
//'
//' # 10000 uniform points from a 2-d arbitrary simplex
//' V = matrix(c(2,3,-1,7,0,0),ncol = 2, nrow = 3, byrow = TRUE)
//' P = Vpolytope$new(V)
//' points = sample_points(P, N = 10000, exact = TRUE)
//' @export
// [[Rcpp::export]]
Rcpp::NumericMatrix sample_points(Rcpp::Nullable<Rcpp::Reference> P = R_NilValue,
                                  Rcpp::Nullable<unsigned int> N = R_NilValue,
                                  Rcpp::Nullable<std::string> distribution = R_NilValue,
                                  Rcpp::Nullable<std::string> WalkType = R_NilValue,
                                  Rcpp::Nullable<unsigned int> walk_step = R_NilValue,
                                  Rcpp::Nullable<bool> exact = R_NilValue,
                                  Rcpp::Nullable<std::string> body = R_NilValue,
                                  Rcpp::Nullable<Rcpp::List> Parameters = R_NilValue,
                                  Rcpp::Nullable<Rcpp::NumericVector> InnerPoint = R_NilValue){

    typedef double NT;
    typedef Cartesian<NT>    Kernel;
    typedef typename Kernel::Point    Point;
    typedef boost::mt19937 RNGType;
    typedef HPolytope <Point> Hpolytope;
    typedef VPolytope <Point, RNGType> Vpolytope;
    typedef Zonotope <Point> zonotope;
    typedef IntersectionOfVpoly<Vpolytope> InterVP;
    typedef Eigen::Matrix<NT,Eigen::Dynamic,1> VT;
    typedef Eigen::Matrix<NT,Eigen::Dynamic,Eigen::Dynamic> MT;

    Hpolytope HP;
    Vpolytope VP;
    zonotope ZP;
    InterVP VPcVP;

    int type, dim, numpoints;
    NT radius = 1.0, delta = -1.0;
    bool set_mean_point = false, cdhr = true, rdhr = false, ball_walk = false, gaussian = false;
    std::list<Point> randPoints;
    std::pair<Point, NT> InnerBall;

    if (!N.isNotNull()) {
        numpoints = 100;
    } else {
        numpoints = Rcpp::as<unsigned int>(N);
    }

    if (exact.isNotNull()) {
        if (P.isNotNull()) {
            type = Rcpp::as<Rcpp::Reference>(P).field("type");
            dim = Rcpp::as<Rcpp::Reference>(P).field("dimension");
            if (Rcpp::as<bool>(exact) && type==2) {
                if (Rcpp::as<MT>(Rcpp::as<Rcpp::Reference>(P).field("V")).rows() == Rcpp::as<MT>(Rcpp::as<Rcpp::Reference>(P).field("V")).cols()+1) {
                    Vpolytope VP;
                    VP.init(dim,
                            Rcpp::as<MT>(Rcpp::as<Rcpp::Reference>(P).field("V")),
                            VT::Ones(Rcpp::as<MT>(Rcpp::as<Rcpp::Reference>(P).field("V")).rows()));
                    Sam_arb_simplex(VP, numpoints, randPoints);
                } else {
                    throw Rcpp::exception("Not a simplex!");
                }
            } else if (Rcpp::as<bool>(exact) && type!=2) {
                throw Rcpp::exception("Not a simplex in V-representation!");
            }
        } else {
            if (!body.isNotNull()) {

                throw Rcpp::exception("Wrong input!");

            } else if (!Rcpp::as<Rcpp::List>(Parameters).containsElementNamed("dimension")) {

                throw Rcpp::exception("Wrong input!");

            }
            dim = Rcpp::as<NT>(Rcpp::as<Rcpp::List>(Parameters)["dimension"]);
            if (Rcpp::as<Rcpp::List>(Parameters).containsElementNamed("radius")) {

                radius = Rcpp::as<NT>(Rcpp::as<Rcpp::List>(Parameters)["radius"]);

            }
            if (Rcpp::as<std::string>(body).compare(std::string("hypersphere"))==0) {

                for (unsigned int k = 0; k < numpoints; ++k) {
                    randPoints.push_back(get_point_on_Dsphere<RNGType , Point >(dim, radius));
                }

            } else if (Rcpp::as<std::string>(body).compare(std::string("ball"))==0) {

                for (unsigned int k = 0; k < numpoints; ++k) {
                    randPoints.push_back(get_point_in_Dsphere<RNGType , Point >(dim, radius));
                }

            } else if (Rcpp::as<std::string>(body).compare(std::string("unit simplex"))==0) {

                Sam_Unit<NT, RNGType >(dim, numpoints, randPoints);

            } else if (Rcpp::as<std::string>(body).compare(std::string("canonical simplex"))==0) {

                Sam_Canon_Unit<NT, RNGType >(dim, numpoints, randPoints);

            } else {

                throw Rcpp::exception("Wrong input!");

            }
        }
    } else if (P.isNotNull()) {

        type = Rcpp::as<Rcpp::Reference>(P).field("type");
        dim = Rcpp::as<Rcpp::Reference>(P).field("dimension");
        unsigned int walkL = 10+dim/10;
        Point MeanPoint;
        if (InnerPoint.isNotNull()) {
            if (Rcpp::as<Rcpp::NumericVector>(InnerPoint).size()!=dim) {
                Rcpp::warning("Internal Point has to lie in the same dimension as the polytope P");
            } else {
                set_mean_point = true;
                MeanPoint = Point( dim , Rcpp::as<std::vector<NT> >(InnerPoint).begin(),
                        Rcpp::as<std::vector<NT> >(InnerPoint).end() );
            }
        }
        if(walk_step.isNotNull()) walkL = Rcpp::as<unsigned int>(walk_step);

        NT a = 0.5;

        if (Rcpp::as<Rcpp::List>(Parameters).containsElementNamed("variance"))
            a = 1.0 / (2.0 * Rcpp::as<NT>(Rcpp::as<Rcpp::List>(Parameters)["variance"]));

        if(!WalkType.isNotNull() || Rcpp::as<std::string>(WalkType).compare(std::string("CDHR"))==0){
            cdhr = true;
            rdhr = false;
            ball_walk = false;
        } else if (Rcpp::as<std::string>(WalkType).compare(std::string("RDHR"))==0) {
            cdhr = false;
            rdhr = true;
            ball_walk = false;
        } else if (Rcpp::as<std::string>(WalkType).compare(std::string("BW"))==0) {
            if (Rcpp::as<Rcpp::List>(Parameters).containsElementNamed("BW_rad")) {
                delta = Rcpp::as<NT>(Rcpp::as<Rcpp::List>(Parameters)["BW_rad"]);
            }
            cdhr = false;
            rdhr = false;
            ball_walk = true;
        } else {
            throw Rcpp::exception("Unknown walk type!");
        }

        if (distribution.isNotNull()) {
            if (Rcpp::as<std::string>(distribution).compare(std::string("gaussian"))==0) {
                gaussian = true;
            } else if(Rcpp::as<std::string>(distribution).compare(std::string("uniform"))!=0) {
                throw Rcpp::exception("Wrong distribution!");
            }
        }
        bool rand_only=false,
                NN=false,
                birk=false,
                verbose=false;

        unsigned seed = std::chrono::system_clock::now().time_since_epoch().count();
        // the random engine with this seed
        RNGType rng(seed);
        boost::random::uniform_real_distribution<>(urdist);
        boost::random::uniform_real_distribution<> urdist1(-1,1);

        switch(type) {
            case 1: {
                // Hpolytope
                HP.init(dim, Rcpp::as<MT>(Rcpp::as<Rcpp::Reference>(P).field("A")),
                        Rcpp::as<VT>(Rcpp::as<Rcpp::Reference>(P).field("b")));

                if (!set_mean_point || ball_walk) {
                    InnerBall = HP.ComputeInnerBall();
                    if (!set_mean_point) MeanPoint = InnerBall.first;
                }
                break;
            }
            case 2: {
                // Vpolytope
                VP.init(dim, Rcpp::as<MT>(Rcpp::as<Rcpp::Reference>(P).field("V")),
                        VT::Ones(Rcpp::as<MT>(Rcpp::as<Rcpp::Reference>(P).field("V")).rows()));

                if (!set_mean_point || ball_walk) {
                    InnerBall = VP.ComputeInnerBall();
                    if (!set_mean_point) MeanPoint = InnerBall.first;
                }
                break;
            }
            case 3: {
                // Zonotope
                ZP.init(dim, Rcpp::as<MT>(Rcpp::as<Rcpp::Reference>(P).field("G")),
                        VT::Ones(Rcpp::as<MT>(Rcpp::as<Rcpp::Reference>(P).field("G")).rows()));

                if (!set_mean_point || ball_walk) {
                    InnerBall = ZP.ComputeInnerBall();
                    if (!set_mean_point) MeanPoint = InnerBall.first;
                }
                break;
            }
            case 4: {
                // Intersection of two V-polytopes
                Vpolytope VP1;
                Vpolytope VP2;
                VP1.init(dim, Rcpp::as<MT>(Rcpp::as<Rcpp::Reference>(P).field("V1")),
                         VT::Ones(Rcpp::as<MT>(Rcpp::as<Rcpp::Reference>(P).field("V1")).rows()));
                VP2.init(dim, Rcpp::as<MT>(Rcpp::as<Rcpp::Reference>(P).field("V2")),
                         VT::Ones(Rcpp::as<MT>(Rcpp::as<Rcpp::Reference>(P).field("V2")).rows()));
                VPcVP.init(VP1, VP2);

                if (!set_mean_point || ball_walk) {
                    InnerBall = VP.ComputeInnerBall();
                    if (!set_mean_point) MeanPoint = InnerBall.first;
                }
                break;
            }
        }

        if (ball_walk) {
            if (gaussian) {
                delta = 4.0 * InnerBall.second / std::sqrt(std::max(NT(1.0), a) * NT(dim));
            } else {
                delta = 4.0 * InnerBall.second / std::sqrt(NT(dim));
            }
        }
        vars<NT, RNGType> var1(1,dim,walkL,1,0.0,0.0,0,0.0,0,InnerBall.second,rng,urdist,urdist1,
                               delta,verbose,rand_only,false,NN,birk,ball_walk,cdhr,rdhr);
        vars_g<NT, RNGType> var2(dim, walkL, 0, 0, 1, 0, InnerBall.second, rng, 0, 0, 0, delta, false, verbose,
                                 rand_only, false, NN, birk, ball_walk, cdhr, rdhr);

        switch (type) {
            case 1: {
                sampling_only<Point>(randPoints, HP, walkL, numpoints, gaussian,
                                     a, MeanPoint, var1, var2);
                break;
            }
            case 2: {
                sampling_only<Point>(randPoints, VP, walkL, numpoints, gaussian,
                                     a, MeanPoint, var1, var2);
                break;
            }
            case 3: {
                sampling_only<Point>(randPoints, ZP, walkL, numpoints, gaussian,
                                     a, MeanPoint, var1, var2);
                break;
            }
            case 4: {
                sampling_only<Point>(randPoints, VPcVP, walkL, numpoints, gaussian,
                                     a, MeanPoint, var1, var2);
                break;
            }
        }

    } else {

        throw Rcpp::exception("Wrong input!");

    }

    Rcpp::NumericMatrix PointSet(dim,numpoints);
    typename std::list<Point>::iterator rpit=randPoints.begin();
    typename std::vector<NT>::iterator qit;
    unsigned int j = 0, i;
    for ( ; rpit!=randPoints.end(); rpit++, j++) {
        qit = (*rpit).iter_begin(); i=0;
        for ( ; qit!=(*rpit).iter_end(); qit++, i++){
            PointSet(i,j)=*qit;
        }
    }
    return PointSet;

}