Skip to content
Snippets Groups Projects
Commit af5fc7ce authored by Gareth Tribello's avatar Gareth Tribello
Browse files

Added multicolvar actions to calculate the angles/torsions between bonds and the cell axes

parent c02f3ea0
No related branches found
No related tags found
No related merge requests found
include ../../scripts/test.make
#! FIELDS time dx.lessthan dy.between
0.000000 1.4487 0.7038
0.050000 1.4590 0.7042
0.100000 1.4215 0.9236
0.150000 1.3815 0.9333
0.200000 1.3622 0.9406
type=driver
# this is to test a different name
arg="--plumed plumed.dat --trajectory-stride 10 --timestep 0.005 --ixyz trajectory.xyz --dump-forces forces --dump-forces-fmt=%8.4f"
extra_files="../../trajectories/trajectory.xyz"
#! FIELDS time parameter dx.lessthan dy.between
0.000000 0 -0.0531 0.0021
0.000000 1 -0.0009 0.1631
0.000000 2 0.0602 0.0000
0.000000 3 0.0531 -0.0021
0.000000 4 0.0009 -0.1631
0.000000 5 -0.0602 0.0000
0.000000 6 -0.8154 -0.1562
0.000000 7 0.0183 0.0068
0.000000 8 -0.0178 0.0000
0.000000 9 0.8154 0.1562
0.000000 10 -0.0183 -0.0068
0.000000 11 0.0178 0.0000
0.000000 12 -0.0806 -0.0038
0.000000 13 -0.0002 0.1547
0.000000 14 0.0563 0.0000
0.000000 15 -0.6937 -0.1331
0.000000 16 0.0156 0.0038
0.000000 17 -0.0159 0.0000
0.000000 18 0.6294 0.1308
0.000000 19 -0.0159 0.1306
0.000000 20 0.0650 0.0000
0.050000 0 -0.0357 0.0012
0.050000 1 -0.0004 0.1507
0.050000 2 0.0442 0.0000
0.050000 3 0.0357 -0.0012
0.050000 4 0.0004 -0.1507
0.050000 5 -0.0442 0.0000
0.050000 6 -0.8378 -0.1609
0.050000 7 0.0227 0.0086
0.050000 8 -0.0223 0.0000
0.050000 9 0.8378 0.1609
0.050000 10 -0.0227 -0.0086
0.050000 11 0.0223 0.0000
0.050000 12 -0.0735 -0.0059
0.050000 13 0.0006 0.1549
0.050000 14 0.0443 0.0000
0.050000 15 -0.6936 -0.1333
0.050000 16 0.0188 0.0059
0.050000 17 -0.0188 0.0000
0.050000 18 0.6526 0.1320
0.050000 19 -0.0188 0.1178
0.050000 20 0.0547 0.0000
0.100000 0 -0.0278 0.0011
0.100000 1 -0.0004 0.1401
0.100000 2 0.0361 0.0000
0.100000 3 0.0278 -0.0011
0.100000 4 0.0004 -0.1401
0.100000 5 -0.0361 0.0000
0.100000 6 -0.8591 0.2132
0.100000 7 -0.0017 0.0009
0.100000 8 0.0018 0.0000
0.100000 9 0.8591 -0.2132
0.100000 10 0.0017 -0.0009
0.100000 11 -0.0018 0.0000
0.100000 12 -0.0280 0.0005
0.100000 13 -0.0004 0.1545
0.100000 14 0.0398 0.0000
0.100000 15 -0.6276 0.1558
0.100000 16 -0.0012 -0.0005
0.100000 17 0.0010 0.0000
0.100000 18 0.6675 -0.1706
0.100000 19 0.0010 0.1184
0.100000 20 0.0292 0.0000
0.150000 0 -0.0233 -0.0006
0.150000 1 0.0002 0.1311
0.150000 2 0.0311 0.0000
0.150000 3 0.0233 0.0006
0.150000 4 -0.0002 -0.1311
0.150000 5 -0.0311 0.0000
0.150000 6 -0.8413 0.2435
0.150000 7 -0.0244 0.0193
0.150000 8 0.0321 0.0000
0.150000 9 0.8413 -0.2435
0.150000 10 0.0244 -0.0193
0.150000 11 -0.0321 0.0000
0.150000 12 0.0146 -0.0128
0.150000 13 0.0014 0.1538
0.150000 14 0.0351 0.0000
0.150000 15 -0.5308 0.1536
0.150000 16 -0.0154 0.0128
0.150000 17 0.0204 0.0000
0.150000 18 0.6788 -0.2029
0.150000 19 0.0204 0.1000
0.150000 20 0.0008 0.0000
0.200000 0 -0.0182 -0.0004
0.200000 1 0.0001 0.1273
0.200000 2 0.0256 0.0000
0.200000 3 0.0182 0.0004
0.200000 4 -0.0001 -0.1273
0.200000 5 -0.0256 0.0000
0.200000 6 -0.7791 0.2441
0.200000 7 -0.0321 0.0307
0.200000 8 0.0460 0.0000
0.200000 9 0.7791 -0.2441
0.200000 10 0.0321 -0.0307
0.200000 11 -0.0460 0.0000
0.200000 12 0.0386 -0.0195
0.200000 13 0.0026 0.1524
0.200000 14 0.0276 0.0000
0.200000 15 -0.4828 0.1513
0.200000 16 -0.0199 0.0195
0.200000 17 0.0286 0.0000
0.200000 18 0.6769 -0.2173
0.200000 19 0.0286 0.0826
0.200000 20 -0.0187 0.0000
dx: XANGLES ATOMS1=1,2 ATOMS2=5,4 LESS_THAN={RATIONAL R_0=1.5}
#dxn: XANGLES ATOMS1=1,2 ATOMS2=5,4 LESS_THAN={RATIONAL R_0=1.5} NUMERICAL_DERIVATIVES
dy: ZYTORSIONS ATOMS1=1,2 ATOMS2=5,4 BETWEEN={GAUSSIAN LOWER=-3 UPPER=1.5 SMEAR=0.5}
#dyn: ZYTORSIONS ATOMS1=1,2 ATOMS2=5,4 BETWEEN={GAUSSIAN LOWER=-3 UPPER=1.5 SMEAR=0.5} NUMERICAL_DERIVATIVES
PRINT ARG=dx.*,dy.* FILE=colvar FMT=%8.4f
DUMPDERIVATIVES ARG=dx.*,dy.* FILE=deriv FMT=%8.4f
/* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Copyright (c) 2012-2016 The plumed team
(see the PEOPLE file at the root of the distribution for a list of names)
See http://www.plumed.org for more information.
This file is part of plumed, version 2.
plumed is free software: you can redistribute it and/or modify
it under the terms of the GNU Lesser General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
plumed is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License
along with plumed. If not, see <http://www.gnu.org/licenses/>.
+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
#include "MultiColvar.h"
#include "core/ActionRegister.h"
#include "tools/Angle.h"
#include "tools/SwitchingFunction.h"
#include <string>
#include <cmath>
using namespace std;
namespace PLMD{
namespace multicolvar{
//+PLUMEDOC MCOLVAR XANGLES
/*
Calculate the angles between the vector connecting two atoms and the x axis.
\par Examples
The following input tells plumed to calculate the angles between the x-axis and the vector connecting atom 3 to atom 5 and between the x-axis
and the vector connecting atom 1 to atom 2. The minimum of these two quantities is then
\verbatim
XANLGES ATOMS1=3,5 ATOMS2=1,2 MIN={BETA=0.1} LABEL=d1
PRINT ARG=d1.min
\endverbatim
(See also \ref PRINT).
*/
//+ENDPLUMEDOC
//+PLUMEDOC MCOLVAR YANGLES
/*
Calculate the angles between the vector connecting two atoms and the y axis.
\par Examples
The following input tells plumed to calculate the angles between the y-axis and the vector connecting atom 3 to atom 5 and between the y-axis
and the vector connecting atom 1 to atom 2. The minimum of these two quantities is then
\verbatim
YANLGES ATOMS1=3,5 ATOMS2=1,2 MIN={BETA=0.1} LABEL=d1
PRINT ARG=d1.min
\endverbatim
(See also \ref PRINT).
*/
//+ENDPLUMEDOC
//+PLUMEDOC MCOLVAR ZANGLES
/*
Calculate the angles between the vector connecting two atoms and the z axis.
\par Examples
The following input tells plumed to calculate the angles between the z-axis and the vector connecting atom 3 to atom 5 and between the z-axis
and the vector connecting atom 1 to atom 2. The minimum of these two quantities is then
\verbatim
ZANLGES ATOMS1=3,5 ATOMS2=1,2 MIN={BETA=0.1} LABEL=d1
PRINT ARG=d1.min
\endverbatim
(See also \ref PRINT).
*/
//+ENDPLUMEDOC
class XAngles : public MultiColvar {
private:
bool use_sf;
unsigned myc;
SwitchingFunction sf1;
public:
static void registerKeywords( Keywords& keys );
explicit XAngles(const ActionOptions&);
// active methods:
virtual double compute( const unsigned& tindex, AtomValuePack& myatoms ) const ;
double calculateWeight( const unsigned& taskCode, const double& weight, AtomValuePack& ) const ;
/// Returns the number of coordinates of the field
bool isPeriodic(){ return false; }
};
PLUMED_REGISTER_ACTION(XAngles,"XANGLES")
PLUMED_REGISTER_ACTION(XAngles,"YANGLES")
PLUMED_REGISTER_ACTION(XAngles,"ZANGLES")
void XAngles::registerKeywords( Keywords& keys ){
MultiColvar::registerKeywords( keys );
keys.use("ATOMS"); keys.use("MAX"); keys.use("ALT_MIN");
keys.use("MEAN"); keys.use("MIN"); keys.use("LESS_THAN");
keys.use("LOWEST"); keys.use("HIGHEST");
keys.use("MORE_THAN"); keys.use("BETWEEN"); keys.use("HISTOGRAM"); keys.use("MOMENTS");
keys.add("atoms-1","GROUP","Calculate the distance between each distinct pair of atoms in the group");
keys.add("atoms-2","GROUPA","Calculate the distances between all the atoms in GROUPA and all "
"the atoms in GROUPB. This must be used in conjuction with GROUPB.");
keys.add("atoms-2","GROUPB","Calculate the distances between all the atoms in GROUPA and all the atoms "
"in GROUPB. This must be used in conjuction with GROUPA.");
keys.add("optional","SWITCH","A switching function that ensures that only angles are only computed when atoms are within "
"are within a certain fixed cutoff. The following provides information on the \\ref switchingfunction that are available.");
}
XAngles::XAngles(const ActionOptions&ao):
PLUMED_MULTICOLVAR_INIT(ao),
use_sf(false)
{
if( getName().find("X")!=std::string::npos) myc=0;
else if( getName().find("Y")!=std::string::npos) myc=1;
else if( getName().find("Z")!=std::string::npos) myc=2;
else plumed_error();
// Read in switching function
std::string sfinput, errors; parse("SWITCH",sfinput);
if( sfinput.length()>0 ){
use_sf=true; weightHasDerivatives=true;
sf1.set(sfinput,errors);
if( errors.length()!=0 ) error("problem reading SWITCH keyword : " + errors );
log.printf(" only calculating angles for atoms separated by less than %s\n", sf1.description().c_str() );
setLinkCellCutoff( sf1.get_dmax() );
}
// Read in the atoms
std::vector<AtomNumber> all_atoms;
readTwoGroups( "GROUP", "GROUPA", "GROUPB", all_atoms );
int natoms=2; readAtoms( natoms, all_atoms );
// And check everything has been read in correctly
checkRead();
}
double XAngles::calculateWeight( const unsigned& taskCode, const double& weight, AtomValuePack& myatoms ) const {
if(!use_sf) return 1.0;
Vector distance=getSeparation( myatoms.getPosition(0), myatoms.getPosition(1) );
double dw, w = sf1.calculateSqr( distance.modulo2(), dw );
addAtomDerivatives( 0, 0, (-dw)*distance, myatoms );
addAtomDerivatives( 0, 1, (+dw)*distance, myatoms );
myatoms.addBoxDerivatives( 0, (-dw)*Tensor(distance,distance) );
return w;
}
double XAngles::compute( const unsigned& tindex, AtomValuePack& myatoms ) const {
Vector ddij, ddik, axis, distance; axis.zero(); axis[myc]=1;
distance=getSeparation( myatoms.getPosition(0), myatoms.getPosition(1) );
PLMD::Angle a; double angle=a.compute( distance, axis, ddij, ddik );
addAtomDerivatives( 1, 0, -ddij, myatoms );
addAtomDerivatives( 1, 1, ddij, myatoms );
myatoms.addBoxDerivatives( 1, -Tensor( distance,ddij ) );
return angle;
}
}
}
/* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
Copyright (c) 2012-2016 The plumed team
(see the PEOPLE file at the root of the distribution for a list of names)
See http://www.plumed.org for more information.
This file is part of plumed, version 2.
plumed is free software: you can redistribute it and/or modify
it under the terms of the GNU Lesser General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
plumed is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License
along with plumed. If not, see <http://www.gnu.org/licenses/>.
+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++ */
#include "MultiColvar.h"
#include "core/ActionRegister.h"
#include "tools/Torsion.h"
#include "tools/SwitchingFunction.h"
#include <string>
#include <cmath>
using namespace std;
namespace PLMD{
namespace multicolvar{
//+PLUMEDOC MCOLVAR XYTORSIONS
/*
Calculate the torsional angle around the x axis from the positive y direction.
\par Examples
The following input tells plumed to calculate the angle around the x direction between the positive y-axis and the vector connecting atom 3 to atom 5 and
the angle around the x direction between the positive y axis and the vector connecting atom 1 to atom 2. The minimum of these two quantities is then output
\verbatim
XYTORSIONS ATOMS1=3,5 ATOMS2=1,2 MIN={BETA=0.1} LABEL=d1
PRINT ARG=d1.min
\endverbatim
(See also \ref PRINT).
*/
//+ENDPLUMEDOC
//+PLUMEDOC MCOLVAR XZTORSIONS
/*
Calculate the torsional angle around the x axis from the positive z direction.
\par Examples
The following input tells plumed to calculate the angle around the x direction between the positive z-axis and the vector connecting atom 3 to atom 5 and
the angle around the x direction between the positive z direction and the vector connecting atom 1 to atom 2. The minimum of these two quantities is then output
\verbatim
XZTORSIONS ATOMS1=3,5 ATOMS2=1,2 MIN={BETA=0.1} LABEL=d1
PRINT ARG=d1.min
\endverbatim
(See also \ref PRINT).
*/
//+ENDPLUMEDOC
//+PLUMEDOC MCOLVAR YXTORSIONS
/*
Calculate the torsional angle around the y axis from the positive x direction.
\par Examples
The following input tells plumed to calculate the angle around the y direction between the positive x-direction and the vector connecting atom 3 to atom 5 and
the angle around the y direction between the positive x axis and the vector connecting atom 1 to atom 2. The minimum of these two quantities is then output
\verbatim
YXTORSIONS ATOMS1=3,5 ATOMS2=1,2 MIN={BETA=0.1} LABEL=d1
PRINT ARG=d1.min
\endverbatim
(See also \ref PRINT).
*/
//+ENDPLUMEDOC
//+PLUMEDOC MCOLVAR YZTORSIONS
/*
Calculate the torsional angle around the y axis from the positive z direction.
\par Examples
The following input tells plumed to calculate the angle around the y direction between the positive z-direction and the vector connecting atom 3 to atom 5 and
the angle around the y direction between the positive z direction and the vector connecting atom 1 to atom 2. The minimum of these two quantities is then output
\verbatim
YZTORSIONS ATOMS1=3,5 ATOMS2=1,2 MIN={BETA=0.1} LABEL=d1
PRINT ARG=d1.min
\endverbatim
(See also \ref PRINT).
*/
//+ENDPLUMEDOC
//+PLUMEDOC MCOLVAR ZXTORSIONS
/*
Calculate the torsional angle around the z axis from the positive x direction.
\par Examples
The following input tells plumed to calculate the angle around the z direction between the positive x-direction and the vector connecting atom 3 to atom 5 and
the angle around the z direction between the positive x-direction and the vector connecting atom 1 to atom 2. The minimum of these two quantities is then output
\verbatim
ZXTORSIONS ATOMS1=3,5 ATOMS2=1,2 MIN={BETA=0.1} LABEL=d1
PRINT ARG=d1.min
\endverbatim
(See also \ref PRINT).
*/
//+ENDPLUMEDOC
//+PLUMEDOC MCOLVAR ZYTORSIONS
/*
Calculate the torsional angle around the z axis from the positive y direction.
\par Examples
The following input tells plumed to calculate the angle around the z direction between the positive y-axis and the vector connecting atom 3 to atom 5 and
the angle around the z direction between the positive y axis and the vector connecting atom 1 to atom 2. The minimum of these two quantities is then output
\verbatim
ZYTORSIONS ATOMS1=3,5 ATOMS2=1,2 MIN={BETA=0.1} LABEL=d1
PRINT ARG=d1.min
\endverbatim
(See also \ref PRINT).
*/
//+ENDPLUMEDOC
class XYTorsion : public MultiColvar {
private:
bool use_sf;
unsigned myc1, myc2;
SwitchingFunction sf1;
public:
static void registerKeywords( Keywords& keys );
explicit XYTorsion(const ActionOptions&);
// active methods:
virtual double compute( const unsigned& tindex, AtomValuePack& myatoms ) const ;
double calculateWeight( const unsigned& taskCode, const double& weight, AtomValuePack& ) const ;
/// Returns the number of coordinates of the field
bool isPeriodic(){ return false; }
};
PLUMED_REGISTER_ACTION(XYTorsion,"XYTORSIONS")
PLUMED_REGISTER_ACTION(XYTorsion,"XZTORSIONS")
PLUMED_REGISTER_ACTION(XYTorsion,"YXTORSIONS")
PLUMED_REGISTER_ACTION(XYTorsion,"YZTORSIONS")
PLUMED_REGISTER_ACTION(XYTorsion,"ZXTORSIONS")
PLUMED_REGISTER_ACTION(XYTorsion,"ZYTORSIONS")
void XYTorsion::registerKeywords( Keywords& keys ){
MultiColvar::registerKeywords( keys );
keys.use("ATOMS"); keys.use("MAX"); keys.use("ALT_MIN");
keys.use("MEAN"); keys.use("MIN");
keys.use("LOWEST"); keys.use("HIGHEST");
keys.use("BETWEEN"); keys.use("HISTOGRAM"); keys.use("MOMENTS");
keys.add("atoms-1","GROUP","Calculate the distance between each distinct pair of atoms in the group");
keys.add("atoms-2","GROUPA","Calculate the distances between all the atoms in GROUPA and all "
"the atoms in GROUPB. This must be used in conjuction with GROUPB.");
keys.add("atoms-2","GROUPB","Calculate the distances between all the atoms in GROUPA and all the atoms "
"in GROUPB. This must be used in conjuction with GROUPA.");
keys.add("optional","SWITCH","A switching function that ensures that only angles are only computed when atoms are within "
"are within a certain fixed cutoff. The following provides information on the \\ref switchingfunction that are available.");
}
XYTorsion::XYTorsion(const ActionOptions&ao):
PLUMED_MULTICOLVAR_INIT(ao),
use_sf(false)
{
if( getName().find("XY")!=std::string::npos){ myc1=0; myc2=1; }
else if( getName().find("XZ")!=std::string::npos){ myc1=0; myc2=2; }
else if( getName().find("YX")!=std::string::npos){ myc1=1; myc2=0; }
else if( getName().find("YZ")!=std::string::npos){ myc1=1; myc2=2; }
else if( getName().find("ZX")!=std::string::npos){ myc1=2; myc2=0; }
else if( getName().find("ZY")!=std::string::npos){ myc1=2; myc2=1; }
else plumed_error();
// Read in switching function
std::string sfinput, errors; parse("SWITCH",sfinput);
if( sfinput.length()>0 ){
use_sf=true; weightHasDerivatives=true;
sf1.set(sfinput,errors);
if( errors.length()!=0 ) error("problem reading SWITCH keyword : " + errors );
log.printf(" only calculating angles for atoms separated by less than %s\n", sf1.description().c_str() );
setLinkCellCutoff( sf1.get_dmax() );
}
// Read in the atoms
std::vector<AtomNumber> all_atoms;
readTwoGroups( "GROUP", "GROUPA", "GROUPB", all_atoms );
int natoms=2; readAtoms( natoms, all_atoms );
// And check everything has been read in correctly
checkRead();
}
double XYTorsion::calculateWeight( const unsigned& taskCode, const double& weight, AtomValuePack& myatoms ) const {
if(!use_sf) return 1.0;
Vector distance=getSeparation( myatoms.getPosition(0), myatoms.getPosition(1) );
double dw, w = sf1.calculateSqr( distance.modulo2(), dw );
addAtomDerivatives( 0, 0, (-dw)*distance, myatoms );
addAtomDerivatives( 0, 1, (+dw)*distance, myatoms );
myatoms.addBoxDerivatives( 0, (-dw)*Tensor(distance,distance) );
return w;
}
double XYTorsion::compute( const unsigned& tindex, AtomValuePack& myatoms ) const {
Vector dd0, dd1, dd2, axis, rot, distance;
axis.zero(); rot.zero(); rot[myc1]=1; axis[myc2]=1;
distance=getSeparation( myatoms.getPosition(0), myatoms.getPosition(1) );
PLMD::Torsion t; double torsion=t.compute( distance, rot, axis, dd0, dd1, dd2 );
addAtomDerivatives( 1, 0, -dd0, myatoms );
addAtomDerivatives( 1, 1, dd0, myatoms );
myatoms.addBoxDerivatives( 1, -extProduct(distance,dd0) );
return torsion;
}
}
}
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment