CASM
AClustersApproachtoStatisticalMechanics
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules
ConfigIOStrain.cc
Go to the documentation of this file.
1 #include <functional>
2 #include "casm/external/boost.hh"
3 #include "casm/clex/ConfigIO.hh"
6 
7 namespace CASM {
8 
9  namespace ConfigIO {
10  bool RelaxationStrain::parse_args(const std::string &args) {
11 
12  std::string tmetric_name = "GL";
13  std::string index_expr;
14  auto it = args.cbegin(), it_end = args.cend();
15  while(it != it_end) {
16  if(std::isdigit(*it) || (*it) == ':') {
17  auto it2 = it;
18  while(it2 != it_end && (std::isdigit(*it2) || isspace(*it2) || (*it2) == ':'))
19  ++it2;
20  index_expr = std::string(it, it2);
21  it = it2;
22  }
23  if(it != it_end && std::isalpha(*it)) {
24  auto it2 = it;
25  while(it2 != it_end && std::isalpha(*it2))
26  ++it2;
27  tmetric_name = std::string(it, it2);
28  it = it2;
29  }
30  while(it != it_end && !std::isalnum(*it) && *it != ':')
31  ++it;
32  }
33  if(m_metric_name.size() > 0 && tmetric_name != m_metric_name) {
34  return false;
35  }
36  m_metric_name = tmetric_name;
37  if(index_expr.size() > 0) {
38  _parse_index_expression(index_expr);
39  }
40 
41  return true;
42  }
43 
44  //****************************************************************************************
45 
46  void RelaxationStrain::init(const Configuration &_tmplt) const {
47  if(m_metric_name.size() == 0)
48  m_metric_name = "GL";
50  if(_index_rules().size() > 0)
51  return;
52 
53  for(Index i = 0; i < 6; i++)
54  _add_rule(std::vector<Index>({i}));
55 
56  }
57  //****************************************************************************************
58 
59  bool RelaxationStrain::validate(const Configuration &_config) const {
60  return _config.calc_properties().contains("relaxation_deformation");
61  }
62 
63  //****************************************************************************************
64 
65  std::vector<std::string> RelaxationStrain::col_header(const Configuration &_tmplt) const {
66  std::vector<std::string> col;
67  auto it(_index_rules().cbegin()), end_it(_index_rules().cend());
68  // Index s = max(8 - int(name().size()), 0);
69  for(; it != end_it; ++it) {
70  std::stringstream t_ss;
71  t_ss << " " << name() << '(' << m_metric_name << ',' << (*it)[0] << ')';
72  col.push_back(t_ss.str());
73  }
74  return col;
75  }
76 
77 
78  //****************************************************************************************
79 
80  std::string RelaxationStrain::short_header(const Configuration &_tmplt) const {
81  return name() + "(" + m_metric_name + ")";
82  }
83 
84  //****************************************************************************************
86  return m_straincalc.unrolled_strain_metric(_config.calc_properties()["relaxation_deformation"].get<Eigen::Matrix3d>());
87  }
88 
89  //****************************************************************************************
90 
91  bool DoFStrain::parse_args(const std::string &args) {
92 
93  std::string tmetric_name = "GL";
94  std::string index_expr;
95  auto it = args.cbegin(), it_end = args.cend();
96  while(it != it_end) {
97  if(std::isdigit(*it) || (*it) == ':') {
98  auto it2 = it;
99  while(it2 != it_end && (std::isdigit(*it2) || isspace(*it2) || (*it2) == ':'))
100  ++it2;
101  index_expr = std::string(it, it2);
102  it = it2;
103  }
104  if(it != it_end && std::isalpha(*it)) {
105  auto it2 = it;
106  while(it2 != it_end && std::isalpha(*it2))
107  ++it2;
108  tmetric_name = std::string(it, it2);
109  it = it2;
110  }
111  while(it != it_end && !std::isalnum(*it) && *it != ':')
112  ++it;
113  }
114  if(m_metric_name.size() > 0 && tmetric_name != m_metric_name) {
115  return false;
116  }
117  m_metric_name = tmetric_name;
118  if(index_expr.size() > 0) {
119  _parse_index_expression(index_expr);
120  }
121 
122  return true;
123  }
124 
125  //****************************************************************************************
126 
127  void DoFStrain::init(const Configuration &_tmplt) const {
128  if(m_metric_name.size() == 0)
129  m_metric_name = "GL";
131  if(_index_rules().size() > 0)
132  return;
133 
134  for(Index i = 0; i < 6; i++)
135  _add_rule(std::vector<Index>({i}));
136 
137  }
138 
139  //****************************************************************************************
140 
141  std::vector<std::string> DoFStrain::col_header(const Configuration &_tmplt) const {
142  std::vector<std::string> col;
143  auto it(_index_rules().cbegin()), end_it(_index_rules().cend());
144  // Index s = max(8 - int(name().size()), 0);
145  for(; it != end_it; ++it) {
146  std::stringstream t_ss;
147  t_ss << name() << '(' << m_metric_name << ',' << (*it)[0] << ')';
148  col.push_back(t_ss.str());
149  }
150  return col;
151  }
152 
153 
154  //****************************************************************************************
155 
156  std::string DoFStrain::short_header(const Configuration &_tmplt) const {
157  return name() + "(" + m_metric_name + ")";
158  }
159 
160  //****************************************************************************************
163  }
164 
165  }
166 
167 }
void set_mode(const std::string &mode_name)
void init(const Configuration &_tmplt) const override
const std::string & name() const
Returns a name for the formatter, which becomes the tag used for parsing.
Eigen::VectorXd evaluate(const Configuration &_config) const override
void _add_rule(const std::vector< Index > &new_rule) const
void init(const Configuration &_tmplt) const override
Main CASM namespace.
Definition: complete.cpp:8
Eigen::VectorXd evaluate(const Configuration &_config) const override
void _parse_index_expression(const std::string &_expr)
T get(Args...args) const
Get data from json, using one of several alternatives.
Definition: jsonParser.hh:729
bool parse_args(const std::string &args) override
Default implementation calls _parse_index_expression.
std::string short_header(const Configuration &_config) const override
std::string short_header(const Configuration &_config) const override
const Properties & calc_properties() const
std::vector< std::string > col_header(const Configuration &_config) const override
EigenIndex Index
For long integer indexing:
StrainConverter m_straincalc
Eigen::VectorXd VectorXd
std::vector< std::string > col_header(const Configuration &_config) const override
Eigen::Matrix3d Matrix3d
bool validate(const Configuration &_config) const override
bool contains(const std::string &name) const
Return true if JSON object contains 'name'.
Definition: jsonParser.cc:500
const IndexContainer & _index_rules() const
const Eigen::Matrix3d & deformation() const
Applied strain.
bool parse_args(const std::string &args) override
Default implementation calls _parse_index_expression.
VectorXd unrolled_strain_metric(const Matrix3d &F) const
A Configuration represents the values of all degrees of freedom in a Supercell.