Line data Source code
1 : /* +++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
2 : Copyright (c) 2012-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 "PathMSDBase.h"
23 : #include "Colvar.h"
24 : #include "ActionRegister.h"
25 : #include "core/PlumedMain.h"
26 : #include "core/Atoms.h"
27 : #include "tools/PDB.h"
28 : #include "tools/RMSD.h"
29 : #include "tools/Tools.h"
30 :
31 : namespace PLMD {
32 : namespace colvar {
33 :
34 27 : void PathMSDBase::registerKeywords(Keywords& keys) {
35 27 : Colvar::registerKeywords(keys);
36 54 : keys.add("compulsory","LAMBDA","the lambda parameter is needed for smoothing, is in the units of plumed");
37 54 : keys.add("compulsory","REFERENCE","the pdb is needed to provide the various milestones");
38 54 : keys.add("optional","NEIGH_SIZE","size of the neighbor list");
39 54 : keys.add("optional","NEIGH_STRIDE","how often the neighbor list needs to be calculated in time units");
40 54 : keys.add("optional", "EPSILON", "(default=-1) the maximum distance between the close and the current structure, the positive value turn on the close structure method");
41 54 : keys.add("optional", "LOG_CLOSE", "(default=0) value 1 enables logging regarding the close structure");
42 54 : keys.add("optional", "DEBUG_CLOSE", "(default=0) value 1 enables extensive debugging info regarding the close structure, the simulation will run much slower");
43 27 : }
44 :
45 25 : PathMSDBase::PathMSDBase(const ActionOptions&ao):
46 : PLUMED_COLVAR_INIT(ao),
47 25 : nopbc(false),
48 25 : neigh_size(-1),
49 25 : neigh_stride(-1),
50 25 : epsilonClose(-1),
51 25 : debugClose(0),
52 25 : logClose(0),
53 25 : computeRefClose(false),
54 25 : nframes(0)
55 : {
56 25 : parse("LAMBDA",lambda);
57 25 : parse("NEIGH_SIZE",neigh_size);
58 25 : parse("NEIGH_STRIDE",neigh_stride);
59 25 : parse("REFERENCE",reference);
60 25 : parse("EPSILON", epsilonClose);
61 25 : parse("LOG_CLOSE", logClose);
62 25 : parse("DEBUG_CLOSE", debugClose);
63 25 : parseFlag("NOPBC",nopbc);
64 :
65 : // open the file
66 25 : FILE* fp=fopen(reference.c_str(),"r");
67 : std::vector<AtomNumber> aaa;
68 25 : if (fp!=NULL)
69 : {
70 25 : log<<"Opening reference file "<<reference.c_str()<<"\n";
71 : bool do_read=true;
72 : unsigned nat=0;
73 1107 : while (do_read) {
74 1107 : PDB mypdb;
75 1107 : RMSD mymsd;
76 2214 : do_read=mypdb.readFromFilepointer(fp,plumed.getAtoms().usingNaturalUnits(),0.1/atoms.getUnits().getLength());
77 1107 : if(do_read) {
78 1082 : nframes++;
79 1082 : if(mypdb.getAtomNumbers().size()==0) error("number of atoms in a frame should be more than zero");
80 1082 : if(nat==0) nat=mypdb.getAtomNumbers().size();
81 1082 : if(nat!=mypdb.getAtomNumbers().size()) error("frames should have the same number of atoms");
82 1082 : if(aaa.empty()) {
83 25 : aaa=mypdb.getAtomNumbers();
84 25 : log.printf(" found %zu atoms in input \n",aaa.size());
85 25 : log.printf(" with indices : ");
86 346 : for(unsigned i=0; i<aaa.size(); ++i) {
87 321 : if(i%25==0) log<<"\n";
88 321 : log.printf("%d ",aaa[i].serial());
89 : }
90 25 : log.printf("\n");
91 : }
92 2164 : if(aaa!=mypdb.getAtomNumbers()) error("frames should contain same atoms in same order");
93 1082 : log<<"Found PDB: "<<nframes<<" containing "<<mypdb.getAtomNumbers().size()<<" atoms\n";
94 1082 : pdbv.push_back(mypdb);
95 1082 : derivs_s.resize(mypdb.getAtomNumbers().size());
96 1082 : derivs_z.resize(mypdb.getAtomNumbers().size());
97 1082 : mymsd.set(mypdb,"OPTIMAL");
98 1082 : msdv.push_back(mymsd); // the vector that stores the frames
99 : } else {break ;}
100 1107 : }
101 25 : fclose (fp);
102 25 : log<<"Found TOTAL "<<nframes<< " PDB in the file "<<reference.c_str()<<" \n";
103 25 : if(nframes==0) error("at least one frame expected");
104 : //set up rmsdRefClose, initialize it to the first structure loaded from reference file
105 25 : rmsdPosClose.set(pdbv[0], "OPTIMAL");
106 25 : firstPosClose = true;
107 : }
108 25 : if(neigh_stride>0 || neigh_size>0) {
109 14 : if(neigh_size>int(nframes)) {
110 0 : log.printf(" List size required ( %d ) is too large: resizing to the maximum number of frames required: %u \n",neigh_size,nframes);
111 0 : neigh_size=nframes;
112 : }
113 14 : log.printf(" Neighbor list enabled: \n");
114 14 : log.printf(" size : %d elements\n",neigh_size);
115 14 : log.printf(" stride : %d timesteps \n",neigh_stride);
116 : } else {
117 11 : log.printf(" Neighbor list NOT enabled \n");
118 : }
119 25 : if (epsilonClose > 0) {
120 2 : log.printf(" Computing with the close structure, epsilon = %lf\n", epsilonClose);
121 6 : log << " Bibliography " << plumed.cite("Pazurikova J, Krenek A, Spiwok V, Simkova M J. Chem. Phys. 146, 115101 (2017)") << "\n";
122 : }
123 : else {
124 23 : debugClose = 0;
125 23 : logClose = 0;
126 : }
127 25 : if (debugClose)
128 2 : log.printf(" Extensive debug info regarding close structure turned on\n");
129 :
130 25 : rotationRefClose.resize(nframes);
131 25 : savedIndices = std::vector<unsigned>(nframes);
132 :
133 25 : if(nopbc) log.printf(" without periodic boundary conditions\n");
134 24 : else log.printf(" using periodic boundary conditions\n");
135 :
136 25 : }
137 :
138 25 : PathMSDBase::~PathMSDBase() {
139 75 : }
140 :
141 11179 : void PathMSDBase::calculate() {
142 :
143 11179 : if(neigh_size>0 && getExchangeStep()) error("Neighbor lists for this collective variable are not compatible with replica exchange, sorry for that!");
144 :
145 : //log.printf("NOW CALCULATE! \n");
146 :
147 11179 : if(!nopbc) makeWhole();
148 :
149 :
150 : // resize the list to full
151 11179 : if(imgVec.empty()) { // this is the signal that means: recalculate all
152 7164 : imgVec.resize(nframes);
153 : #pragma omp simd
154 7164 : for(unsigned i=0; i<nframes; i++) {
155 300920 : imgVec[i].property=indexvec[i];
156 300920 : imgVec[i].index=i;
157 : }
158 : }
159 :
160 : // THIS IS THE HEAVY PART (RMSD STUFF)
161 11179 : unsigned stride=comm.Get_size();
162 11179 : unsigned rank=comm.Get_rank();
163 11179 : size_t nat=pdbv[0].size();
164 11179 : plumed_assert(nat>0);
165 11179 : plumed_assert(nframes>0);
166 11179 : plumed_assert(imgVec.size()>0);
167 :
168 11179 : std::vector<Tensor> tmp_rotationRefClose(nframes);
169 :
170 11179 : if (epsilonClose > 0) {
171 : //compute rmsd between positions and close structure, save rotation matrix, drotation_drr01
172 1092 : double posclose = rmsdPosClose.calc_Rot_DRotDRr01(getPositions(), rotationPosClose, drotationPosCloseDrr01, true);
173 : //if we compute for the first time or the existing close structure is too far from current structure
174 1092 : if (firstPosClose || (posclose > epsilonClose)) {
175 : //set the current structure as close one for a few next steps
176 16 : if (logClose)
177 16 : log << "PLUMED_CLOSE: new close structure, rmsd pos close " << posclose << "\n";
178 16 : rmsdPosClose.clear();
179 16 : rmsdPosClose.setReference(getPositions());
180 : //as this is a new close structure, we need to save the rotation matrices fitted to the reference structures
181 : // and we need to accurately recalculate for all reference structures
182 16 : computeRefClose = true;
183 16 : imgVec.resize(nframes);
184 688 : for(unsigned i=0; i<nframes; i++) {
185 672 : imgVec[i].property=indexvec[i];
186 672 : imgVec[i].index=i;
187 : }
188 16 : firstPosClose = false;
189 16 : }
190 : else {
191 : //the current structure is pretty close to the close structure, so we use saved rotation matrices to decrease the complexity of rmsd comuptation
192 1076 : if (debugClose)
193 1076 : log << "PLUMED-CLOSE: old close structure, rmsd pos close " << posclose << "\n";
194 1076 : computeRefClose = false;
195 : }
196 : }
197 :
198 11179 : std::vector<double> tmp_distances(imgVec.size(),0.0);
199 : std::vector<Vector> tmp_derivs;
200 : // this array is a merge of all tmp_derivs, so as to allow a single comm.Sum below
201 11179 : std::vector<Vector> tmp_derivs2(imgVec.size()*nat);
202 :
203 : // if imgVec.size() is less than nframes, it means that only some msd will be calculated
204 11179 : if (epsilonClose > 0) {
205 1092 : if (computeRefClose) {
206 : //recompute rotation matrices accurately
207 688 : for(unsigned i=rank; i<imgVec.size(); i+=stride) {
208 672 : tmp_distances[i] = msdv[imgVec[i].index].calc_Rot(getPositions(), tmp_derivs, tmp_rotationRefClose[imgVec[i].index], true);
209 672 : plumed_assert(tmp_derivs.size()==nat);
210 : #pragma omp simd
211 8736 : for(unsigned j=0; j<nat; j++) tmp_derivs2[i*nat+j]=tmp_derivs[j];
212 : }
213 : }
214 : else {
215 : //approximate distance with saved rotation matrices
216 46268 : for(unsigned i=rank; i<imgVec.size(); i+=stride) {
217 45192 : tmp_distances[i] = msdv[imgVec[i].index].calculateWithCloseStructure(getPositions(), tmp_derivs, rotationPosClose, rotationRefClose[imgVec[i].index], drotationPosCloseDrr01, true);
218 45192 : plumed_assert(tmp_derivs.size()==nat);
219 : #pragma omp simd
220 587496 : for(unsigned j=0; j<nat; j++) tmp_derivs2[i*nat+j]=tmp_derivs[j];
221 45192 : if (debugClose) {
222 45192 : double withclose = tmp_distances[i];
223 45192 : RMSD opt;
224 45192 : opt.setType("OPTIMAL");
225 90384 : opt.setReference(msdv[imgVec[i].index].getReference());
226 : std::vector<Vector> ders;
227 45192 : double withoutclose = opt.calculate(getPositions(), ders, true);
228 45192 : float difference = std::abs(withoutclose-withclose);
229 45192 : log<<"PLUMED-CLOSE: difference original "<<withoutclose;
230 45192 : log<<" - with close "<<withclose<<" = "<<difference<<", step "<<getStep()<<", i "<<i<<" imgVec[i].index "<<imgVec[i].index<<"\n";
231 45192 : }
232 : }
233 : }
234 : }
235 : else {
236 : // store temporary local results
237 297781 : for(unsigned i=rank; i<imgVec.size(); i+=stride) {
238 287694 : tmp_distances[i]=msdv[imgVec[i].index].calculate(getPositions(),tmp_derivs,true);
239 287694 : plumed_assert(tmp_derivs.size()==nat);
240 : #pragma omp simd
241 3729822 : for(unsigned j=0; j<nat; j++) tmp_derivs2[i*nat+j]=tmp_derivs[j];
242 : }
243 : }
244 :
245 : // reduce over all processors
246 11179 : comm.Sum(tmp_distances);
247 11179 : comm.Sum(tmp_derivs2);
248 11179 : if (epsilonClose > 0 && computeRefClose) {
249 16 : comm.Sum(tmp_rotationRefClose);
250 688 : for (unsigned i=0; i<nframes; i++) {
251 672 : rotationRefClose[i] = tmp_rotationRefClose[i];
252 : }
253 : }
254 : // assign imgVec[i].distance and imgVec[i].distder
255 482329 : for(size_t i=0; i<imgVec.size(); i++) {
256 471150 : imgVec[i].distance=tmp_distances[i];
257 471150 : imgVec[i].distder.assign(&tmp_derivs2[i*nat],nat+&tmp_derivs2[i*nat]);
258 : }
259 :
260 : // END OF THE HEAVY PART
261 :
262 : std::vector<Value*> val_s_path;
263 11179 : if(labels.size()>0) {
264 18018 : for(unsigned i=0; i<labels.size(); i++) { val_s_path.push_back(getPntrToComponent(labels[i].c_str()));}
265 : } else {
266 5173 : val_s_path.push_back(getPntrToComponent("sss"));
267 : }
268 22358 : Value* val_z_path=getPntrToComponent("zzz");
269 :
270 28364 : std::vector<double> s_path(val_s_path.size()); for(unsigned i=0; i<s_path.size(); i++)s_path[i]=0.;
271 : double partition=0.;
272 : double tmp;
273 :
274 : // clean vector
275 156302 : for(unsigned i=0; i< derivs_z.size(); i++) {derivs_z[i].zero();}
276 :
277 482329 : for(auto & it : imgVec) {
278 471150 : it.similarity=std::exp(-lambda*(it.distance));
279 1194552 : for(unsigned i=0; i<s_path.size(); i++) {
280 723402 : s_path[i]+=(it.property[i])*it.similarity;
281 : }
282 471150 : partition+=it.similarity;
283 : }
284 28364 : for(unsigned i=0; i<s_path.size(); i++) { s_path[i]/=partition; val_s_path[i]->set(s_path[i]) ;}
285 11179 : val_z_path->set(-(1./lambda)*std::log(partition));
286 28364 : for(unsigned j=0; j<s_path.size(); j++) {
287 : // clean up
288 : #pragma omp simd
289 240386 : for(unsigned i=0; i< derivs_s.size(); i++) {derivs_s[i].zero();}
290 : // do the derivative
291 740587 : for(const auto & it : imgVec) {
292 723402 : double expval=it.similarity;
293 723402 : tmp=lambda*expval*(s_path[j]-it.property[j])/partition;
294 : #pragma omp simd
295 10117428 : for(unsigned i=0; i< derivs_s.size(); i++) { derivs_s[i]+=tmp*it.distder[i] ;}
296 723402 : if(j==0) {
297 : #pragma omp simd
298 6585900 : for(unsigned i=0; i< derivs_z.size(); i++) { derivs_z[i]+=it.distder[i]*expval/partition;}
299 : }
300 : }
301 240386 : for(unsigned i=0; i< derivs_s.size(); i++) {
302 223201 : setAtomsDerivatives (val_s_path[j],i,derivs_s[i]);
303 223201 : if(j==0) {setAtomsDerivatives (val_z_path,i,derivs_z[i]);}
304 : }
305 : }
306 28364 : for(unsigned i=0; i<val_s_path.size(); ++i) setBoxDerivativesNoPbc(val_s_path[i]);
307 11179 : setBoxDerivativesNoPbc(val_z_path);
308 : //
309 : // here set next round neighbors
310 : //
311 11179 : if (neigh_size>0) {
312 : //if( int(getStep())%int(neigh_stride/getTimeStep())==0 ){
313 : // enforce consistency: the stride is in time steps
314 7153 : if( int(getStep())%int(neigh_stride)==0 ) {
315 :
316 : // next round do it all:empty the vector
317 7153 : imgVec.clear();
318 : }
319 : // time to analyze the results:
320 7153 : if(imgVec.size()==nframes) {
321 : //sort by msd
322 0 : sort(imgVec.begin(), imgVec.end(), imgOrderByDist());
323 : //resize
324 0 : imgVec.resize(neigh_size);
325 : }
326 : }
327 : //log.printf("CALCULATION DONE! \n");
328 11179 : }
329 :
330 : }
331 :
332 : }
|