occ
Loading...
Searching...
No Matches
occ::mults::CrystalEnergy Class Reference

Crystal energy evaluator for rigid molecule assemblies. More...

#include <crystal_energy.h>

Classes

struct  MoleculeCache
 Cached per-molecule data to avoid redundant recomputation. More...
 
struct  MoleculeGeometry
 Cached molecule geometry (atom positions in body frame). More...
 

Public Member Functions

 CrystalEnergy (CrystalEnergySetup setup)
 Construct crystal energy evaluator.
 
 CrystalEnergy (const crystal::Crystal &crystal, std::vector< MultipoleSource > multipoles, double cutoff_radius=20.0, ForceFieldType ff=ForceFieldType::BuckinghamDE, bool use_cartesian=true, bool use_ewald=true, double ewald_accuracy=1e-6, double ewald_eta=0.0, int ewald_kmax=0, bool defer_setup=false)
 Legacy constructor from Crystal + MultipleSources.
 
 ~CrystalEnergy ()
 Destructor (defined in .cpp for unique_ptr with forward-declared type).
 
 CrystalEnergy (CrystalEnergy &&) noexcept
 Move operations (defined in .cpp where EwaldLatticeCache is complete).
 
CrystalEnergyoperator= (CrystalEnergy &&) noexcept
 
CrystalEnergyResult compute (const std::vector< MoleculeState > &molecules)
 Compute energy and gradient for given molecular states.
 
CrystalEnergyResultWithHessian compute_with_hessian (const std::vector< MoleculeState > &molecules)
 Compute energy, gradient, and Hessian for given molecular states.
 
bool can_compute_exact_hessian () const
 Returns true when compute_with_hessian provides an exact Hessian for the current configured energy model.
 
double compute_energy (const std::vector< MoleculeState > &molecules)
 Compute energy only (faster, for line search).
 
std::vector< MoleculeStateinitial_states () const
 Get initial states from crystal geometry.
 
void update_lattice (const crystal::Crystal &strained_crystal, std::vector< MoleculeState > new_states)
 Update the crystal lattice for strained evaluation.
 
int num_molecules () const
 Number of molecules in asymmetric unit.
 
const std::vector< MoleculeGeometry > & molecule_geometry () const
 Access molecule geometry (body-frame atom positions).
 
const std::vector< NeighborPair > & neighbor_pairs () const
 Get neighbor list.
 
ForceFieldParamsforce_field ()
 Access force field parameters.
 
const ForceFieldParamsforce_field () const
 
void set_typed_aniso_params (const std::map< std::pair< int, int >, AnisotropicRepulsionParams > &params)
 Set/query anisotropic repulsion parameters (typed).
 
bool has_aniso_params (int type1, int type2) const
 
AnisotropicRepulsionParams get_aniso_params (int type1, int type2) const
 
bool has_any_aniso_params () const
 
void set_buckingham_params (int Z1, int Z2, const BuckinghamParams &params)
 Get/set Buckingham parameters for atom pair.
 
void set_typed_buckingham_params (int type1, int type2, const BuckinghamParams &params)
 
void set_typed_buckingham_params (const std::map< std::pair< int, int >, BuckinghamParams > &params)
 
void clear_typed_buckingham_params ()
 
void set_short_range_type_labels (const std::map< int, std::string > &labels)
 
BuckinghamParams get_buckingham_params (int Z1, int Z2) const
 
bool has_buckingham_params (int Z1, int Z2) const
 
bool uses_williams_atom_typing () const
 
bool uses_short_range_typing () const
 
bool has_typed_buckingham_params (int type1, int type2) const
 
BuckinghamParams get_buckingham_params_for_types (int type1, int type2) const
 
std::string short_range_type_name (int type_code) const
 
const crystal::Crystalcrystal () const
 Get the underlying crystal.
 
double cutoff_radius () const
 
void set_cutoff_radius (double cutoff)
 
void update_neighbors ()
 Update neighbor list (called if molecules move significantly).
 
void update_neighbors (const std::vector< MoleculeState > &states)
 Update neighbor list using current molecule states.
 
void build_neighbor_list_from_positions (const std::vector< Vec3 > &mol_coms, bool force_com_cutoff=false, const std::vector< MoleculeState > *orientation_states=nullptr)
 Build neighbor list from explicit molecule COM positions.
 
void set_neighbor_list (const std::vector< NeighborPair > &neighbors)
 Set neighbor list directly (e.g. to reuse a list from another CrystalEnergy).
 
const std::vector< NeighborPair > & neighbor_list () const
 Get the neighbor list.
 
void set_molecule_geometry (std::vector< MoleculeGeometry > geometry)
 Set molecule geometry directly (bypasses Crystal's molecule detection).
 
void set_initial_states (std::vector< MoleculeState > states)
 Set initial molecule states directly.
 
void set_ewald_dipole (bool enable)
 Enable/disable dipole Ewald (charge-dipole + dipole-dipole).
 
bool use_ewald () const
 
void set_max_interaction_order (int max_order)
 Set maximum interaction order for electrostatics (rankA + rankB <= max).
 
int max_interaction_order () const
 
void set_elec_site_cutoff (double cutoff)
 Set per-site cutoff for electrostatic interactions (Angstrom).
 
double elec_site_cutoff () const
 
void set_use_com_elec_gate (bool enable)
 Enable COM-based gate for electrostatic interactions.
 
bool use_com_elec_gate () const
 
void set_buckingham_site_cutoff (double cutoff)
 Set a separate per-site cutoff for Buckingham (default: same as neighbor cutoff).
 
double buckingham_site_cutoff () const
 
void set_electrostatic_taper (double r_on, double r_off, int order=3)
 Set DMACRYS-style electrostatic spline taper.
 
void clear_electrostatic_taper ()
 
const CutoffSplineelectrostatic_taper () const
 
void set_electrostatic_taper_in_hessian (bool enable)
 When false, the electrostatic taper is applied to energy and gradient but NOT to the pair Hessian.
 
bool electrostatic_taper_in_hessian () const
 
void set_short_range_taper (double r_on, double r_off, int order=3)
 Set DMACRYS-style short-range spline taper.
 
void clear_short_range_taper ()
 
const CutoffSplineshort_range_taper () const
 
std::vector< std::vector< bool > > compute_buckingham_site_masks (const std::vector< MoleculeState > &states) const
 Compute which atom-atom pairs are within the Buckingham cutoff for each neighbor pair at the given molecular states.
 
void set_fixed_site_masks (std::vector< std::vector< bool > > masks)
 Set frozen site-pair masks.
 
void clear_fixed_site_masks ()
 Clear frozen site-pair masks (revert to distance-based cutoff).
 
int max_multipole_rank () const
 Get maximum multipole rank across all sites.
 
size_t num_neighbor_pairs () const
 Get number of neighbor pairs.
 
std::vector< PairEnergyDebugdebug_pair_energies (const std::vector< MoleculeState > &molecules)
 Compute per-pair debug breakdown at current states.
 
std::vector< int > neighbor_shell_histogram () const
 Neighbor shell histogram for quick sanity checks.
 
size_t num_sites () const
 Get total number of multipole sites.
 

Static Public Member Functions

static const char * short_range_type_label (int type_code)
 Convert a Williams/NEIGHCRYS type code to a short label (e.g. 512 -> C_W3).
 
static int short_range_type_atomic_number (int type_code)
 Map a Williams/NEIGHCRYS type code to element Z, or 0 if unknown.
 
static std::map< std::pair< int, int >, BuckinghamParamswilliams_de_params ()
 Williams DE Buckingham parameters (built-in).
 

Detailed Description

Crystal energy evaluator for rigid molecule assemblies.

Computes electrostatic energy using multipole interactions and short-range repulsion-dispersion using Buckingham or LJ potentials.

Usage:

Crystal crystal = ...; // Load from CIF
std::vector<MultipoleSource> multipoles = ...; // From DMA
CrystalEnergy energy_calc(crystal, multipoles, 20.0, ForceFieldType::BuckinghamDE);
std::vector<MoleculeState> states = energy_calc.initial_states();
auto result = energy_calc.compute(states);
Crystal energy evaluator for rigid molecule assemblies.
Definition crystal_energy.h:126
const crystal::Crystal & crystal() const
Get the underlying crystal.
Definition crystal_energy.h:247
@ BuckinghamDE
Williams DE Buckingham parameters.

Constructor & Destructor Documentation

◆ CrystalEnergy() [1/3]

occ::mults::CrystalEnergy::CrystalEnergy ( CrystalEnergySetup  setup)
explicit

Construct crystal energy evaluator.

Parameters
crystalCrystal structure (provides geometry and symmetry)
multipolesBody-frame multipoles for each symmetry-unique molecule
cutoff_radiusNeighbor cutoff in Angstrom (default 20.0)
ffForce field for short-range (default BuckinghamDE)
use_cartesianUse Cartesian T-tensor engine (default true) Construct from explicit molecule data. No Crystal molecule assembly. This is the preferred constructor.

◆ CrystalEnergy() [2/3]

occ::mults::CrystalEnergy::CrystalEnergy ( const crystal::Crystal crystal,
std::vector< MultipoleSource multipoles,
double  cutoff_radius = 20.0,
ForceFieldType  ff = ForceFieldType::BuckinghamDE,
bool  use_cartesian = true,
bool  use_ewald = true,
double  ewald_accuracy = 1e-6,
double  ewald_eta = 0.0,
int  ewald_kmax = 0,
bool  defer_setup = false 
)

Legacy constructor from Crystal + MultipleSources.

◆ ~CrystalEnergy()

occ::mults::CrystalEnergy::~CrystalEnergy ( )

Destructor (defined in .cpp for unique_ptr with forward-declared type).

◆ CrystalEnergy() [3/3]

occ::mults::CrystalEnergy::CrystalEnergy ( CrystalEnergy &&  )
noexcept

Move operations (defined in .cpp where EwaldLatticeCache is complete).

Member Function Documentation

◆ buckingham_site_cutoff()

double occ::mults::CrystalEnergy::buckingham_site_cutoff ( ) const
inline

◆ build_neighbor_list_from_positions()

void occ::mults::CrystalEnergy::build_neighbor_list_from_positions ( const std::vector< Vec3 > &  mol_coms,
bool  force_com_cutoff = false,
const std::vector< MoleculeState > *  orientation_states = nullptr 
)

Build neighbor list from explicit molecule COM positions.

Bypasses Crystal's molecule detection (useful when loading external data). When force_com_cutoff=true, uses COM distance for the molecule pair gate even when atom geometry is available (matches DMACRYS TBLCNT behavior).

◆ can_compute_exact_hessian()

bool occ::mults::CrystalEnergy::can_compute_exact_hessian ( ) const

Returns true when compute_with_hessian provides an exact Hessian for the current configured energy model.

◆ clear_electrostatic_taper()

void occ::mults::CrystalEnergy::clear_electrostatic_taper ( )
inline

◆ clear_fixed_site_masks()

void occ::mults::CrystalEnergy::clear_fixed_site_masks ( )
inline

Clear frozen site-pair masks (revert to distance-based cutoff).

◆ clear_short_range_taper()

void occ::mults::CrystalEnergy::clear_short_range_taper ( )
inline

◆ clear_typed_buckingham_params()

void occ::mults::CrystalEnergy::clear_typed_buckingham_params ( )
inline

◆ compute()

CrystalEnergyResult occ::mults::CrystalEnergy::compute ( const std::vector< MoleculeState > &  molecules)

Compute energy and gradient for given molecular states.

◆ compute_buckingham_site_masks()

std::vector< std::vector< bool > > occ::mults::CrystalEnergy::compute_buckingham_site_masks ( const std::vector< MoleculeState > &  states) const

Compute which atom-atom pairs are within the Buckingham cutoff for each neighbor pair at the given molecular states.

Returns a per-neighbor-pair boolean mask (row-major, nA*nB) indicating inclusion.

◆ compute_energy()

double occ::mults::CrystalEnergy::compute_energy ( const std::vector< MoleculeState > &  molecules)

Compute energy only (faster, for line search).

◆ compute_with_hessian()

CrystalEnergyResultWithHessian occ::mults::CrystalEnergy::compute_with_hessian ( const std::vector< MoleculeState > &  molecules)

Compute energy, gradient, and Hessian for given molecular states.

Hessian is assembled analytically from pair terms (short-range full rigid-body and full Cartesian multipole electrostatics). Ewald site-position Hessian terms are included when enabled; exactness for the full model is reported by exact_for_model.

◆ crystal()

const crystal::Crystal & occ::mults::CrystalEnergy::crystal ( ) const
inline

Get the underlying crystal.

◆ cutoff_radius()

double occ::mults::CrystalEnergy::cutoff_radius ( ) const
inline

◆ debug_pair_energies()

std::vector< PairEnergyDebug > occ::mults::CrystalEnergy::debug_pair_energies ( const std::vector< MoleculeState > &  molecules)

Compute per-pair debug breakdown at current states.

◆ elec_site_cutoff()

double occ::mults::CrystalEnergy::elec_site_cutoff ( ) const
inline

◆ electrostatic_taper()

const CutoffSpline & occ::mults::CrystalEnergy::electrostatic_taper ( ) const
inline

◆ electrostatic_taper_in_hessian()

bool occ::mults::CrystalEnergy::electrostatic_taper_in_hessian ( ) const
inline

◆ force_field() [1/2]

ForceFieldParams & occ::mults::CrystalEnergy::force_field ( )
inline

Access force field parameters.

◆ force_field() [2/2]

const ForceFieldParams & occ::mults::CrystalEnergy::force_field ( ) const
inline

◆ get_aniso_params()

AnisotropicRepulsionParams occ::mults::CrystalEnergy::get_aniso_params ( int  type1,
int  type2 
) const
inline

◆ get_buckingham_params()

BuckinghamParams occ::mults::CrystalEnergy::get_buckingham_params ( int  Z1,
int  Z2 
) const
inline

◆ get_buckingham_params_for_types()

BuckinghamParams occ::mults::CrystalEnergy::get_buckingham_params_for_types ( int  type1,
int  type2 
) const
inline

◆ has_aniso_params()

bool occ::mults::CrystalEnergy::has_aniso_params ( int  type1,
int  type2 
) const
inline

◆ has_any_aniso_params()

bool occ::mults::CrystalEnergy::has_any_aniso_params ( ) const
inline

◆ has_buckingham_params()

bool occ::mults::CrystalEnergy::has_buckingham_params ( int  Z1,
int  Z2 
) const
inline

◆ has_typed_buckingham_params()

bool occ::mults::CrystalEnergy::has_typed_buckingham_params ( int  type1,
int  type2 
) const
inline

◆ initial_states()

std::vector< MoleculeState > occ::mults::CrystalEnergy::initial_states ( ) const

Get initial states from crystal geometry.

◆ max_interaction_order()

int occ::mults::CrystalEnergy::max_interaction_order ( ) const
inline

◆ max_multipole_rank()

int occ::mults::CrystalEnergy::max_multipole_rank ( ) const

Get maximum multipole rank across all sites.

◆ molecule_geometry()

const std::vector< MoleculeGeometry > & occ::mults::CrystalEnergy::molecule_geometry ( ) const
inline

Access molecule geometry (body-frame atom positions).

◆ neighbor_list()

const std::vector< NeighborPair > & occ::mults::CrystalEnergy::neighbor_list ( ) const
inline

Get the neighbor list.

◆ neighbor_pairs()

const std::vector< NeighborPair > & occ::mults::CrystalEnergy::neighbor_pairs ( ) const
inline

Get neighbor list.

◆ neighbor_shell_histogram()

std::vector< int > occ::mults::CrystalEnergy::neighbor_shell_histogram ( ) const

Neighbor shell histogram for quick sanity checks.

◆ num_molecules()

int occ::mults::CrystalEnergy::num_molecules ( ) const
inline

Number of molecules in asymmetric unit.

◆ num_neighbor_pairs()

size_t occ::mults::CrystalEnergy::num_neighbor_pairs ( ) const
inline

Get number of neighbor pairs.

◆ num_sites()

size_t occ::mults::CrystalEnergy::num_sites ( ) const

Get total number of multipole sites.

◆ operator=()

CrystalEnergy & occ::mults::CrystalEnergy::operator= ( CrystalEnergy &&  )
noexcept

◆ set_buckingham_params()

void occ::mults::CrystalEnergy::set_buckingham_params ( int  Z1,
int  Z2,
const BuckinghamParams params 
)
inline

Get/set Buckingham parameters for atom pair.

◆ set_buckingham_site_cutoff()

void occ::mults::CrystalEnergy::set_buckingham_site_cutoff ( double  cutoff)
inline

Set a separate per-site cutoff for Buckingham (default: same as neighbor cutoff).

Use a larger value than the neighbor cutoff to reduce derivative discontinuities near the neighbor boundary.

◆ set_cutoff_radius()

void occ::mults::CrystalEnergy::set_cutoff_radius ( double  cutoff)

◆ set_elec_site_cutoff()

void occ::mults::CrystalEnergy::set_elec_site_cutoff ( double  cutoff)
inline

Set per-site cutoff for electrostatic interactions (Angstrom).

Default 0.0 means no per-site cutoff (all site pairs within included molecule pairs are computed). DMACRYS applies RANG2 per-site cutoff in its PAIR module for higher multipole terms.

◆ set_electrostatic_taper()

void occ::mults::CrystalEnergy::set_electrostatic_taper ( double  r_on,
double  r_off,
int  order = 3 
)

Set DMACRYS-style electrostatic spline taper.

f=1 for r<=r_on, spline to 0 over (r_on,r_off], zero beyond r_off. Order must be 3 (cubic) or 5 (quintic).

◆ set_electrostatic_taper_in_hessian()

void occ::mults::CrystalEnergy::set_electrostatic_taper_in_hessian ( bool  enable)
inline

When false, the electrostatic taper is applied to energy and gradient but NOT to the pair Hessian.

This matches DMACRYS behaviour where the multipole SEC array does not include spline chain-rule terms.

◆ set_ewald_dipole()

void occ::mults::CrystalEnergy::set_ewald_dipole ( bool  enable)
inline

Enable/disable dipole Ewald (charge-dipole + dipole-dipole).

Default is true (DMACRYS computes Ewald for qq, qμ, and μμ).

◆ set_fixed_site_masks()

void occ::mults::CrystalEnergy::set_fixed_site_masks ( std::vector< std::vector< bool > >  masks)
inline

Set frozen site-pair masks.

When non-empty, compute_short_range_pair uses only the atom pairs marked true, bypassing the distance cutoff. This ensures a smoother energy surface under cell/molecular perturbations.

◆ set_initial_states()

void occ::mults::CrystalEnergy::set_initial_states ( std::vector< MoleculeState states)

Set initial molecule states directly.

◆ set_max_interaction_order()

void occ::mults::CrystalEnergy::set_max_interaction_order ( int  max_order)
inline

Set maximum interaction order for electrostatics (rankA + rankB <= max).

Default -1 means no truncation (compute all orders up to 2*max_rank). Set to 4 to match DMACRYS truncation.

◆ set_molecule_geometry()

void occ::mults::CrystalEnergy::set_molecule_geometry ( std::vector< MoleculeGeometry geometry)

Set molecule geometry directly (bypasses Crystal's molecule detection).

◆ set_neighbor_list()

void occ::mults::CrystalEnergy::set_neighbor_list ( const std::vector< NeighborPair > &  neighbors)

Set neighbor list directly (e.g. to reuse a list from another CrystalEnergy).

◆ set_short_range_taper()

void occ::mults::CrystalEnergy::set_short_range_taper ( double  r_on,
double  r_off,
int  order = 3 
)

Set DMACRYS-style short-range spline taper.

◆ set_short_range_type_labels()

void occ::mults::CrystalEnergy::set_short_range_type_labels ( const std::map< int, std::string > &  labels)
inline

◆ set_typed_aniso_params()

void occ::mults::CrystalEnergy::set_typed_aniso_params ( const std::map< std::pair< int, int >, AnisotropicRepulsionParams > &  params)
inline

Set/query anisotropic repulsion parameters (typed).

◆ set_typed_buckingham_params() [1/2]

void occ::mults::CrystalEnergy::set_typed_buckingham_params ( const std::map< std::pair< int, int >, BuckinghamParams > &  params)
inline

◆ set_typed_buckingham_params() [2/2]

void occ::mults::CrystalEnergy::set_typed_buckingham_params ( int  type1,
int  type2,
const BuckinghamParams params 
)
inline

◆ set_use_com_elec_gate()

void occ::mults::CrystalEnergy::set_use_com_elec_gate ( bool  enable)
inline

Enable COM-based gate for electrostatic interactions.

When enabled, only molecule pairs with COM distance <= cutoff_radius contribute to electrostatics (matching DMACRYS TBLCNT behavior). Buckingham is unaffected — uses full atom-based neighbor list. Default is true for DMACRYS compatibility.

◆ short_range_taper()

const CutoffSpline & occ::mults::CrystalEnergy::short_range_taper ( ) const
inline

◆ short_range_type_atomic_number()

static int occ::mults::CrystalEnergy::short_range_type_atomic_number ( int  type_code)
inlinestatic

Map a Williams/NEIGHCRYS type code to element Z, or 0 if unknown.

◆ short_range_type_label()

static const char * occ::mults::CrystalEnergy::short_range_type_label ( int  type_code)
inlinestatic

Convert a Williams/NEIGHCRYS type code to a short label (e.g. 512 -> C_W3).

◆ short_range_type_name()

std::string occ::mults::CrystalEnergy::short_range_type_name ( int  type_code) const
inline

◆ update_lattice()

void occ::mults::CrystalEnergy::update_lattice ( const crystal::Crystal strained_crystal,
std::vector< MoleculeState new_states 
)

Update the crystal lattice for strained evaluation.

Keeps: neighbors, force field, geometry (body-frame), multipoles. Updates: crystal (for cell translations + Ewald), initial states.

◆ update_neighbors() [1/2]

void occ::mults::CrystalEnergy::update_neighbors ( )

Update neighbor list (called if molecules move significantly).

◆ update_neighbors() [2/2]

void occ::mults::CrystalEnergy::update_neighbors ( const std::vector< MoleculeState > &  states)

Update neighbor list using current molecule states.

For explicit neighbor mode (all UC molecules), this rebuilds the list from current COM positions and orientations.

◆ use_com_elec_gate()

bool occ::mults::CrystalEnergy::use_com_elec_gate ( ) const
inline

◆ use_ewald()

bool occ::mults::CrystalEnergy::use_ewald ( ) const
inline

◆ uses_short_range_typing()

bool occ::mults::CrystalEnergy::uses_short_range_typing ( ) const
inline

◆ uses_williams_atom_typing()

bool occ::mults::CrystalEnergy::uses_williams_atom_typing ( ) const
inline

◆ williams_de_params()

static std::map< std::pair< int, int >, BuckinghamParams > occ::mults::CrystalEnergy::williams_de_params ( )
inlinestatic

Williams DE Buckingham parameters (built-in).


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