Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2014-2020 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 "core/ActionRegister.h"
23 : #include "core/PlumedMain.h"
24 : #include "tools/Units.h"
25 : #include "tools/Pbc.h"
26 : #include "ActionVolume.h"
27 : #include "VolumeShortcut.h"
28 :
29 : //+PLUMEDOC VOLUMES CAVITY
30 : /*
31 : This quantity can be used to calculate functions of the distribution of collective variables for the atoms that lie in a box defined by the positions of four atoms.
32 :
33 : Each of the base quantities calculated by a multicolvar can can be assigned to a particular point in three
34 : dimensional space. For example, if we have the coordination numbers for all the atoms in the
35 : system each coordination number can be assumed to lie on the position of the central atom.
36 : Because each base quantity can be assigned to a particular point in space we can calculate functions of the
37 : distribution of base quantities in a particular part of the box by using:
38 :
39 : \f[
40 : \overline{s}_{\tau} = \frac{ \sum_i f(s_i) w(u_i,v_i,w_i) }{ \sum_i w(u_i,v_i,w_i) }
41 : \f]
42 :
43 : where the sum is over the collective variables, \f$s_i\f$, each of which can be thought to be at \f$ (u_i,v_i,z_i)\f$.
44 : The function \f$(s_i)\f$ can be any of the usual LESS_THAN, MORE_THAN, WITHIN etc that are used in all other multicolvars.
45 : Notice that here (at variance with what is done in \ref AROUND) we have transformed from the usual \f$(x_i,y_i,z_i)\f$
46 : position to a position in \f$ (u_i,v_i,z_i)\f$. This is done using a rotation matrix as follows:
47 :
48 : \f[
49 : \left(
50 : \begin{matrix}
51 : u_i \\
52 : v_i \\
53 : w_i
54 : \end{matrix}
55 : \right) = \mathbf{R}
56 : \left(
57 : \begin{matrix}
58 : x_i - x_o \\
59 : y_i - y_o \\
60 : z_i - z_o
61 : \end{matrix}
62 : \right)
63 : \f]
64 :
65 : where \f$\mathbf{R}\f$ is a rotation matrix that is calculated by constructing a set of three orthonormal vectors from the
66 : reference positions specified by the user. The first of these unit vectors points from the first reference atom to the second.
67 : The second is then the normal to the plane containing atoms 1,2 and 3 and the the third is the unit vector orthogonal to
68 : these first two vectors. \f$(x_o,y_o,z_o)\f$, meanwhile, specifies the position of the first reference atom.
69 :
70 : In the previous function \f$ w(u_i,v_i,w_i) \f$ measures whether or not the system is in the subregion of interest. It
71 : is equal to:
72 :
73 : \f[
74 : w(u_i,v_i,w_i) = \int_{0}^{u'} \int_{0}^{v'} \int_{0}^{w'} \textrm{d}u\textrm{d}v\textrm{d}w
75 : K\left( \frac{u - u_i}{\sigma} \right)K\left( \frac{v - v_i}{\sigma} \right)K\left( \frac{w - w_i}{\sigma} \right)
76 : \f]
77 :
78 : where \f$K\f$ is one of the kernel functions described on \ref histogrambead and \f$\sigma\f$ is a bandwidth parameter.
79 : The vector connecting atom 1 to atom 4 is used to define the extent of the box in each of the \f$u\f$, \f$v\f$ and \f$w\f$
80 : directions. Essentially the vector connecting atom 1 to atom 4 is projected onto the three unit vectors
81 : described above and the resulting projections determine the \f$u'\f$, \f$v'\f$ and \f$w'\f$ parameters in the above expression.
82 :
83 : \par Examples
84 :
85 : The following commands tell plumed to calculate the number of atoms in an ion channel in a protein.
86 : The extent of the channel is calculated from the positions of atoms 1, 4, 5 and 11. The final value will be labeled cav.
87 :
88 : \plumedfile
89 : d1: DENSITY SPECIES=20-500
90 : CAVITY DATA=d1 ATOMS=1,4,5,11 SIGMA=0.1 LABEL=cav
91 : \endplumedfile
92 :
93 : The following command tells plumed to calculate the coordination numbers (with other water molecules) for the water
94 : molecules in the protein channel described above. The average coordination number and the number of coordination
95 : numbers more than 4 is then calculated. The values of these two quantities are given the labels cav.mean and cav.morethan
96 :
97 : \plumedfile
98 : d1: COORDINATIONNUMBER SPECIES=20-500 R_0=0.1
99 : CAVITY DATA=d1 ATOMS=1,4,5,11 SIGMA=0.1 MEAN MORE_THAN={RATIONAL R_0=4} LABEL=cav
100 : \endplumedfile
101 :
102 : */
103 : //+ENDPLUMEDOC
104 :
105 : //+PLUMEDOC MCOLVAR CAVITY_CALC
106 : /*
107 : Calculate a vector from the input positions with elements equal to one when the positions are in a particular part of the cell and elements equal to zero otherwise
108 :
109 : \par Examples
110 :
111 : */
112 : //+ENDPLUMEDOC
113 :
114 : namespace PLMD {
115 : namespace volumes {
116 :
117 : class VolumeCavity : public ActionVolume {
118 : private:
119 : bool boxout;
120 : OFile boxfile;
121 : double lenunit;
122 : double jacob_det;
123 : double len_bi, len_cross, len_perp, sigma;
124 : Vector origin, bi, cross, perp;
125 : std::vector<Vector> dlbi, dlcross, dlperp;
126 : std::vector<Tensor> dbi, dcross, dperp;
127 : public:
128 : static void registerKeywords( Keywords& keys );
129 : explicit VolumeCavity(const ActionOptions& ao);
130 : ~VolumeCavity();
131 : void setupRegions() override;
132 : void update() override;
133 : double calculateNumberInside( const Vector& cpos, Vector& derivatives, Tensor& vir, std::vector<Vector>& refders ) const override;
134 : };
135 :
136 : PLUMED_REGISTER_ACTION(VolumeCavity,"CAVITY_CALC")
137 : char glob_cavity[] = "CAVITY";
138 : typedef VolumeShortcut<glob_cavity> VolumeCavityShortcut;
139 : PLUMED_REGISTER_ACTION(VolumeCavityShortcut,"CAVITY")
140 :
141 7 : void VolumeCavity::registerKeywords( Keywords& keys ) {
142 7 : ActionVolume::registerKeywords( keys ); keys.setDisplayName("CAVITY");
143 14 : keys.add("atoms","BOX","the positions of four atoms that define spatial extent of the cavity");
144 14 : keys.addFlag("PRINT_BOX",false,"write out the positions of the corners of the box to an xyz file");
145 14 : keys.add("optional","FILE","the file on which to write out the box coordinates");
146 14 : keys.add("optional","UNITS","( default=nm ) the units in which to write out the corners of the box");
147 7 : }
148 :
149 1 : VolumeCavity::VolumeCavity(const ActionOptions& ao):
150 : Action(ao),
151 : ActionVolume(ao),
152 1 : boxout(false),
153 1 : lenunit(1.0),
154 1 : dlbi(4),
155 1 : dlcross(4),
156 1 : dlperp(4),
157 1 : dbi(3),
158 1 : dcross(3),
159 2 : dperp(3)
160 : {
161 2 : std::vector<AtomNumber> atoms; parseAtomList("BOX",atoms);
162 1 : if( atoms.size()!=4 ) error("number of atoms in box should be equal to four");
163 :
164 1 : log.printf(" boundaries for region are calculated based on positions of atoms : ");
165 5 : for(unsigned i=0; i<atoms.size(); ++i) log.printf("%d ",atoms[i].serial() );
166 1 : log.printf("\n"); requestAtoms( atoms );
167 :
168 1 : boxout=false; parseFlag("PRINT_BOX",boxout);
169 1 : if(boxout) {
170 0 : std::string boxfname; parse("FILE",boxfname);
171 0 : if(boxfname.length()==0) error("no name for box file specified");
172 0 : std::string unitname; parse("UNITS",unitname);
173 0 : if ( unitname.length()>0 ) {
174 0 : Units u; u.setLength(unitname);
175 0 : lenunit=getUnits().getLength()/u.getLength();
176 0 : } else {
177 : unitname="nm";
178 : }
179 0 : boxfile.link(*this);
180 0 : boxfile.open( boxfname );
181 0 : log.printf(" printing box coordinates on file named %s in %s \n",boxfname.c_str(), unitname.c_str() );
182 : }
183 :
184 1 : checkRead();
185 1 : }
186 :
187 2 : VolumeCavity::~VolumeCavity() {
188 2 : }
189 :
190 1560 : void VolumeCavity::setupRegions() {
191 : // Make some space for things
192 1560 : Vector d1, d2, d3;
193 :
194 : // Retrieve the sigma value
195 1560 : sigma=getSigma();
196 : // Set the position of the origin
197 1560 : origin=getPosition(0);
198 :
199 : // Get two vectors
200 1560 : d1 = pbcDistance(origin,getPosition(1));
201 1560 : double d1l=d1.modulo();
202 1560 : d2 = pbcDistance(origin,getPosition(2));
203 :
204 : // Find the vector connecting the origin to the top corner of
205 : // the subregion
206 1560 : d3 = pbcDistance(origin,getPosition(3));
207 :
208 : // Create a set of unit vectors
209 1560 : bi = d1 / d1l; len_bi=dotProduct( d3, bi );
210 1560 : cross = crossProduct( d1, d2 ); double crossmod=cross.modulo();
211 1560 : cross = cross / crossmod; len_cross=dotProduct( d3, cross );
212 1560 : perp = crossProduct( cross, bi ); len_perp=dotProduct( d3, perp );
213 :
214 : // Calculate derivatives of box shape with respect to atoms
215 1560 : double d1l3=d1l*d1l*d1l;
216 1560 : dbi[0](0,0) = ( -(d1[1]*d1[1]+d1[2]*d1[2])/d1l3 ); // dx/dx
217 1560 : dbi[0](0,1) = ( d1[0]*d1[1]/d1l3 ); // dx/dy
218 1560 : dbi[0](0,2) = ( d1[0]*d1[2]/d1l3 ); // dx/dz
219 1560 : dbi[0](1,0) = ( d1[1]*d1[0]/d1l3 ); // dy/dx
220 1560 : dbi[0](1,1) = ( -(d1[0]*d1[0]+d1[2]*d1[2])/d1l3 ); // dy/dy
221 1560 : dbi[0](1,2) = ( d1[1]*d1[2]/d1l3 );
222 1560 : dbi[0](2,0) = ( d1[2]*d1[0]/d1l3 );
223 1560 : dbi[0](2,1) = ( d1[2]*d1[1]/d1l3 );
224 1560 : dbi[0](2,2) = ( -(d1[1]*d1[1]+d1[0]*d1[0])/d1l3 );
225 :
226 1560 : dbi[1](0,0) = ( (d1[1]*d1[1]+d1[2]*d1[2])/d1l3 );
227 1560 : dbi[1](0,1) = ( -d1[0]*d1[1]/d1l3 );
228 1560 : dbi[1](0,2) = ( -d1[0]*d1[2]/d1l3 );
229 1560 : dbi[1](1,0) = ( -d1[1]*d1[0]/d1l3 );
230 1560 : dbi[1](1,1) = ( (d1[0]*d1[0]+d1[2]*d1[2])/d1l3 );
231 1560 : dbi[1](1,2) = ( -d1[1]*d1[2]/d1l3 );
232 1560 : dbi[1](2,0) = ( -d1[2]*d1[0]/d1l3 );
233 1560 : dbi[1](2,1) = ( -d1[2]*d1[1]/d1l3 );
234 1560 : dbi[1](2,2) = ( (d1[1]*d1[1]+d1[0]*d1[0])/d1l3 );
235 1560 : dbi[2].zero();
236 :
237 1560 : Tensor tcderiv; double cmod3=crossmod*crossmod*crossmod; Vector ucross=crossmod*cross;
238 1560 : tcderiv.setCol( 0, crossProduct( d1, Vector(-1.0,0.0,0.0) ) + crossProduct( Vector(-1.0,0.0,0.0), d2 ) );
239 1560 : tcderiv.setCol( 1, crossProduct( d1, Vector(0.0,-1.0,0.0) ) + crossProduct( Vector(0.0,-1.0,0.0), d2 ) );
240 1560 : tcderiv.setCol( 2, crossProduct( d1, Vector(0.0,0.0,-1.0) ) + crossProduct( Vector(0.0,0.0,-1.0), d2 ) );
241 1560 : dcross[0](0,0)=( tcderiv(0,0)/crossmod - ucross[0]*(ucross[0]*tcderiv(0,0) + ucross[1]*tcderiv(1,0) + ucross[2]*tcderiv(2,0))/cmod3 ); // dx/dx
242 1560 : dcross[0](0,1)=( tcderiv(0,1)/crossmod - ucross[0]*(ucross[0]*tcderiv(0,1) + ucross[1]*tcderiv(1,1) + ucross[2]*tcderiv(2,1))/cmod3 ); // dx/dy
243 1560 : dcross[0](0,2)=( tcderiv(0,2)/crossmod - ucross[0]*(ucross[0]*tcderiv(0,2) + ucross[1]*tcderiv(1,2) + ucross[2]*tcderiv(2,2))/cmod3 ); // dx/dz
244 1560 : dcross[0](1,0)=( tcderiv(1,0)/crossmod - ucross[1]*(ucross[0]*tcderiv(0,0) + ucross[1]*tcderiv(1,0) + ucross[2]*tcderiv(2,0))/cmod3 ); // dy/dx
245 1560 : dcross[0](1,1)=( tcderiv(1,1)/crossmod - ucross[1]*(ucross[0]*tcderiv(0,1) + ucross[1]*tcderiv(1,1) + ucross[2]*tcderiv(2,1))/cmod3 ); // dy/dy
246 1560 : dcross[0](1,2)=( tcderiv(1,2)/crossmod - ucross[1]*(ucross[0]*tcderiv(0,2) + ucross[1]*tcderiv(1,2) + ucross[2]*tcderiv(2,2))/cmod3 ); // dy/dz
247 1560 : dcross[0](2,0)=( tcderiv(2,0)/crossmod - ucross[2]*(ucross[0]*tcderiv(0,0) + ucross[1]*tcderiv(1,0) + ucross[2]*tcderiv(2,0))/cmod3 ); // dz/dx
248 1560 : dcross[0](2,1)=( tcderiv(2,1)/crossmod - ucross[2]*(ucross[0]*tcderiv(0,1) + ucross[1]*tcderiv(1,1) + ucross[2]*tcderiv(2,1))/cmod3 ); // dz/dy
249 1560 : dcross[0](2,2)=( tcderiv(2,2)/crossmod - ucross[2]*(ucross[0]*tcderiv(0,2) + ucross[1]*tcderiv(1,2) + ucross[2]*tcderiv(2,2))/cmod3 ); // dz/dz
250 :
251 1560 : tcderiv.setCol( 0, crossProduct( Vector(1.0,0.0,0.0), d2 ) );
252 1560 : tcderiv.setCol( 1, crossProduct( Vector(0.0,1.0,0.0), d2 ) );
253 1560 : tcderiv.setCol( 2, crossProduct( Vector(0.0,0.0,1.0), d2 ) );
254 1560 : dcross[1](0,0)=( tcderiv(0,0)/crossmod - ucross[0]*(ucross[0]*tcderiv(0,0) + ucross[1]*tcderiv(1,0) + ucross[2]*tcderiv(2,0))/cmod3 ); // dx/dx
255 1560 : dcross[1](0,1)=( tcderiv(0,1)/crossmod - ucross[0]*(ucross[0]*tcderiv(0,1) + ucross[1]*tcderiv(1,1) + ucross[2]*tcderiv(2,1))/cmod3 ); // dx/dy
256 1560 : dcross[1](0,2)=( tcderiv(0,2)/crossmod - ucross[0]*(ucross[0]*tcderiv(0,2) + ucross[1]*tcderiv(1,2) + ucross[2]*tcderiv(2,2))/cmod3 ); // dx/dz
257 1560 : dcross[1](1,0)=( tcderiv(1,0)/crossmod - ucross[1]*(ucross[0]*tcderiv(0,0) + ucross[1]*tcderiv(1,0) + ucross[2]*tcderiv(2,0))/cmod3 ); // dy/dx
258 1560 : dcross[1](1,1)=( tcderiv(1,1)/crossmod - ucross[1]*(ucross[0]*tcderiv(0,1) + ucross[1]*tcderiv(1,1) + ucross[2]*tcderiv(2,1))/cmod3 ); // dy/dy
259 1560 : dcross[1](1,2)=( tcderiv(1,2)/crossmod - ucross[1]*(ucross[0]*tcderiv(0,2) + ucross[1]*tcderiv(1,2) + ucross[2]*tcderiv(2,2))/cmod3 ); // dy/dz
260 1560 : dcross[1](2,0)=( tcderiv(2,0)/crossmod - ucross[2]*(ucross[0]*tcderiv(0,0) + ucross[1]*tcderiv(1,0) + ucross[2]*tcderiv(2,0))/cmod3 ); // dz/dx
261 1560 : dcross[1](2,1)=( tcderiv(2,1)/crossmod - ucross[2]*(ucross[0]*tcderiv(0,1) + ucross[1]*tcderiv(1,1) + ucross[2]*tcderiv(2,1))/cmod3 ); // dz/dy
262 1560 : dcross[1](2,2)=( tcderiv(2,2)/crossmod - ucross[2]*(ucross[0]*tcderiv(0,2) + ucross[1]*tcderiv(1,2) + ucross[2]*tcderiv(2,2))/cmod3 ); // dz/dz
263 :
264 1560 : tcderiv.setCol( 0, crossProduct( d1, Vector(1.0,0.0,0.0) ) );
265 1560 : tcderiv.setCol( 1, crossProduct( d1, Vector(0.0,1.0,0.0) ) );
266 1560 : tcderiv.setCol( 2, crossProduct( d1, Vector(0.0,0.0,1.0) ) );
267 1560 : dcross[2](0,0)=( tcderiv(0,0)/crossmod - ucross[0]*(ucross[0]*tcderiv(0,0) + ucross[1]*tcderiv(1,0) + ucross[2]*tcderiv(2,0))/cmod3 ); // dx/dx
268 1560 : dcross[2](0,1)=( tcderiv(0,1)/crossmod - ucross[0]*(ucross[0]*tcderiv(0,1) + ucross[1]*tcderiv(1,1) + ucross[2]*tcderiv(2,1))/cmod3 ); // dx/dy
269 1560 : dcross[2](0,2)=( tcderiv(0,2)/crossmod - ucross[0]*(ucross[0]*tcderiv(0,2) + ucross[1]*tcderiv(1,2) + ucross[2]*tcderiv(2,2))/cmod3 ); // dx/dz
270 1560 : dcross[2](1,0)=( tcderiv(1,0)/crossmod - ucross[1]*(ucross[0]*tcderiv(0,0) + ucross[1]*tcderiv(1,0) + ucross[2]*tcderiv(2,0))/cmod3 ); // dy/dx
271 1560 : dcross[2](1,1)=( tcderiv(1,1)/crossmod - ucross[1]*(ucross[0]*tcderiv(0,1) + ucross[1]*tcderiv(1,1) + ucross[2]*tcderiv(2,1))/cmod3 ); // dy/dy
272 1560 : dcross[2](1,2)=( tcderiv(1,2)/crossmod - ucross[1]*(ucross[0]*tcderiv(0,2) + ucross[1]*tcderiv(1,2) + ucross[2]*tcderiv(2,2))/cmod3 ); // dy/dz
273 1560 : dcross[2](2,0)=( tcderiv(2,0)/crossmod - ucross[2]*(ucross[0]*tcderiv(0,0) + ucross[1]*tcderiv(1,0) + ucross[2]*tcderiv(2,0))/cmod3 ); // dz/dx
274 1560 : dcross[2](2,1)=( tcderiv(2,1)/crossmod - ucross[2]*(ucross[0]*tcderiv(0,1) + ucross[1]*tcderiv(1,1) + ucross[2]*tcderiv(2,1))/cmod3 ); // dz/dy
275 1560 : dcross[2](2,2)=( tcderiv(2,2)/crossmod - ucross[2]*(ucross[0]*tcderiv(0,2) + ucross[1]*tcderiv(1,2) + ucross[2]*tcderiv(2,2))/cmod3 ); // dz/dz
276 :
277 1560 : dperp[0].setCol( 0, ( crossProduct( dcross[0].getCol(0), bi ) + crossProduct( cross, dbi[0].getCol(0) ) ) );
278 1560 : dperp[0].setCol( 1, ( crossProduct( dcross[0].getCol(1), bi ) + crossProduct( cross, dbi[0].getCol(1) ) ) );
279 1560 : dperp[0].setCol( 2, ( crossProduct( dcross[0].getCol(2), bi ) + crossProduct( cross, dbi[0].getCol(2) ) ) );
280 :
281 1560 : dperp[1].setCol( 0, ( crossProduct( dcross[1].getCol(0), bi ) + crossProduct( cross, dbi[1].getCol(0) ) ) );
282 1560 : dperp[1].setCol( 1, ( crossProduct( dcross[1].getCol(1), bi ) + crossProduct( cross, dbi[1].getCol(1) ) ) );
283 1560 : dperp[1].setCol( 2, ( crossProduct( dcross[1].getCol(2), bi ) + crossProduct( cross, dbi[1].getCol(2) ) ) );
284 :
285 1560 : dperp[2].setCol( 0, ( crossProduct( dcross[2].getCol(0), bi ) ) );
286 1560 : dperp[2].setCol( 1, ( crossProduct( dcross[2].getCol(1), bi ) ) );
287 1560 : dperp[2].setCol( 2, ( crossProduct( dcross[2].getCol(2), bi ) ) );
288 :
289 : // Ensure that all lengths are positive
290 1560 : if( len_bi<0 ) {
291 0 : bi=-bi; len_bi=-len_bi;
292 0 : for(unsigned i=0; i<3; ++i) dbi[i]*=-1.0;
293 : }
294 1560 : if( len_cross<0 ) {
295 0 : cross=-cross; len_cross=-len_cross;
296 0 : for(unsigned i=0; i<3; ++i) dcross[i]*=-1.0;
297 : }
298 1560 : if( len_perp<0 ) {
299 0 : perp=-perp; len_perp=-len_perp;
300 0 : for(unsigned i=0; i<3; ++i) dperp[i]*=-1.0;
301 : }
302 1560 : if( len_bi<=0 || len_cross<=0 || len_bi<=0 ) plumed_merror("Invalid box coordinates");
303 :
304 : // Now derivatives of lengths
305 1560 : Tensor dd3( Tensor::identity() );
306 1560 : dlbi[0] = matmul(d3,dbi[0]) - matmul(bi,dd3);
307 1560 : dlbi[1] = matmul(d3,dbi[1]);
308 1560 : dlbi[2] = matmul(d3,dbi[2]);
309 1560 : dlbi[3] = matmul(bi,dd3);
310 :
311 1560 : dlcross[0] = matmul(d3,dcross[0]) - matmul(cross,dd3);
312 1560 : dlcross[1] = matmul(d3,dcross[1]);
313 1560 : dlcross[2] = matmul(d3,dcross[2]);
314 1560 : dlcross[3] = matmul(cross,dd3);
315 :
316 1560 : dlperp[0] = matmul(d3,dperp[0]) - matmul(perp,dd3);
317 1560 : dlperp[1] = matmul(d3,dperp[1]);
318 1560 : dlperp[2] = matmul(d3,dperp[2]);
319 1560 : dlperp[3] = matmul(perp,dd3);
320 :
321 : // Need to calculate the jacobian
322 1560 : Tensor jacob;
323 1560 : jacob(0,0)=bi[0]; jacob(1,0)=bi[1]; jacob(2,0)=bi[2];
324 1560 : jacob(0,1)=cross[0]; jacob(1,1)=cross[1]; jacob(2,1)=cross[2];
325 1560 : jacob(0,2)=perp[0]; jacob(1,2)=perp[1]; jacob(2,2)=perp[2];
326 1560 : jacob_det = fabs( jacob.determinant() );
327 1560 : }
328 :
329 60 : void VolumeCavity::update() {
330 60 : if(boxout) {
331 0 : boxfile.printf("%d\n",8);
332 0 : const Tensor & t(getPbc().getBox());
333 0 : if(getPbc().isOrthorombic()) {
334 0 : boxfile.printf(" %f %f %f\n",lenunit*t(0,0),lenunit*t(1,1),lenunit*t(2,2));
335 : } else {
336 0 : boxfile.printf(" %f %f %f %f %f %f %f %f %f\n",
337 0 : lenunit*t(0,0),lenunit*t(0,1),lenunit*t(0,2),
338 0 : lenunit*t(1,0),lenunit*t(1,1),lenunit*t(1,2),
339 0 : lenunit*t(2,0),lenunit*t(2,1),lenunit*t(2,2)
340 : );
341 : }
342 0 : boxfile.printf("AR %f %f %f \n",lenunit*origin[0],lenunit*origin[1],lenunit*origin[2]);
343 0 : Vector ut, vt, wt;
344 0 : ut = origin + len_bi*bi;
345 0 : vt = origin + len_cross*cross;
346 0 : wt = origin + len_perp*perp;
347 0 : boxfile.printf("AR %f %f %f \n",lenunit*(ut[0]), lenunit*(ut[1]), lenunit*(ut[2]) );
348 0 : boxfile.printf("AR %f %f %f \n",lenunit*(vt[0]), lenunit*(vt[1]), lenunit*(vt[2]) );
349 0 : boxfile.printf("AR %f %f %f \n",lenunit*(wt[0]), lenunit*(wt[1]), lenunit*(wt[2]) );
350 0 : boxfile.printf("AR %f %f %f \n",lenunit*(vt[0]+len_bi*bi[0]),
351 0 : lenunit*(vt[1]+len_bi*bi[1]),
352 0 : lenunit*(vt[2]+len_bi*bi[2]) );
353 0 : boxfile.printf("AR %f %f %f \n",lenunit*(ut[0]+len_perp*perp[0]),
354 0 : lenunit*(ut[1]+len_perp*perp[1]),
355 0 : lenunit*(ut[2]+len_perp*perp[2]) );
356 0 : boxfile.printf("AR %f %f %f \n",lenunit*(vt[0]+len_perp*perp[0]),
357 0 : lenunit*(vt[1]+len_perp*perp[1]),
358 0 : lenunit*(vt[2]+len_perp*perp[2]) );
359 0 : boxfile.printf("AR %f %f %f \n",lenunit*(vt[0]+len_perp*perp[0]+len_bi*bi[0]),
360 0 : lenunit*(vt[1]+len_perp*perp[1]+len_bi*bi[1]),
361 0 : lenunit*(vt[2]+len_perp*perp[2]+len_bi*bi[2]) );
362 : }
363 60 : }
364 :
365 1560 : double VolumeCavity::calculateNumberInside( const Vector& cpos, Vector& derivatives, Tensor& vir, std::vector<Vector>& rderiv ) const {
366 : // Setup the histogram bead
367 3120 : HistogramBead bead; bead.isNotPeriodic(); bead.setKernelType( getKernelType() );
368 :
369 : // Calculate distance of atom from origin of new coordinate frame
370 1560 : Vector datom=pbcDistance( origin, cpos );
371 : double ucontr, uder, vcontr, vder, wcontr, wder;
372 :
373 : // Calculate contribution from integral along bi
374 1560 : bead.set( 0, len_bi, sigma );
375 1560 : double upos=dotProduct( datom, bi );
376 1560 : ucontr=bead.calculate( upos, uder );
377 1560 : double udlen=bead.uboundDerivative( upos );
378 1560 : double uder2 = bead.lboundDerivative( upos ) - udlen;
379 :
380 : // Calculate contribution from integral along cross
381 1560 : bead.set( 0, len_cross, sigma );
382 1560 : double vpos=dotProduct( datom, cross );
383 1560 : vcontr=bead.calculate( vpos, vder );
384 1560 : double vdlen=bead.uboundDerivative( vpos );
385 1560 : double vder2 = bead.lboundDerivative( vpos ) - vdlen;
386 :
387 : // Calculate contribution from integral along perp
388 1560 : bead.set( 0, len_perp, sigma );
389 1560 : double wpos=dotProduct( datom, perp );
390 1560 : wcontr=bead.calculate( wpos, wder );
391 1560 : double wdlen=bead.uboundDerivative( wpos );
392 1560 : double wder2 = bead.lboundDerivative( wpos ) - wdlen;
393 :
394 1560 : Vector dfd; dfd[0]=uder*vcontr*wcontr; dfd[1]=ucontr*vder*wcontr; dfd[2]=ucontr*vcontr*wder;
395 1560 : derivatives[0] = (dfd[0]*bi[0]+dfd[1]*cross[0]+dfd[2]*perp[0]);
396 1560 : derivatives[1] = (dfd[0]*bi[1]+dfd[1]*cross[1]+dfd[2]*perp[1]);
397 1560 : derivatives[2] = (dfd[0]*bi[2]+dfd[1]*cross[2]+dfd[2]*perp[2]);
398 1560 : double tot = ucontr*vcontr*wcontr*jacob_det;
399 :
400 : // Add reference atom derivatives
401 1560 : dfd[0]=uder2*vcontr*wcontr; dfd[1]=ucontr*vder2*wcontr; dfd[2]=ucontr*vcontr*wder2;
402 1560 : Vector dfld; dfld[0]=udlen*vcontr*wcontr; dfld[1]=ucontr*vdlen*wcontr; dfld[2]=ucontr*vcontr*wdlen;
403 1560 : rderiv[0] = dfd[0]*matmul(datom,dbi[0]) + dfd[1]*matmul(datom,dcross[0]) + dfd[2]*matmul(datom,dperp[0]) +
404 3120 : dfld[0]*dlbi[0] + dfld[1]*dlcross[0] + dfld[2]*dlperp[0] - derivatives;
405 1560 : rderiv[1] = dfd[0]*matmul(datom,dbi[1]) + dfd[1]*matmul(datom,dcross[1]) + dfd[2]*matmul(datom,dperp[1]) +
406 3120 : dfld[0]*dlbi[1] + dfld[1]*dlcross[1] + dfld[2]*dlperp[1];
407 1560 : rderiv[2] = dfd[0]*matmul(datom,dbi[2]) + dfd[1]*matmul(datom,dcross[2]) + dfd[2]*matmul(datom,dperp[2]) +
408 3120 : dfld[0]*dlbi[2] + dfld[1]*dlcross[2] + dfld[2]*dlperp[2];
409 1560 : rderiv[3] = dfld[0]*dlbi[3] + dfld[1]*dlcross[3] + dfld[2]*dlperp[3];
410 :
411 1560 : vir.zero(); vir-=Tensor( cpos,derivatives );
412 7800 : for(unsigned i=0; i<4; ++i) {
413 6240 : vir -= Tensor( getPosition(i), rderiv[i] );
414 : }
415 :
416 1560 : return tot;
417 : }
418 :
419 : }
420 : }
|