Descriptors¶
-
class
Descriptor¶ Subclassed by B2, B2_Simple, B3, FourBody, ThreeBody, ThreeBodyWide, TwoBody
Public Functions
-
Descriptor()¶
-
virtual DescriptorValues
compute_struc(Structure &structure) = 0¶
-
virtual
~Descriptor()¶
-
void
write_to_file(std::ofstream &coeff_file, int coeff_size)¶
-
virtual nlohmann::json
return_json() = 0¶
Public Members
-
std::string
descriptor_name¶
-
-
class
DescriptorValues¶ Public Functions
-
DescriptorValues()¶
Public Members
-
int
n_descriptors¶
-
int
n_types¶
-
int
n_atoms¶
-
double
volume¶
-
std::vector<Eigen::MatrixXd>
descriptors¶
-
std::vector<Eigen::MatrixXd>
descriptor_force_dervs¶
-
std::vector<Eigen::MatrixXd>
neighbor_coordinates¶
-
std::vector<Eigen::VectorXd>
descriptor_norms¶
-
std::vector<Eigen::VectorXd>
descriptor_force_dots¶
-
std::vector<Eigen::VectorXd>
cutoff_values¶
-
std::vector<Eigen::VectorXd>
cutoff_dervs¶
-
std::vector<Eigen::VectorXi>
neighbor_counts¶
-
std::vector<Eigen::VectorXi>
cumulative_neighbor_counts¶
-
std::vector<Eigen::VectorXi>
atom_indices¶
-
std::vector<Eigen::VectorXi>
neighbor_indices¶
-
int
n_clusters= 0¶
-
std::vector<int>
n_clusters_by_type¶
-
std::vector<int>
cumulative_type_count¶
-
std::vector<int>
n_neighbors_by_type¶
-
-
class
ClusterDescriptor¶ Public Functions
-
ClusterDescriptor()¶
-
ClusterDescriptor(const DescriptorValues &structure)¶
-
ClusterDescriptor(const DescriptorValues &structure, const std::vector<std::vector<int>> &clusters)¶
-
ClusterDescriptor(const DescriptorValues &structure, const std::vector<int> &clusters)¶
-
void
initialize_cluster(int n_types, int n_descriptors)¶
-
void
add_clusters_by_type(const DescriptorValues &structure, const std::vector<std::vector<int>> &clusters)¶
-
void
add_clusters(const DescriptorValues &structure, const std::vector<int> &clusters)¶
-
void
add_all_clusters(const DescriptorValues &structure)¶
-
-
class
B2: public Descriptor¶ Subclassed by B2_Norm
Public Functions
-
B2()¶
-
B2(const std::string &radial_basis, const std::string &cutoff_function, const std::vector<double> &radial_hyps, const std::vector<double> &cutoff_hyps, const std::vector<int> &descriptor_settings)¶
-
B2(const std::string &radial_basis, const std::string &cutoff_function, const std::vector<double> &radial_hyps, const std::vector<double> &cutoff_hyps, const std::vector<int> &descriptor_settings, const Eigen::MatrixXd &cutoffs)¶ Construct the B2 descriptor with distinct cutoffs for each pair of species.
-
DescriptorValues
compute_struc(Structure &structure)¶
-
void
write_to_file(std::ofstream &coeff_file, int coeff_size)¶
-
nlohmann::json
return_json()¶
Public Members
-
std::function<void(std::vector<double>&, std::vector<double>&, double, int, std::vector<double>)>
radial_pointer¶
-
std::function<void(std::vector<double>&, double, double, std::vector<double>)>
cutoff_pointer¶
-
std::string
radial_basis¶
-
std::string
cutoff_function¶
-
std::vector<double>
radial_hyps¶
-
std::vector<double>
cutoff_hyps¶
-
std::vector<int>
descriptor_settings¶
-
std::string
descriptor_name= "B2"¶
-
Eigen::MatrixXd
cutoffs¶ Matrix of cutoff values, with element (i, j) corresponding to the cutoff assigned to the species pair (i, j), where i is the central species and j is the environment species.
-