LCOV - code coverage report
Current view: top level - symfunc - CylindricalHarmonic.cpp (source / functions) Hit Total Coverage
Test: plumed test coverage Lines: 24 24 100.0 %
Date: 2024-10-18 14:00:25 Functions: 4 4 100.0 %

          Line data    Source code
       1             : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
       2             :    Copyright (c) 2012-2017 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 "function/FunctionTemplateBase.h"
      23             : #include "function/FunctionShortcut.h"
      24             : #include "function/FunctionOfMatrix.h"
      25             : #include "core/ActionRegister.h"
      26             : 
      27             : #include <complex>
      28             : 
      29             : namespace PLMD {
      30             : namespace symfunc {
      31             : 
      32             : //+PLUMEDOC MCOLVAR CYLINDRICAL_HARMONIC
      33             : /*
      34             : Calculate the cylindrical harmonic function
      35             : 
      36             : \par Examples
      37             : 
      38             : 
      39             : */
      40             : //+ENDPLUMEDOC
      41             : 
      42             : //+PLUMEDOC MCOLVAR CYLINDRICAL_HARMONIC_MATRIX
      43             : /*
      44             : Calculate the cylindrical harmonic function from the elements in two input matrices
      45             : 
      46             : \par Examples
      47             : 
      48             : 
      49             : */
      50             : //+ENDPLUMEDOC
      51             : 
      52             : 
      53          12 : class CylindricalHarmonic : public function::FunctionTemplateBase {
      54             : private:
      55             :   int tmom;
      56             : public:
      57             :   void registerKeywords( Keywords& keys ) override;
      58             :   void read( ActionWithArguments* action ) override;
      59             :   void setPeriodicityForOutputs( ActionWithValue* action ) override;
      60             :   void calc( const ActionWithArguments* action, const std::vector<double>& args, std::vector<double>& vals, Matrix<double>& derivatives ) const override;
      61             : };
      62             : 
      63             : typedef function::FunctionShortcut<CylindricalHarmonic> CyHarmShortcut;
      64             : PLUMED_REGISTER_ACTION(CyHarmShortcut,"CYLINDRICAL_HARMONIC")
      65             : typedef function::FunctionOfMatrix<CylindricalHarmonic> MatrixCyHarm;
      66             : PLUMED_REGISTER_ACTION(MatrixCyHarm,"CYLINDRICAL_HARMONIC_MATRIX")
      67             : 
      68           9 : void CylindricalHarmonic::registerKeywords( Keywords& keys ) {
      69          18 :   keys.add("compulsory","DEGREE","the value of the n parameter in the equation above");
      70          18 :   keys.addOutputComponent("rm","default","the real part of the cylindrical harmonic");
      71          18 :   keys.addOutputComponent("im","default","the imaginary part of the cylindrical harmonic");
      72           9 : }
      73             : 
      74           2 : void CylindricalHarmonic::read( ActionWithArguments* action ) {
      75           2 :   parse(action,"DEGREE",tmom);
      76           2 :   action->log.printf("  calculating %dth order cylindrical harmonic with %s and %s as input \n", tmom, action->getPntrToArgument(0)->getName().c_str(), action->getPntrToArgument(1)->getName().c_str() );
      77           2 :   if( action->getNumberOfArguments()==3 ) action->log.printf("  multiplying cylindrical harmonic by weight from %s \n", action->getPntrToArgument(2)->getName().c_str() );
      78           2 : }
      79             : 
      80           2 : void CylindricalHarmonic::setPeriodicityForOutputs( ActionWithValue* action ) {
      81           4 :   action->componentIsNotPeriodic("rm"); action->componentIsNotPeriodic("im");
      82           2 : }
      83             : 
      84      943610 : void CylindricalHarmonic::calc( const ActionWithArguments* action, const std::vector<double>& args, std::vector<double>& vals, Matrix<double>& derivatives ) const {
      85      943610 :   double dlen2 = args[0]*args[0] + args[1]*args[1]; double dlen = sqrt( dlen2 ); double dlen3 = dlen2*dlen;
      86      943610 :   std::complex<double> com1( args[0]/dlen,args[1]/dlen ); double weight=1; if( args.size()==3 ) weight=args[2];
      87      943610 :   std::complex<double> ppp = pow( com1, tmom-1 ), ii( 0, 1 ); double real_z = real( ppp*com1 ), imag_z = imag( ppp*com1 );
      88      943610 :   std::complex<double> dp_x = static_cast<double>(tmom)*ppp*( (1.0/dlen)-(args[0]*args[0])/dlen3-ii*(args[0]*args[1])/dlen3 );
      89      943610 :   std::complex<double> dp_y = static_cast<double>(tmom)*ppp*( ii*(1.0/dlen)-(args[0]*args[1])/dlen3-ii*(args[1]*args[1])/dlen3 );
      90      943610 :   vals[0] = weight*real_z; derivatives(0,0) = weight*real(dp_x); derivatives(0,1) = weight*real(dp_y);
      91      943610 :   vals[1] = weight*imag_z; derivatives(1,0) = weight*imag(dp_x); derivatives(1,1) = weight*imag(dp_y);
      92      943610 :   if( args.size()==3 ) { derivatives(0,2) = real_z; derivatives(1,2) = imag_z; }
      93      943610 : }
      94             : 
      95             : }
      96             : }
      97             : 

Generated by: LCOV version 1.16