#pragma once#include<Eigen/Eigen>#include<iostream>#include<fstream>#include<vector>#include<string>#include<cstdio>#include"itkParticlePositionReader.h"#include"itkParticlePositionWriter.h"#include"Shapeworks.h"#include"ParticleSystem.h"namespaceshapeworks{classProject;classParticleShapeStatistics{public:constexprstaticintVDimension=3;ParticleShapeStatistics(){};ParticleShapeStatistics(std::shared_ptr<Project>project);~ParticleShapeStatistics(){};intDoPCA(std::vector<std::vector<Point>>global_pts,intdomainsPerShape=1);intDoPCA(ParticleSystemparticleSystem,intdomainsPerShape=1);itkStaticConstMacro(Dimension,unsignedint,VDimension);intImportPoints(std::vector<Eigen::VectorXd>points,std::vector<int>group_ids);intReadPointFiles(conststd::string&s);intReloadPointFiles();intWriteCSVFile(conststd::string&s);intWriteCSVFile2(conststd::string&s);intComputeModes();intPrincipalComponentProjections();intFisherLinearDiscriminant(unsignedintnumModes);intSampleSize()const{returnm_numSamples;}intGroup1SampleSize()const{returnm_numSamples1;}intGroup2SampleSize()const{returnm_numSamples2;}intNumberOfDimensions()const{returnm_numDimensions;}intGroupID(unsignedinti)const{returnm_groupIDs[i];}conststd::vector<int>&GroupID()const{returnm_groupIDs;}constEigen::MatrixXd&Eigenvectors()const{returnm_eigenvectors;}conststd::vector<double>&Eigenvalues()const{returnm_eigenvalues;}constEigen::VectorXd&Mean()const{returnm_mean;}constEigen::VectorXd&Group1Mean()const{returnm_mean1;}constEigen::VectorXd&Group2Mean()const{returnm_mean2;}constEigen::VectorXd&NormalizedGroupDifference()const{returnm_groupdiffnorm;}constEigen::VectorXd&GroupDifference()const{returnm_groupdiff;}intComputeMedianShape(constintID);doubleL1Norm(unsignedinta,unsignedintb);Eigen::MatrixXd&PCALoadings(){returnm_principals;}constEigen::VectorXd&FishersLDA()const{returnm_fishersLD;}constEigen::MatrixXd&ShapeMatrix()const{returnm_shapes;}constEigen::MatrixXd&RecenteredShape()const{returnm_pointsMinusMean;}conststd::vector<double>&PercentVarByMode()const{returnm_percentVarByMode;}intSimpleLinearRegression(conststd::vector<double>&y,conststd::vector<double>&x,double&a,double&b)const;Eigen::VectorXdget_compactness(std::function<void(float)>progress_callback);Eigen::VectorXdget_specificity(std::function<void(float)>progress_callback);Eigen::VectorXdget_generalization(std::function<void(float)>progress_callback);Eigen::MatrixXdget_group1_matrix();Eigen::MatrixXdget_group2_matrix();private:voidcompute_good_bad_points();unsignedintm_numSamples1;unsignedintm_numSamples2;unsignedintm_numSamples;unsignedintm_domainsPerShape;unsignedintm_numDimensions;std::vector<int>m_groupIDs;Eigen::MatrixXdm_eigenvectors;std::vector<double>m_eigenvalues;Eigen::VectorXdm_mean;Eigen::VectorXdm_mean1;Eigen::VectorXdm_mean2;Eigen::MatrixXdm_pointsMinusMean;Eigen::MatrixXdm_shapes;Eigen::MatrixXdm_projectedPMM1;Eigen::MatrixXdm_projectedPMM2;Eigen::VectorXdm_projectedMean1;Eigen::VectorXdm_projectedMean2;std::vector<double>m_fishersProjection;std::vector<double>m_percentVarByMode;Eigen::VectorXdm_fishersLD;Eigen::MatrixXdm_principals;Eigen::VectorXdm_groupdiff;Eigen::VectorXdm_groupdiffnorm;// used to keep the points' files that needs to be reloaded when new updates come in.std::vector<std::string>m_pointsfiles;Eigen::MatrixXdm_Matrix;Eigen::MatrixXdm_group_1_matrix;Eigen::MatrixXdm_group_2_matrix;// 0 = bad, 1 = goodstd::vector<bool>m_goodPoints;std::vector<Eigen::VectorXd>points_;};}// shapeworks