CASM
AClustersApproachtoStatisticalMechanics
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules
CASM::NiggliRep Class Reference

#include <Niggli.hh>

Detailed Description

Returns a values of A, B, C, ksi, eta and zeta. If the initialization cell has vector lenghts a, b and c, with angles alpha, beta, gamma, then the niggli values are

A=a*a B=b*b C=c*c ksi=2bc*cos(alpha) eta=2ac*cos(beta) zeta=2ab*cos(gamma)

Put another way, if the initialization lattice has vectors a, b and c then

A=aa B=bb C=cc ksi=2*bc eta=2*ac zeta=2*ab

See also
I. Krivy and B. Gruber, Acta Cryst. (1976). A32, 297. [doi:10.1107/S0567739476000636] R. W. Grosse-Kunstleve, N. K. Sauter and P. D. Adams, Acta Cryst. (2004). A60, 1. [doi:10.1107/S010876730302186X]

Definition at line 46 of file Niggli.hh.

Public Member Functions

 NiggliRep (const Lattice &init_lat)
 
 NiggliRep (const Eigen::Matrix3d &init_lat_col_mat)
 
double A () const
 Square of lattice length a. More...
 
double B () const
 Square of lattice length b. More...
 
double C () const
 Square of lattice length c. More...
 
double ksi () const
 2bc*cos(alpha) More...
 
double eta () const
 2ac*cos(beta) More...
 
double zeta () const
 2ab*cos(gamma) More...
 
const Eigen::Matrix3d & metrical_matrix () const
 
bool meets_criteria_1 (double compare_tol) const
 A<=B OR (A==B, |ksi| <= |eta|) More...
 
bool meets_criteria_2 (double compare_tol) const
 B<=C OR (B==C, |eta| <= |zeta|) More...
 
bool meets_criteria_3 (double compare_tol) const
 For type I: ksi>0 && eta>0 && zeta>0 (all angles < 90) More...
 
bool meets_criteria_4 (double compare_tol) const
 For type II: ksi<=0 && eta<=0 && zeta<=0 (all angles >= 90) More...
 
bool meets_criteria_5 (double compare_tol) const
 |ksi|<=B OR (ksi==B, zeta<=2*eta) OR (ksi==-B, zeta==0) More...
 
bool meets_criteria_6 (double compare_tol) const
 |eta|<=A OR (eta==A, zeta<=2*ksi) OR (eta==-A, zeta==0) More...
 
bool meets_criteria_7 (double compare_tol) const
 |zeta|<=A OR (zeta==A, eta<=2*ksi) OR (zeta==-A, eta==0) More...
 
bool meets_criteria_8 (double compare_tol) const
 
bool is_niggli (double compare_tol) const
 True if all conditions are true, and either 4 OR 3 is false. More...
 
bool is_niggli_type1 (double compare_tol) const
 True if all conditions except 4 are true. More...
 
bool is_niggli_type2 (double compare_tol) const
 True if all conditions except 3 are true. More...
 
void debug_criteria (double compare_tol) const
 

Private Attributes

const Eigen::Matrix3d m_metrical_matrix
 Transpose of initialization lattice dotted with itself. More...
 

Constructor & Destructor Documentation

CASM::NiggliRep::NiggliRep ( const Lattice init_lat)

Definition at line 11 of file Niggli.cc.

CASM::NiggliRep::NiggliRep ( const Eigen::Matrix3d &  init_lat_col_mat)

Definition at line 7 of file Niggli.cc.

Member Function Documentation

double CASM::NiggliRep::A ( ) const

Square of lattice length a.

Definition at line 15 of file Niggli.cc.

double CASM::NiggliRep::B ( ) const

Square of lattice length b.

Definition at line 19 of file Niggli.cc.

double CASM::NiggliRep::C ( ) const

Square of lattice length c.

Definition at line 23 of file Niggli.cc.

void CASM::NiggliRep::debug_criteria ( double  compare_tol) const

Definition at line 202 of file Niggli.cc.

double CASM::NiggliRep::eta ( ) const

2ac*cos(beta)

Definition at line 31 of file Niggli.cc.

bool CASM::NiggliRep::is_niggli ( double  compare_tol) const

True if all conditions are true, and either 4 OR 3 is false.

Definition at line 186 of file Niggli.cc.

bool CASM::NiggliRep::is_niggli_type1 ( double  compare_tol) const

True if all conditions except 4 are true.

Definition at line 178 of file Niggli.cc.

bool CASM::NiggliRep::is_niggli_type2 ( double  compare_tol) const

True if all conditions except 3 are true.

Definition at line 182 of file Niggli.cc.

double CASM::NiggliRep::ksi ( ) const

2bc*cos(alpha)

Definition at line 27 of file Niggli.cc.

bool CASM::NiggliRep::meets_criteria_1 ( double  compare_tol) const

A<=B OR (A==B, |ksi| <= |eta|)

Definition at line 43 of file Niggli.cc.

bool CASM::NiggliRep::meets_criteria_2 ( double  compare_tol) const

B<=C OR (B==C, |eta| <= |zeta|)

Definition at line 59 of file Niggli.cc.

bool CASM::NiggliRep::meets_criteria_3 ( double  compare_tol) const

For type I: ksi>0 && eta>0 && zeta>0 (all angles < 90)

Definition at line 75 of file Niggli.cc.

bool CASM::NiggliRep::meets_criteria_4 ( double  compare_tol) const

For type II: ksi<=0 && eta<=0 && zeta<=0 (all angles >= 90)

Definition at line 86 of file Niggli.cc.

bool CASM::NiggliRep::meets_criteria_5 ( double  compare_tol) const

|ksi|<=B OR (ksi==B, zeta<=2*eta) OR (ksi==-B, zeta==0)

Definition at line 97 of file Niggli.cc.

bool CASM::NiggliRep::meets_criteria_6 ( double  compare_tol) const

|eta|<=A OR (eta==A, zeta<=2*ksi) OR (eta==-A, zeta==0)

Definition at line 118 of file Niggli.cc.

bool CASM::NiggliRep::meets_criteria_7 ( double  compare_tol) const

|zeta|<=A OR (zeta==A, eta<=2*ksi) OR (zeta==-A, eta==0)

Definition at line 139 of file Niggli.cc.

bool CASM::NiggliRep::meets_criteria_8 ( double  compare_tol) const

ksi+eta+zeta+A+B<0 OR (ksi+eta+zeta+A+B==0, 2*A+2*eta+zeta>0) C<=A+B+C+ksi+eta+zeta OR (C==A+B+C+ksi+eta+zeta, 2*A+2*eta+zeta<=0)

Definition at line 160 of file Niggli.cc.

const Eigen::Matrix3d & CASM::NiggliRep::metrical_matrix ( ) const

Definition at line 39 of file Niggli.cc.

double CASM::NiggliRep::zeta ( ) const

2ab*cos(gamma)

Definition at line 35 of file Niggli.cc.

Member Data Documentation

const Eigen::Matrix3d CASM::NiggliRep::m_metrical_matrix
private

Transpose of initialization lattice dotted with itself.

Definition at line 112 of file Niggli.hh.


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