Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2015-2019 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 :
28 47 : void HistogramOnGrid::registerKeywords( Keywords& keys ) {
29 47 : GridVessel::registerKeywords( keys );
30 188 : keys.add("compulsory","KERNEL","the type of kernel to use");
31 188 : keys.add("compulsory","BANDWIDTH","the bandwidths");
32 188 : keys.add("compulsory","CONCENTRATION","the concentration parameter for Von Mises-Fisher distributions");
33 47 : }
34 :
35 36 : HistogramOnGrid::HistogramOnGrid( const vesselbase::VesselOptions& da ):
36 : GridVessel(da),
37 : neigh_tot(0),
38 : addOneKernelAtATime(false),
39 36 : bandwidths(dimension),
40 144 : discrete(false)
41 : {
42 72 : if( getType()=="flat" ) {
43 68 : parse("KERNEL",kerneltype);
44 68 : if( kerneltype=="discrete" || kerneltype=="DISCRETE" ) {
45 4 : discrete=true; setNoDerivatives();
46 : } else {
47 60 : parseVector("BANDWIDTH",bandwidths);
48 : }
49 : } else {
50 4 : parse("CONCENTRATION",von_misses_concentration);
51 2 : von_misses_norm = von_misses_concentration / ( 4*pi*sinh( von_misses_concentration ) );
52 : }
53 36 : }
54 :
55 3 : double HistogramOnGrid::getFibonacciCutoff() const {
56 3 : return std::log( epsilon / von_misses_norm ) / von_misses_concentration;
57 : }
58 :
59 18 : bool HistogramOnGrid::noDiscreteKernels() const {
60 18 : return !discrete;
61 : }
62 :
63 41 : 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 41 : GridVessel::setBounds( smin, smax, nbins, spacing );
66 41 : if( !discrete ) {
67 37 : std::vector<double> point(dimension,0);
68 74 : KernelFunctions kernel( point, bandwidths, kerneltype, false, 1.0, true ); neigh_tot=1;
69 74 : nneigh=kernel.getSupport( dx ); std::vector<double> support( kernel.getContinuousSupport() );
70 155 : for(unsigned i=0; i<dimension; ++i) {
71 178 : if( pbc[i] && 2*support[i]>getGridExtent(i) ) error("bandwidth is too large for periodic grid");
72 59 : neigh_tot *= (2*nneigh[i]+1);
73 : }
74 : }
75 41 : }
76 :
77 21475 : KernelFunctions* HistogramOnGrid::getKernelAndNeighbors( std::vector<double>& point, unsigned& num_neigh, std::vector<unsigned>& neighbors ) const {
78 21475 : if( discrete ) {
79 10 : plumed_assert( getType()=="flat" );
80 15 : num_neigh=1; for(unsigned i=0; i<dimension; ++i) point[i] += 0.5*dx[i];
81 5 : neighbors[0] = getIndex( point ); return NULL;
82 42940 : } else if( getType()=="flat" ) {
83 21413 : KernelFunctions* kernel = new KernelFunctions( point, bandwidths, kerneltype, false, 1.0, true );
84 64239 : getNeighbors( kernel->getCenter(), nneigh, num_neigh, neighbors );
85 21413 : return kernel;
86 114 : } else if( getType()=="fibonacci" ) {
87 57 : getNeighbors( point, nneigh, num_neigh, neighbors );
88 57 : return NULL;
89 : } else {
90 0 : plumed_error();
91 : }
92 : return NULL;
93 : }
94 :
95 36118 : std::vector<Value*> HistogramOnGrid::getVectorOfValues() const {
96 : std::vector<Value*> vv;
97 242634 : for(unsigned i=0; i<dimension; ++i) {
98 206516 : vv.push_back(new Value());
99 316440 : if( pbc[i] ) vv[i]->setDomain( str_min[i], str_max[i] );
100 48296 : else vv[i]->setNotPeriodic();
101 : }
102 36118 : return vv;
103 : }
104 :
105 36118 : void HistogramOnGrid::calculate( const unsigned& current, MultiValue& myvals, std::vector<double>& buffer, std::vector<unsigned>& der_list ) const {
106 36118 : if( addOneKernelAtATime ) {
107 : plumed_dbg_assert( myvals.getNumberOfValues()==2 && !wasforced );
108 14713 : std::vector<double> der( dimension );
109 93491 : for(unsigned i=0; i<dimension; ++i) der[i]=myvals.getDerivative( 1, i );
110 58852 : accumulate( getAction()->getPositionInCurrentTaskList(current), myvals.get(0), myvals.get(1), der, buffer );
111 : } else {
112 : plumed_dbg_assert( myvals.getNumberOfValues()==dimension+2 );
113 64215 : std::vector<double> point( dimension ); double weight=myvals.get(0)*myvals.get( 1+dimension );
114 149143 : for(unsigned i=0; i<dimension; ++i) point[i]=myvals.get( 1+i );
115 :
116 : // Get the kernel
117 : unsigned num_neigh; std::vector<unsigned> neighbors;
118 21405 : std::vector<double> der( dimension );
119 21405 : KernelFunctions* kernel=getKernelAndNeighbors( point, num_neigh, neighbors );
120 :
121 42867 : if( !kernel && getType()=="flat" ) {
122 0 : plumed_dbg_assert( num_neigh==1 ); der.resize(0);
123 0 : accumulate( neighbors[0], weight, 1.0, der, buffer );
124 : } else {
125 : double totwforce=0.0;
126 21405 : std::vector<double> intforce( 2*dimension, 0.0 );
127 21405 : std::vector<Value*> vv( getVectorOfValues() );
128 :
129 21405 : double newval; std::vector<unsigned> tindices( dimension ); std::vector<double> xx( dimension );
130 41141097 : for(unsigned i=0; i<num_neigh; ++i) {
131 41119692 : unsigned ineigh=neighbors[i];
132 20559846 : if( inactive( ineigh ) ) continue ;
133 11847474 : getGridPointCoordinates( ineigh, tindices, xx );
134 11847474 : if( kernel ) {
135 153875677 : for(unsigned j=0; j<dimension; ++j) vv[j]->set(xx[j]);
136 11842701 : newval = kernel->evaluate( vv, der, true );
137 : } else {
138 : // Evalulate dot product
139 62049 : double dot=0; for(unsigned j=0; j<dimension; ++j) { dot+=xx[j]*point[j]; der[j]=xx[j]; }
140 : // Von misses distribution for concentration parameter
141 4773 : newval = von_misses_norm*exp( von_misses_concentration*dot );
142 : // And final derivatives
143 33411 : for(unsigned j=0; j<dimension; ++j) der[j] *= von_misses_concentration*newval;
144 : }
145 11847474 : accumulate( ineigh, weight, newval, der, buffer );
146 11847474 : if( wasForced() ) {
147 5744 : accumulateForce( ineigh, weight, der, intforce );
148 17232 : totwforce += myvals.get( 1+dimension )*newval*forces[ineigh];
149 : }
150 : }
151 21405 : if( wasForced() ) {
152 : // Minus sign for kernel here as we are taking derivative with respect to position of center of
153 : // kernel NOT derivative wrt to grid point
154 110 : double pref = 1; if( kernel ) pref = -1;
155 110 : unsigned nder = getAction()->getNumberOfDerivatives();
156 220 : unsigned gridbuf = getNumberOfBufferPoints()*getNumberOfQuantities();
157 490 : for(unsigned j=0; j<dimension; ++j) {
158 67690 : for(unsigned k=0; k<myvals.getNumberActive(); ++k) {
159 : unsigned kder=myvals.getActiveIndex(k);
160 101250 : buffer[ bufstart + gridbuf + kder ] += pref*intforce[j]*myvals.getDerivative( j+1, kder );
161 : }
162 : }
163 : // Accumulate the sum of all the weights
164 220 : buffer[ bufstart + gridbuf + nder ] += myvals.get(0);
165 : // Add the derivatives of the weights into the force -- this is separate loop as weights of all parts are considered together
166 64010 : for(unsigned k=0; k<myvals.getNumberActive(); ++k) {
167 : unsigned kder=myvals.getActiveIndex(k);
168 63900 : buffer[ bufstart + gridbuf + kder ] += totwforce*myvals.getDerivative( 0, kder );
169 63900 : buffer[ bufstart + gridbuf + nder + 1 + kder ] += myvals.getDerivative( 0, kder );
170 : }
171 : }
172 21405 : if( kernel ) delete kernel;
173 149143 : for(unsigned i=0; i<dimension; ++i) delete vv[i];
174 : }
175 : }
176 36118 : }
177 :
178 37076 : void HistogramOnGrid::accumulate( const unsigned& ipoint, const double& weight, const double& dens, const std::vector<double>& der, std::vector<double>& buffer ) const {
179 74152 : buffer[bufstart+nper*ipoint] += weight*dens;
180 211130 : if( der.size()>0 ) for(unsigned j=0; j<dimension; ++j) buffer[bufstart+nper*ipoint + 1 + j] += weight*der[j];
181 37076 : }
182 :
183 5744 : void HistogramOnGrid::accumulateForce( const unsigned& ipoint, const double& weight, const std::vector<double>& der, std::vector<double>& intforce ) const {
184 74838 : for(unsigned j=0; j<der.size(); ++j) intforce[j] += forces[ipoint]*weight*der[j];
185 5744 : }
186 :
187 20 : void HistogramOnGrid::getFinalForces( const std::vector<double>& buffer, std::vector<double>& finalForces ) {
188 24 : if( finalForces.size()!=getAction()->getNumberOfDerivatives() ) finalForces.resize( getAction()->getNumberOfDerivatives() );
189 : // And the final force
190 20 : unsigned nder = getAction()->getNumberOfDerivatives();
191 : // Derivatives due to normalization
192 40 : unsigned gridbuf = getNumberOfBufferPoints()*getNumberOfQuantities();
193 11425 : for(unsigned i=0; i<finalForces.size(); ++i) finalForces[i] = buffer[ bufstart + gridbuf + i ];
194 : // Derivatives due to normalization
195 20 : if( !noAverage() ) {
196 15 : unsigned wderstart = bufstart + gridbuf + nder; double pref=0;
197 7225 : for(unsigned ipoint=0; ipoint<getNumberOfPoints(); ++ipoint) {
198 14420 : pref += forces[ipoint]*buffer[ bufstart + ipoint*nper ] / buffer[wderstart];
199 : }
200 14250 : for(unsigned j=0; j<finalForces.size(); ++j) finalForces[j] -= pref*buffer[ wderstart + 1 + j ];
201 : }
202 20 : }
203 :
204 143 : void HistogramOnGrid::finish( const std::vector<double>& buffer ) {
205 143 : if( addOneKernelAtATime ) {
206 29491 : for(unsigned i=0; i<getAction()->getCurrentNumberOfActiveTasks(); ++i) {
207 177019 : for(unsigned j=0; j<nper; ++j) addDataElement( nper*getAction()->getActiveTask(i)+j, buffer[bufstart+i*nper+j] );
208 : }
209 : } else {
210 78 : GridVessel::finish( buffer );
211 : }
212 143 : }
213 :
214 : }
215 4839 : }
|