Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2013-2018 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 :
29 18 : void Steinhardt::registerKeywords( Keywords& keys ) {
30 18 : VectorMultiColvar::registerKeywords( keys );
31 18 : keys.add("compulsory","NN","12","The n parameter of the switching function ");
32 18 : keys.add("compulsory","MM","0","The m parameter of the switching function; 0 implies 2*NN");
33 18 : keys.add("compulsory","D_0","0.0","The d_0 parameter of the switching function");
34 18 : 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 swiching function defined above. "
36 : "The following provides information on the \\ref switchingfunction that are available. "
37 18 : "When this keyword is present you no longer need the NN, MM, D_0 and R_0 keywords.");
38 18 : keys.use("SPECIES"); keys.use("SPECIESA"); keys.use("SPECIESB");
39 18 : keys.use("MEAN"); keys.use("LESS_THAN"); keys.use("MORE_THAN"); keys.use("VMEAN");
40 18 : keys.use("BETWEEN"); keys.use("HISTOGRAM"); keys.use("MOMENTS"); keys.use("MIN"); keys.use("ALT_MIN");
41 18 : keys.use("LOWEST"); keys.use("HIGHEST");
42 18 : }
43 :
44 15 : Steinhardt::Steinhardt( const ActionOptions& ao ):
45 : Action(ao),
46 : VectorMultiColvar(ao),
47 15 : tmom(0)
48 : {
49 : // Read in the switching function
50 30 : std::string sw, errors; parse("SWITCH",sw);
51 15 : if(sw.length()>0) {
52 3 : switchingFunction.set(sw,errors);
53 : } else {
54 12 : double r_0=-1.0, d_0; int nn, mm;
55 12 : parse("NN",nn); parse("MM",mm);
56 12 : parse("R_0",r_0); parse("D_0",d_0);
57 12 : if( r_0<0.0 ) error("you must set a value for R_0");
58 12 : switchingFunction.set(nn,mm,r_0,d_0);
59 : }
60 15 : log.printf(" Steinhardt parameter of central atom and those within %s\n",( switchingFunction.description() ).c_str() );
61 15 : 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 15 : setLinkCellCutoff( switchingFunction.get_dmax() );
64 15 : rcut = switchingFunction.get_dmax(); rcut2 = rcut*rcut;
65 30 : std::vector<AtomNumber> all_atoms; setupMultiColvarBase( all_atoms );
66 15 : }
67 :
68 15 : void Steinhardt::setAngularMomentum( const unsigned& ang ) {
69 15 : tmom=ang; setVectorDimensionality( 2*(2*ang + 1) );
70 15 : }
71 :
72 66492 : void Steinhardt::calculateVector( multicolvar::AtomValuePack& myatoms ) const {
73 : double dfunc, dpoly_ass, md, tq6, itq6, real_z, imag_z;
74 66492 : Vector dz, myrealvec, myimagvec, real_dz, imag_dz;
75 : // The square root of -1
76 66499 : std::complex<double> ii( 0.0, 1.0 ), dp_x, dp_y, dp_z;
77 :
78 66500 : unsigned ncomp=2*tmom+1;
79 66500 : double sw, poly_ass, d2, dlen; std::complex<double> powered;
80 4513452 : for(unsigned i=1; i<myatoms.getNumberOfAtoms(); ++i) {
81 4447198 : Vector& distance=myatoms.getPosition(i); // getSeparation( myatoms.getPosition(0), myatoms.getPosition(i) );
82 21268504 : if ( (d2=distance[0]*distance[0])<rcut2 &&
83 14425051 : (d2+=distance[1]*distance[1])<rcut2 &&
84 17744070 : (d2+=distance[2]*distance[2])<rcut2 &&
85 2352349 : d2>epsilon ) {
86 :
87 2352351 : dlen = sqrt(d2);
88 2352351 : sw = switchingFunction.calculate( dlen, dfunc );
89 2352355 : accumulateSymmetryFunction( -1, i, sw, (+dfunc)*distance, (-dfunc)*Tensor( distance,distance ), myatoms );
90 2352345 : double dlen3 = d2*dlen;
91 : // Do stuff for m=0
92 2352345 : poly_ass=deriv_poly( 0, distance[2]/dlen, dpoly_ass );
93 : // Derivatives of z/r wrt x, y, z
94 2352347 : dz = -( distance[2] / dlen3 )*distance; dz[2] += (1.0 / dlen);
95 : // Derivative wrt to the vector connecting the two atoms
96 2352348 : myrealvec = (+sw)*dpoly_ass*dz + poly_ass*(+dfunc)*distance;
97 : // Accumulate the derivatives
98 2352345 : accumulateSymmetryFunction( 2 + tmom, i, sw*poly_ass, myrealvec, Tensor( -myrealvec,distance ), myatoms );
99 :
100 : // The complex number of which we have to take powers
101 2352346 : std::complex<double> com1( distance[0]/dlen,distance[1]/dlen );
102 :
103 : // Do stuff for all other m values
104 15044375 : for(unsigned m=1; m<=tmom; ++m) {
105 : // Calculate Legendre Polynomial
106 12692031 : poly_ass=deriv_poly( m, distance[2]/dlen, dpoly_ass );
107 : // Calculate powe of complex number
108 12692053 : powered=pow(com1,m-1); md=static_cast<double>(m);
109 : // Real and imaginary parts of z
110 12692048 : real_z = real(com1*powered); imag_z = imag(com1*powered );
111 :
112 : // Calculate steinhardt parameter
113 12692079 : tq6=poly_ass*real_z; // Real part of steinhardt parameter
114 12692079 : itq6=poly_ass*imag_z; // Imaginary part of steinhardt parameter
115 :
116 : // Derivatives wrt ( x/r + iy )^m
117 12692079 : dp_x = md*powered*( (1.0/dlen)-(distance[0]*distance[0])/dlen3-ii*(distance[0]*distance[1])/dlen3 );
118 12692008 : dp_y = md*powered*( ii*(1.0/dlen)-(distance[0]*distance[1])/dlen3-ii*(distance[1]*distance[1])/dlen3 );
119 12692066 : dp_z = md*powered*( -(distance[0]*distance[2])/dlen3-ii*(distance[1]*distance[2])/dlen3 );
120 :
121 : // Derivatives of real and imaginary parts of above
122 12692013 : real_dz[0] = real( dp_x ); real_dz[1] = real( dp_y ); real_dz[2] = real( dp_z );
123 12692024 : imag_dz[0] = imag( dp_x ); imag_dz[1] = imag( dp_y ); imag_dz[2] = imag( dp_z );
124 :
125 : // Complete derivative of steinhardt parameter
126 12692044 : myrealvec = (+sw)*dpoly_ass*real_z*dz + (+dfunc)*distance*tq6 + (+sw)*poly_ass*real_dz;
127 12692080 : myimagvec = (+sw)*dpoly_ass*imag_z*dz + (+dfunc)*distance*itq6 + (+sw)*poly_ass*imag_dz;
128 :
129 : // Real part
130 12692068 : accumulateSymmetryFunction( 2 + tmom + m, i, sw*tq6, myrealvec, Tensor( -myrealvec,distance ), myatoms );
131 : // Imaginary part
132 12692163 : accumulateSymmetryFunction( 2+ncomp+tmom+m, i, sw*itq6, myimagvec, Tensor( -myimagvec,distance ), myatoms );
133 : // Store -m part of vector
134 12692137 : double pref=pow(-1.0,m);
135 : // -m part of vector is just +m part multiplied by (-1.0)**m and multiplied by complex
136 : // conjugate of Legendre polynomial
137 : // Real part
138 12692137 : accumulateSymmetryFunction( 2+tmom-m, i, pref*sw*tq6, pref*myrealvec, pref*Tensor( -myrealvec,distance ), myatoms );
139 : // Imaginary part
140 12692030 : accumulateSymmetryFunction( 2+ncomp+tmom-m, i, -pref*sw*itq6, -pref*myimagvec, pref*Tensor( myimagvec,distance ), myatoms );
141 : }
142 : }
143 : }
144 :
145 : // Normalize
146 66499 : updateActiveAtoms( myatoms );
147 66500 : for(unsigned i=0; i<getNumberOfComponentsInVector(); ++i) myatoms.getUnderlyingMultiValue().quotientRule( 2+i, 2+i );
148 66500 : }
149 :
150 15044041 : double Steinhardt::deriv_poly( const unsigned& m, const double& val, double& df ) const {
151 15044041 : double fact=1.0;
152 15044041 : for(unsigned j=1; j<=m; ++j) fact=fact*j;
153 15044041 : double res=coeff_poly[m]*fact;
154 :
155 15044429 : double pow=1.0, xi=val, dxi=1.0; df=0.0;
156 56621348 : for(int i=m+1; i<=tmom; ++i) {
157 41577063 : double fact=1.0;
158 41577063 : for(unsigned j=i-m+1; j<=i; ++j) fact=fact*j;
159 41577063 : res=res+coeff_poly[i]*fact*xi;
160 41575289 : df = df + pow*coeff_poly[i]*fact*dxi;
161 41576919 : xi=xi*val; dxi=dxi*val; pow+=1.0;
162 : }
163 15044285 : df = df*normaliz[m];
164 15044195 : return normaliz[m]*res;
165 : }
166 :
167 : }
168 2523 : }
|