CASM
AClustersApproachtoStatisticalMechanics
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules
GrandCanonical.hh
Go to the documentation of this file.
1 #ifndef CASM_GrandCanonical_HH
2 #define CASM_GrandCanonical_HH
3 
4 #include "casm/clex/Clex.hh"
12 
13 
14 namespace CASM {
15 
31  class GrandCanonical : public MonteCarlo {
32 
33  public:
34 
35  static const Monte::ENSEMBLE ensemble;
39 
40 
42  GrandCanonical(PrimClex &primclex, const SettingsType &settings, Log &_log);
43 
44 
46  Index steps_per_pass() const;
47 
48 
50  const CondType &conditions() const;
51 
53  void set_conditions(const CondType &new_conditions);
54 
56  void set_configdof(const ConfigDoF &configdof, const std::string &msg = "");
57 
59  std::pair<ConfigDoF, std::string> set_state(
60  const GrandCanonicalConditions &new_conditions,
62 
64  void set_state(const CondType &new_conditions,
65  const ConfigDoF &configdof,
66  const std::string &msg = "");
67 
69  const EventType &propose();
70 
72  bool check(const EventType &event);
73 
75  void accept(const EventType &event);
76 
78  void reject(const EventType &event);
79 
80  void check_corr() {
81  std::cout << "corr:" << std::endl;
82  std::cout << correlations_vec(_configdof(), supercell(), _clexulator()) << std::endl;
83  std::cout << "OK corr" << std::endl;
84  }
85 
87  void write_results(Index cond_index) const;
88 
89 
91  double lte_grand_canonical_free_energy() const;
92 
94  const double &formation_energy() const {
95  return *m_formation_energy;
96  }
97 
99  const double &potential_energy() const {
100  return *m_potential_energy;
101  }
102 
104  const Eigen::VectorXd &corr() const {
105  return *m_corr;
106  }
107 
109  const Eigen::VectorXd &comp_n() const {
110  return *m_comp_n;
111  }
112 
114  double potential_energy(const Configuration &config) const;
115 
116 
117  private:
118 
120  double &_formation_energy() {
121  return *m_formation_energy;
122  }
123 
125  double &_potential_energy() {
126  return *m_potential_energy;
127  }
128 
131  return *m_corr;
132  }
133 
136  return *m_comp_n;
137  }
138 
141  }
142 
143  const ECIContainer &_eci() const {
144  return m_formation_energy_clex.eci();
145  }
146 
148  void _set_dCorr(GrandCanonicalEvent &event,
149  Index mutating_site,
150  int sublat,
151  int current_occupant,
152  int new_occupant,
153  bool use_deltas,
154  bool all_correlations) const;
155 
158  std::string title,
159  std::string colheader,
160  bool all_correlations) const;
161 
164  Index mutating_site,
165  int sublat,
166  int current_occupant,
167  int new_occupant) const;
168 
170  void _update_properties();
171 
173  ConfigDoF _default_motif() const;
174 
176  std::pair<ConfigDoF, std::string> _auto_motif(const GrandCanonicalConditions &cond) const;
177 
179  std::pair<ConfigDoF, std::string> _restricted_auto_motif(const GrandCanonicalConditions &cond) const;
180 
182  ConfigDoF _configname_motif(const std::string &configname) const;
183 
184 
187 
190 
193 
196 
198  EventType m_event;
199 
202 
203 
204  // ---- Pointers to properties for faster access
205 
208 
211 
214 
217 
218  };
219 
220 }
221 
222 #endif
223 
224 
225 
226 
void set_conditions(const CondType &new_conditions)
Set conditions and clear previously collected data.
Eigen::VectorXd & _corr()
Correlations, normalized per primitive cell.
ConfigDoF & _configdof() const
Access current microstate.
Definition: MonteCarlo.hh:204
Data structure for storing information regarding a proposed grand canonical Monte Carlo event...
const SiteExchanger m_site_swaps
Keeps track of what sites can change to what.
GrandCanonicalConditions CondType
Clexulator & clexulator(Log &status_log=null_log()) const
Definition: Clex.cc:24
const Configuration & config() const
const Access current microstate
Definition: MonteCarlo.hh:88
const Eigen::VectorXd & comp_n() const
Number of atoms of each type, normalized per primitive cell.
const ECIContainer & _eci() const
bool check(const EventType &event)
Based on a random number, decide if the change in energy from the proposed event is low enough to be ...
ENSEMBLE
Monte Carlo ensemble type.
const ConfigDoF & configdof() const
const Access current microstate
Definition: MonteCarlo.hh:93
double * m_formation_energy
Formation energy, normalized per primitive cell.
EventType m_event
Event to propose, check, accept/reject:
Main CASM namespace.
Definition: complete.cpp:8
void _print_correlations(const Eigen::VectorXd &corr, std::string title, std::string colheader, bool all_correlations) const
Print correlations to _log()
void reject(const EventType &event)
Nothing needs to be done to reject a GrandCanonicalEvent.
bool m_all_correlations
If true, calculate all correlations; if false, calculate correlations with non-zero eci...
Eigen::VectorXd & _comp_n()
Number of atoms of each type, normalized per primitive cell.
void write_results(Index cond_index) const
Write results to files.
Eigen::VectorXd * m_corr
Correlations, normalized per primitive cell.
void _update_properties()
Calculate properties given current conditions.
void set_configdof(const ConfigDoF &configdof, const std::string &msg="")
Set configdof and clear previously collected data.
Log & _log() const
Definition: MonteCarlo.hh:208
double lte_grand_canonical_free_energy() const
Calculate the single spin flip low temperature expansion of the grand canonical potential.
double & _formation_energy()
Formation energy, normalized per primitive cell.
const double & potential_energy() const
Potential energy, normalized per primitive cell.
double & _potential_energy()
Potential energy, normalized per primitive cell.
GrandCanonicalEvent EventType
const Supercell & supercell() const
const Access the Supercell that *this is based on
Definition: MonteCarlo.hh:73
const EventType & propose()
Propose a new event, calculate delta properties, and return reference to it.
EigenIndex Index
For long integer indexing:
Evaluates correlations.
Definition: Clexulator.hh:240
A container class for the different degrees of freedom a Configuration might have.
Definition: ConfigDoF.hh:27
const Eigen::VectorXd & corr() const
Correlations, normalized per primitive cell.
const PrimClex & primclex() const
const Access the PrimClex that *this is based on
Definition: MonteCarlo.hh:68
std::pair< ConfigDoF, std::string > _restricted_auto_motif(const GrandCanonicalConditions &cond) const
Generate minimum potential energy ConfigDoF for this supercell.
Eigen::VectorXd VectorXd
PrimClex is the top-level data structure for a CASM project.
Definition: PrimClex.hh:52
bool m_use_deltas
If the supercell is large enough, calculate delta correlations directly.
Clexulator & _clexulator() const
ConfigIO::GenericConfigFormatter< std::string > configname()
Constructs DataFormmaterDictionary containing all Configuration DatumFormatters.
Definition: ConfigIO.cc:340
Eigen::VectorXd * m_comp_n
Number of atoms of each type, normalized per primitive cell.
Interface base class for all types of Monte Carlo simulations (not meant to be used polymorphically) ...
Definition: MonteCarlo.hh:32
std::pair< ConfigDoF, std::string > _auto_motif(const GrandCanonicalConditions &cond) const
Generate minimum potential energy ConfigDoF.
std::pair< ConfigDoF, std::string > set_state(const GrandCanonicalConditions &new_conditions, const GrandCanonicalSettings &settings)
Set configdof and conditions and clear previously collected data.
void accept(const EventType &event)
Accept proposed event. Change configuration accordingly and update energies etc.
const double & formation_energy() const
Formation energy, normalized per primitive cell.
ConfigDoF _default_motif() const
Generate supercell filling ConfigDoF from default configuration.
const CondType & conditions() const
Return current conditions.
static const Monte::ENSEMBLE ensemble
void _update_deltas(GrandCanonicalEvent &event, Index mutating_site, int sublat, int current_occupant, int new_occupant) const
Calculate delta properties for an event and update the event with those properties.
void _set_dCorr(GrandCanonicalEvent &event, Index mutating_site, int sublat, int current_occupant, int new_occupant, bool use_deltas, bool all_correlations) const
Calculate delta correlations for an event.
const ECIContainer & eci() const
Definition: Clex.cc:32
GrandCanonicalSettings SettingsType
Definition: Log.hh:9
Index steps_per_pass() const
Return number of steps per pass. Equals number of sites with variable occupation. ...
double * m_potential_energy
Potential energy, normalized per primitive cell.
ConfigDoF _configname_motif(const std::string &configname) const
Generate supercell filling ConfigDoF from configuration.
Clex m_formation_energy_clex
Holds Clexulator and ECI references.
Data structure used for cluster expansions.
Definition: Clex.hh:24
Eigen::VectorXd correlations_vec(const ConfigDoF &configdof, const Supercell &scel, Clexulator &clexulator)
Returns correlations using 'clexulator'. Supercell needs a correctly populated neighbor list...
Definition: ConfigDoF.cc:242
const MonteSettings & settings() const
const Access settings used for construction
Definition: MonteCarlo.hh:63
GrandCanonical(PrimClex &primclex, const SettingsType &settings, Log &_log)
Constructs a GrandCanonical object and prepares it for running based on MonteSettings.
GrandCanonicalConditions m_condition
Conditions (T, mu). Initially determined by m_settings, but can be changed halfway through the run...
A Configuration represents the values of all degrees of freedom in a Supercell.
A sparse container of ECI values and their corresponding orbit indices.
Definition: ECIContainer.hh:12