SSAGES  0.9.3
Software Suite for Advanced General Ensemble Simulations
Public Member Functions | Static Public Member Functions | Private Attributes | List of all members
SSAGES::ANNCV Class Reference

ANN (artifical neural network) collective variables. More...

#include <ANNCV.h>

Inheritance diagram for SSAGES::ANNCV:
Inheritance graph
[legend]

Public Member Functions

 ANNCV (Label atomids, double scaling_factor, std::vector< unsigned int > num_nodes, std::string coeff_file, std::vector< std::string > activations, int out_index)
 Constructor. More...
 
void Initialize (const Snapshot &snapshot) override
 Initialize necessary variables. More...
 
std::vector< Vectorforward_prop (Vector &input_vec)
 
std::vector< Vectorback_prop (std::vector< Vector > &output_of_layers)
 
void Evaluate (const Snapshot &snapshot) override
 Evaluate the CV. More...
 
- Public Member Functions inherited from SSAGES::CollectiveVariable
 CollectiveVariable ()
 Constructor.
 
virtual ~CollectiveVariable ()
 Destructor.
 
virtual void Initialize (const class Snapshot &)
 Initialize CV. More...
 
virtual void Evaluate (const class Snapshot &)=0
 Evaluate CV. More...
 
double GetValue () const
 Get current value of the CV. More...
 
virtual double GetMinimumImage (double) const
 Returns the minimum image of a CV based on the input location. More...
 
virtual double GetPeriodicValue (double location) const
 Apply periodic boundaries to a given value. More...
 
const std::vector< Vector3 > & GetGradient () const
 Get current gradient of the CV. More...
 
const Matrix3GetBoxGradient () const
 Get gradient contribution to box. More...
 
const std::array< double, 2 > & GetBoundaries ()
 Get CV boundaries. More...
 
virtual double GetDifference (double location) const
 

Static Public Member Functions

static ANNCVBuild (const Json::Value &json, const std::string &path)
 Set up collective variable. More...
 
- Static Public Member Functions inherited from SSAGES::CollectiveVariable
static CollectiveVariableBuildCV (const Json::Value &json, const std::string &path)
 Set up collective variable. More...
 

Private Attributes

Label atomids_
 
double scaling_factor_
 
std::vector< unsigned int > num_nodes_
 
std::vector< Eigen::MatrixXd > weight_coeff_
 
std::vector< Vectorbias_
 
std::vector< std::string > activations_
 
int out_index_
 

Additional Inherited Members

- Protected Attributes inherited from SSAGES::CollectiveVariable
std::vector< Vector3grad_
 Gradient vector dCv/dxi.
 
Matrix3 boxgrad_
 Gradient w.r.t box vectors dCv/dHij.
 
double val_
 Current value of CV.
 
std::array< double, 2 > bounds_
 Bounds on CV.
 

Detailed Description

ANN (artifical neural network) collective variables.

This CV takes scaled (specified by "scaling_factor") Cartesian coordinates of a group of atoms (specified by "atomids") as inputs to a neural network (its number of nodes, connection weights, and activation functions are specified by "num_nodes", "coeff_file", "activations", respectively), computes one component (specified by "out_index") of the final neural network outputs as the CV value.

Definition at line 42 of file ANNCV.h.

Constructor & Destructor Documentation

◆ ANNCV()

SSAGES::ANNCV::ANNCV ( Label  atomids,
double  scaling_factor,
std::vector< unsigned int >  num_nodes,
std::string  coeff_file,
std::vector< std::string >  activations,
int  out_index 
)
inline

Constructor.

Parameters
atomidsatom IDs used as inputs to ANN
scaling_factorscaling factor for input coordinates (input coordinates will be divided by this factor)
num_nodesnumbers of nodes for each layer of ANN
coeff_filefile storing weights and bias of ANN
activationsactivations functions (as strings) of ANN
out_indexindex of output component of ANN
Todo:
Bounds needs to be an input and periodic boundary conditions

Definition at line 66 of file ANNCV.h.

73  :
74  atomids_(atomids), scaling_factor_(scaling_factor), num_nodes_(num_nodes),
75  activations_(activations), out_index_(out_index)
76  {
77  // read coefficients from file
78  std::ifstream my_f(coeff_file);
79  std::string temp_vec;
80  int layer_index = 0;
81  if (num_nodes_[0] != atomids_.size() * 3) {
82  throw BuildException({
83  "WARNING: input dim should be " + std::to_string(atomids_.size() * 3) + " found: "
84  + std::to_string(num_nodes_[0])
85  });
86  }
87  while (std::getline(my_f, temp_vec) )
88  {
89  std::istringstream ss(temp_vec);
90  std::string token;
91  std::vector<double> temp_weight, temp_bias;
92  while (std::getline(ss, token, ',')) // coefficients are separated by comma
93  {
94  temp_weight.push_back(stod(token));
95  }
96  if (temp_weight.size() != num_nodes_[layer_index] * num_nodes_[layer_index + 1]) {
97  throw BuildException({
98  "WARNING: layer weight size = " + std::to_string(temp_weight.size()) + " expected: "
99  + std::to_string(num_nodes_[layer_index] * num_nodes_[layer_index + 1])
100  });
101  }
102  Eigen::Map<Eigen::MatrixXd> temp_weight_v(
103  &temp_weight[0], num_nodes_[layer_index], num_nodes_[layer_index + 1]
104  );
105  std::getline(my_f, temp_vec);
106  std::istringstream ss2(temp_vec);
107  while (std::getline(ss2, token, ','))
108  {
109  temp_bias.push_back(stod(token));
110  }
111  if (temp_bias.size() != num_nodes_[layer_index + 1]) {
112  throw BuildException({
113  "WARNING: layer bias size = " + std::to_string(temp_bias.size())
114  + " expected: " + std::to_string(num_nodes_[layer_index + 1])
115  });
116  }
117  Vector temp_bias_v = Vector::Map(temp_bias.data(), temp_bias.size());
118  weight_coeff_.push_back(temp_weight_v);
119  bias_.push_back(temp_bias_v);
120  layer_index ++;
121  }
122  }
Eigen::VectorXd Vector
Arbitrary length vector.
Definition: types.h:30

Referenced by Build().

Here is the caller graph for this function:

Member Function Documentation

◆ Build()

static ANNCV* SSAGES::ANNCV::Build ( const Json::Value &  json,
const std::string &  path 
)
inlinestatic

Set up collective variable.

Parameters
jsonJSON input value.
pathPath for JSON path specification.
Returns
Pointer to the CV built by this function. nullptr in case of unknown error.

Builds a CV from a JSON node. Returns a pointer to the built cv. If an unknown error is encountered, this function will return a nullptr, but generally it will throw a BuildException on failure.

Warning
Object lifetime is the caller's responsibility.

Definition at line 244 of file ANNCV.h.

245  {
246  Json::ObjectRequirement validator;
247  Json::Value schema;
248  Json::CharReaderBuilder rbuilder;
249  Json::CharReader* reader = rbuilder.newCharReader();
250 
251  reader->parse(JsonSchema::ANNCV.c_str(),
252  JsonSchema::ANNCV.c_str() + JsonSchema::ANNCV.size(),
253  &schema, NULL);
254  validator.Parse(schema, path);
255 
256  // Validate inputs.
257  validator.Validate(json, path);
258  if(validator.HasErrors())
259  throw BuildException(validator.GetErrors());
260 
261  std::vector<int> atomids;
262  for(auto& s : json["atom_ids"])
263  atomids.push_back(s.asInt());
264  double scaling_factor = json["scaling_factor"].asDouble();
265  std::vector<unsigned int> num_nodes;
266  for (auto &s : json["num_nodes"])
267  num_nodes.push_back(s.asUInt());
268  std::vector<std::string> activations;
269  std::string coeff_file = json["coeff_file"].asString();
270  for (auto &s : json["activations"]) {
271  auto name = s.asString();
272  if (name == "Tanh" || name == "ReLU" || name == "Linear")
273  activations.push_back(name);
274  else
275  throw std::runtime_error("invalid activation function " + name + " provided");
276  }
277  int out_index = json["index"].asInt();
278  return new ANNCV(atomids, scaling_factor, num_nodes, coeff_file, activations, out_index);
279  }
Requirements on an object.
virtual void Parse(Value json, const std::string &path) override
Parse JSON value to generate Requirement(s).
virtual void Validate(const Value &json, const std::string &path) override
Validate JSON value.
std::vector< std::string > GetErrors()
Get list of error messages.
Definition: Requirement.h:92
bool HasErrors()
Check if errors have occured.
Definition: Requirement.h:86
ANNCV(Label atomids, double scaling_factor, std::vector< unsigned int > num_nodes, std::string coeff_file, std::vector< std::string > activations, int out_index)
Constructor.
Definition: ANNCV.h:66

References ANNCV(), Json::Requirement::GetErrors(), Json::Requirement::HasErrors(), Json::ObjectRequirement::Parse(), and Json::ObjectRequirement::Validate().

Referenced by SSAGES::CollectiveVariable::BuildCV().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ Evaluate()

void SSAGES::ANNCV::Evaluate ( const Snapshot snapshot)
inlineoverride

Evaluate the CV.

Parameters
snapshotCurrent simulation snapshot.

Definition at line 194 of file ANNCV.h.

195  {
196  // Get data from snapshot.
197  auto n = snapshot.GetNumAtoms();
198  const auto& pos = snapshot.GetPositions();
199  auto& comm = snapshot.GetCommunicator();
200 
201  // Initialize gradient.
202  std::fill(grad_.begin(), grad_.end(), Vector3{0,0,0});
203  grad_.resize(n, Vector3{0,0,0});
204 
205  // Vector3 xi{0, 0, 0}, xj{0, 0, 0}, xk{0, 0, 0};
206  std::vector<Vector3> positions(atomids_.size(), Vector3({0,0,0}));
207  Label local_idx;
208  snapshot.GetLocalIndices(atomids_, &local_idx);
209  for (size_t ii = 0; ii < atomids_.size(); ii ++) {
210  if (local_idx[ii] != -1) {
211  positions[ii] = pos[local_idx[ii]];
212  }
213  }
214  // By performing a reduce, we actually collect all. This can
215  // be converted to a more intelligent allgater on rank then bcast.
216  MPI_Allreduce(MPI_IN_PLACE, positions.data(), positions.size(), MPI_DOUBLE, MPI_SUM, comm);
217  auto com = snapshot.CenterOfMass(atomids_, false); // center of mass coordinates (not weighted with mass)
218  // remove translation degree of freedom
219  for (auto& temp_pos: positions) {
220  temp_pos = temp_pos - com;
221  }
222  Vector input_vec(positions.size() * 3); // flatten input vector
223  for (size_t ii = 0; ii < positions.size(); ii ++) {
224  for (size_t jj = 0; jj < 3; jj ++) {
225  input_vec[ii * 3 + jj] = positions[ii][jj];
226  }
227  }
228  input_vec = input_vec / scaling_factor_;
229  auto output_of_layers = forward_prop(input_vec);
230  val_ = output_of_layers[output_of_layers.size() - 1][out_index_];
231  auto deriv_back = back_prop(output_of_layers);
232  // subtract mean from deriv_back[0]
233  int num_atoms = deriv_back[0].size() / 3;
234  double average[3] = {0};
235  for (int kk = 0; kk < 3; kk ++) {
236  for (int ii = 0; ii < num_atoms; ii ++) {
237  average[kk] += deriv_back[0][ii * 3 + kk];
238  }
239  average[kk] /= num_atoms;
240  }
241  }
std::vector< Vector3 > grad_
Gradient vector dCv/dxi.
double val_
Current value of CV.
Eigen::Vector3d Vector3
Three-dimensional vector.
Definition: types.h:33
std::vector< int > Label
List of integers.
Definition: types.h:48

References SSAGES::Snapshot::CenterOfMass(), SSAGES::Snapshot::GetCommunicator(), SSAGES::Snapshot::GetLocalIndices(), SSAGES::Snapshot::GetNumAtoms(), SSAGES::Snapshot::GetPositions(), SSAGES::CollectiveVariable::grad_, and SSAGES::CollectiveVariable::val_.

Here is the call graph for this function:

◆ Initialize()

void SSAGES::ANNCV::Initialize ( const Snapshot snapshot)
inlineoverride

Initialize necessary variables.

Parameters
snapshotCurrent simulation snapshot.

Definition at line 128 of file ANNCV.h.

129  {
130  using std::to_string;
131 
132  std::vector<int> found;
133  snapshot.GetLocalIndices(atomids_, &found);
134  size_t nfound = found.size();
135  MPI_Allreduce(MPI_IN_PLACE, &nfound, 1, MPI_INT, MPI_SUM, snapshot.GetCommunicator());
136 
137  if(nfound != atomids_.size())
138  throw BuildException({
139  "ANNCV: Expected to find " +
140  to_string(atomids_.size()) +
141  " atoms, but only found " +
142  to_string(nfound) + "."
143  });
144  }

References SSAGES::Snapshot::GetCommunicator(), and SSAGES::Snapshot::GetLocalIndices().

Here is the call graph for this function:

The documentation for this class was generated from the following file: