CASM  1.1.0
A Clusters Approach to Statistical Mechanics
MonteIO_impl.hh
Go to the documentation of this file.
1 #ifndef CASM_MonteIO_impl_HH
2 #define CASM_MonteIO_impl_HH
3 
4 #include <boost/filesystem.hpp>
5 #include <boost/filesystem/fstream.hpp>
6 
7 #include "casm/casm_io/Log.hh"
9 #include "casm/clex/PrimClex.hh"
15 
16 namespace CASM {
17 namespace Monte {
18 
19 // --- Template definitions ---------------------
20 
23 template <typename MonteType>
25  auto evaluator = [=](const ConstMonteCarloPtr &mc) {
26  ConstMonteCarloPtr ptr = mc;
27  return static_cast<const MonteType *>(ptr)->conditions().temperature();
28  };
29  return GenericDatumFormatter<double, ConstMonteCarloPtr>("T", "Temperature",
30  evaluator);
31 }
32 
35 template <typename MonteType>
37  auto evaluator = [=](const ConstMonteCarloPtr &mc) {
38  ConstMonteCarloPtr ptr = mc;
39  return static_cast<const MonteType *>(ptr)->conditions().beta();
40  };
42  evaluator);
43 }
44 
47 template <typename MonteType>
49 MonteCarloParamChemPotFormatter(const MonteType &mc, int index) {
50  auto evaluator = [=](const ConstMonteCarloPtr &mc) {
51  ConstMonteCarloPtr ptr = mc;
52  return static_cast<const MonteType *>(ptr)->conditions().param_chem_pot()(
53  index);
54  };
55  std::string header = std::string("param_chem_pot(") +
56  CompositionConverter::comp_var(index) + ")";
58  evaluator);
59 }
60 
63 template <typename MonteType>
65  const MonteType &mc, int index) {
66  auto evaluator = [=](const ConstMonteCarloPtr &mc) {
67  ConstMonteCarloPtr ptr = mc;
68  return static_cast<const MonteType *>(ptr)->conditions().chem_pot(index);
69  };
70  std::string header = std::string("chem_pot(") +
71  struc_molecule_name(mc.primclex().prim())[index] + ")";
73  evaluator);
74 }
75 
78 template <typename MonteType>
80  const MonteType &mc, int index) {
81  auto evaluator = [=](const ConstMonteCarloPtr &mc) {
82  ConstMonteCarloPtr ptr = mc;
83  return static_cast<const MonteType *>(ptr)
84  ->conditions()
85  .param_composition()(index);
86  };
87  std::string header =
88  std::string("comp(") + CompositionConverter::comp_var(index) + ")";
90  evaluator);
91 }
92 
95 template <typename MonteType>
97  const MonteType &mc, int index) {
98  auto evaluator = [=](const ConstMonteCarloPtr &mc) {
99  ConstMonteCarloPtr ptr = mc;
100  return static_cast<const MonteType *>(ptr)->conditions().mol_composition()(
101  index);
102  };
103  std::string header = std::string("comp_n(") +
104  xtal::struc_molecule_name(mc.primclex().prim())[index] +
105  ")";
107  evaluator);
108 }
109 
112 template <typename MonteType>
114  const MonteType &mc, int index) {
115  auto evaluator = [=](const ConstMonteCarloPtr &mc) {
116  const MonteType *mc_ptr = static_cast<const MonteType *>(mc);
117  Eigen::VectorXd comp_n = mc_ptr->comp_n();
118  Eigen::VectorXd site_frac = comp_n / mc->primclex().prim().basis().size();
119  return site_frac(index);
120  };
121  std::string header = std::string("site_frac(") +
122  struc_molecule_name(mc.primclex().prim())[index] + ")";
124  evaluator);
125 }
126 
129 template <typename MonteType>
131  const MonteType &mc, int index) {
132  auto evaluator = [=](const ConstMonteCarloPtr &mc) {
133  const MonteType *mc_ptr = static_cast<const MonteType *>(mc);
134  Eigen::VectorXd comp_n = mc_ptr->comp_n();
135  if (mc->primclex().vacancy_allowed()) {
136  comp_n(mc->primclex().vacancy_index()) = 0.0;
137  }
138  auto atom_frac = comp_n / comp_n.sum();
139  return atom_frac[index];
140  };
141  std::string header = std::string("atom_frac(") +
142  struc_molecule_name(mc.primclex().prim())[index] + ")";
144  evaluator);
145 }
146 
151 template <typename MonteType>
154  auto evaluator = [=](const ConstMonteCarloPtr &mc) {
155  CovEvaluator cov_evaluator("potential_energy", "potential_energy");
156  auto N = mc->supercell().volume();
157  ConstMonteCarloPtr ptr = mc;
158  auto T = static_cast<const MonteType *>(ptr)->conditions().temperature();
159  return cov_evaluator(mc) * N / (KB * T * T);
160  };
161 
162  auto validator = [=](const ConstMonteCarloPtr &mc) {
163  return mc->is_equilibrated().first;
164  };
165 
166  std::string header = std::string("heat_capacity");
167 
169  header, header, evaluator, validator);
170 }
171 
178 template <typename MonteType>
180  std::string comp_var_i, std::string comp_var_j) {
181  std::string name_i = "comp(" + comp_var_i + ")";
182  std::string name_j = "comp(" + comp_var_j + ")";
183 
184  auto evaluator = [=](const ConstMonteCarloPtr &mc) {
185  CovEvaluator cov_evaluator(name_i, name_j);
186  auto N = mc->supercell().volume();
187  ConstMonteCarloPtr ptr = mc;
188  auto T = static_cast<const MonteType *>(ptr)->conditions().temperature();
189  return cov_evaluator(mc) * N / (KB * T);
190  };
191 
192  auto validator = [=](const ConstMonteCarloPtr &mc) {
193  return mc->is_equilibrated().first;
194  };
195 
196  std::string header =
197  std::string("susc_x(" + comp_var_i + "," + comp_var_j + ")");
198 
200  header, header, evaluator, validator);
201 }
202 
209 template <typename MonteType>
211  std::string species_i, std::string species_j) {
212  std::string name_i = "comp_n(" + species_i + ")";
213  std::string name_j = "comp_n(" + species_j + ")";
214 
215  auto evaluator = [=](const ConstMonteCarloPtr &mc) {
216  CovEvaluator cov_evaluator(name_i, name_j);
217  auto N = mc->supercell().volume();
218  ConstMonteCarloPtr ptr = mc;
219  auto T = static_cast<const MonteType *>(ptr)->conditions().temperature();
220  return cov_evaluator(mc) * N / (KB * T);
221  };
222 
223  auto validator = [=](const ConstMonteCarloPtr &mc) {
224  return mc->is_equilibrated().first;
225  };
226 
227  std::string header =
228  std::string("susc_n(" + species_i + "," + species_j + ")");
229 
231  header, header, evaluator, validator);
232 }
233 
240 template <typename MonteType>
242 MonteCarloThermoChemSuscXFormatter(std::string comp_var_i) {
243  std::string name_i = "comp(" + comp_var_i + ")";
244 
245  auto evaluator = [=](const ConstMonteCarloPtr &mc) {
246  CovEvaluator cov_evaluator("potential_energy", name_i);
247  auto N = mc->supercell().volume();
248  ConstMonteCarloPtr ptr = mc;
249  auto T = static_cast<const MonteType *>(ptr)->conditions().temperature();
250  return cov_evaluator(mc) * N / (KB * T);
251  };
252 
253  auto validator = [=](const ConstMonteCarloPtr &mc) {
254  return mc->is_equilibrated().first;
255  };
256 
257  std::string header = std::string("susc_x(S," + comp_var_i + ")");
258 
260  header, header, evaluator, validator);
261 }
262 
264 template <typename MonteType>
266 MonteCarloThermoChemSuscNFormatter(std::string species_i) {
267  std::string name_i = "comp_n(" + species_i + ")";
268 
269  auto evaluator = [=](const ConstMonteCarloPtr &mc) {
270  CovEvaluator cov_evaluator("potential_energy", name_i);
271  auto N = mc->supercell().volume();
272  ConstMonteCarloPtr ptr = mc;
273  auto T = static_cast<const MonteType *>(ptr)->conditions().temperature();
274  return cov_evaluator(mc) * N / (KB * T);
275  };
276 
277  auto validator = [=](const ConstMonteCarloPtr &mc) {
278  return mc->is_equilibrated().first;
279  };
280 
281  std::string header = std::string("susc_n(S," + species_i + ")");
282 
284  header, header, evaluator, validator);
285 }
286 
289 template <typename MonteType>
290 void write_results(const MonteSettings &settings, const MonteType &mc,
291  Log &_log) {
292  try {
293  fs::create_directories(settings.output_directory());
295  auto formatter = make_results_formatter(mc);
296 
297  // write csv path results
298  if (settings.write_csv()) {
299  _log << "write: " << dir.results_csv() << "\n";
300  fs::path file = dir.results_csv();
301  fs::ofstream sout;
302 
303  if (!fs::exists(file)) {
304  sout.open(file);
305  formatter.print_header(&mc, sout);
306  } else {
307  sout.open(file, std::ios::app);
308  }
309 
310  formatter.print(&mc, sout);
311 
312  sout.close();
313  }
314 
315  // write json path results
316  if (settings.write_json()) {
317  _log << "write: " << dir.results_json() << "\n";
318  fs::path file = dir.results_json();
319 
320  jsonParser results;
321  if (fs::exists(file)) {
322  results.read(file);
323  } else {
324  results = jsonParser::object();
325  }
326 
327  formatter.to_json_arrays(&mc, results);
328  results.write(file);
329  }
330  } catch (...) {
331  std::cerr << "ERROR writing results" << std::endl;
332  throw;
333  }
334 }
335 
337 template <typename MonteType>
338 void write_conditions_json(const MonteSettings &settings, const MonteType &mc,
339  size_type cond_index, Log &_log) {
340  try {
342  fs::create_directories(dir.conditions_dir(cond_index));
343  jsonParser json;
344  to_json(mc.conditions(), json);
345  _log << "write: " << dir.conditions_json(cond_index) << "\n";
346  json.write(dir.conditions_json(cond_index));
347  } catch (...) {
348  std::cerr << "ERROR writing conditions.json" << std::endl;
349  throw;
350  }
351 }
352 
353 } // namespace Monte
354 } // namespace CASM
355 
356 #endif
static std::string comp_var(size_type i)
Composition variable names: "a", "b", ...
A DatumFormatter that returns a value of specified type, via functions that may be specified at runti...
Definition: Log.hh:48
Interface base class for all types of Monte Carlo simulations (not meant to be used polymorphically)
Definition: MonteCarlo.hh:44
Settings for Monte Carlo calculations.
const fs::path output_directory() const
Directory where output should go.
bool write_csv() const
Write csv versions of files? (csv is the default format if no 'output_format' given)
bool write_json() const
Write json versions of files?
static jsonParser object()
Returns an empty json object.
Definition: jsonParser.hh:395
bool read(std::istream &stream)
Reads json from the stream.
Definition: jsonParser.cc:168
void write(const std::string &file_name, unsigned int indent=2, unsigned int prec=12) const
Write json to file.
Definition: jsonParser.cc:196
Eigen::VectorXd comp_n(const Configuration &config)
Returns the composition, as number of each species per unit cell.
Eigen::VectorXd site_frac(const Configuration &config)
Returns the composition as site fraction, in the order of Structure::get_struc_molecule.
std::vector< std::string > struc_molecule_name(BasicStructure const &_struc)
Returns an Array of each possible Molecule in this Structure.
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloSiteFracFormatter(const MonteType &mc, int index)
Print site_frac(x) for any class MonteType with valid 'Eigen::VectorXd MonteType::comp_n'.
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloSuscXFormatter(std::string comp_var_i, std::string comp_var_j)
Print parametric susceptibility, 'susc_x(a,b)'.
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloParamChemPotFormatter(const MonteType &mc, int index)
Print param_chem_pot(x)
Definition: MonteIO_impl.hh:49
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloChemPotFormatter(const MonteType &mc, int index)
Print chem_pot(N)
Definition: MonteIO_impl.hh:64
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloTFormatter()
Print Temperature.
Definition: MonteIO_impl.hh:24
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloAtomFracFormatter(const MonteType &mc, int index)
Print atom_frac(X) for any class MonteType with valid 'Eigen::VectorXd MonteType::comp_n'.
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloThermoChemSuscNFormatter(std::string species_i)
Print thermo-chemical susceptibility, 'susc_n(S,A)'.
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloSuscNFormatter(std::string species_i, std::string species_j)
Print susceptibility, 'susc_n(A,B)'.
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloBetaFormatter()
Print Beta.
Definition: MonteIO_impl.hh:36
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloThermoChemSuscXFormatter(std::string comp_var_i)
Print parametric thermo-chemical susceptibility, 'susc_x(S,a)'.
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloCompNFormatter(const MonteType &mc, int index)
Print comp_n(N)
Definition: MonteIO_impl.hh:96
void write_results(const MonteSettings &settings, const MonteType &mc, Log &_log)
Will create new file or append to existing file results of the latest run.
void write_conditions_json(const MonteSettings &settings, const MonteType &mc, size_type cond_index, Log &_log)
Write conditions to conditions.cond_index directory.
jsonParser & to_json(const CanonicalConditions &conditions, jsonParser &json)
Store CanonicalConditions in JSON format.
Definition: CanonicalIO.cc:102
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloHeatCapacityFormatter()
Print heat capacity, 'heat_capacity'.
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloCompFormatter(const MonteType &mc, int index)
Print comp(x)
Definition: MonteIO_impl.hh:79
DataFormatter< ConstMonteCarloPtr > make_results_formatter(const Canonical &mc)
Make a LTE results formatter.
Definition: CanonicalIO.cc:31
Main CASM namespace.
Definition: APICommand.hh:8
const double KB
Definition: definitions.hh:33
Eigen::VectorXd VectorXd
DirectoryStructure const & dir
Definition: settings.cc:136
Functor to help evaluate covariance.
Definition: MonteIO.hh:144