1 /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 Copyright (c) 2013-2021 The plumed team
3 (see the PEOPLE file at the root of the distribution for a list of names)
4
5 See http://www.plumed.org for more information.
6
7 This file is part of plumed, version 2.
8
9 plumed is free software: you can redistribute it and/or modify
10 it under the terms of the GNU Lesser General Public License as published by
11 the Free Software Foundation, either version 3 of the License, or
12 (at your option) any later version.
13
14 plumed is distributed in the hope that it will be useful,
15 but WITHOUT ANY WARRANTY; without even the implied warranty of
16 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
17 GNU Lesser General Public License for more details.
18
19 You should have received a copy of the GNU Lesser General Public License
20 along with plumed. If not, see <http://www.gnu.org/licenses/>.
21 +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
22 #include "Steinhardt.h"
23 #include "core/PlumedMain.h"
24 #include <complex>
25
26 namespace PLMD {
27 namespace crystallization {
28
registerKeywords(Keywords & keys)29 void Steinhardt::registerKeywords( Keywords& keys ) {
30 VectorMultiColvar::registerKeywords( keys );
31 keys.add("compulsory","NN","12","The n parameter of the switching function ");
32 keys.add("compulsory","MM","0","The m parameter of the switching function; 0 implies 2*NN");
33 keys.add("compulsory","D_0","0.0","The d_0 parameter of the switching function");
34 keys.add("compulsory","R_0","The r_0 parameter of the switching function");
35 keys.add("optional","SWITCH","This keyword is used if you want to employ an alternative to the continuous switching function defined above. "
36 "The following provides information on the \\ref switchingfunction that are available. "
37 "When this keyword is present you no longer need the NN, MM, D_0 and R_0 keywords.");
38 keys.use("SPECIES"); keys.use("SPECIESA"); keys.use("SPECIESB");
39 keys.use("MEAN"); keys.use("LESS_THAN"); keys.use("MORE_THAN"); keys.use("VMEAN");
40 keys.use("BETWEEN"); keys.use("HISTOGRAM"); keys.use("MOMENTS"); keys.use("MIN"); keys.use("ALT_MIN");
41 keys.use("LOWEST"); keys.use("HIGHEST");
42 }
43
Steinhardt(const ActionOptions & ao)44 Steinhardt::Steinhardt( const ActionOptions& ao ):
45 Action(ao),
46 VectorMultiColvar(ao),
47 tmom(0)
48 {
49 // Read in the switching function
50 std::string sw, errors; parse("SWITCH",sw);
51 if(sw.length()>0) {
52 switchingFunction.set(sw,errors);
53 } else {
54 double r_0=-1.0, d_0; int nn, mm;
55 parse("NN",nn); parse("MM",mm);
56 parse("R_0",r_0); parse("D_0",d_0);
57 if( r_0<0.0 ) error("you must set a value for R_0");
58 switchingFunction.set(nn,mm,r_0,d_0);
59 }
60 log.printf(" Steinhardt parameter of central atom and those within %s\n",( switchingFunction.description() ).c_str() );
61 log<<" Bibliography "<<plumed.cite("Tribello, Giberti, Sosso, Salvalaglio and Parrinello, J. Chem. Theory Comput. 13, 1317 (2017)")<<"\n";
62 // Set the link cell cutoff
63 setLinkCellCutoff( switchingFunction.get_dmax() );
64 rcut = switchingFunction.get_dmax(); rcut2 = rcut*rcut;
65 std::vector<AtomNumber> all_atoms; setupMultiColvarBase( all_atoms );
66 }
67
setAngularMomentum(const unsigned & ang)68 void Steinhardt::setAngularMomentum( const unsigned& ang ) {
69 tmom=ang; setVectorDimensionality( 2*(2*ang + 1) );
70 }
71
calculateVector(multicolvar::AtomValuePack & myatoms) const72 void Steinhardt::calculateVector( multicolvar::AtomValuePack& myatoms ) const {
73 double dfunc, dpoly_ass, md, tq6, itq6, real_z, imag_z;
74 Vector dz, myrealvec, myimagvec, real_dz, imag_dz;
75 // The square root of -1
76 std::complex<double> ii( 0.0, 1.0 ), dp_x, dp_y, dp_z;
77
78 unsigned ncomp=2*tmom+1;
79 double sw, poly_ass, dlen; std::complex<double> powered;
80 for(unsigned i=1; i<myatoms.getNumberOfAtoms(); ++i) {
81 Vector& distance=myatoms.getPosition(i); // getSeparation( myatoms.getPosition(0), myatoms.getPosition(i) );
82 double d2;
83 if ( (d2=distance[0]*distance[0])<rcut2 &&
84 (d2+=distance[1]*distance[1])<rcut2 &&
85 (d2+=distance[2]*distance[2])<rcut2 &&
86 d2>epsilon ) {
87
88 dlen = sqrt(d2);
89 sw = switchingFunction.calculate( dlen, dfunc );
90 accumulateSymmetryFunction( -1, i, sw, (+dfunc)*distance, (-dfunc)*Tensor( distance,distance ), myatoms );
91 double dlen3 = d2*dlen;
92 // Do stuff for m=0
93 poly_ass=deriv_poly( 0, distance[2]/dlen, dpoly_ass );
94 // Derivatives of z/r wrt x, y, z
95 dz = -( distance[2] / dlen3 )*distance; dz[2] += (1.0 / dlen);
96 // Derivative wrt to the vector connecting the two atoms
97 myrealvec = (+sw)*dpoly_ass*dz + poly_ass*(+dfunc)*distance;
98 // Accumulate the derivatives
99 accumulateSymmetryFunction( 2 + tmom, i, sw*poly_ass, myrealvec, Tensor( -myrealvec,distance ), myatoms );
100
101 // The complex number of which we have to take powers
102 std::complex<double> com1( distance[0]/dlen,distance[1]/dlen );
103 powered = std::complex<double>(1.0,0.0);
104
105 // Do stuff for all other m values
106 for(unsigned m=1; m<=tmom; ++m) {
107 // Calculate Legendre Polynomial
108 poly_ass=deriv_poly( m, distance[2]/dlen, dpoly_ass );
109 // Calculate power of complex number
110 // if(std::abs(com1)>epsilon) powered=pow(com1,m-1);
111 // else if(m==1) powered=std::complex<double>(1.,0);
112 // else powered = std::complex<double>(0.,0.);
113 // Real and imaginary parts of z
114 real_z = real(com1*powered); imag_z = imag(com1*powered );
115
116 // Calculate steinhardt parameter
117 tq6=poly_ass*real_z; // Real part of steinhardt parameter
118 itq6=poly_ass*imag_z; // Imaginary part of steinhardt parameter
119
120 // Derivatives wrt ( x/r + iy )^m
121 md=static_cast<double>(m);
122 dp_x = md*powered*( (1.0/dlen)-(distance[0]*distance[0])/dlen3-ii*(distance[0]*distance[1])/dlen3 );
123 dp_y = md*powered*( ii*(1.0/dlen)-(distance[0]*distance[1])/dlen3-ii*(distance[1]*distance[1])/dlen3 );
124 dp_z = md*powered*( -(distance[0]*distance[2])/dlen3-ii*(distance[1]*distance[2])/dlen3 );
125
126 // Derivatives of real and imaginary parts of above
127 real_dz[0] = real( dp_x ); real_dz[1] = real( dp_y ); real_dz[2] = real( dp_z );
128 imag_dz[0] = imag( dp_x ); imag_dz[1] = imag( dp_y ); imag_dz[2] = imag( dp_z );
129
130 // Complete derivative of steinhardt parameter
131 myrealvec = (+sw)*dpoly_ass*real_z*dz + (+dfunc)*distance*tq6 + (+sw)*poly_ass*real_dz;
132 myimagvec = (+sw)*dpoly_ass*imag_z*dz + (+dfunc)*distance*itq6 + (+sw)*poly_ass*imag_dz;
133
134 // Real part
135 accumulateSymmetryFunction( 2 + tmom + m, i, sw*tq6, myrealvec, Tensor( -myrealvec,distance ), myatoms );
136 // Imaginary part
137 accumulateSymmetryFunction( 2+ncomp+tmom+m, i, sw*itq6, myimagvec, Tensor( -myimagvec,distance ), myatoms );
138 // Store -m part of vector
139 double pref=pow(-1.0,m);
140 // -m part of vector is just +m part multiplied by (-1.0)**m and multiplied by complex
141 // conjugate of Legendre polynomial
142 // Real part
143 accumulateSymmetryFunction( 2+tmom-m, i, pref*sw*tq6, pref*myrealvec, pref*Tensor( -myrealvec,distance ), myatoms );
144 // Imaginary part
145 accumulateSymmetryFunction( 2+ncomp+tmom-m, i, -pref*sw*itq6, -pref*myimagvec, pref*Tensor( myimagvec,distance ), myatoms );
146 // Calculate next power of complex number
147 powered *= com1;
148 }
149 }
150 }
151
152 // Normalize
153 updateActiveAtoms( myatoms );
154 for(unsigned i=0; i<getNumberOfComponentsInVector(); ++i) myatoms.getUnderlyingMultiValue().quotientRule( 2+i, 2+i );
155 }
156
deriv_poly(const unsigned & m,const double & val,double & df) const157 double Steinhardt::deriv_poly( const unsigned& m, const double& val, double& df ) const {
158 double fact=1.0;
159 for(unsigned j=1; j<=m; ++j) fact=fact*j;
160 double res=coeff_poly[m]*fact;
161
162 double pow=1.0, xi=val, dxi=1.0; df=0.0;
163 for(int i=m+1; i<=tmom; ++i) {
164 double fact=1.0;
165 for(unsigned j=i-m+1; j<=i; ++j) fact=fact*j;
166 res=res+coeff_poly[i]*fact*xi;
167 df = df + pow*coeff_poly[i]*fact*dxi;
168 xi=xi*val; dxi=dxi*val; pow+=1.0;
169 }
170 df = df*normaliz[m];
171 return normaliz[m]*res;
172 }
173
174 }
175 }
176