1 #ifndef DISTRIBUTEDCOLLECTION_H_
2 #define DISTRIBUTEDCOLLECTION_H_
12 #include <parcer/Communicator.h>
17 namespace SamplingAlgorithms {
29 virtual void Add(std::shared_ptr<SamplingState> newSamp)
override;
35 std::shared_ptr<SamplingState>
LocalAt(
unsigned i);
41 const std::shared_ptr<SamplingState>
LocalAt(
unsigned i)
const;
47 std::shared_ptr<SamplingState>
GlobalAt(
unsigned i);
53 const std::shared_ptr<SamplingState>
GlobalAt(
unsigned i)
const;
59 virtual std::shared_ptr<SamplingState>
at(
unsigned i)
override;
65 virtual const std::shared_ptr<SamplingState>
at(
unsigned i)
const override;
84 virtual unsigned int size()
const override;
109 virtual Eigen::VectorXd
CentralMoment(
unsigned order,
int blockDim=-1)
const override;
116 Eigen::VectorXd
LocalMean(
int blockDim=-1)
const;
123 Eigen::VectorXd
GlobalMean(
int blockDim=-1)
const;
131 virtual Eigen::VectorXd
Mean(
int blockDim=-1)
const override;
153 virtual Eigen::VectorXd
Variance(
int blockDim=-1)
const override;
175 virtual Eigen::MatrixXd
Covariance(
int blockDim=-1)
const override;
182 Eigen::VectorXd
LocalESS(
int blockDim=-1)
const;
189 Eigen::VectorXd
GlobalESS(
int blockDim=-1)
const;
197 virtual Eigen::VectorXd
ESS(
int blockDim=-1)
const override;
218 virtual Eigen::MatrixXd
AsMatrix(
int blockDim=-1)
const override;
239 virtual Eigen::VectorXd
Weights()
const override;
246 virtual void WriteToFile(std::string
const& filename, std::string
const& dataset =
"/")
const override;
248 Eigen::VectorXd
GlobalExpectedValue(std::shared_ptr<muq::Modeling::ModPiece>
const& f, std::vector<std::string>
const& metains = std::vector<std::string>())
const;
250 Eigen::VectorXd
LocalExpectedValue(std::shared_ptr<muq::Modeling::ModPiece>
const& f, std::vector<std::string>
const& metains = std::vector<std::string>())
const;
252 virtual Eigen::VectorXd
ExpectedValue(std::shared_ptr<muq::Modeling::ModPiece>
const& f, std::vector<std::string>
const& metains = std::vector<std::string>())
const override;
256 template <
class scalar,
int rows,
int cols,
int options,
int maxRows,
int maxCols>
257 Eigen::Matrix<scalar, rows, cols, options, maxRows, maxCols>
GlobalEigenMean(Eigen::Matrix<scalar, rows, cols, options, maxRows, maxCols>
const& local)
const {
259 typedef Eigen::Matrix<scalar, rows, cols, options, maxRows, maxCols> MatType;
260 MatType global = MatType::Zero(local.rows(), local.cols());
262 for(
unsigned int i=0; i<
comm->GetSize(); ++i ) {
263 MatType l(local.rows(), local.cols());
264 if(
comm->GetRank()==i ) { l = local; }
267 global += l/(double)
comm->GetSize();
277 std::shared_ptr<parcer::Communicator>
comm;
std::shared_ptr< parcer::Communicator > comm
The communicator for this collection.
Eigen::VectorXd LocalExpectedValue(std::shared_ptr< muq::Modeling::ModPiece > const &f, std::vector< std::string > const &metains=std::vector< std::string >()) const
std::shared_ptr< SampleCollection > collection
The local sample collection (stored on this processor)
virtual Eigen::VectorXd ESS(int blockDim=-1) const override
Compute the ESS using all of the samples.
Eigen::Matrix< scalar, rows, cols, options, maxRows, maxCols > GlobalEigenMean(Eigen::Matrix< scalar, rows, cols, options, maxRows, maxCols > const &local) const
virtual Eigen::VectorXd Weights() const override
Return all of the weights as a vector.
Eigen::MatrixXd GlobalCovariance(int blockDim=-1) const
Compute the covariance using all of the samples.
virtual void Add(std::shared_ptr< SamplingState > newSamp) override
Add a sample that is stored on this processor.
Eigen::VectorXd GlobalMean(int blockDim=-1) const
Compute the mean using all of the samples.
Eigen::VectorXd LocalESS(int blockDim=-1) const
Compute the ESS using only local samples.
virtual void WriteToFile(std::string const &filename, std::string const &dataset="/") const override
virtual std::shared_ptr< SamplingState > at(unsigned i) override
Get the local state at the index.
virtual Eigen::VectorXd ExpectedValue(std::shared_ptr< muq::Modeling::ModPiece > const &f, std::vector< std::string > const &metains=std::vector< std::string >()) const override
Eigen::MatrixXd AsGlobalMatrix(int blockDim=-1) const
Return all of the samples as a matrix.
DistributedCollection(std::shared_ptr< SampleCollection > collection, std::shared_ptr< parcer::Communicator > comm)
Eigen::VectorXd LocalMean(int blockDim=-1) const
Compute the mean using only local samples.
virtual Eigen::VectorXd CentralMoment(unsigned order, int blockDim=-1) const override
Computes the componentwise central moments (e.g., variance, skewness, kurtosis, etc....
virtual ~DistributedCollection()=default
std::shared_ptr< SamplingState > GlobalAt(unsigned i)
Get the global state at the index.
virtual Eigen::VectorXd Mean(int blockDim=-1) const override
Compute the mean using all of the samples.
Eigen::VectorXd GlobalCentralMoment(unsigned order, int blockDim=-1) const
Computes the componentwise central moments (e.g., variance, skewness, kurtosis, etc....
unsigned int GlobalSize() const
The total number of samples.
virtual unsigned int size() const override
The total number of samples.
Eigen::VectorXd LocalCentralMoment(unsigned order, int blockDim=-1) const
Computes the componentwise central moments (e.g., variance, skewness, kurtosis, etc....
std::shared_ptr< SamplingState > LocalAt(unsigned i)
Get the local state at the index.
Eigen::VectorXd GlobalWeights() const
Return all of the weights as a vector.
virtual Eigen::MatrixXd AsMatrix(int blockDim=-1) const override
Return all of the samples as a matrix.
unsigned int LocalSize() const
The number of samples stored locally.
virtual Eigen::MatrixXd Covariance(int blockDim=-1) const override
Compute the covariance using all of the samples.
virtual Eigen::VectorXd Variance(int blockDim=-1) const override
Compute the variance using all of the samples.
Eigen::VectorXd LocalWeights() const
Return all of the weights on this processor as a vector.
Eigen::VectorXd LocalVariance(int blockDim=-1) const
Compute the variance using only local samples.
Eigen::VectorXd GlobalVariance(int blockDim=-1) const
Compute the variance using all of the samples.
Eigen::VectorXd GlobalESS(int blockDim=-1) const
Compute the ESS using all of the samples.
Eigen::MatrixXd AsLocalMatrix(int blockDim=-1) const
Return all of the samples on the local processor as a matrix.
Eigen::VectorXd GlobalExpectedValue(std::shared_ptr< muq::Modeling::ModPiece > const &f, std::vector< std::string > const &metains=std::vector< std::string >()) const
Eigen::MatrixXd LocalCovariance(int blockDim=-1) const
Compute the covariance using only local samples.
A class to hold and analyze a collection of SamplingState objects.