Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ 2 : Copyright (c) 2015-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 "core/ActionRegister.h" 23 : #include "tools/Pbc.h" 24 : #include "tools/SwitchingFunction.h" 25 : #include "ActionVolume.h" 26 : 27 : //+PLUMEDOC VOLUMES INSPHERE 28 : /* 29 : This quantity can be used to calculate functions of the distribution of collective variables for the atoms that lie in a particular, user-specified part of of the cell. 30 : 31 : Each of the base quantities calculated by a multicolvar can can be assigned to a particular point in three 32 : dimensional space. For example, if we have the coordination numbers for all the atoms in the 33 : system each coordination number can be assumed to lie on the position of the central atom. 34 : Because each base quantity can be assigned to a particular point in space we can calculate functions of the 35 : distribution of base quantities in a particular part of the box by using: 36 : 37 : \f[ 38 : \overline{s}_{\tau} = \frac{ \sum_i f(s_i) \sigma(r) }{ \sum_i \sigma(r) } 39 : \f] 40 : 41 : where the sum is over the collective variables, \f$s_i\f$, each of which can be thought to be at \f$ (x_i,y_i,z_i)\f$. 42 : The function \f$\sigma\f$ is a \ref switchingfunction that acts on the distance between the point at which the 43 : collective is located \f$(x_i,y_i,z_i)\f$ and the position of the atom that was specified using the ORIGIN keyword. 44 : In other words: 45 : \f[ 46 : r = sqrt{ ( x_i - x_0)^2 + ( y_i - y_0)^2 + ( z_i - z_0)^2} 47 : \f] 48 : In short this function, \f$\sigma(r_{xy})\f$, measures whether or not the CV is within a sphere that is 49 : centered on the position of the atom specified using the keyword ORIGIN. 50 : 51 : 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. 52 : 53 : When INCYLINDER is used with the \ref DENSITY action the number of atoms in the specified region is calculated 54 : 55 : \par Examples 56 : 57 : The input below can be use to calculate the average coordination numbers for those atoms that are within a sphere 58 : of radius 1.5 nm that is centered on the position of atom 101. 59 : 60 : \plumedfile 61 : c1: COORDINATIONNUMBER SPECIES=1-100 SWITCH={RATIONAL R_0=0.1} 62 : d2: INSPHERE ATOM=101 DATA=c1 RADIUS={TANH R_0=1.5} MEAN 63 : PRINT ARG=d2.* FILE=colvar 64 : \endplumedfile 65 : 66 : */ 67 : //+ENDPLUMEDOC 68 : 69 : namespace PLMD { 70 : namespace multicolvar { 71 : 72 : class VolumeInSphere : public ActionVolume { 73 : private: 74 : Vector origin; 75 : SwitchingFunction switchingFunction; 76 : public: 77 : static void registerKeywords( Keywords& keys ); 78 : explicit VolumeInSphere(const ActionOptions& ao); 79 : void setupRegions() override; 80 : double calculateNumberInside( const Vector& cpos, Vector& derivatives, Tensor& vir, std::vector<Vector>& refders ) const override; 81 : }; 82 : 83 10425 : PLUMED_REGISTER_ACTION(VolumeInSphere,"INSPHERE") 84 : 85 4 : void VolumeInSphere::registerKeywords( Keywords& keys ) { 86 4 : ActionVolume::registerKeywords( keys ); 87 8 : keys.add("atoms","ATOM","the atom whose vicinity we are interested in examining"); 88 8 : keys.add("compulsory","RADIUS","the switching function that tells us the extent of the spherical region of interest"); 89 4 : keys.remove("SIGMA"); 90 4 : } 91 : 92 3 : VolumeInSphere::VolumeInSphere(const ActionOptions& ao): 93 : Action(ao), 94 3 : ActionVolume(ao) 95 : { 96 : std::vector<AtomNumber> atom; 97 6 : parseAtomList("ATOM",atom); 98 3 : if( atom.size()!=1 ) error("should only be one atom specified"); 99 3 : log.printf(" center of sphere is at position of atom : %d\n",atom[0].serial() ); 100 : 101 6 : std::string sw, errors; parse("RADIUS",sw); 102 3 : if(sw.length()==0) error("missing RADIUS keyword"); 103 3 : switchingFunction.set(sw,errors); 104 3 : if( errors.length()!=0 ) error("problem reading RADIUS keyword : " + errors ); 105 3 : log.printf(" radius of sphere is given by %s \n", ( switchingFunction.description() ).c_str() ); 106 : 107 3 : checkRead(); requestAtoms(atom); 108 3 : } 109 : 110 82 : void VolumeInSphere::setupRegions() { } 111 : 112 155474 : double VolumeInSphere::calculateNumberInside( const Vector& cpos, Vector& derivatives, Tensor& vir, std::vector<Vector>& refders ) const { 113 : // Calculate position of atom wrt to origin 114 155474 : Vector fpos=pbcDistance( getPosition(0), cpos ); 115 155474 : double dfunc, value = switchingFunction.calculateSqr( fpos.modulo2(), dfunc ); 116 155474 : derivatives.zero(); derivatives = dfunc*fpos; refders[0] = -derivatives; 117 : // Add a virial contribution 118 155474 : vir -= Tensor(fpos,derivatives); 119 155474 : return value; 120 : } 121 : 122 : } 123 : }