Latin Hypercube Samples (lhs)  1.0
R, C++, and Rcpp code to generate Latin hypercube samples
LHSCommonDefines.h
Go to the documentation of this file.
1 
21 #ifndef LHSCOMMONDEFINES_H
22 #define LHSCOMMONDEFINES_H
23 
24 #include <cstdlib>
25 #include <cmath>
26 #include <exception>
27 #include <vector>
28 #include <algorithm>
29 #include <functional>
30 #include <numeric>
31 #include <cfloat>
32 #include <climits>
33 #include <cstdio>
34 #include <iostream>
35 #include "matrix.h"
36 #include "order.h"
37 #include "CRandom.h"
38 
39 #ifdef RCOMPILE
40 #include <Rcpp.h>
41 #define PRINT_MACRO Rcpp::Rcout
42 #define ERROR_MACRO Rcpp::Rcerr
43 #else // RCOMPILE
44 
45 #define PRINT_MACRO std::cout
46 
47 #define ERROR_MACRO std::cerr
48 #endif // RCOMPILE
49 
51 #define PRINT_RESULT 0
52 
54 #define START_RNG Rcpp::RNGScope * tempRNG = new Rcpp::RNGScope(); // instantiate a pointer so that the destructor is not implicitly called
55 
56 #define END_RNG delete tempRNG; // explicitly release the RNG state to avoid memory corruption
57 
61 namespace lhslib
62 {
71  void improvedLHS(int n, int k, int dup, bclib::matrix<int> & result,
72  bclib::CRandom<double> & oRandom);
81  void maximinLHS(int n, int k, int dup, bclib::matrix<int> & result,
82  bclib::CRandom<double> & oRandom);
94  void optimumLHS(int n, int k, int maxSweeps, double eps,
95  bclib::matrix<int> & outlhs, int optimalityRecordLength,
96  bclib::CRandom<double> & oRandom, bool bVerbose);
107  void optSeededLHS(int n, int k, int maxSweeps, double eps,
108  bclib::matrix<double> & oldHypercube, int optimalityRecordLength, bool bVerbose);
109 
114  typedef bclib::matrix<int>::size_type msize_type;
119  typedef std::vector<int>::size_type vsize_type;
120 
129  void randomLHS(int n, int k, bool bPreserveDraw, bclib::matrix<double> & result, bclib::CRandom<double> & oRandom);
130 
138  void randomLHS(int n, int k, bclib::matrix<int> & result, bclib::CRandom<double> & oRandom);
139 
152  void geneticLHS(int n, int k, int pop, int gen, double pMut, const std::string & criterium,
153  bool bVerbose, bclib::matrix<double> & result, bclib::CRandom<double> & oRandom);
154 }
155 
156 #endif /* LHSCOMMONDEFINES_H */
CRandom.h
lhslib::optimumLHS
void optimumLHS(int n, int k, int maxSweeps, double eps, bclib::matrix< int > &outlhs, int optimalityRecordLength, bclib::CRandom< double > &oRandom, bool bVerbose)
Definition: optimumLHS.cpp:49
lhslib
Definition: geneticLHS.cpp:25
lhslib::improvedLHS
void improvedLHS(int n, int k, int dup, bclib::matrix< int > &result, bclib::CRandom< double > &oRandom)
Definition: improvedLHS.cpp:45
lhslib::maximinLHS
void maximinLHS(int n, int k, int dup, bclib::matrix< int > &result, bclib::CRandom< double > &oRandom)
Definition: maximinLHS.cpp:40
lhslib::msize_type
bclib::matrix< int >::size_type msize_type
Definition: LHSCommonDefines.h:114
lhslib::vsize_type
std::vector< int >::size_type vsize_type
Definition: LHSCommonDefines.h:119
order.h
lhslib::randomLHS
void randomLHS(int n, int k, bool bPreserveDraw, bclib::matrix< double > &result, bclib::CRandom< double > &oRandom)
Definition: randomLHS.cpp:44
lhslib::geneticLHS
void geneticLHS(int n, int k, int pop, int gen, double pMut, const std::string &criterium, bool bVerbose, bclib::matrix< double > &result, bclib::CRandom< double > &oRandom)
Definition: geneticLHS.cpp:27
bclib::CRandom< double >
lhslib::optSeededLHS
void optSeededLHS(int n, int k, int maxSweeps, double eps, bclib::matrix< double > &oldHypercube, int optimalityRecordLength, bool bVerbose)
Definition: optSeededLHS.cpp:53