1 /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2    Copyright (c) 2015-2021 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 "HistogramOnGrid.h"
23 #include "tools/KernelFunctions.h"
24 
25 namespace PLMD {
26 namespace gridtools {
27 
registerKeywords(Keywords & keys)28 void HistogramOnGrid::registerKeywords( Keywords& keys ) {
29   GridVessel::registerKeywords( keys );
30   keys.add("compulsory","KERNEL","the type of kernel to use");
31   keys.add("compulsory","BANDWIDTH","the bandwidths");
32   keys.add("compulsory","CONCENTRATION","the concentration parameter for Von Mises-Fisher distributions");
33 }
34 
HistogramOnGrid(const vesselbase::VesselOptions & da)35 HistogramOnGrid::HistogramOnGrid( const vesselbase::VesselOptions& da ):
36   GridVessel(da),
37   neigh_tot(0),
38   addOneKernelAtATime(false),
39   bandwidths(dimension),
40   discrete(false)
41 {
42   if( getType()=="flat" ) {
43     parse("KERNEL",kerneltype);
44     if( kerneltype=="discrete" || kerneltype=="DISCRETE" ) {
45       discrete=true; setNoDerivatives();
46     } else {
47       parseVector("BANDWIDTH",bandwidths);
48     }
49   } else {
50     parse("CONCENTRATION",von_misses_concentration);
51     von_misses_norm = von_misses_concentration / ( 4*pi*sinh( von_misses_concentration ) );
52   }
53 }
54 
getFibonacciCutoff() const55 double HistogramOnGrid::getFibonacciCutoff() const {
56   return std::log( epsilon / von_misses_norm ) / von_misses_concentration;
57 }
58 
noDiscreteKernels() const59 bool HistogramOnGrid::noDiscreteKernels() const {
60   return !discrete;
61 }
62 
setBounds(const std::vector<std::string> & smin,const std::vector<std::string> & smax,const std::vector<unsigned> & nbins,const std::vector<double> & spacing)63 void HistogramOnGrid::setBounds( const std::vector<std::string>& smin, const std::vector<std::string>& smax,
64                                  const std::vector<unsigned>& nbins, const std::vector<double>& spacing ) {
65   GridVessel::setBounds( smin, smax, nbins, spacing );
66   if( !discrete ) {
67     std::vector<double> point(dimension,0);
68     KernelFunctions kernel( point, bandwidths, kerneltype, "DIAGONAL", 1.0 ); neigh_tot=1;
69     nneigh=kernel.getSupport( dx ); std::vector<double> support( kernel.getContinuousSupport() );
70     for(unsigned i=0; i<dimension; ++i) {
71       if( pbc[i] && 2*support[i]>getGridExtent(i) ) error("bandwidth is too large for periodic grid");
72       neigh_tot *= (2*nneigh[i]+1);
73     }
74   }
75 }
76 
getKernelAndNeighbors(std::vector<double> & point,unsigned & num_neigh,std::vector<unsigned> & neighbors) const77 std::unique_ptr<KernelFunctions> HistogramOnGrid::getKernelAndNeighbors( std::vector<double>& point, unsigned& num_neigh, std::vector<unsigned>& neighbors ) const {
78   if( discrete ) {
79     plumed_assert( getType()=="flat" );
80     num_neigh=1; for(unsigned i=0; i<dimension; ++i) point[i] += 0.5*dx[i];
81     neighbors[0] = getIndex( point ); return NULL;
82   } else if( getType()=="flat" ) {
83     std::unique_ptr<KernelFunctions> kernel(new KernelFunctions( point, bandwidths, kerneltype, "DIAGONAL", 1.0 ));
84 // GB: Now values is destroyed when exiting this function.
85 // I think before there was a leak
86     std::vector<std::unique_ptr<Value>> values=getVectorOfValues();
87     kernel->normalize( Tools::unique2raw(values) ); getNeighbors( kernel->getCenter(), nneigh, num_neigh, neighbors );
88     return kernel;
89   } else if( getType()=="fibonacci" ) {
90     getNeighbors( point, nneigh, num_neigh, neighbors );
91     return NULL;
92   } else {
93     plumed_error();
94   }
95   return NULL;
96 }
97 
getVectorOfValues() const98 std::vector<std::unique_ptr<Value>> HistogramOnGrid::getVectorOfValues() const {
99   std::vector<std::unique_ptr<Value>> vv;
100   for(unsigned i=0; i<dimension; ++i) {
101     vv.emplace_back(new Value());
102     if( pbc[i] ) vv[i]->setDomain( str_min[i], str_max[i] );
103     else vv[i]->setNotPeriodic();
104   }
105   return vv;
106 }
107 
calculate(const unsigned & current,MultiValue & myvals,std::vector<double> & buffer,std::vector<unsigned> & der_list) const108 void HistogramOnGrid::calculate( const unsigned& current, MultiValue& myvals, std::vector<double>& buffer, std::vector<unsigned>& der_list ) const {
109   if( addOneKernelAtATime ) {
110     plumed_dbg_assert( myvals.getNumberOfValues()==2 && !wasforced );
111     std::vector<double> der( dimension );
112     for(unsigned i=0; i<dimension; ++i) der[i]=myvals.getDerivative( 1, i );
113     accumulate( getAction()->getPositionInCurrentTaskList(current), myvals.get(0), myvals.get(1), der, buffer );
114   } else {
115     plumed_dbg_assert( myvals.getNumberOfValues()==dimension+2 );
116     std::vector<double> point( dimension ); double weight=myvals.get(0)*myvals.get( 1+dimension );
117     for(unsigned i=0; i<dimension; ++i) point[i]=myvals.get( 1+i );
118     // Get the kernel
119     unsigned num_neigh; std::vector<unsigned> neighbors(1);
120     std::vector<double> der( dimension );
121     std::unique_ptr<KernelFunctions> kernel=getKernelAndNeighbors( point, num_neigh, neighbors );
122 
123     if( !kernel && getType()=="flat" ) {
124       plumed_dbg_assert( num_neigh==1 ); der.resize(0);
125       accumulate( neighbors[0], weight, 1.0, der, buffer );
126     } else {
127       double totwforce=0.0;
128       std::vector<double> intforce( 2*dimension, 0.0 );
129       std::vector<std::unique_ptr<Value>> vv( getVectorOfValues() );
130 
131       double newval; std::vector<unsigned> tindices( dimension ); std::vector<double> xx( dimension );
132       for(unsigned i=0; i<num_neigh; ++i) {
133         unsigned ineigh=neighbors[i];
134         if( inactive( ineigh ) ) continue ;
135         getGridPointCoordinates( ineigh, tindices, xx );
136         if( kernel ) {
137           for(unsigned j=0; j<dimension; ++j) vv[j]->set(xx[j]);
138           newval = kernel->evaluate( Tools::unique2raw(vv), der, true );
139         } else {
140           // Evalulate dot product
141           double dot=0; for(unsigned j=0; j<dimension; ++j) { dot+=xx[j]*point[j]; der[j]=xx[j]; }
142           // Von misses distribution for concentration parameter
143           newval = von_misses_norm*exp( von_misses_concentration*dot );
144           // And final derivatives
145           for(unsigned j=0; j<dimension; ++j) der[j] *= von_misses_concentration*newval;
146         }
147         accumulate( ineigh, weight, newval, der, buffer );
148         if( wasForced() ) {
149           accumulateForce( ineigh, weight, der, intforce );
150           totwforce += myvals.get( 1+dimension )*newval*forces[ineigh];
151         }
152       }
153       if( wasForced() ) {
154         // Minus sign for kernel here as we are taking derivative with respect to position of center of
155         // kernel NOT derivative wrt to grid point
156         double pref = 1; if( kernel ) pref = -1;
157         unsigned nder = getAction()->getNumberOfDerivatives();
158         unsigned gridbuf = getNumberOfBufferPoints()*getNumberOfQuantities();
159         for(unsigned j=0; j<dimension; ++j) {
160           for(unsigned k=0; k<myvals.getNumberActive(); ++k) {
161             unsigned kder=myvals.getActiveIndex(k);
162             buffer[ bufstart + gridbuf + kder ] += pref*intforce[j]*myvals.getDerivative( j+1, kder );
163           }
164         }
165         // Accumulate the sum of all the weights
166         buffer[ bufstart + gridbuf + nder ] += myvals.get(0);
167         // Add the derivatives of the weights into the force -- this is separate loop as weights of all parts are considered together
168         for(unsigned k=0; k<myvals.getNumberActive(); ++k) {
169           unsigned kder=myvals.getActiveIndex(k);
170           buffer[ bufstart + gridbuf + kder ] += totwforce*myvals.getDerivative( 0, kder );
171           buffer[ bufstart + gridbuf + nder + 1 + kder ] += myvals.getDerivative( 0, kder );
172         }
173       }
174     }
175   }
176 }
177 
accumulate(const unsigned & ipoint,const double & weight,const double & dens,const std::vector<double> & der,std::vector<double> & buffer) const178 void HistogramOnGrid::accumulate( const unsigned& ipoint, const double& weight, const double& dens, const std::vector<double>& der, std::vector<double>& buffer ) const {
179   buffer[bufstart+nper*ipoint] += weight*dens;
180   if( der.size()>0 ) for(unsigned j=0; j<dimension; ++j) buffer[bufstart+nper*ipoint + 1 + j] += weight*der[j];
181 }
182 
accumulateForce(const unsigned & ipoint,const double & weight,const std::vector<double> & der,std::vector<double> & intforce) const183 void HistogramOnGrid::accumulateForce( const unsigned& ipoint, const double& weight, const std::vector<double>& der, std::vector<double>& intforce ) const {
184   for(unsigned j=0; j<der.size(); ++j) intforce[j] += forces[ipoint]*weight*der[j];
185 }
186 
getFinalForces(const std::vector<double> & buffer,std::vector<double> & finalForces)187 void HistogramOnGrid::getFinalForces( const std::vector<double>& buffer, std::vector<double>& finalForces ) {
188   if( finalForces.size()!=getAction()->getNumberOfDerivatives() ) finalForces.resize( getAction()->getNumberOfDerivatives() );
189   // And the final force
190   unsigned nder = getAction()->getNumberOfDerivatives();
191   // Derivatives due to normalization
192   unsigned gridbuf = getNumberOfBufferPoints()*getNumberOfQuantities();
193   for(unsigned i=0; i<finalForces.size(); ++i) finalForces[i] = buffer[ bufstart + gridbuf + i ];
194   // Derivatives due to normalization
195   if( !noAverage() ) {
196     unsigned wderstart = bufstart + gridbuf + nder; double pref=0;
197     for(unsigned ipoint=0; ipoint<getNumberOfPoints(); ++ipoint) {
198       pref += forces[ipoint]*buffer[ bufstart + ipoint*nper ] / buffer[wderstart];
199     }
200     for(unsigned j=0; j<finalForces.size(); ++j) finalForces[j] -= pref*buffer[ wderstart + 1 + j ];
201   }
202 }
203 
finish(const std::vector<double> & buffer)204 void HistogramOnGrid::finish( const std::vector<double>& buffer ) {
205   if( addOneKernelAtATime ) {
206     for(unsigned i=0; i<getAction()->getCurrentNumberOfActiveTasks(); ++i) {
207       for(unsigned j=0; j<nper; ++j) addDataElement( nper*getAction()->getActiveTask(i)+j, buffer[bufstart+i*nper+j] );
208     }
209   } else {
210     GridVessel::finish( buffer );
211   }
212 }
213 
214 }
215 }
216