CASM
AClustersApproachtoStatisticalMechanics
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules
MonteIO_impl.hh
Go to the documentation of this file.
1 #ifndef CASM_MonteIO_impl_HH
2 #define CASM_MonteIO_impl_HH
3 
5 
6 namespace CASM {
7 
8  // --- Template definitions ---------------------
9 
11  template<typename MonteType>
13  auto evaluator = [ = ](const ConstMonteCarloPtr & mc) {
14  ConstMonteCarloPtr ptr = mc;
15  return static_cast<const MonteType *>(ptr)->conditions().temperature();
16  };
17  return GenericDatumFormatter<double, ConstMonteCarloPtr>("T", "Temperature", evaluator);
18  }
19 
21  template<typename MonteType>
23  auto evaluator = [ = ](const ConstMonteCarloPtr & mc) {
24  ConstMonteCarloPtr ptr = mc;
25  return static_cast<const MonteType *>(ptr)->conditions().beta();
26  };
27  return GenericDatumFormatter<double, ConstMonteCarloPtr>("Beta", "Beta", evaluator);
28  }
29 
31  template<typename MonteType>
33  auto evaluator = [ = ](const ConstMonteCarloPtr & mc) {
34  ConstMonteCarloPtr ptr = mc;
35  return static_cast<const MonteType *>(ptr)->conditions().param_chem_pot()(index);
36  };
37  std::string header = std::string("param_chem_pot(") + CompositionConverter::comp_var(index) + ")";
38  return GenericDatumFormatter<double, ConstMonteCarloPtr>(header, header, evaluator);
39  }
40 
42  template<typename MonteType>
44  auto evaluator = [ = ](const ConstMonteCarloPtr & mc) {
45  ConstMonteCarloPtr ptr = mc;
46  return static_cast<const MonteType *>(ptr)->conditions().chem_pot(index);
47  };
48  std::string header = std::string("chem_pot(") + mc.primclex().get_prim().get_struc_molecule_name()[index] + ")";
49  return GenericDatumFormatter<double, ConstMonteCarloPtr>(header, header, evaluator);
50  }
51 
53  template<typename MonteType>
55  auto evaluator = [ = ](const ConstMonteCarloPtr & mc) {
56  ConstMonteCarloPtr ptr = mc;
57  return static_cast<const MonteType *>(ptr)->conditions().param_composition()(index);
58  };
59  std::string header = std::string("comp(") + CompositionConverter::comp_var(index) + ")";
60  return GenericDatumFormatter<double, ConstMonteCarloPtr>(header, header, evaluator);
61  }
62 
64  template<typename MonteType>
66  auto evaluator = [ = ](const ConstMonteCarloPtr & mc) {
67  ConstMonteCarloPtr ptr = mc;
68  return static_cast<const MonteType *>(ptr)->conditions().mol_composition()(index);
69  };
70  std::string header = std::string("comp_n(") + mc.primclex().get_prim().get_struc_molecule_name()[index] + ")";
71  return GenericDatumFormatter<double, ConstMonteCarloPtr>(header, header, evaluator);
72  }
73 
75  template<typename MonteType>
77  auto evaluator = [ = ](const ConstMonteCarloPtr & mc) {
78  const MonteType *mc_ptr = static_cast<const MonteType *>(mc);
79  Eigen::VectorXd comp_n = mc_ptr->comp_n();
80  Eigen::VectorXd site_frac = comp_n / mc->primclex().get_prim().basis.size();
81  return site_frac(index);
82  };
83  std::string header = std::string("site_frac(") + mc.primclex().get_prim().get_struc_molecule_name()[index] + ")";
84  return GenericDatumFormatter<double, ConstMonteCarloPtr>(header, header, evaluator);
85  }
86 
88  template<typename MonteType>
90  auto evaluator = [ = ](const ConstMonteCarloPtr & mc) {
91  const MonteType *mc_ptr = static_cast<const MonteType *>(mc);
92  Eigen::VectorXd comp_n = mc_ptr->comp_n();
93  if(mc->primclex().vacancy_allowed()) {
94  comp_n(mc->primclex().vacancy_index()) = 0.0;
95  }
96  auto atom_frac = comp_n / comp_n.sum();
97  return atom_frac[index];
98  };
99  std::string header = std::string("atom_frac(") + mc.primclex().get_prim().get_struc_molecule_name()[index] + ")";
100  return GenericDatumFormatter<double, ConstMonteCarloPtr>(header, header, evaluator);
101  }
102 
106  template<typename MonteType>
108 
109  auto evaluator = [ = ](const ConstMonteCarloPtr & mc) {
110  CovEvaluator cov_evaluator("potential_energy", "potential_energy");
111  auto N = mc->supercell().volume();
112  ConstMonteCarloPtr ptr = mc;
113  auto T = static_cast<const MonteType *>(ptr)->conditions().temperature();
114  return cov_evaluator(mc) * N / (KB * T * T);
115  };
116 
117  auto validator = [ = ](const ConstMonteCarloPtr & mc) {
118  return mc->is_equilibrated().first;
119  };
120 
121  std::string header = std::string("heat_capacity");
122 
123  return GenericDatumFormatter<double, ConstMonteCarloPtr>(header, header, evaluator, validator);
124  }
125 
132  template<typename MonteType>
133  GenericDatumFormatter<double, ConstMonteCarloPtr> MonteCarloSuscXFormatter(std::string comp_var_i, std::string comp_var_j) {
134 
135  std::string name_i = "comp(" + comp_var_i + ")";
136  std::string name_j = "comp(" + comp_var_j + ")";
137 
138  auto evaluator = [ = ](const ConstMonteCarloPtr & mc) {
139  CovEvaluator cov_evaluator(name_i, name_j);
140  auto N = mc->supercell().volume();
141  ConstMonteCarloPtr ptr = mc;
142  auto T = static_cast<const MonteType *>(ptr)->conditions().temperature();
143  return cov_evaluator(mc) * N / (KB * T);
144  };
145 
146  auto validator = [ = ](const ConstMonteCarloPtr & mc) {
147  return mc->is_equilibrated().first;
148  };
149 
150  std::string header = std::string("susc_x(" + comp_var_i + "," + comp_var_j + ")");
151 
152  return GenericDatumFormatter<double, ConstMonteCarloPtr>(header, header, evaluator, validator);
153  }
154 
161  template<typename MonteType>
162  GenericDatumFormatter<double, ConstMonteCarloPtr> MonteCarloSuscNFormatter(std::string species_i, std::string species_j) {
163 
164  std::string name_i = "comp_n(" + species_i + ")";
165  std::string name_j = "comp_n(" + species_j + ")";
166 
167  auto evaluator = [ = ](const ConstMonteCarloPtr & mc) {
168  CovEvaluator cov_evaluator(name_i, name_j);
169  auto N = mc->supercell().volume();
170  ConstMonteCarloPtr ptr = mc;
171  auto T = static_cast<const MonteType *>(ptr)->conditions().temperature();
172  return cov_evaluator(mc) * N / (KB * T);
173  };
174 
175  auto validator = [ = ](const ConstMonteCarloPtr & mc) {
176  return mc->is_equilibrated().first;
177  };
178 
179  std::string header = std::string("susc_n(" + species_i + "," + species_j + ")");
180 
181  return GenericDatumFormatter<double, ConstMonteCarloPtr>(header, header, evaluator, validator);
182  }
183 
189  template<typename MonteType>
190  GenericDatumFormatter<double, ConstMonteCarloPtr>
191  MonteCarloThermoChemSuscXFormatter(std::string comp_var_i) {
192 
193  std::string name_i = "comp(" + comp_var_i + ")";
194 
195  auto evaluator = [ = ](const ConstMonteCarloPtr & mc) {
196  CovEvaluator cov_evaluator("potential_energy", name_i);
197  auto N = mc->supercell().volume();
198  ConstMonteCarloPtr ptr = mc;
199  auto T = static_cast<const MonteType *>(ptr)->conditions().temperature();
200  return cov_evaluator(mc) * N / (KB * T);
201  };
202 
203  auto validator = [ = ](const ConstMonteCarloPtr & mc) {
204  return mc->is_equilibrated().first;
205  };
206 
207  std::string header = std::string("susc_x(S," + comp_var_i + ")");
208 
209  return GenericDatumFormatter<double, ConstMonteCarloPtr>(header, header, evaluator, validator);
210  }
211 
213  template<typename MonteType>
214  GenericDatumFormatter<double, ConstMonteCarloPtr>
215  MonteCarloThermoChemSuscNFormatter(std::string species_i) {
216 
217  std::string name_i = "comp_n(" + species_i + ")";
218 
219  auto evaluator = [ = ](const ConstMonteCarloPtr & mc) {
220  CovEvaluator cov_evaluator("potential_energy", name_i);
221  auto N = mc->supercell().volume();
222  ConstMonteCarloPtr ptr = mc;
223  auto T = static_cast<const MonteType *>(ptr)->conditions().temperature();
224  return cov_evaluator(mc) * N / (KB * T);
225  };
226 
227  auto validator = [ = ](const ConstMonteCarloPtr & mc) {
228  return mc->is_equilibrated().first;
229  };
230 
231  std::string header = std::string("susc_n(S," + species_i + ")");
232 
233  return GenericDatumFormatter<double, ConstMonteCarloPtr>(header, header, evaluator, validator);
234  }
235 
237  template<typename MonteType>
238  void write_results(const MonteSettings &settings, const MonteType &mc, Log &_log) {
239  try {
240 
241  fs::create_directories(settings.output_directory());
243  auto formatter = make_results_formatter(mc);
244 
245  // write csv path results
246  if(settings.write_csv()) {
247  _log << "write: " << dir.results_csv() << "\n";
248  fs::path file = dir.results_csv();
249  fs::ofstream sout;
250 
251  if(!fs::exists(file)) {
252  sout.open(file);
253  formatter.print_header(&mc, sout);
254  }
255  else {
256  sout.open(file, std::ios::app);
257  }
258 
259  formatter.print(&mc, sout);
260 
261  sout.close();
262  }
263 
264  // write json path results
265  if(settings.write_json()) {
266  _log << "write: " << dir.results_json() << "\n";
267  fs::path file = dir.results_json();
268 
269  jsonParser results;
270  if(fs::exists(file)) {
271  results.read(file);
272  }
273  else {
274  results = jsonParser::object();
275  }
276 
277  formatter.to_json_arrays(&mc, results);
278  results.write(file);
279  }
280  }
281  catch(...) {
282  std::cerr << "ERROR writing results" << std::endl;
283  throw;
284  }
285 
286  }
287 
289  template<typename MonteType>
290  void write_conditions_json(const MonteSettings &settings, const MonteType &mc, Index cond_index, Log &_log) {
291  try {
293  fs::create_directories(dir.conditions_dir(cond_index));
294  jsonParser json;
295  to_json(mc.conditions(), json);
296  _log << "write: " << dir.conditions_json(cond_index) << "\n";
297  json.write(dir.conditions_json(cond_index));
298  }
299  catch(...) {
300  std::cerr << "ERROR writing conditions.json" << std::endl;
301  throw;
302  }
303  }
304 
305 }
306 
307 #endif
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloThermoChemSuscXFormatter(std::string comp_var_i)
Print parametric thermo-chemical susceptibility, 'susc_x(S,a)'.
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.
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloAtomFracFormatter(const MonteType &mc, int index)
Print atom_frac(X) for any class MonteType with valid 'Eigen::VectorXd MonteType::comp_n'.
Definition: MonteIO_impl.hh:89
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloParamChemPotFormatter(const MonteType &mc, int index)
Print param_chem_pot(x)
Definition: MonteIO_impl.hh:32
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloBetaFormatter()
Print Beta.
Definition: MonteIO_impl.hh:22
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloCompFormatter(const MonteType &mc, int index)
Print comp(x)
Definition: MonteIO_impl.hh:54
void write_conditions_json(const MonteSettings &settings, const MonteType &mc, Index cond_index, Log &_log)
Write conditions to conditions.cond_index directory.
void write(const std::string &file_name, unsigned int indent=2, unsigned int prec=12) const
Write json to file.
Definition: jsonParser.cc:191
jsonParser & to_json(const ClexDescription &desc, jsonParser &json)
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloCompNFormatter(const MonteType &mc, int index)
Print comp_n(N)
Definition: MonteIO_impl.hh:65
bool read(std::istream &stream)
Reads json from the stream.
Definition: jsonParser.cc:165
Main CASM namespace.
Definition: complete.cpp:8
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloHeatCapacityFormatter()
Print heat capacity, 'heat_capacity'.
const fs::path output_directory() const
Directory where output should go.
std::vector< std::string > get_struc_molecule_name() const
Returns an Array of each possible Molecule in this Structure.
Definition: Structure.cc:166
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloSuscXFormatter(std::string comp_var_i, std::string comp_var_j)
Print parametric susceptibility, 'susc_x(a,b)'.
Eigen::VectorXd comp_n(const ConfigDoF &configdof, const Supercell &scel)
Returns comp_n, the number of each molecule per primitive cell, ordered as Structure::get_struc_molec...
Definition: ConfigDoF.cc:312
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloChemPotFormatter(const MonteType &mc, int index)
Print chem_pot(N)
Definition: MonteIO_impl.hh:43
static std::string comp_var(size_type i)
Composition variable names: "a", "b", ...
EigenIndex Index
For long integer indexing:
const PrimClex & primclex() const
const Access the PrimClex that *this is based on
Definition: MonteCarlo.hh:68
Eigen::VectorXd VectorXd
const double KB
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloSiteFracFormatter(const MonteType &mc, int index)
Print site_frac(x) for any class MonteType with valid 'Eigen::VectorXd MonteType::comp_n'.
Definition: MonteIO_impl.hh:76
Interface base class for all types of Monte Carlo simulations (not meant to be used polymorphically) ...
Definition: MonteCarlo.hh:32
Eigen::VectorXd site_frac(const Configuration &config)
Returns the composition as site fraction, in the order of Structure::get_struc_molecule.
DirectoryStructure & dir
Definition: settings.cc:102
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloSuscNFormatter(std::string species_i, std::string species_j)
Print susceptibility, 'susc_n(A,B)'.
Functor to help evaluate covariance.
Definition: MonteIO.hh:144
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?
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloTFormatter()
Print Temperature.
Definition: MonteIO_impl.hh:12
GenericDatumFormatter< double, ConstMonteCarloPtr > MonteCarloThermoChemSuscNFormatter(std::string species_i)
Print thermo-chemical susceptibility, 'susc_n(S,A)'.
Settings for Monte Carlo calculations.
Definition: Log.hh:9
A DatumFormatter that returns a value of specified type, via functions that may be specified at runti...
static jsonParser object()
Returns an empty json object.
Definition: jsonParser.hh:329
const Structure & get_prim() const
const Access to primitive Structure
Definition: PrimClex.cc:260
DataFormatter< ConstMonteCarloPtr > make_results_formatter(const GrandCanonical &mc)
Make a LTE results formatter.