Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2013-2023 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 30 : void Steinhardt::registerKeywords( Keywords& keys ) {
30 30 : VectorMultiColvar::registerKeywords( keys );
31 60 : keys.add("compulsory","NN","12","The n parameter of the switching function ");
32 60 : keys.add("compulsory","MM","0","The m parameter of the switching function; 0 implies 2*NN");
33 60 : keys.add("compulsory","D_0","0.0","The d_0 parameter of the switching function");
34 60 : keys.add("compulsory","R_0","The r_0 parameter of the switching function");
35 60 : 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 30 : keys.use("SPECIES");
39 30 : keys.use("SPECIESA");
40 30 : keys.use("SPECIESB");
41 30 : keys.use("MEAN");
42 30 : keys.use("LESS_THAN");
43 30 : keys.use("MORE_THAN");
44 30 : keys.use("VMEAN");
45 30 : keys.use("BETWEEN");
46 30 : keys.use("HISTOGRAM");
47 30 : keys.use("MOMENTS");
48 30 : keys.use("MIN");
49 30 : keys.use("ALT_MIN");
50 30 : keys.use("LOWEST");
51 30 : keys.use("HIGHEST");
52 30 : }
53 :
54 18 : Steinhardt::Steinhardt( const ActionOptions& ao ):
55 : Action(ao),
56 : VectorMultiColvar(ao),
57 18 : tmom(0) {
58 : // Read in the switching function
59 : std::string sw, errors;
60 36 : parse("SWITCH",sw);
61 18 : if(sw.length()>0) {
62 6 : switchingFunction.set(sw,errors);
63 : } else {
64 12 : double r_0=-1.0, d_0;
65 : int nn, mm;
66 12 : parse("NN",nn);
67 12 : parse("MM",mm);
68 12 : parse("R_0",r_0);
69 12 : parse("D_0",d_0);
70 12 : if( r_0<0.0 ) {
71 0 : error("you must set a value for R_0");
72 : }
73 12 : switchingFunction.set(nn,mm,r_0,d_0);
74 : }
75 18 : log.printf(" Steinhardt parameter of central atom and those within %s\n",( switchingFunction.description() ).c_str() );
76 36 : log<<" Bibliography "<<plumed.cite("Tribello, Giberti, Sosso, Salvalaglio and Parrinello, J. Chem. Theory Comput. 13, 1317 (2017)")<<"\n";
77 : // Set the link cell cutoff
78 18 : setLinkCellCutoff( switchingFunction.get_dmax() );
79 18 : rcut = switchingFunction.get_dmax();
80 18 : rcut2 = rcut*rcut;
81 : std::vector<AtomNumber> all_atoms;
82 18 : setupMultiColvarBase( all_atoms );
83 18 : }
84 :
85 18 : void Steinhardt::setAngularMomentum( const unsigned& ang ) {
86 18 : tmom=ang;
87 18 : setVectorDimensionality( 2*(2*ang + 1) );
88 18 : }
89 :
90 76317 : void Steinhardt::calculateVector( multicolvar::AtomValuePack& myatoms ) const {
91 : double dfunc, dpoly_ass, md, tq6, itq6, real_z, imag_z;
92 76317 : Vector dz, myrealvec, myimagvec, real_dz, imag_dz;
93 : // The square root of -1
94 : std::complex<double> ii( 0.0, 1.0 ), dp_x, dp_y, dp_z;
95 :
96 76317 : unsigned ncomp=2*tmom+1;
97 : double sw, poly_ass, dlen;
98 : std::complex<double> powered;
99 5297031 : for(unsigned i=1; i<myatoms.getNumberOfAtoms(); ++i) {
100 : Vector& distance=myatoms.getPosition(i); // getSeparation( myatoms.getPosition(0), myatoms.getPosition(i) );
101 : double d2;
102 9591512 : if ( (d2=distance[0]*distance[0])<rcut2 &&
103 4370798 : (d2+=distance[1]*distance[1])<rcut2 &&
104 8662510 : (d2+=distance[2]*distance[2])<rcut2 &&
105 : d2>epsilon ) {
106 :
107 2437950 : dlen = sqrt(d2);
108 2437950 : sw = switchingFunction.calculate( dlen, dfunc );
109 2437950 : accumulateSymmetryFunction( -1, i, sw, (+dfunc)*distance, (-dfunc)*Tensor( distance,distance ), myatoms );
110 2437950 : double dlen3 = d2*dlen;
111 : // Do stuff for m=0
112 2437950 : poly_ass=deriv_poly( 0, distance[2]/dlen, dpoly_ass );
113 : // Derivatives of z/r wrt x, y, z
114 2437950 : dz = -( distance[2] / dlen3 )*distance;
115 2437950 : dz[2] += (1.0 / dlen);
116 : // Derivative wrt to the vector connecting the two atoms
117 2437950 : myrealvec = (+sw)*dpoly_ass*dz + poly_ass*(+dfunc)*distance;
118 : // Accumulate the derivatives
119 2437950 : accumulateSymmetryFunction( 2 + tmom, i, sw*poly_ass, myrealvec, Tensor( -myrealvec,distance ), myatoms );
120 :
121 : // The complex number of which we have to take powers
122 2437950 : std::complex<double> com1( distance[0]/dlen,distance[1]/dlen );
123 : powered = std::complex<double>(1.0,0.0);
124 :
125 : // Do stuff for all other m values
126 15643770 : for(unsigned m=1; m<=tmom; ++m) {
127 : // Calculate Legendre Polynomial
128 13205820 : poly_ass=deriv_poly( m, distance[2]/dlen, dpoly_ass );
129 : // Calculate power of complex number
130 : // if(std::abs(com1)>epsilon) powered=pow(com1,m-1);
131 : // else if(m==1) powered=std::complex<double>(1.,0);
132 : // else powered = std::complex<double>(0.,0.);
133 : // Real and imaginary parts of z
134 : real_z = real(com1*powered);
135 : imag_z = imag(com1*powered );
136 :
137 : // Calculate steinhardt parameter
138 13205820 : tq6=poly_ass*real_z; // Real part of steinhardt parameter
139 13205820 : itq6=poly_ass*imag_z; // Imaginary part of steinhardt parameter
140 :
141 : // Derivatives wrt ( x/r + iy )^m
142 13205820 : md=static_cast<double>(m);
143 13205820 : dp_x = md*powered*( (1.0/dlen)-(distance[0]*distance[0])/dlen3-ii*(distance[0]*distance[1])/dlen3 );
144 13205820 : dp_y = md*powered*( ii*(1.0/dlen)-(distance[0]*distance[1])/dlen3-ii*(distance[1]*distance[1])/dlen3 );
145 13205820 : dp_z = md*powered*( -(distance[0]*distance[2])/dlen3-ii*(distance[1]*distance[2])/dlen3 );
146 :
147 : // Derivatives of real and imaginary parts of above
148 13205820 : real_dz[0] = real( dp_x );
149 13205820 : real_dz[1] = real( dp_y );
150 13205820 : real_dz[2] = real( dp_z );
151 13205820 : imag_dz[0] = imag( dp_x );
152 13205820 : imag_dz[1] = imag( dp_y );
153 13205820 : imag_dz[2] = imag( dp_z );
154 :
155 : // Complete derivative of steinhardt parameter
156 13205820 : myrealvec = (+sw)*dpoly_ass*real_z*dz + (+dfunc)*distance*tq6 + (+sw)*poly_ass*real_dz;
157 13205820 : myimagvec = (+sw)*dpoly_ass*imag_z*dz + (+dfunc)*distance*itq6 + (+sw)*poly_ass*imag_dz;
158 :
159 : // Real part
160 13205820 : accumulateSymmetryFunction( 2 + tmom + m, i, sw*tq6, myrealvec, Tensor( -myrealvec,distance ), myatoms );
161 : // Imaginary part
162 13205820 : accumulateSymmetryFunction( 2+ncomp+tmom+m, i, sw*itq6, myimagvec, Tensor( -myimagvec,distance ), myatoms );
163 : // Store -m part of vector
164 13205820 : double pref=pow(-1.0,m);
165 : // -m part of vector is just +m part multiplied by (-1.0)**m and multiplied by complex
166 : // conjugate of Legendre polynomial
167 : // Real part
168 13205820 : accumulateSymmetryFunction( 2+tmom-m, i, pref*sw*tq6, pref*myrealvec, pref*Tensor( -myrealvec,distance ), myatoms );
169 : // Imaginary part
170 13205820 : accumulateSymmetryFunction( 2+ncomp+tmom-m, i, -pref*sw*itq6, -pref*myimagvec, pref*Tensor( myimagvec,distance ), myatoms );
171 : // Calculate next power of complex number
172 : powered *= com1;
173 : }
174 : }
175 : }
176 :
177 : // Normalize
178 76317 : updateActiveAtoms( myatoms );
179 1929999 : for(unsigned i=0; i<getNumberOfComponentsInVector(); ++i) {
180 1853682 : myatoms.getUnderlyingMultiValue().quotientRule( 2+i, 2+i );
181 : }
182 76317 : }
183 :
184 15643770 : double Steinhardt::deriv_poly( const unsigned& m, const double& val, double& df ) const {
185 : double fact=1.0;
186 59020380 : for(unsigned j=1; j<=m; ++j) {
187 43376610 : fact=fact*j;
188 : }
189 15643770 : double res=coeff_poly[m]*fact;
190 :
191 15643770 : double pow=1.0, xi=val, dxi=1.0;
192 15643770 : df=0.0;
193 59020380 : for(int i=m+1; i<=tmom; ++i) {
194 : double fact=1.0;
195 110931360 : for(unsigned j=i-m+1; j<=i; ++j) {
196 67554750 : fact=fact*j;
197 : }
198 43376610 : res=res+coeff_poly[i]*fact*xi;
199 43376610 : df = df + pow*coeff_poly[i]*fact*dxi;
200 43376610 : xi=xi*val;
201 43376610 : dxi=dxi*val;
202 43376610 : pow+=1.0;
203 : }
204 15643770 : df = df*normaliz[m];
205 15643770 : return normaliz[m]*res;
206 : }
207 :
208 : }
209 : }
|