Add function to get residue start and end
[gromacs.git] / src / gromacs / essentialdynamics / edsam.cpp
blobaebd672980e33fa4eadd1843f523a45086755c60
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
7 * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 #include "gmxpre.h"
40 #include "edsam.h"
42 #include <cmath>
43 #include <cstdio>
44 #include <cstring>
45 #include <ctime>
47 #include <memory>
49 #include "gromacs/commandline/filenm.h"
50 #include "gromacs/domdec/domdec_struct.h"
51 #include "gromacs/fileio/gmxfio.h"
52 #include "gromacs/fileio/xvgr.h"
53 #include "gromacs/gmxlib/network.h"
54 #include "gromacs/gmxlib/nrnb.h"
55 #include "gromacs/linearalgebra/nrjac.h"
56 #include "gromacs/math/functions.h"
57 #include "gromacs/math/utilities.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/math/vectypes.h"
60 #include "gromacs/mdlib/broadcaststructs.h"
61 #include "gromacs/mdlib/constr.h"
62 #include "gromacs/mdlib/groupcoord.h"
63 #include "gromacs/mdlib/stat.h"
64 #include "gromacs/mdlib/update.h"
65 #include "gromacs/mdrunutility/handlerestart.h"
66 #include "gromacs/mdtypes/commrec.h"
67 #include "gromacs/mdtypes/edsamhistory.h"
68 #include "gromacs/mdtypes/inputrec.h"
69 #include "gromacs/mdtypes/md_enums.h"
70 #include "gromacs/mdtypes/observableshistory.h"
71 #include "gromacs/mdtypes/state.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/topology/mtop_lookup.h"
74 #include "gromacs/topology/topology.h"
75 #include "gromacs/utility/cstringutil.h"
76 #include "gromacs/utility/fatalerror.h"
77 #include "gromacs/utility/gmxassert.h"
78 #include "gromacs/utility/logger.h"
79 #include "gromacs/utility/smalloc.h"
80 #include "gromacs/utility/strconvert.h"
82 namespace
84 /*! \brief Identifies the type of ED: none, normal ED, flooding. */
85 enum class EssentialDynamicsType
87 //! No essential dynamics
88 None,
89 //! Essential dynamics sampling
90 EDSampling,
91 //! Essential dynamics flooding
92 Flooding
95 /*! \brief Identify on which structure to operate. */
96 enum class EssentialDynamicsStructure
98 //! Reference structure
99 Reference,
100 //! Average Structure
101 Average,
102 //! Origin Structure
103 Origin,
104 //! Target Structure
105 Target
108 /*! \brief Essential dynamics vector.
109 * TODO: split into working data and input data
110 * NOTE: called eigvec, because traditionally eigenvectors from PCA analysis
111 * were used as essential dynamics vectors, however, vectors used for ED need
112 * not have any special properties
114 struct t_eigvec
116 //! nr of eigenvectors
117 int neig = 0;
118 //! index nrs of eigenvectors
119 int* ieig = nullptr;
120 //! stepsizes (per eigenvector)
121 real* stpsz = nullptr;
122 //! eigenvector components
123 rvec** vec = nullptr;
124 //! instantaneous x projections
125 real* xproj = nullptr;
126 //! instantaneous f projections
127 real* fproj = nullptr;
128 //! instantaneous radius
129 real radius = 0.;
130 //! starting or target projections
131 real* refproj = nullptr;
134 /*! \brief Essential dynamics vectors per method implementation.
136 struct t_edvecs
138 //! only monitored, no constraints
139 t_eigvec mon = {};
140 //! fixed linear constraints
141 t_eigvec linfix = {};
142 //! acceptance linear constraints
143 t_eigvec linacc = {};
144 //! fixed radial constraints (exp)
145 t_eigvec radfix = {};
146 //! acceptance radial constraints (exp)
147 t_eigvec radacc = {};
148 //! acceptance rad. contraction constr.
149 t_eigvec radcon = {};
152 /*! \brief Essential dynamics flooding parameters and work data.
153 * TODO: split into working data and input parameters
154 * NOTE: The implementation here follows:
155 * O.E. Lange, L.V. Schafer, and H. Grubmuller,
156 * “Flooding in GROMACS: Accelerated barrier crossings in molecular dynamics,”
157 * J. Comp. Chem., 27 1693–1702 (2006)
159 struct t_edflood
161 /*! \brief Target destabilisation free energy.
162 * Controls the flooding potential strength.
163 * Linked to the expected speed-up of mean first passage time out of flooded minimum */
164 real deltaF0 = 0;
165 //! Do not calculate a flooding potential, instead flood with a constant force
166 bool bConstForce = false;
167 //! Relaxation time scale for slowest flooded degree of freedom
168 real tau = 0;
169 //! Current estimated destabilisation free energy
170 real deltaF = 0;
171 //! Flooding energy, acting as a proportionality factor for the flooding potential
172 real Efl = 0;
173 //! Boltzmann constant times temperature, provided by user
174 real kT = 0;
175 //! The flooding potential
176 real Vfl = 0;
177 //! Integration time step
178 real dt = 0;
179 //! Inital flooding strenth
180 real constEfl = 0;
181 //! Empirical global scaling parameter of the essential dynamics vectors.
182 real alpha2 = 0;
183 //! The forces from flooding in atom coordinate space (in contrast to projections onto essential dynamics vectors)
184 rvec* forces_cartesian = nullptr;
185 //! The vectors along which to flood
186 t_eigvec vecs = {};
187 //! Use flooding for harmonic restraint on the eigenvector
188 bool bHarmonic = false;
189 //! The initial reference projection of the flooding vectors. Only with harmonic restraint.
190 real* initialReferenceProjection = nullptr;
191 //! The current reference projection is the initialReferenceProjection + step * slope. Only with harmonic restraint.
192 real* referenceProjectionSlope = nullptr;
194 } // namespace
197 /* This type is for the average, reference, target, and origin structure */
198 struct gmx_edx
200 int nr = 0; /* number of atoms this structure contains */
201 int nr_loc = 0; /* number of atoms on local node */
202 int* anrs = nullptr; /* atom index numbers */
203 int* anrs_loc = nullptr; /* local atom index numbers */
204 int nalloc_loc = 0; /* allocation size of anrs_loc */
205 int* c_ind = nullptr; /* at which position of the whole anrs
206 * array is a local atom?, i.e.
207 * c_ind[0...nr_loc-1] gives the atom index
208 * with respect to the collective
209 * anrs[0...nr-1] array */
210 rvec* x = nullptr; /* positions for this structure */
211 rvec* x_old = nullptr; /* Last positions which have the correct PBC
212 representation of the ED group. In
213 combination with keeping track of the
214 shift vectors, the ED group can always
215 be made whole */
216 real* m = nullptr; /* masses */
217 real mtot = 0.; /* total mass (only used in sref) */
218 real* sqrtm = nullptr; /* sqrt of the masses used for mass-
219 * weighting of analysis (only used in sav) */
223 typedef struct edpar
225 int nini = 0; /* total Nr of atoms */
226 gmx_bool fitmas = false; /* true if trans fit with cm */
227 gmx_bool pcamas = false; /* true if mass-weighted PCA */
228 int presteps = 0; /* number of steps to run without any
229 * perturbations ... just monitoring */
230 int outfrq = 0; /* freq (in steps) of writing to edo */
231 int maxedsteps = 0; /* max nr of steps per cycle */
233 /* all gmx_edx datasets are copied to all nodes in the parallel case */
234 gmx_edx sref = {}; /* reference positions, to these fitting
235 * will be done */
236 gmx_bool bRefEqAv = false; /* If true, reference & average indices
237 * are the same. Used for optimization */
238 gmx_edx sav = {}; /* average positions */
239 gmx_edx star = {}; /* target positions */
240 gmx_edx sori = {}; /* origin positions */
242 t_edvecs vecs = {}; /* eigenvectors */
243 real slope = 0; /* minimal slope in acceptance radexp */
245 t_edflood flood = {}; /* parameters especially for flooding */
246 struct t_ed_buffer* buf = nullptr; /* handle to local buffers */
247 } t_edpar;
250 struct gmx_edsam
252 ~gmx_edsam();
253 //! The type of ED
254 EssentialDynamicsType eEDtype = EssentialDynamicsType::None;
255 //! output file pointer
256 FILE* edo = nullptr;
257 std::vector<t_edpar> edpar;
258 gmx_bool bFirst = false;
260 gmx_edsam::~gmx_edsam()
262 if (edo)
264 /* edo is opened sometimes with xvgropen, sometimes with
265 * gmx_fio_fopen, so we use the least common denominator for
266 * closing. */
267 gmx_fio_fclose(edo);
271 struct t_do_edsam
273 matrix old_rotmat;
274 real oldrad;
275 rvec old_transvec, older_transvec, transvec_compact;
276 rvec* xcoll; /* Positions from all nodes, this is the
277 collective set we work on.
278 These are the positions of atoms with
279 average structure indices */
280 rvec* xc_ref; /* same but with reference structure indices */
281 ivec* shifts_xcoll; /* Shifts for xcoll */
282 ivec* extra_shifts_xcoll; /* xcoll shift changes since last NS step */
283 ivec* shifts_xc_ref; /* Shifts for xc_ref */
284 ivec* extra_shifts_xc_ref; /* xc_ref shift changes since last NS step */
285 gmx_bool bUpdateShifts; /* TRUE in NS steps to indicate that the
286 ED shifts for this ED group need to
287 be updated */
291 /* definition of ED buffer structure */
292 struct t_ed_buffer
294 struct t_fit_to_ref* fit_to_ref;
295 struct t_do_edfit* do_edfit;
296 struct t_do_edsam* do_edsam;
297 struct t_do_radcon* do_radcon;
300 namespace gmx
302 class EssentialDynamics::Impl
304 public:
305 // TODO: move all fields from gmx_edsam here and remove gmx_edsam
306 gmx_edsam essentialDynamics_;
308 EssentialDynamics::EssentialDynamics() : impl_(new Impl) {}
309 EssentialDynamics::~EssentialDynamics() = default;
311 gmx_edsam* EssentialDynamics::getLegacyED()
313 return &impl_->essentialDynamics_;
315 } // namespace gmx
317 /* Function declarations */
318 static void fit_to_reference(rvec* xcoll, rvec transvec, matrix rotmat, t_edpar* edi);
319 static void translate_and_rotate(rvec* x, int nat, rvec transvec, matrix rotmat);
320 static real rmsd_from_structure(rvec* x, struct gmx_edx* s);
321 namespace
323 /*! \brief Read in the essential dynamics input file and return its contents.
324 * \param[in] fn the file name of the edi file to be read
325 * \param[in] nr_mdatoms the number of atoms in the simulation
326 * \returns A vector of containing the essentiail dyanmics parameters.
327 * NOTE: edi files may that it may contain several ED data sets from concatenated edi files.
328 * The standard case would be a single ED data set, though. */
329 std::vector<t_edpar> read_edi_file(const char* fn, int nr_mdatoms);
330 } // namespace
331 static void crosscheck_edi_file_vs_checkpoint(const gmx_edsam& ed, edsamhistory_t* EDstate);
332 static void init_edsamstate(const gmx_edsam& ed, edsamhistory_t* EDstate);
333 static void write_edo_legend(gmx_edsam* ed, int nED, const gmx_output_env_t* oenv);
334 /* End function declarations */
336 /* Define formats for the column width in the output file */
337 const char EDcol_sfmt[] = "%17s";
338 const char EDcol_efmt[] = "%17.5e";
339 const char EDcol_ffmt[] = "%17f";
341 /* Define a formats for reading, otherwise cppcheck complains for scanf without width limits */
342 const char max_ev_fmt_d[] = "%7d"; /* Number of eigenvectors. 9,999,999 should be enough */
343 const char max_ev_fmt_lf[] = "%12lf";
344 const char max_ev_fmt_dlf[] = "%7d%12lf";
345 const char max_ev_fmt_dlflflf[] = "%7d%12lf%12lf%12lf";
346 const char max_ev_fmt_lelele[] = "%12le%12le%12le";
348 namespace
350 /*!\brief Determines whether to perform essential dynamics constraints other than flooding.
351 * \param[in] edi the essential dynamics parameters
352 * \returns true if essential dyanmics constraints need to be performed
354 bool bNeedDoEdsam(const t_edpar& edi)
356 return (edi.vecs.mon.neig != 0) || (edi.vecs.linfix.neig != 0) || (edi.vecs.linacc.neig != 0)
357 || (edi.vecs.radfix.neig != 0) || (edi.vecs.radacc.neig != 0) || (edi.vecs.radcon.neig != 0);
359 } // namespace
362 /* Multiple ED groups will be labeled with letters instead of numbers
363 * to avoid confusion with eigenvector indices */
364 static char get_EDgroupChar(int nr_edi, int nED)
366 if (nED == 1)
368 return ' ';
371 /* nr_edi = 1 -> A
372 * nr_edi = 2 -> B ...
374 return 'A' + nr_edi - 1;
377 namespace
379 /*! \brief The mass-weighted inner product of two coordinate vectors.
380 * Does not subtract average positions, projection on single eigenvector is returned
381 * used by: do_linfix, do_linacc, do_radfix, do_radacc, do_radcon
382 * Average position is subtracted in ed_apply_constraints prior to calling projectx
383 * \param[in] edi Essential dynamics parameters
384 * \param[in] xcoll vector of atom coordinates
385 * \param[in] vec vector of coordinates to project onto
386 * \return mass-weighted projection.
388 real projectx(const t_edpar& edi, rvec* xcoll, rvec* vec)
390 int i;
391 real proj = 0.0;
394 for (i = 0; i < edi.sav.nr; i++)
396 proj += edi.sav.sqrtm[i] * iprod(vec[i], xcoll[i]);
399 return proj;
401 /*!\brief Project coordinates onto vector after substracting average position.
402 * projection is stored in vec->refproj which is used for radacc, radfix,
403 * radcon and center of flooding potential.
404 * Average positions are first substracted from x, then added back again.
405 * \param[in] edi essential dynamics parameters with average position
406 * \param[in] x Coordinates to be projected
407 * \param[out] vec eigenvector, radius and refproj are overwritten here
409 void rad_project(const t_edpar& edi, rvec* x, t_eigvec* vec)
411 int i;
412 real rad = 0.0;
414 /* Subtract average positions */
415 for (i = 0; i < edi.sav.nr; i++)
417 rvec_dec(x[i], edi.sav.x[i]);
420 for (i = 0; i < vec->neig; i++)
422 vec->refproj[i] = projectx(edi, x, vec->vec[i]);
423 rad += gmx::square((vec->refproj[i] - vec->xproj[i]));
425 vec->radius = sqrt(rad);
427 /* Add average positions */
428 for (i = 0; i < edi.sav.nr; i++)
430 rvec_inc(x[i], edi.sav.x[i]);
434 /*!\brief Projects coordinates onto eigenvectors and stores result in vec->xproj.
435 * Mass-weighting is applied. Subtracts average positions prior to projection and add
436 * them afterwards to retain the unchanged vector.
437 * \param[in] x The coordinates to project to an eigenvector
438 * \param[in,out] vec The eigenvectors
439 * \param[in] edi essential dynamics parameters holding average structure and masses
441 void project_to_eigvectors(rvec* x, t_eigvec* vec, const t_edpar& edi)
443 if (!vec->neig)
445 return;
448 /* Subtract average positions */
449 for (int i = 0; i < edi.sav.nr; i++)
451 rvec_dec(x[i], edi.sav.x[i]);
454 for (int i = 0; i < vec->neig; i++)
456 vec->xproj[i] = projectx(edi, x, vec->vec[i]);
459 /* Add average positions */
460 for (int i = 0; i < edi.sav.nr; i++)
462 rvec_inc(x[i], edi.sav.x[i]);
465 } // namespace
467 /* Project vector x onto all edi->vecs (mon, linfix,...) */
468 static void project(rvec* x, /* positions to project */
469 t_edpar* edi) /* edi data set */
471 /* It is not more work to subtract the average position in every
472 * subroutine again, because these routines are rarely used simultaneously */
473 project_to_eigvectors(x, &edi->vecs.mon, *edi);
474 project_to_eigvectors(x, &edi->vecs.linfix, *edi);
475 project_to_eigvectors(x, &edi->vecs.linacc, *edi);
476 project_to_eigvectors(x, &edi->vecs.radfix, *edi);
477 project_to_eigvectors(x, &edi->vecs.radacc, *edi);
478 project_to_eigvectors(x, &edi->vecs.radcon, *edi);
481 namespace
483 /*!\brief Evaluates the distance from reference to current eigenvector projection.
484 * \param[in] vec eigenvector
485 * \returns distance
487 real calc_radius(const t_eigvec& vec)
489 real rad = 0.0;
491 for (int i = 0; i < vec.neig; i++)
493 rad += gmx::square((vec.refproj[i] - vec.xproj[i]));
496 return rad = sqrt(rad);
498 } // namespace
500 struct t_do_edfit
502 double** omega;
503 double** om;
506 static void do_edfit(int natoms, rvec* xp, rvec* x, matrix R, t_edpar* edi)
508 /* this is a copy of do_fit with some modifications */
509 int c, r, n, j, i, irot;
510 double d[6], xnr, xpc;
511 matrix vh, vk, u;
512 int index;
513 real max_d;
515 struct t_do_edfit* loc;
516 gmx_bool bFirst;
518 if (edi->buf->do_edfit != nullptr)
520 bFirst = FALSE;
522 else
524 bFirst = TRUE;
525 snew(edi->buf->do_edfit, 1);
527 loc = edi->buf->do_edfit;
529 if (bFirst)
531 snew(loc->omega, 2 * DIM);
532 snew(loc->om, 2 * DIM);
533 for (i = 0; i < 2 * DIM; i++)
535 snew(loc->omega[i], 2 * DIM);
536 snew(loc->om[i], 2 * DIM);
540 for (i = 0; (i < 6); i++)
542 d[i] = 0;
543 for (j = 0; (j < 6); j++)
545 loc->omega[i][j] = 0;
546 loc->om[i][j] = 0;
550 /* calculate the matrix U */
551 clear_mat(u);
552 for (n = 0; (n < natoms); n++)
554 for (c = 0; (c < DIM); c++)
556 xpc = xp[n][c];
557 for (r = 0; (r < DIM); r++)
559 xnr = x[n][r];
560 u[c][r] += xnr * xpc;
565 /* construct loc->omega */
566 /* loc->omega is symmetric -> loc->omega==loc->omega' */
567 for (r = 0; (r < 6); r++)
569 for (c = 0; (c <= r); c++)
571 if ((r >= 3) && (c < 3))
573 loc->omega[r][c] = u[r - 3][c];
574 loc->omega[c][r] = u[r - 3][c];
576 else
578 loc->omega[r][c] = 0;
579 loc->omega[c][r] = 0;
584 /* determine h and k */
585 jacobi(loc->omega, 6, d, loc->om, &irot);
587 if (irot == 0)
589 fprintf(stderr, "IROT=0\n");
592 index = 0; /* For the compiler only */
594 for (j = 0; (j < 3); j++)
596 max_d = -1000;
597 for (i = 0; (i < 6); i++)
599 if (d[i] > max_d)
601 max_d = d[i];
602 index = i;
605 d[index] = -10000;
606 for (i = 0; (i < 3); i++)
608 vh[j][i] = M_SQRT2 * loc->om[i][index];
609 vk[j][i] = M_SQRT2 * loc->om[i + DIM][index];
613 /* determine R */
614 for (c = 0; (c < 3); c++)
616 for (r = 0; (r < 3); r++)
618 R[c][r] = vk[0][r] * vh[0][c] + vk[1][r] * vh[1][c] + vk[2][r] * vh[2][c];
621 if (det(R) < 0)
623 for (c = 0; (c < 3); c++)
625 for (r = 0; (r < 3); r++)
627 R[c][r] = vk[0][r] * vh[0][c] + vk[1][r] * vh[1][c] - vk[2][r] * vh[2][c];
634 static void rmfit(int nat, rvec* xcoll, const rvec transvec, matrix rotmat)
636 rvec vec;
637 matrix tmat;
640 /* Remove rotation.
641 * The inverse rotation is described by the transposed rotation matrix */
642 transpose(rotmat, tmat);
643 rotate_x(xcoll, nat, tmat);
645 /* Remove translation */
646 vec[XX] = -transvec[XX];
647 vec[YY] = -transvec[YY];
648 vec[ZZ] = -transvec[ZZ];
649 translate_x(xcoll, nat, vec);
653 /**********************************************************************************
654 ******************** FLOODING ****************************************************
655 **********************************************************************************
657 The flooding ability was added later to edsam. Many of the edsam functionality could be reused for that purpose.
658 The flooding covariance matrix, i.e. the selected eigenvectors and their corresponding eigenvalues are
659 read as 7th Component Group. The eigenvalues are coded into the stepsize parameter (as used by -linfix or -linacc).
661 do_md clls right in the beginning the function init_edsam, which reads the edi file, saves all the necessary information in
662 the edi structure and calls init_flood, to initialise some extra fields in the edi->flood structure.
664 since the flooding acts on forces do_flood is called from the function force() (force.c), while the other
665 edsam functionality is hooked into md via the update() (update.c) function acting as constraint on positions.
667 do_flood makes a copy of the positions,
668 fits them, projects them computes flooding_energy, and flooding forces. The forces are computed in the
669 space of the eigenvectors and are then blown up to the full cartesian space and rotated back to remove the
670 fit. Then do_flood adds these forces to the forcefield-forces
671 (given as parameter) and updates the adaptive flooding parameters Efl and deltaF.
673 To center the flooding potential at a different location one can use the -ori option in make_edi. The ori
674 structure is projected to the system of eigenvectors and then this position in the subspace is used as
675 center of the flooding potential. If the option is not used, the center will be zero in the subspace,
676 i.e. the average structure as given in the make_edi file.
678 To use the flooding potential as restraint, make_edi has the option -restrain, which leads to inverted
679 signs of alpha2 and Efl, such that the sign in the exponential of Vfl is not inverted but the sign of
680 Vfl is inverted. Vfl = Efl * exp (- .../Efl/alpha2*x^2...) With tau>0 the negative Efl will grow slowly
681 so that the restraint is switched off slowly. When Efl==0 and inverted flooding is ON is reached no
682 further adaption is applied, Efl will stay constant at zero.
684 To use restraints with harmonic potentials switch -restrain and -harmonic. Then the eigenvalues are
685 used as spring constants for the harmonic potential.
686 Note that eq3 in the flooding paper (J. Comp. Chem. 2006, 27, 1693-1702) defines the parameter lambda \
687 as the inverse of the spring constant, whereas the implementation uses lambda as the spring constant.
689 To use more than one flooding matrix just concatenate several .edi files (cat flood1.edi flood2.edi > flood_all.edi)
690 the routine read_edi_file reads all of theses flooding files.
691 The structure t_edi is now organized as a list of t_edis and the function do_flood cycles through the list
692 calling the do_single_flood() routine for every single entry. Since every state variables have been kept in one
693 edi there is no interdependence whatsoever. The forces are added together.
695 To write energies into the .edr file, call the function
696 get_flood_enx_names(char**, int *nnames) to get the Header (Vfl1 Vfl2... Vfln)
697 and call
698 get_flood_energies(real Vfl[],int nnames);
700 TODO:
701 - one could program the whole thing such that Efl, Vfl and deltaF is written to the .edr file. -- i dont know how to do that, yet.
703 Maybe one should give a range of atoms for which to remove motion, so that motion is removed with
704 two edsam files from two peptide chains
707 // TODO split this into multiple files
708 namespace
710 /*!\brief Output flooding simulation settings and results to file.
711 * \param[in] edi Essential dynamics input parameters
712 * \param[in] fp output file
713 * \param[in] rmsd rmsd to reference structure
715 void write_edo_flood(const t_edpar& edi, FILE* fp, real rmsd)
717 /* Output how well we fit to the reference structure */
718 fprintf(fp, EDcol_ffmt, rmsd);
720 for (int i = 0; i < edi.flood.vecs.neig; i++)
722 fprintf(fp, EDcol_efmt, edi.flood.vecs.xproj[i]);
724 /* Check whether the reference projection changes with time (this can happen
725 * in case flooding is used as harmonic restraint). If so, output the
726 * current reference projection */
727 if (edi.flood.bHarmonic && edi.flood.referenceProjectionSlope[i] != 0.0)
729 fprintf(fp, EDcol_efmt, edi.flood.vecs.refproj[i]);
732 /* Output Efl if we are doing adaptive flooding */
733 if (0 != edi.flood.tau)
735 fprintf(fp, EDcol_efmt, edi.flood.Efl);
737 fprintf(fp, EDcol_efmt, edi.flood.Vfl);
739 /* Output deltaF if we are doing adaptive flooding */
740 if (0 != edi.flood.tau)
742 fprintf(fp, EDcol_efmt, edi.flood.deltaF);
744 fprintf(fp, EDcol_efmt, edi.flood.vecs.fproj[i]);
747 } // namespace
750 /* From flood.xproj compute the Vfl(x) at this point */
751 static real flood_energy(t_edpar* edi, int64_t step)
753 /* compute flooding energy Vfl
754 Vfl = Efl * exp( - \frac {kT} {2Efl alpha^2} * sum_i { \lambda_i c_i^2 } )
755 \lambda_i is the reciprocal eigenvalue 1/\sigma_i
756 it is already computed by make_edi and stored in stpsz[i]
757 bHarmonic:
758 Vfl = - Efl * 1/2(sum _i {\frac 1{\lambda_i} c_i^2})
760 real sum;
761 real Vfl;
762 int i;
765 /* Each time this routine is called (i.e. each time step), we add a small
766 * value to the reference projection. This way a harmonic restraint towards
767 * a moving reference is realized. If no value for the additive constant
768 * is provided in the edi file, the reference will not change. */
769 if (edi->flood.bHarmonic)
771 for (i = 0; i < edi->flood.vecs.neig; i++)
773 edi->flood.vecs.refproj[i] = edi->flood.initialReferenceProjection[i]
774 + step * edi->flood.referenceProjectionSlope[i];
778 sum = 0.0;
779 /* Compute sum which will be the exponent of the exponential */
780 for (i = 0; i < edi->flood.vecs.neig; i++)
782 /* stpsz stores the reciprocal eigenvalue 1/sigma_i */
783 sum += edi->flood.vecs.stpsz[i] * (edi->flood.vecs.xproj[i] - edi->flood.vecs.refproj[i])
784 * (edi->flood.vecs.xproj[i] - edi->flood.vecs.refproj[i]);
787 /* Compute the Gauss function*/
788 if (edi->flood.bHarmonic)
790 Vfl = -0.5 * edi->flood.Efl * sum; /* minus sign because Efl is negative, if restrain is on. */
792 else
794 Vfl = edi->flood.Efl != 0
795 ? edi->flood.Efl
796 * std::exp(-edi->flood.kT / 2 / edi->flood.Efl / edi->flood.alpha2 * sum)
797 : 0;
800 return Vfl;
804 /* From the position and from Vfl compute forces in subspace -> store in edi->vec.flood.fproj */
805 static void flood_forces(t_edpar* edi)
807 /* compute the forces in the subspace of the flooding eigenvectors
808 * by the formula F_i= V_{fl}(c) * ( \frac {kT} {E_{fl}} \lambda_i c_i */
810 int i;
811 real energy = edi->flood.Vfl;
814 if (edi->flood.bHarmonic)
816 for (i = 0; i < edi->flood.vecs.neig; i++)
818 edi->flood.vecs.fproj[i] = edi->flood.Efl * edi->flood.vecs.stpsz[i]
819 * (edi->flood.vecs.xproj[i] - edi->flood.vecs.refproj[i]);
822 else
824 for (i = 0; i < edi->flood.vecs.neig; i++)
826 /* if Efl is zero the forces are zero if not use the formula */
827 edi->flood.vecs.fproj[i] =
828 edi->flood.Efl != 0
829 ? edi->flood.kT / edi->flood.Efl / edi->flood.alpha2 * energy
830 * edi->flood.vecs.stpsz[i]
831 * (edi->flood.vecs.xproj[i] - edi->flood.vecs.refproj[i])
832 : 0;
837 namespace
839 /*!\brief Raise forces from subspace into cartesian space.
840 * This function lifts the forces from the subspace to the cartesian space
841 * all the values not contained in the subspace are assumed to be zero and then
842 * a coordinate transformation from eigenvector to cartesian vectors is performed
843 * The nonexistent values don't have to be set to zero explicitly, they would occur
844 * as zero valued summands, hence we just stop to compute this part of the sum.
845 * For every atom we add all the contributions to this atom from all the different eigenvectors.
846 * NOTE: one could add directly to the forcefield forces, would mean we wouldn't have to clear the
847 * field forces_cart prior the computation, but we compute the forces separately
848 * to have them accessible for diagnostics
850 * \param[in] edi Essential dynamics input parameters
851 * \param[out] forces_cart The cartesian forces
854 void flood_blowup(const t_edpar& edi, rvec* forces_cart)
856 const real* forces_sub = edi.flood.vecs.fproj;
857 /* Calculate the cartesian forces for the local atoms */
859 /* Clear forces first */
860 for (int j = 0; j < edi.sav.nr_loc; j++)
862 clear_rvec(forces_cart[j]);
865 /* Now compute atomwise */
866 for (int j = 0; j < edi.sav.nr_loc; j++)
868 /* Compute forces_cart[edi.sav.anrs[j]] */
869 for (int eig = 0; eig < edi.flood.vecs.neig; eig++)
871 rvec addedForce;
872 /* Force vector is force * eigenvector (compute only atom j) */
873 svmul(forces_sub[eig], edi.flood.vecs.vec[eig][edi.sav.c_ind[j]], addedForce);
874 /* Add this vector to the cartesian forces */
875 rvec_inc(forces_cart[j], addedForce);
880 } // namespace
883 /* Update the values of Efl, deltaF depending on tau and Vfl */
884 static void update_adaption(t_edpar* edi)
886 /* this function updates the parameter Efl and deltaF according to the rules given in
887 * 'predicting unimolecular chemical reactions: chemical flooding' M Mueller et al,
888 * J. chem Phys. */
890 if ((edi->flood.tau < 0 ? -edi->flood.tau : edi->flood.tau) > 0.00000001)
892 edi->flood.Efl = edi->flood.Efl
893 + edi->flood.dt / edi->flood.tau * (edi->flood.deltaF0 - edi->flood.deltaF);
894 /* check if restrain (inverted flooding) -> don't let EFL become positive */
895 if (edi->flood.alpha2 < 0 && edi->flood.Efl > -0.00000001)
897 edi->flood.Efl = 0;
900 edi->flood.deltaF = (1 - edi->flood.dt / edi->flood.tau) * edi->flood.deltaF
901 + edi->flood.dt / edi->flood.tau * edi->flood.Vfl;
906 static void do_single_flood(FILE* edo,
907 const rvec x[],
908 rvec force[],
909 t_edpar* edi,
910 int64_t step,
911 const matrix box,
912 const t_commrec* cr,
913 gmx_bool bNS) /* Are we in a neighbor searching step? */
915 int i;
916 matrix rotmat; /* rotation matrix */
917 matrix tmat; /* inverse rotation */
918 rvec transvec; /* translation vector */
919 real rmsdev;
920 struct t_do_edsam* buf;
923 buf = edi->buf->do_edsam;
926 /* Broadcast the positions of the AVERAGE structure such that they are known on
927 * every processor. Each node contributes its local positions x and stores them in
928 * the collective ED array buf->xcoll */
929 communicate_group_positions(cr, buf->xcoll, buf->shifts_xcoll, buf->extra_shifts_xcoll, bNS, x,
930 edi->sav.nr, edi->sav.nr_loc, edi->sav.anrs_loc, edi->sav.c_ind,
931 edi->sav.x_old, box);
933 /* Only assembly REFERENCE positions if their indices differ from the average ones */
934 if (!edi->bRefEqAv)
936 communicate_group_positions(cr, buf->xc_ref, buf->shifts_xc_ref, buf->extra_shifts_xc_ref,
937 bNS, x, edi->sref.nr, edi->sref.nr_loc, edi->sref.anrs_loc,
938 edi->sref.c_ind, edi->sref.x_old, box);
941 /* If bUpdateShifts was TRUE, the shifts have just been updated in get_positions.
942 * We do not need to update the shifts until the next NS step */
943 buf->bUpdateShifts = FALSE;
945 /* Now all nodes have all of the ED/flooding positions in edi->sav->xcoll,
946 * as well as the indices in edi->sav.anrs */
948 /* Fit the reference indices to the reference structure */
949 if (edi->bRefEqAv)
951 fit_to_reference(buf->xcoll, transvec, rotmat, edi);
953 else
955 fit_to_reference(buf->xc_ref, transvec, rotmat, edi);
958 /* Now apply the translation and rotation to the ED structure */
959 translate_and_rotate(buf->xcoll, edi->sav.nr, transvec, rotmat);
961 /* Project fitted structure onto supbspace -> store in edi->flood.vecs.xproj */
962 project_to_eigvectors(buf->xcoll, &edi->flood.vecs, *edi);
964 if (!edi->flood.bConstForce)
966 /* Compute Vfl(x) from flood.xproj */
967 edi->flood.Vfl = flood_energy(edi, step);
969 update_adaption(edi);
971 /* Compute the flooding forces */
972 flood_forces(edi);
975 /* Translate them into cartesian positions */
976 flood_blowup(*edi, edi->flood.forces_cartesian);
978 /* Rotate forces back so that they correspond to the given structure and not to the fitted one */
979 /* Each node rotates back its local forces */
980 transpose(rotmat, tmat);
981 rotate_x(edi->flood.forces_cartesian, edi->sav.nr_loc, tmat);
983 /* Finally add forces to the main force variable */
984 for (i = 0; i < edi->sav.nr_loc; i++)
986 rvec_inc(force[edi->sav.anrs_loc[i]], edi->flood.forces_cartesian[i]);
989 /* Output is written by the master process */
990 if (do_per_step(step, edi->outfrq) && MASTER(cr))
992 /* Output how well we fit to the reference */
993 if (edi->bRefEqAv)
995 /* Indices of reference and average structures are identical,
996 * thus we can calculate the rmsd to SREF using xcoll */
997 rmsdev = rmsd_from_structure(buf->xcoll, &edi->sref);
999 else
1001 /* We have to translate & rotate the reference atoms first */
1002 translate_and_rotate(buf->xc_ref, edi->sref.nr, transvec, rotmat);
1003 rmsdev = rmsd_from_structure(buf->xc_ref, &edi->sref);
1006 write_edo_flood(*edi, edo, rmsdev);
1011 /* Main flooding routine, called from do_force */
1012 extern void do_flood(const t_commrec* cr,
1013 const t_inputrec* ir,
1014 const rvec x[],
1015 rvec force[],
1016 gmx_edsam* ed,
1017 const matrix box,
1018 int64_t step,
1019 gmx_bool bNS)
1021 /* Write time to edo, when required. Output the time anyhow since we need
1022 * it in the output file for ED constraints. */
1023 if (MASTER(cr) && do_per_step(step, ed->edpar.begin()->outfrq))
1025 fprintf(ed->edo, "\n%12f", ir->init_t + step * ir->delta_t);
1028 if (ed->eEDtype != EssentialDynamicsType::Flooding)
1030 return;
1033 for (auto& edi : ed->edpar)
1035 /* Call flooding for one matrix */
1036 if (edi.flood.vecs.neig)
1038 do_single_flood(ed->edo, x, force, &edi, step, box, cr, bNS);
1044 /* Called by init_edi, configure some flooding related variables and structures,
1045 * print headers to output files */
1046 static void init_flood(t_edpar* edi, gmx_edsam* ed, real dt)
1048 int i;
1051 edi->flood.Efl = edi->flood.constEfl;
1052 edi->flood.Vfl = 0;
1053 edi->flood.dt = dt;
1055 if (edi->flood.vecs.neig)
1057 /* If in any of the ED groups we find a flooding vector, flooding is turned on */
1058 ed->eEDtype = EssentialDynamicsType::Flooding;
1060 fprintf(stderr, "ED: Flooding %d eigenvector%s.\n", edi->flood.vecs.neig,
1061 edi->flood.vecs.neig > 1 ? "s" : "");
1063 if (edi->flood.bConstForce)
1065 /* We have used stpsz as a vehicle to carry the fproj values for constant
1066 * force flooding. Now we copy that to flood.vecs.fproj. Note that
1067 * in const force flooding, fproj is never changed. */
1068 for (i = 0; i < edi->flood.vecs.neig; i++)
1070 edi->flood.vecs.fproj[i] = edi->flood.vecs.stpsz[i];
1072 fprintf(stderr, "ED: applying on eigenvector %d a constant force of %g\n",
1073 edi->flood.vecs.ieig[i], edi->flood.vecs.fproj[i]);
1080 #ifdef DEBUGHELPERS
1081 /*********** Energy book keeping ******/
1082 static void get_flood_enx_names(t_edpar* edi, char** names, int* nnames) /* get header of energies */
1084 t_edpar* actual;
1085 int count;
1086 char buf[STRLEN];
1087 actual = edi;
1088 count = 1;
1089 while (actual)
1091 srenew(names, count);
1092 sprintf(buf, "Vfl_%d", count);
1093 names[count - 1] = gmx_strdup(buf);
1094 actual = actual->next_edi;
1095 count++;
1097 *nnames = count - 1;
1101 static void get_flood_energies(t_edpar* edi, real Vfl[], int nnames)
1103 /*fl has to be big enough to capture nnames-many entries*/
1104 t_edpar* actual;
1105 int count;
1108 actual = edi;
1109 count = 1;
1110 while (actual)
1112 Vfl[count - 1] = actual->flood.Vfl;
1113 actual = actual->next_edi;
1114 count++;
1116 if (nnames != count - 1)
1118 gmx_fatal(FARGS, "Number of energies is not consistent with t_edi structure");
1121 /************* END of FLOODING IMPLEMENTATION ****************************/
1122 #endif
1125 /* This function opens the ED input and output files, reads in all datasets it finds
1126 * in the input file, and cross-checks whether the .edi file information is consistent
1127 * with the essential dynamics data found in the checkpoint file (if present).
1128 * gmx make_edi can be used to create an .edi input file.
1130 static std::unique_ptr<gmx::EssentialDynamics> ed_open(int natoms,
1131 ObservablesHistory* oh,
1132 const char* ediFileName,
1133 const char* edoFileName,
1134 const gmx::StartingBehavior startingBehavior,
1135 const gmx_output_env_t* oenv,
1136 const t_commrec* cr)
1138 auto edHandle = std::make_unique<gmx::EssentialDynamics>();
1139 auto ed = edHandle->getLegacyED();
1140 /* We want to perform ED (this switch might later be upgraded to EssentialDynamicsType::Flooding) */
1141 ed->eEDtype = EssentialDynamicsType::EDSampling;
1143 if (MASTER(cr))
1145 // If we start from a checkpoint file, we already have an edsamHistory struct
1146 if (oh->edsamHistory == nullptr)
1148 oh->edsamHistory = std::make_unique<edsamhistory_t>(edsamhistory_t{});
1150 edsamhistory_t* EDstate = oh->edsamHistory.get();
1152 /* Read the edi input file: */
1153 ed->edpar = read_edi_file(ediFileName, natoms);
1155 /* Make sure the checkpoint was produced in a run using this .edi file */
1156 if (EDstate->bFromCpt)
1158 crosscheck_edi_file_vs_checkpoint(*ed, EDstate);
1160 else
1162 EDstate->nED = ed->edpar.size();
1164 init_edsamstate(*ed, EDstate);
1166 /* The master opens the ED output file */
1167 if (startingBehavior == gmx::StartingBehavior::RestartWithAppending)
1169 ed->edo = gmx_fio_fopen(edoFileName, "a+");
1171 else
1173 ed->edo = xvgropen(edoFileName, "Essential dynamics / flooding output", "Time (ps)",
1174 "RMSDs (nm), projections on EVs (nm), ...", oenv);
1176 /* Make a descriptive legend */
1177 write_edo_legend(ed, EDstate->nED, oenv);
1180 return edHandle;
1184 /* Broadcasts the structure data */
1185 static void bc_ed_positions(const t_commrec* cr, struct gmx_edx* s, EssentialDynamicsStructure stype)
1187 snew_bc(MASTER(cr), s->anrs, s->nr); /* Index numbers */
1188 snew_bc(MASTER(cr), s->x, s->nr); /* Positions */
1189 nblock_bc(cr->mpi_comm_mygroup, s->nr, s->anrs);
1190 nblock_bc(cr->mpi_comm_mygroup, s->nr, s->x);
1192 /* For the average & reference structures we need an array for the collective indices,
1193 * and we need to broadcast the masses as well */
1194 if (stype == EssentialDynamicsStructure::Average || stype == EssentialDynamicsStructure::Reference)
1196 /* We need these additional variables in the parallel case: */
1197 snew(s->c_ind, s->nr); /* Collective indices */
1198 /* Local atom indices get assigned in dd_make_local_group_indices.
1199 * There, also memory is allocated */
1200 s->nalloc_loc = 0; /* allocation size of s->anrs_loc */
1201 snew_bc(MASTER(cr), s->x_old, s->nr); /* To be able to always make the ED molecule whole, ... */
1202 nblock_bc(cr->mpi_comm_mygroup, s->nr,
1203 s->x_old); /* ... keep track of shift changes with the help of old coords */
1206 /* broadcast masses for the reference structure (for mass-weighted fitting) */
1207 if (stype == EssentialDynamicsStructure::Reference)
1209 snew_bc(MASTER(cr), s->m, s->nr);
1210 nblock_bc(cr->mpi_comm_mygroup, s->nr, s->m);
1213 /* For the average structure we might need the masses for mass-weighting */
1214 if (stype == EssentialDynamicsStructure::Average)
1216 snew_bc(MASTER(cr), s->sqrtm, s->nr);
1217 nblock_bc(cr->mpi_comm_mygroup, s->nr, s->sqrtm);
1218 snew_bc(MASTER(cr), s->m, s->nr);
1219 nblock_bc(cr->mpi_comm_mygroup, s->nr, s->m);
1224 /* Broadcasts the eigenvector data */
1225 static void bc_ed_vecs(const t_commrec* cr, t_eigvec* ev, int length)
1227 int i;
1229 snew_bc(MASTER(cr), ev->ieig, ev->neig); /* index numbers of eigenvector */
1230 snew_bc(MASTER(cr), ev->stpsz, ev->neig); /* stepsizes per eigenvector */
1231 snew_bc(MASTER(cr), ev->xproj, ev->neig); /* instantaneous x projection */
1232 snew_bc(MASTER(cr), ev->fproj, ev->neig); /* instantaneous f projection */
1233 snew_bc(MASTER(cr), ev->refproj, ev->neig); /* starting or target projection */
1235 nblock_bc(cr->mpi_comm_mygroup, ev->neig, ev->ieig);
1236 nblock_bc(cr->mpi_comm_mygroup, ev->neig, ev->stpsz);
1237 nblock_bc(cr->mpi_comm_mygroup, ev->neig, ev->xproj);
1238 nblock_bc(cr->mpi_comm_mygroup, ev->neig, ev->fproj);
1239 nblock_bc(cr->mpi_comm_mygroup, ev->neig, ev->refproj);
1241 snew_bc(MASTER(cr), ev->vec, ev->neig); /* Eigenvector components */
1242 for (i = 0; i < ev->neig; i++)
1244 snew_bc(MASTER(cr), ev->vec[i], length);
1245 nblock_bc(cr->mpi_comm_mygroup, length, ev->vec[i]);
1250 /* Broadcasts the ED / flooding data to other nodes
1251 * and allocates memory where needed */
1252 static void broadcast_ed_data(const t_commrec* cr, gmx_edsam* ed)
1254 /* Master lets the other nodes know if its ED only or also flooding */
1255 gmx_bcast(sizeof(ed->eEDtype), &(ed->eEDtype), cr->mpi_comm_mygroup);
1257 int numedis = ed->edpar.size();
1258 /* First let everybody know how many ED data sets to expect */
1259 gmx_bcast(sizeof(numedis), &numedis, cr->mpi_comm_mygroup);
1260 nblock_abc(MASTER(cr), cr->mpi_comm_mygroup, numedis, &(ed->edpar));
1261 /* Now transfer the ED data set(s) */
1262 for (auto& edi : ed->edpar)
1264 /* Broadcast a single ED data set */
1265 block_bc(cr->mpi_comm_mygroup, edi);
1267 /* Broadcast positions */
1268 bc_ed_positions(cr, &(edi.sref),
1269 EssentialDynamicsStructure::Reference); /* reference positions (don't broadcast masses) */
1270 bc_ed_positions(cr, &(edi.sav),
1271 EssentialDynamicsStructure::Average); /* average positions (do broadcast masses as well) */
1272 bc_ed_positions(cr, &(edi.star), EssentialDynamicsStructure::Target); /* target positions */
1273 bc_ed_positions(cr, &(edi.sori), EssentialDynamicsStructure::Origin); /* origin positions */
1275 /* Broadcast eigenvectors */
1276 bc_ed_vecs(cr, &edi.vecs.mon, edi.sav.nr);
1277 bc_ed_vecs(cr, &edi.vecs.linfix, edi.sav.nr);
1278 bc_ed_vecs(cr, &edi.vecs.linacc, edi.sav.nr);
1279 bc_ed_vecs(cr, &edi.vecs.radfix, edi.sav.nr);
1280 bc_ed_vecs(cr, &edi.vecs.radacc, edi.sav.nr);
1281 bc_ed_vecs(cr, &edi.vecs.radcon, edi.sav.nr);
1282 /* Broadcast flooding eigenvectors and, if needed, values for the moving reference */
1283 bc_ed_vecs(cr, &edi.flood.vecs, edi.sav.nr);
1285 /* For harmonic restraints the reference projections can change with time */
1286 if (edi.flood.bHarmonic)
1288 snew_bc(MASTER(cr), edi.flood.initialReferenceProjection, edi.flood.vecs.neig);
1289 snew_bc(MASTER(cr), edi.flood.referenceProjectionSlope, edi.flood.vecs.neig);
1290 nblock_bc(cr->mpi_comm_mygroup, edi.flood.vecs.neig, edi.flood.initialReferenceProjection);
1291 nblock_bc(cr->mpi_comm_mygroup, edi.flood.vecs.neig, edi.flood.referenceProjectionSlope);
1297 /* init-routine called for every *.edi-cycle, initialises t_edpar structure */
1298 static void init_edi(const gmx_mtop_t* mtop, t_edpar* edi)
1300 int i;
1301 real totalmass = 0.0;
1302 rvec com;
1304 /* NOTE Init_edi is executed on the master process only
1305 * The initialized data sets are then transmitted to the
1306 * other nodes in broadcast_ed_data */
1308 /* evaluate masses (reference structure) */
1309 snew(edi->sref.m, edi->sref.nr);
1310 int molb = 0;
1311 for (i = 0; i < edi->sref.nr; i++)
1313 if (edi->fitmas)
1315 edi->sref.m[i] = mtopGetAtomMass(mtop, edi->sref.anrs[i], &molb);
1317 else
1319 edi->sref.m[i] = 1.0;
1322 /* Check that every m > 0. Bad things will happen otherwise. */
1323 if (edi->sref.m[i] <= 0.0)
1325 gmx_fatal(FARGS,
1326 "Reference structure atom %d (sam.edi index %d) has a mass of %g.\n"
1327 "For a mass-weighted fit, all reference structure atoms need to have a mass "
1328 ">0.\n"
1329 "Either make the covariance analysis non-mass-weighted, or exclude massless\n"
1330 "atoms from the reference structure by creating a proper index group.\n",
1331 i, edi->sref.anrs[i] + 1, edi->sref.m[i]);
1334 totalmass += edi->sref.m[i];
1336 edi->sref.mtot = totalmass;
1338 /* Masses m and sqrt(m) for the average structure. Note that m
1339 * is needed if forces have to be evaluated in do_edsam */
1340 snew(edi->sav.sqrtm, edi->sav.nr);
1341 snew(edi->sav.m, edi->sav.nr);
1342 for (i = 0; i < edi->sav.nr; i++)
1344 edi->sav.m[i] = mtopGetAtomMass(mtop, edi->sav.anrs[i], &molb);
1345 if (edi->pcamas)
1347 edi->sav.sqrtm[i] = sqrt(edi->sav.m[i]);
1349 else
1351 edi->sav.sqrtm[i] = 1.0;
1354 /* Check that every m > 0. Bad things will happen otherwise. */
1355 if (edi->sav.sqrtm[i] <= 0.0)
1357 gmx_fatal(FARGS,
1358 "Average structure atom %d (sam.edi index %d) has a mass of %g.\n"
1359 "For ED with mass-weighting, all average structure atoms need to have a mass "
1360 ">0.\n"
1361 "Either make the covariance analysis non-mass-weighted, or exclude massless\n"
1362 "atoms from the average structure by creating a proper index group.\n",
1363 i, edi->sav.anrs[i] + 1, edi->sav.m[i]);
1367 /* put reference structure in origin */
1368 get_center(edi->sref.x, edi->sref.m, edi->sref.nr, com);
1369 com[XX] = -com[XX];
1370 com[YY] = -com[YY];
1371 com[ZZ] = -com[ZZ];
1372 translate_x(edi->sref.x, edi->sref.nr, com);
1374 /* Init ED buffer */
1375 snew(edi->buf, 1);
1379 static void check(const char* line, const char* label)
1381 if (!strstr(line, label))
1383 gmx_fatal(FARGS,
1384 "Could not find input parameter %s at expected position in edsam input-file "
1385 "(.edi)\nline read instead is %s",
1386 label, line);
1391 static int read_checked_edint(FILE* file, const char* label)
1393 char line[STRLEN + 1];
1394 int idum;
1396 fgets2(line, STRLEN, file);
1397 check(line, label);
1398 fgets2(line, STRLEN, file);
1399 sscanf(line, max_ev_fmt_d, &idum);
1400 return idum;
1403 static bool read_checked_edbool(FILE* file, const char* label)
1405 return read_checked_edint(file, label) != 0;
1408 static int read_edint(FILE* file, bool* bEOF)
1410 char line[STRLEN + 1];
1411 int idum;
1412 char* eof;
1415 eof = fgets2(line, STRLEN, file);
1416 if (eof == nullptr)
1418 *bEOF = TRUE;
1419 return -1;
1421 eof = fgets2(line, STRLEN, file);
1422 if (eof == nullptr)
1424 *bEOF = TRUE;
1425 return -1;
1427 sscanf(line, max_ev_fmt_d, &idum);
1428 *bEOF = FALSE;
1429 return idum;
1433 static real read_checked_edreal(FILE* file, const char* label)
1435 char line[STRLEN + 1];
1436 double rdum;
1439 fgets2(line, STRLEN, file);
1440 check(line, label);
1441 fgets2(line, STRLEN, file);
1442 sscanf(line, max_ev_fmt_lf, &rdum);
1443 return static_cast<real>(rdum); /* always read as double and convert to single */
1447 static void read_edx(FILE* file, int number, int* anrs, rvec* x)
1449 int i, j;
1450 char line[STRLEN + 1];
1451 double d[3];
1454 for (i = 0; i < number; i++)
1456 fgets2(line, STRLEN, file);
1457 sscanf(line, max_ev_fmt_dlflflf, &anrs[i], &d[0], &d[1], &d[2]);
1458 anrs[i]--; /* we are reading FORTRAN indices */
1459 for (j = 0; j < 3; j++)
1461 x[i][j] = d[j]; /* always read as double and convert to single */
1466 namespace
1468 /*!\brief Read eigenvector coordinates for multiple eigenvectors from a file.
1469 * \param[in] in the file to read from
1470 * \param[in] numAtoms the number of atoms/coordinates for the eigenvector
1471 * \param[out] vec the eigen-vectors
1472 * \param[in] nEig the number of eigenvectors
1474 void scan_edvec(FILE* in, int numAtoms, rvec*** vec, int nEig)
1476 snew(*vec, nEig);
1477 for (int iEigenvector = 0; iEigenvector < nEig; iEigenvector++)
1479 snew((*vec)[iEigenvector], numAtoms);
1480 for (int iAtom = 0; iAtom < numAtoms; iAtom++)
1482 char line[STRLEN + 1];
1483 double x, y, z;
1484 fgets2(line, STRLEN, in);
1485 sscanf(line, max_ev_fmt_lelele, &x, &y, &z);
1486 (*vec)[iEigenvector][iAtom][XX] = x;
1487 (*vec)[iEigenvector][iAtom][YY] = y;
1488 (*vec)[iEigenvector][iAtom][ZZ] = z;
1492 /*!\brief Read an essential dynamics eigenvector and corresponding step size.
1493 * \param[in] in input file
1494 * \param[in] nrAtoms number of atoms/coordinates
1495 * \param[out] tvec the eigenvector
1497 void read_edvec(FILE* in, int nrAtoms, t_eigvec* tvec)
1499 tvec->neig = read_checked_edint(in, "NUMBER OF EIGENVECTORS");
1500 if (tvec->neig <= 0)
1502 return;
1505 snew(tvec->ieig, tvec->neig);
1506 snew(tvec->stpsz, tvec->neig);
1507 for (int i = 0; i < tvec->neig; i++)
1509 char line[STRLEN + 1];
1510 fgets2(line, STRLEN, in);
1511 int numEigenVectors;
1512 double stepSize;
1513 const int nscan = sscanf(line, max_ev_fmt_dlf, &numEigenVectors, &stepSize);
1514 if (nscan != 2)
1516 gmx_fatal(FARGS, "Expected 2 values for flooding vec: <nr> <stpsz>\n");
1518 tvec->ieig[i] = numEigenVectors;
1519 tvec->stpsz[i] = stepSize;
1520 } /* end of loop over eigenvectors */
1522 scan_edvec(in, nrAtoms, &tvec->vec, tvec->neig);
1524 /*!\brief Read an essential dynamics eigenvector and intial reference projections and slopes if available.
1526 * Eigenvector and its intial reference and slope are stored on the
1527 * same line in the .edi format. To avoid re-winding the .edi file,
1528 * the reference eigen-vector and reference data are read in one go.
1530 * \param[in] in input file
1531 * \param[in] nrAtoms number of atoms/coordinates
1532 * \param[out] tvec the eigenvector
1533 * \param[out] initialReferenceProjection The projection onto eigenvectors as read from file.
1534 * \param[out] referenceProjectionSlope The slope of the reference projections.
1536 bool readEdVecWithReferenceProjection(FILE* in,
1537 int nrAtoms,
1538 t_eigvec* tvec,
1539 real** initialReferenceProjection,
1540 real** referenceProjectionSlope)
1542 bool providesReference = false;
1543 tvec->neig = read_checked_edint(in, "NUMBER OF EIGENVECTORS");
1544 if (tvec->neig <= 0)
1546 return false;
1549 snew(tvec->ieig, tvec->neig);
1550 snew(tvec->stpsz, tvec->neig);
1551 snew(*initialReferenceProjection, tvec->neig);
1552 snew(*referenceProjectionSlope, tvec->neig);
1553 for (int i = 0; i < tvec->neig; i++)
1555 char line[STRLEN + 1];
1556 fgets2(line, STRLEN, in);
1557 int numEigenVectors;
1558 double stepSize = 0.;
1559 double referenceProjection = 0.;
1560 double referenceSlope = 0.;
1561 const int nscan = sscanf(line, max_ev_fmt_dlflflf, &numEigenVectors, &stepSize,
1562 &referenceProjection, &referenceSlope);
1563 /* Zero out values which were not scanned */
1564 switch (nscan)
1566 case 4:
1567 /* Every 4 values read, including reference position */
1568 providesReference = true;
1569 break;
1570 case 3:
1571 /* A reference position is provided */
1572 providesReference = true;
1573 /* No value for slope, set to 0 */
1574 referenceSlope = 0.0;
1575 break;
1576 case 2:
1577 /* No values for reference projection and slope, set to 0 */
1578 referenceProjection = 0.0;
1579 referenceSlope = 0.0;
1580 break;
1581 default:
1582 gmx_fatal(FARGS,
1583 "Expected 2 - 4 (not %d) values for flooding vec: <nr> <spring const> "
1584 "<refproj> <refproj-slope>\n",
1585 nscan);
1587 (*initialReferenceProjection)[i] = referenceProjection;
1588 (*referenceProjectionSlope)[i] = referenceSlope;
1590 tvec->ieig[i] = numEigenVectors;
1591 tvec->stpsz[i] = stepSize;
1592 } /* end of loop over eigenvectors */
1594 scan_edvec(in, nrAtoms, &tvec->vec, tvec->neig);
1595 return providesReference;
1598 /*!\brief Allocate working memory for the eigenvectors.
1599 * \param[in,out] tvec eigenvector for which memory will be allocated
1601 void setup_edvec(t_eigvec* tvec)
1603 snew(tvec->xproj, tvec->neig);
1604 snew(tvec->fproj, tvec->neig);
1605 snew(tvec->refproj, tvec->neig);
1607 } // namespace
1610 /* Check if the same atom indices are used for reference and average positions */
1611 static gmx_bool check_if_same(struct gmx_edx sref, struct gmx_edx sav)
1613 int i;
1616 /* If the number of atoms differs between the two structures,
1617 * they cannot be identical */
1618 if (sref.nr != sav.nr)
1620 return FALSE;
1623 /* Now that we know that both stuctures have the same number of atoms,
1624 * check if also the indices are identical */
1625 for (i = 0; i < sav.nr; i++)
1627 if (sref.anrs[i] != sav.anrs[i])
1629 return FALSE;
1632 fprintf(stderr,
1633 "ED: Note: Reference and average structure are composed of the same atom indices.\n");
1635 return TRUE;
1638 namespace
1641 //! Oldest (lowest) magic number indicating unsupported essential dynamics input
1642 constexpr int c_oldestUnsupportedMagicNumber = 666;
1643 //! Oldest (lowest) magic number indicating supported essential dynamics input
1644 constexpr int c_oldestSupportedMagicNumber = 669;
1645 //! Latest (highest) magic number indicating supported essential dynamics input
1646 constexpr int c_latestSupportedMagicNumber = 670;
1648 /*!\brief Set up essential dynamics work parameters.
1649 * \param[in] edi Essential dynamics working structure.
1651 void setup_edi(t_edpar* edi)
1653 snew(edi->sref.x_old, edi->sref.nr);
1654 edi->sref.sqrtm = nullptr;
1655 snew(edi->sav.x_old, edi->sav.nr);
1656 if (edi->star.nr > 0)
1658 edi->star.sqrtm = nullptr;
1660 if (edi->sori.nr > 0)
1662 edi->sori.sqrtm = nullptr;
1664 setup_edvec(&edi->vecs.linacc);
1665 setup_edvec(&edi->vecs.mon);
1666 setup_edvec(&edi->vecs.linfix);
1667 setup_edvec(&edi->vecs.linacc);
1668 setup_edvec(&edi->vecs.radfix);
1669 setup_edvec(&edi->vecs.radacc);
1670 setup_edvec(&edi->vecs.radcon);
1671 setup_edvec(&edi->flood.vecs);
1674 /*!\brief Check if edi file version supports CONST_FORCE_FLOODING.
1675 * \param[in] magicNumber indicates the essential dynamics input file version
1676 * \returns true if CONST_FORCE_FLOODING is supported
1678 bool ediFileHasConstForceFlooding(int magicNumber)
1680 return magicNumber > c_oldestSupportedMagicNumber;
1683 /*!\brief Reports reading success of the essential dynamics input file magic number.
1684 * \param[in] in pointer to input file
1685 * \param[out] magicNumber determines the edi file version
1686 * \returns true if setting file version from input file was successful.
1688 bool ediFileHasValidDataSet(FILE* in, int* magicNumber)
1690 /* the edi file is not free format, so expect problems if the input is corrupt. */
1691 bool reachedEndOfFile = true;
1692 /* check the magic number */
1693 *magicNumber = read_edint(in, &reachedEndOfFile);
1694 /* Check whether we have reached the end of the input file */
1695 return !reachedEndOfFile;
1698 /*!\brief Trigger fatal error if magic number is unsupported.
1699 * \param[in] magicNumber A magic number read from the edi file
1700 * \param[in] fn name of input file for error reporting
1702 void exitWhenMagicNumberIsInvalid(int magicNumber, const char* fn)
1705 if (magicNumber < c_oldestSupportedMagicNumber || magicNumber > c_latestSupportedMagicNumber)
1707 if (magicNumber >= c_oldestUnsupportedMagicNumber && magicNumber < c_oldestSupportedMagicNumber)
1709 gmx_fatal(FARGS,
1710 "Wrong magic number: Use newest version of make_edi to produce edi file");
1712 else
1714 gmx_fatal(FARGS, "Wrong magic number %d in %s", magicNumber, fn);
1719 /*!\brief reads an essential dynamics input file into a essential dynamics data structure.
1721 * \param[in] in input file
1722 * \param[in] nr_mdatoms the number of md atoms, used for consistency checking
1723 * \param[in] hasConstForceFlooding determines if field "CONST_FORCE_FLOODING" is available in input file
1724 * \param[in] fn the file name of the input file for error reporting
1725 * \returns edi essential dynamics parameters
1727 t_edpar read_edi(FILE* in, int nr_mdatoms, bool hasConstForceFlooding, const char* fn)
1729 t_edpar edi;
1730 bool bEOF;
1731 /* check the number of atoms */
1732 edi.nini = read_edint(in, &bEOF);
1733 if (edi.nini != nr_mdatoms)
1735 gmx_fatal(FARGS, "Nr of atoms in %s (%d) does not match nr of md atoms (%d)", fn, edi.nini,
1736 nr_mdatoms);
1739 /* Done checking. For the rest we blindly trust the input */
1740 edi.fitmas = read_checked_edbool(in, "FITMAS");
1741 edi.pcamas = read_checked_edbool(in, "ANALYSIS_MAS");
1742 edi.outfrq = read_checked_edint(in, "OUTFRQ");
1743 edi.maxedsteps = read_checked_edint(in, "MAXLEN");
1744 edi.slope = read_checked_edreal(in, "SLOPECRIT");
1746 edi.presteps = read_checked_edint(in, "PRESTEPS");
1747 edi.flood.deltaF0 = read_checked_edreal(in, "DELTA_F0");
1748 edi.flood.deltaF = read_checked_edreal(in, "INIT_DELTA_F");
1749 edi.flood.tau = read_checked_edreal(in, "TAU");
1750 edi.flood.constEfl = read_checked_edreal(in, "EFL_NULL");
1751 edi.flood.alpha2 = read_checked_edreal(in, "ALPHA2");
1752 edi.flood.kT = read_checked_edreal(in, "KT");
1753 edi.flood.bHarmonic = read_checked_edbool(in, "HARMONIC");
1754 if (hasConstForceFlooding)
1756 edi.flood.bConstForce = read_checked_edbool(in, "CONST_FORCE_FLOODING");
1758 else
1760 edi.flood.bConstForce = FALSE;
1762 edi.sref.nr = read_checked_edint(in, "NREF");
1764 /* allocate space for reference positions and read them */
1765 snew(edi.sref.anrs, edi.sref.nr);
1766 snew(edi.sref.x, edi.sref.nr);
1767 read_edx(in, edi.sref.nr, edi.sref.anrs, edi.sref.x);
1769 /* average positions. they define which atoms will be used for ED sampling */
1770 edi.sav.nr = read_checked_edint(in, "NAV");
1771 snew(edi.sav.anrs, edi.sav.nr);
1772 snew(edi.sav.x, edi.sav.nr);
1773 read_edx(in, edi.sav.nr, edi.sav.anrs, edi.sav.x);
1775 /* Check if the same atom indices are used for reference and average positions */
1776 edi.bRefEqAv = check_if_same(edi.sref, edi.sav);
1778 /* eigenvectors */
1780 read_edvec(in, edi.sav.nr, &edi.vecs.mon);
1781 read_edvec(in, edi.sav.nr, &edi.vecs.linfix);
1782 read_edvec(in, edi.sav.nr, &edi.vecs.linacc);
1783 read_edvec(in, edi.sav.nr, &edi.vecs.radfix);
1784 read_edvec(in, edi.sav.nr, &edi.vecs.radacc);
1785 read_edvec(in, edi.sav.nr, &edi.vecs.radcon);
1787 /* Was a specific reference point for the flooding/umbrella potential provided in the edi file? */
1788 bool bHaveReference = false;
1789 if (edi.flood.bHarmonic)
1791 bHaveReference = readEdVecWithReferenceProjection(in, edi.sav.nr, &edi.flood.vecs,
1792 &edi.flood.initialReferenceProjection,
1793 &edi.flood.referenceProjectionSlope);
1795 else
1797 read_edvec(in, edi.sav.nr, &edi.flood.vecs);
1800 /* target positions */
1801 edi.star.nr = read_edint(in, &bEOF);
1802 if (edi.star.nr > 0)
1804 snew(edi.star.anrs, edi.star.nr);
1805 snew(edi.star.x, edi.star.nr);
1806 read_edx(in, edi.star.nr, edi.star.anrs, edi.star.x);
1809 /* positions defining origin of expansion circle */
1810 edi.sori.nr = read_edint(in, &bEOF);
1811 if (edi.sori.nr > 0)
1813 if (bHaveReference)
1815 /* Both an -ori structure and a at least one manual reference point have been
1816 * specified. That's ambiguous and probably not intentional. */
1817 gmx_fatal(FARGS,
1818 "ED: An origin structure has been provided and a at least one (moving) "
1819 "reference\n"
1820 " point was manually specified in the edi file. That is ambiguous. "
1821 "Aborting.\n");
1823 snew(edi.sori.anrs, edi.sori.nr);
1824 snew(edi.sori.x, edi.sori.nr);
1825 read_edx(in, edi.sori.nr, edi.sori.anrs, edi.sori.x);
1827 return edi;
1830 std::vector<t_edpar> read_edi_file(const char* fn, int nr_mdatoms)
1832 std::vector<t_edpar> essentialDynamicsParameters;
1833 FILE* in;
1834 /* This routine is executed on the master only */
1836 /* Open the .edi parameter input file */
1837 in = gmx_fio_fopen(fn, "r");
1838 fprintf(stderr, "ED: Reading edi file %s\n", fn);
1840 /* Now read a sequence of ED input parameter sets from the edi file */
1841 int ediFileMagicNumber;
1842 while (ediFileHasValidDataSet(in, &ediFileMagicNumber))
1844 exitWhenMagicNumberIsInvalid(ediFileMagicNumber, fn);
1845 bool hasConstForceFlooding = ediFileHasConstForceFlooding(ediFileMagicNumber);
1846 auto edi = read_edi(in, nr_mdatoms, hasConstForceFlooding, fn);
1847 setup_edi(&edi);
1848 essentialDynamicsParameters.emplace_back(edi);
1850 /* Since we arrived within this while loop we know that there is still another data set to be read in */
1851 /* We need to allocate space for the data: */
1853 if (essentialDynamicsParameters.empty())
1855 gmx_fatal(FARGS, "No complete ED data set found in edi file %s.", fn);
1858 fprintf(stderr, "ED: Found %zu ED group%s.\n", essentialDynamicsParameters.size(),
1859 essentialDynamicsParameters.size() > 1 ? "s" : "");
1861 /* Close the .edi file again */
1862 gmx_fio_fclose(in);
1864 return essentialDynamicsParameters;
1866 } // anonymous namespace
1869 struct t_fit_to_ref
1871 rvec* xcopy; /* Working copy of the positions in fit_to_reference */
1874 /* Fit the current positions to the reference positions
1875 * Do not actually do the fit, just return rotation and translation.
1876 * Note that the COM of the reference structure was already put into
1877 * the origin by init_edi. */
1878 static void fit_to_reference(rvec* xcoll, /* The positions to be fitted */
1879 rvec transvec, /* The translation vector */
1880 matrix rotmat, /* The rotation matrix */
1881 t_edpar* edi) /* Just needed for do_edfit */
1883 rvec com; /* center of mass */
1884 int i;
1885 struct t_fit_to_ref* loc;
1888 /* Allocate memory the first time this routine is called for each edi group */
1889 if (nullptr == edi->buf->fit_to_ref)
1891 snew(edi->buf->fit_to_ref, 1);
1892 snew(edi->buf->fit_to_ref->xcopy, edi->sref.nr);
1894 loc = edi->buf->fit_to_ref;
1896 /* We do not touch the original positions but work on a copy. */
1897 for (i = 0; i < edi->sref.nr; i++)
1899 copy_rvec(xcoll[i], loc->xcopy[i]);
1902 /* Calculate the center of mass */
1903 get_center(loc->xcopy, edi->sref.m, edi->sref.nr, com);
1905 transvec[XX] = -com[XX];
1906 transvec[YY] = -com[YY];
1907 transvec[ZZ] = -com[ZZ];
1909 /* Subtract the center of mass from the copy */
1910 translate_x(loc->xcopy, edi->sref.nr, transvec);
1912 /* Determine the rotation matrix */
1913 do_edfit(edi->sref.nr, edi->sref.x, loc->xcopy, rotmat, edi);
1917 static void translate_and_rotate(rvec* x, /* The positions to be translated and rotated */
1918 int nat, /* How many positions are there? */
1919 rvec transvec, /* The translation vector */
1920 matrix rotmat) /* The rotation matrix */
1922 /* Translation */
1923 translate_x(x, nat, transvec);
1925 /* Rotation */
1926 rotate_x(x, nat, rotmat);
1930 /* Gets the rms deviation of the positions to the structure s */
1931 /* fit_to_structure has to be called before calling this routine! */
1932 static real rmsd_from_structure(rvec* x, /* The positions under consideration */
1933 struct gmx_edx* s) /* The structure from which the rmsd shall be computed */
1935 real rmsd = 0.0;
1936 int i;
1939 for (i = 0; i < s->nr; i++)
1941 rmsd += distance2(s->x[i], x[i]);
1944 rmsd /= static_cast<real>(s->nr);
1945 rmsd = sqrt(rmsd);
1947 return rmsd;
1951 void dd_make_local_ed_indices(gmx_domdec_t* dd, struct gmx_edsam* ed)
1953 if (ed->eEDtype != EssentialDynamicsType::None)
1955 /* Loop over ED groups */
1957 for (auto& edi : ed->edpar)
1959 /* Local atoms of the reference structure (for fitting), need only be assembled
1960 * if their indices differ from the average ones */
1961 if (!edi.bRefEqAv)
1963 dd_make_local_group_indices(dd->ga2la, edi.sref.nr, edi.sref.anrs, &edi.sref.nr_loc,
1964 &edi.sref.anrs_loc, &edi.sref.nalloc_loc, edi.sref.c_ind);
1967 /* Local atoms of the average structure (on these ED will be performed) */
1968 dd_make_local_group_indices(dd->ga2la, edi.sav.nr, edi.sav.anrs, &edi.sav.nr_loc,
1969 &edi.sav.anrs_loc, &edi.sav.nalloc_loc, edi.sav.c_ind);
1971 /* Indicate that the ED shift vectors for this structure need to be updated
1972 * at the next call to communicate_group_positions, since obviously we are in a NS step */
1973 edi.buf->do_edsam->bUpdateShifts = TRUE;
1979 static inline void ed_unshift_single_coord(const matrix box, const rvec x, const ivec is, rvec xu)
1981 int tx, ty, tz;
1984 tx = is[XX];
1985 ty = is[YY];
1986 tz = is[ZZ];
1988 if (TRICLINIC(box))
1990 xu[XX] = x[XX] - tx * box[XX][XX] - ty * box[YY][XX] - tz * box[ZZ][XX];
1991 xu[YY] = x[YY] - ty * box[YY][YY] - tz * box[ZZ][YY];
1992 xu[ZZ] = x[ZZ] - tz * box[ZZ][ZZ];
1994 else
1996 xu[XX] = x[XX] - tx * box[XX][XX];
1997 xu[YY] = x[YY] - ty * box[YY][YY];
1998 xu[ZZ] = x[ZZ] - tz * box[ZZ][ZZ];
2002 namespace
2004 /*!\brief Apply fixed linear constraints to essential dynamics variable.
2005 * \param[in,out] xcoll The collected coordinates.
2006 * \param[in] edi the essential dynamics parameters
2007 * \param[in] step the current simulation step
2009 void do_linfix(rvec* xcoll, const t_edpar& edi, int64_t step)
2011 /* loop over linfix vectors */
2012 for (int i = 0; i < edi.vecs.linfix.neig; i++)
2014 /* calculate the projection */
2015 real proj = projectx(edi, xcoll, edi.vecs.linfix.vec[i]);
2017 /* calculate the correction */
2018 real preFactor = edi.vecs.linfix.refproj[i] + step * edi.vecs.linfix.stpsz[i] - proj;
2020 /* apply the correction */
2021 preFactor /= edi.sav.sqrtm[i];
2022 for (int j = 0; j < edi.sav.nr; j++)
2024 rvec differenceVector;
2025 svmul(preFactor, edi.vecs.linfix.vec[i][j], differenceVector);
2026 rvec_inc(xcoll[j], differenceVector);
2031 /*!\brief Apply acceptance linear constraints to essential dynamics variable.
2032 * \param[in,out] xcoll The collected coordinates.
2033 * \param[in] edi the essential dynamics parameters
2035 void do_linacc(rvec* xcoll, t_edpar* edi)
2037 /* loop over linacc vectors */
2038 for (int i = 0; i < edi->vecs.linacc.neig; i++)
2040 /* calculate the projection */
2041 real proj = projectx(*edi, xcoll, edi->vecs.linacc.vec[i]);
2043 /* calculate the correction */
2044 real preFactor = 0.0;
2045 if (edi->vecs.linacc.stpsz[i] > 0.0)
2047 if ((proj - edi->vecs.linacc.refproj[i]) < 0.0)
2049 preFactor = edi->vecs.linacc.refproj[i] - proj;
2052 if (edi->vecs.linacc.stpsz[i] < 0.0)
2054 if ((proj - edi->vecs.linacc.refproj[i]) > 0.0)
2056 preFactor = edi->vecs.linacc.refproj[i] - proj;
2060 /* apply the correction */
2061 preFactor /= edi->sav.sqrtm[i];
2062 for (int j = 0; j < edi->sav.nr; j++)
2064 rvec differenceVector;
2065 svmul(preFactor, edi->vecs.linacc.vec[i][j], differenceVector);
2066 rvec_inc(xcoll[j], differenceVector);
2069 /* new positions will act as reference */
2070 edi->vecs.linacc.refproj[i] = proj + preFactor;
2073 } // namespace
2076 static void do_radfix(rvec* xcoll, t_edpar* edi)
2078 int i, j;
2079 real *proj, rad = 0.0, ratio;
2080 rvec vec_dum;
2083 if (edi->vecs.radfix.neig == 0)
2085 return;
2088 snew(proj, edi->vecs.radfix.neig);
2090 /* loop over radfix vectors */
2091 for (i = 0; i < edi->vecs.radfix.neig; i++)
2093 /* calculate the projections, radius */
2094 proj[i] = projectx(*edi, xcoll, edi->vecs.radfix.vec[i]);
2095 rad += gmx::square(proj[i] - edi->vecs.radfix.refproj[i]);
2098 rad = sqrt(rad);
2099 ratio = (edi->vecs.radfix.stpsz[0] + edi->vecs.radfix.radius) / rad - 1.0;
2100 edi->vecs.radfix.radius += edi->vecs.radfix.stpsz[0];
2102 /* loop over radfix vectors */
2103 for (i = 0; i < edi->vecs.radfix.neig; i++)
2105 proj[i] -= edi->vecs.radfix.refproj[i];
2107 /* apply the correction */
2108 proj[i] /= edi->sav.sqrtm[i];
2109 proj[i] *= ratio;
2110 for (j = 0; j < edi->sav.nr; j++)
2112 svmul(proj[i], edi->vecs.radfix.vec[i][j], vec_dum);
2113 rvec_inc(xcoll[j], vec_dum);
2117 sfree(proj);
2121 static void do_radacc(rvec* xcoll, t_edpar* edi)
2123 int i, j;
2124 real *proj, rad = 0.0, ratio = 0.0;
2125 rvec vec_dum;
2128 if (edi->vecs.radacc.neig == 0)
2130 return;
2133 snew(proj, edi->vecs.radacc.neig);
2135 /* loop over radacc vectors */
2136 for (i = 0; i < edi->vecs.radacc.neig; i++)
2138 /* calculate the projections, radius */
2139 proj[i] = projectx(*edi, xcoll, edi->vecs.radacc.vec[i]);
2140 rad += gmx::square(proj[i] - edi->vecs.radacc.refproj[i]);
2142 rad = sqrt(rad);
2144 /* only correct when radius decreased */
2145 if (rad < edi->vecs.radacc.radius)
2147 ratio = edi->vecs.radacc.radius / rad - 1.0;
2149 else
2151 edi->vecs.radacc.radius = rad;
2154 /* loop over radacc vectors */
2155 for (i = 0; i < edi->vecs.radacc.neig; i++)
2157 proj[i] -= edi->vecs.radacc.refproj[i];
2159 /* apply the correction */
2160 proj[i] /= edi->sav.sqrtm[i];
2161 proj[i] *= ratio;
2162 for (j = 0; j < edi->sav.nr; j++)
2164 svmul(proj[i], edi->vecs.radacc.vec[i][j], vec_dum);
2165 rvec_inc(xcoll[j], vec_dum);
2168 sfree(proj);
2172 struct t_do_radcon
2174 real* proj;
2177 static void do_radcon(rvec* xcoll, t_edpar* edi)
2179 int i, j;
2180 real rad = 0.0, ratio = 0.0;
2181 struct t_do_radcon* loc;
2182 gmx_bool bFirst;
2183 rvec vec_dum;
2186 if (edi->buf->do_radcon != nullptr)
2188 bFirst = FALSE;
2190 else
2192 bFirst = TRUE;
2193 snew(edi->buf->do_radcon, 1);
2195 loc = edi->buf->do_radcon;
2197 if (edi->vecs.radcon.neig == 0)
2199 return;
2202 if (bFirst)
2204 snew(loc->proj, edi->vecs.radcon.neig);
2207 /* loop over radcon vectors */
2208 for (i = 0; i < edi->vecs.radcon.neig; i++)
2210 /* calculate the projections, radius */
2211 loc->proj[i] = projectx(*edi, xcoll, edi->vecs.radcon.vec[i]);
2212 rad += gmx::square(loc->proj[i] - edi->vecs.radcon.refproj[i]);
2214 rad = sqrt(rad);
2215 /* only correct when radius increased */
2216 if (rad > edi->vecs.radcon.radius)
2218 ratio = edi->vecs.radcon.radius / rad - 1.0;
2220 /* loop over radcon vectors */
2221 for (i = 0; i < edi->vecs.radcon.neig; i++)
2223 /* apply the correction */
2224 loc->proj[i] -= edi->vecs.radcon.refproj[i];
2225 loc->proj[i] /= edi->sav.sqrtm[i];
2226 loc->proj[i] *= ratio;
2228 for (j = 0; j < edi->sav.nr; j++)
2230 svmul(loc->proj[i], edi->vecs.radcon.vec[i][j], vec_dum);
2231 rvec_inc(xcoll[j], vec_dum);
2235 else
2237 edi->vecs.radcon.radius = rad;
2242 static void ed_apply_constraints(rvec* xcoll, t_edpar* edi, int64_t step)
2244 int i;
2247 /* subtract the average positions */
2248 for (i = 0; i < edi->sav.nr; i++)
2250 rvec_dec(xcoll[i], edi->sav.x[i]);
2253 /* apply the constraints */
2254 if (step >= 0)
2256 do_linfix(xcoll, *edi, step);
2258 do_linacc(xcoll, edi);
2259 if (step >= 0)
2261 do_radfix(xcoll, edi);
2263 do_radacc(xcoll, edi);
2264 do_radcon(xcoll, edi);
2266 /* add back the average positions */
2267 for (i = 0; i < edi->sav.nr; i++)
2269 rvec_inc(xcoll[i], edi->sav.x[i]);
2274 namespace
2276 /*!\brief Write out the projections onto the eigenvectors.
2277 * The order of output corresponds to ed_output_legend().
2278 * \param[in] edi The essential dyanmics parameters
2279 * \param[in] fp The output file
2280 * \param[in] rmsd the rmsd to the reference structure
2282 void write_edo(const t_edpar& edi, FILE* fp, real rmsd)
2284 /* Output how well we fit to the reference structure */
2285 fprintf(fp, EDcol_ffmt, rmsd);
2287 for (int i = 0; i < edi.vecs.mon.neig; i++)
2289 fprintf(fp, EDcol_efmt, edi.vecs.mon.xproj[i]);
2292 for (int i = 0; i < edi.vecs.linfix.neig; i++)
2294 fprintf(fp, EDcol_efmt, edi.vecs.linfix.xproj[i]);
2297 for (int i = 0; i < edi.vecs.linacc.neig; i++)
2299 fprintf(fp, EDcol_efmt, edi.vecs.linacc.xproj[i]);
2302 for (int i = 0; i < edi.vecs.radfix.neig; i++)
2304 fprintf(fp, EDcol_efmt, edi.vecs.radfix.xproj[i]);
2306 if (edi.vecs.radfix.neig)
2308 fprintf(fp, EDcol_ffmt, calc_radius(edi.vecs.radfix)); /* fixed increment radius */
2311 for (int i = 0; i < edi.vecs.radacc.neig; i++)
2313 fprintf(fp, EDcol_efmt, edi.vecs.radacc.xproj[i]);
2315 if (edi.vecs.radacc.neig)
2317 fprintf(fp, EDcol_ffmt, calc_radius(edi.vecs.radacc)); /* acceptance radius */
2320 for (int i = 0; i < edi.vecs.radcon.neig; i++)
2322 fprintf(fp, EDcol_efmt, edi.vecs.radcon.xproj[i]);
2324 if (edi.vecs.radcon.neig)
2326 fprintf(fp, EDcol_ffmt, calc_radius(edi.vecs.radcon)); /* contracting radius */
2329 } // namespace
2332 /* Returns if any constraints are switched on */
2333 static bool ed_constraints(EssentialDynamicsType edtype, const t_edpar& edi)
2335 if (edtype == EssentialDynamicsType::EDSampling || edtype == EssentialDynamicsType::Flooding)
2337 return ((edi.vecs.linfix.neig != 0) || (edi.vecs.linacc.neig != 0) || (edi.vecs.radfix.neig != 0)
2338 || (edi.vecs.radacc.neig != 0) || (edi.vecs.radcon.neig != 0));
2340 return false;
2344 /* Copies reference projection 'refproj' to fixed 'initialReferenceProjection' variable for
2345 * flooding/ umbrella sampling simulations. */
2346 static void copyEvecReference(t_eigvec* floodvecs, real* initialReferenceProjection)
2348 int i;
2351 if (nullptr == initialReferenceProjection)
2353 snew(initialReferenceProjection, floodvecs->neig);
2356 for (i = 0; i < floodvecs->neig; i++)
2358 initialReferenceProjection[i] = floodvecs->refproj[i];
2363 /* Call on MASTER only. Check whether the essential dynamics / flooding
2364 * groups of the checkpoint file are consistent with the provided .edi file. */
2365 static void crosscheck_edi_file_vs_checkpoint(const gmx_edsam& ed, edsamhistory_t* EDstate)
2367 if (nullptr == EDstate->nref || nullptr == EDstate->nav)
2369 gmx_fatal(FARGS,
2370 "Essential dynamics and flooding can only be switched on (or off) at the\n"
2371 "start of a new simulation. If a simulation runs with/without ED constraints,\n"
2372 "it must also continue with/without ED constraints when checkpointing.\n"
2373 "To switch on (or off) ED constraints, please prepare a new .tpr to start\n"
2374 "from without a checkpoint.\n");
2377 for (size_t edinum = 0; edinum < ed.edpar.size(); ++edinum)
2379 /* Check number of atoms in the reference and average structures */
2380 if (EDstate->nref[edinum] != ed.edpar[edinum].sref.nr)
2382 gmx_fatal(FARGS,
2383 "The number of reference structure atoms in ED group %c is\n"
2384 "not the same in .cpt (NREF=%d) and .edi (NREF=%d) files!\n",
2385 get_EDgroupChar(edinum + 1, 0), EDstate->nref[edinum], ed.edpar[edinum].sref.nr);
2387 if (EDstate->nav[edinum] != ed.edpar[edinum].sav.nr)
2389 gmx_fatal(FARGS,
2390 "The number of average structure atoms in ED group %c is\n"
2391 "not the same in .cpt (NREF=%d) and .edi (NREF=%d) files!\n",
2392 get_EDgroupChar(edinum + 1, 0), EDstate->nav[edinum], ed.edpar[edinum].sav.nr);
2396 if (gmx::ssize(ed.edpar) != EDstate->nED)
2398 gmx_fatal(FARGS,
2399 "The number of essential dynamics / flooding groups is not consistent.\n"
2400 "There are %d ED groups in the .cpt file, but %zu in the .edi file!\n"
2401 "Are you sure this is the correct .edi file?\n",
2402 EDstate->nED, ed.edpar.size());
2407 /* The edsamstate struct stores the information we need to make the ED group
2408 * whole again after restarts from a checkpoint file. Here we do the following:
2409 * a) If we did not start from .cpt, we prepare the struct for proper .cpt writing,
2410 * b) if we did start from .cpt, we copy over the last whole structures from .cpt,
2411 * c) in any case, for subsequent checkpoint writing, we set the pointers in
2412 * edsamstate to the x_old arrays, which contain the correct PBC representation of
2413 * all ED structures at the last time step. */
2414 static void init_edsamstate(const gmx_edsam& ed, edsamhistory_t* EDstate)
2416 snew(EDstate->old_sref_p, EDstate->nED);
2417 snew(EDstate->old_sav_p, EDstate->nED);
2419 /* If we did not read in a .cpt file, these arrays are not yet allocated */
2420 if (!EDstate->bFromCpt)
2422 snew(EDstate->nref, EDstate->nED);
2423 snew(EDstate->nav, EDstate->nED);
2426 /* Loop over all ED/flooding data sets (usually only one, though) */
2427 for (int nr_edi = 0; nr_edi < EDstate->nED; nr_edi++)
2429 const auto& edi = ed.edpar[nr_edi];
2430 /* We always need the last reference and average positions such that
2431 * in the next time step we can make the ED group whole again
2432 * if the atoms do not have the correct PBC representation */
2433 if (EDstate->bFromCpt)
2435 /* Copy the last whole positions of reference and average group from .cpt */
2436 for (int i = 0; i < edi.sref.nr; i++)
2438 copy_rvec(EDstate->old_sref[nr_edi][i], edi.sref.x_old[i]);
2440 for (int i = 0; i < edi.sav.nr; i++)
2442 copy_rvec(EDstate->old_sav[nr_edi][i], edi.sav.x_old[i]);
2445 else
2447 EDstate->nref[nr_edi] = edi.sref.nr;
2448 EDstate->nav[nr_edi] = edi.sav.nr;
2451 /* For subsequent checkpoint writing, set the edsamstate pointers to the edi arrays: */
2452 EDstate->old_sref_p[nr_edi] = edi.sref.x_old;
2453 EDstate->old_sav_p[nr_edi] = edi.sav.x_old;
2458 /* Adds 'buf' to 'str' */
2459 static void add_to_string(char** str, const char* buf)
2461 int len;
2464 len = strlen(*str) + strlen(buf) + 1;
2465 srenew(*str, len);
2466 strcat(*str, buf);
2470 static void add_to_string_aligned(char** str, const char* buf)
2472 char buf_aligned[STRLEN];
2474 sprintf(buf_aligned, EDcol_sfmt, buf);
2475 add_to_string(str, buf_aligned);
2479 static void nice_legend(const char*** setname,
2480 int* nsets,
2481 char** LegendStr,
2482 const char* value,
2483 const char* unit,
2484 char EDgroupchar)
2486 auto tmp = gmx::formatString("%c %s", EDgroupchar, value);
2487 add_to_string_aligned(LegendStr, tmp.c_str());
2488 tmp += gmx::formatString(" (%s)", unit);
2489 (*setname)[*nsets] = gmx_strdup(tmp.c_str());
2490 (*nsets)++;
2494 static void nice_legend_evec(const char*** setname,
2495 int* nsets,
2496 char** LegendStr,
2497 t_eigvec* evec,
2498 char EDgroupChar,
2499 const char* EDtype)
2501 int i;
2502 char tmp[STRLEN];
2505 for (i = 0; i < evec->neig; i++)
2507 sprintf(tmp, "EV%dprj%s", evec->ieig[i], EDtype);
2508 nice_legend(setname, nsets, LegendStr, tmp, "nm", EDgroupChar);
2513 /* Makes a legend for the xvg output file. Call on MASTER only! */
2514 static void write_edo_legend(gmx_edsam* ed, int nED, const gmx_output_env_t* oenv)
2516 int i;
2517 int nr_edi, nsets, n_flood, n_edsam;
2518 const char** setname;
2519 char buf[STRLEN];
2520 char* LegendStr = nullptr;
2523 auto edi = ed->edpar.begin();
2525 fprintf(ed->edo, "# Output will be written every %d step%s\n", edi->outfrq,
2526 edi->outfrq != 1 ? "s" : "");
2528 for (nr_edi = 1; nr_edi <= nED; nr_edi++)
2530 fprintf(ed->edo, "#\n");
2531 fprintf(ed->edo, "# Summary of applied con/restraints for the ED group %c\n",
2532 get_EDgroupChar(nr_edi, nED));
2533 fprintf(ed->edo, "# Atoms in average structure: %d\n", edi->sav.nr);
2534 fprintf(ed->edo, "# monitor : %d vec%s\n", edi->vecs.mon.neig,
2535 edi->vecs.mon.neig != 1 ? "s" : "");
2536 fprintf(ed->edo, "# LINFIX : %d vec%s\n", edi->vecs.linfix.neig,
2537 edi->vecs.linfix.neig != 1 ? "s" : "");
2538 fprintf(ed->edo, "# LINACC : %d vec%s\n", edi->vecs.linacc.neig,
2539 edi->vecs.linacc.neig != 1 ? "s" : "");
2540 fprintf(ed->edo, "# RADFIX : %d vec%s\n", edi->vecs.radfix.neig,
2541 edi->vecs.radfix.neig != 1 ? "s" : "");
2542 fprintf(ed->edo, "# RADACC : %d vec%s\n", edi->vecs.radacc.neig,
2543 edi->vecs.radacc.neig != 1 ? "s" : "");
2544 fprintf(ed->edo, "# RADCON : %d vec%s\n", edi->vecs.radcon.neig,
2545 edi->vecs.radcon.neig != 1 ? "s" : "");
2546 fprintf(ed->edo, "# FLOODING : %d vec%s ", edi->flood.vecs.neig,
2547 edi->flood.vecs.neig != 1 ? "s" : "");
2549 if (edi->flood.vecs.neig)
2551 /* If in any of the groups we find a flooding vector, flooding is turned on */
2552 ed->eEDtype = EssentialDynamicsType::Flooding;
2554 /* Print what flavor of flooding we will do */
2555 if (0 == edi->flood.tau) /* constant flooding strength */
2557 fprintf(ed->edo, "Efl_null = %g", edi->flood.constEfl);
2558 if (edi->flood.bHarmonic)
2560 fprintf(ed->edo, ", harmonic");
2563 else /* adaptive flooding */
2565 fprintf(ed->edo, ", adaptive");
2568 fprintf(ed->edo, "\n");
2570 ++edi;
2573 /* Print a nice legend */
2574 snew(LegendStr, 1);
2575 LegendStr[0] = '\0';
2576 sprintf(buf, "# %6s", "time");
2577 add_to_string(&LegendStr, buf);
2579 /* Calculate the maximum number of columns we could end up with */
2580 edi = ed->edpar.begin();
2581 nsets = 0;
2582 for (nr_edi = 1; nr_edi <= nED; nr_edi++)
2584 nsets += 5 + edi->vecs.mon.neig + edi->vecs.linfix.neig + edi->vecs.linacc.neig
2585 + edi->vecs.radfix.neig + edi->vecs.radacc.neig + edi->vecs.radcon.neig
2586 + 6 * edi->flood.vecs.neig;
2587 ++edi;
2589 snew(setname, nsets);
2591 /* In the mdrun time step in a first function call (do_flood()) the flooding
2592 * forces are calculated and in a second function call (do_edsam()) the
2593 * ED constraints. To get a corresponding legend, we need to loop twice
2594 * over the edi groups and output first the flooding, then the ED part */
2596 /* The flooding-related legend entries, if flooding is done */
2597 nsets = 0;
2598 if (EssentialDynamicsType::Flooding == ed->eEDtype)
2600 edi = ed->edpar.begin();
2601 for (nr_edi = 1; nr_edi <= nED; nr_edi++)
2603 /* Always write out the projection on the flooding EVs. Of course, this can also
2604 * be achieved with the monitoring option in do_edsam() (if switched on by the
2605 * user), but in that case the positions need to be communicated in do_edsam(),
2606 * which is not necessary when doing flooding only. */
2607 nice_legend(&setname, &nsets, &LegendStr, "RMSD to ref", "nm", get_EDgroupChar(nr_edi, nED));
2609 for (i = 0; i < edi->flood.vecs.neig; i++)
2611 sprintf(buf, "EV%dprjFLOOD", edi->flood.vecs.ieig[i]);
2612 nice_legend(&setname, &nsets, &LegendStr, buf, "nm", get_EDgroupChar(nr_edi, nED));
2614 /* Output the current reference projection if it changes with time;
2615 * this can happen when flooding is used as harmonic restraint */
2616 if (edi->flood.bHarmonic && edi->flood.referenceProjectionSlope[i] != 0.0)
2618 sprintf(buf, "EV%d ref.prj.", edi->flood.vecs.ieig[i]);
2619 nice_legend(&setname, &nsets, &LegendStr, buf, "nm", get_EDgroupChar(nr_edi, nED));
2622 /* For flooding we also output Efl, Vfl, deltaF, and the flooding forces */
2623 if (0 != edi->flood.tau) /* only output Efl for adaptive flooding (constant otherwise) */
2625 sprintf(buf, "EV%d-Efl", edi->flood.vecs.ieig[i]);
2626 nice_legend(&setname, &nsets, &LegendStr, buf, "kJ/mol", get_EDgroupChar(nr_edi, nED));
2629 sprintf(buf, "EV%d-Vfl", edi->flood.vecs.ieig[i]);
2630 nice_legend(&setname, &nsets, &LegendStr, buf, "kJ/mol", get_EDgroupChar(nr_edi, nED));
2632 if (0 != edi->flood.tau) /* only output deltaF for adaptive flooding (zero otherwise) */
2634 sprintf(buf, "EV%d-deltaF", edi->flood.vecs.ieig[i]);
2635 nice_legend(&setname, &nsets, &LegendStr, buf, "kJ/mol", get_EDgroupChar(nr_edi, nED));
2638 sprintf(buf, "EV%d-FLforces", edi->flood.vecs.ieig[i]);
2639 nice_legend(&setname, &nsets, &LegendStr, buf, "kJ/mol/nm", get_EDgroupChar(nr_edi, nED));
2642 ++edi;
2643 } /* End of flooding-related legend entries */
2645 n_flood = nsets;
2647 /* Now the ED-related entries, if essential dynamics is done */
2648 edi = ed->edpar.begin();
2649 for (nr_edi = 1; nr_edi <= nED; nr_edi++)
2651 if (bNeedDoEdsam(*edi)) /* Only print ED legend if at least one ED option is on */
2653 nice_legend(&setname, &nsets, &LegendStr, "RMSD to ref", "nm", get_EDgroupChar(nr_edi, nED));
2655 /* Essential dynamics, projections on eigenvectors */
2656 nice_legend_evec(&setname, &nsets, &LegendStr, &edi->vecs.mon,
2657 get_EDgroupChar(nr_edi, nED), "MON");
2658 nice_legend_evec(&setname, &nsets, &LegendStr, &edi->vecs.linfix,
2659 get_EDgroupChar(nr_edi, nED), "LINFIX");
2660 nice_legend_evec(&setname, &nsets, &LegendStr, &edi->vecs.linacc,
2661 get_EDgroupChar(nr_edi, nED), "LINACC");
2662 nice_legend_evec(&setname, &nsets, &LegendStr, &edi->vecs.radfix,
2663 get_EDgroupChar(nr_edi, nED), "RADFIX");
2664 if (edi->vecs.radfix.neig)
2666 nice_legend(&setname, &nsets, &LegendStr, "RADFIX radius", "nm",
2667 get_EDgroupChar(nr_edi, nED));
2669 nice_legend_evec(&setname, &nsets, &LegendStr, &edi->vecs.radacc,
2670 get_EDgroupChar(nr_edi, nED), "RADACC");
2671 if (edi->vecs.radacc.neig)
2673 nice_legend(&setname, &nsets, &LegendStr, "RADACC radius", "nm",
2674 get_EDgroupChar(nr_edi, nED));
2676 nice_legend_evec(&setname, &nsets, &LegendStr, &edi->vecs.radcon,
2677 get_EDgroupChar(nr_edi, nED), "RADCON");
2678 if (edi->vecs.radcon.neig)
2680 nice_legend(&setname, &nsets, &LegendStr, "RADCON radius", "nm",
2681 get_EDgroupChar(nr_edi, nED));
2684 ++edi;
2685 } /* end of 'pure' essential dynamics legend entries */
2686 n_edsam = nsets - n_flood;
2688 xvgr_legend(ed->edo, nsets, setname, oenv);
2689 sfree(setname);
2691 fprintf(ed->edo,
2692 "#\n"
2693 "# Legend for %d column%s of flooding plus %d column%s of essential dynamics data:\n",
2694 n_flood, 1 == n_flood ? "" : "s", n_edsam, 1 == n_edsam ? "" : "s");
2695 fprintf(ed->edo, "%s", LegendStr);
2696 sfree(LegendStr);
2698 fflush(ed->edo);
2702 /* Init routine for ED and flooding. Calls init_edi in a loop for every .edi-cycle
2703 * contained in the input file, creates a NULL terminated list of t_edpar structures */
2704 std::unique_ptr<gmx::EssentialDynamics> init_edsam(const gmx::MDLogger& mdlog,
2705 const char* ediFileName,
2706 const char* edoFileName,
2707 const gmx_mtop_t* mtop,
2708 const t_inputrec* ir,
2709 const t_commrec* cr,
2710 gmx::Constraints* constr,
2711 const t_state* globalState,
2712 ObservablesHistory* oh,
2713 const gmx_output_env_t* oenv,
2714 const gmx::StartingBehavior startingBehavior)
2716 int i, avindex;
2717 rvec* x_pbc = nullptr; /* positions of the whole MD system with pbc removed */
2718 rvec * xfit = nullptr, *xstart = nullptr; /* dummy arrays to determine initial RMSDs */
2719 rvec fit_transvec; /* translation ... */
2720 matrix fit_rotmat; /* ... and rotation from fit to reference structure */
2721 rvec* ref_x_old = nullptr; /* helper pointer */
2724 if (MASTER(cr))
2726 fprintf(stderr, "ED: Initializing essential dynamics constraints.\n");
2728 if (oh->edsamHistory != nullptr && oh->edsamHistory->bFromCpt && !gmx_fexist(ediFileName))
2730 gmx_fatal(FARGS,
2731 "The checkpoint file you provided is from an essential dynamics or flooding\n"
2732 "simulation. Please also set the .edi file on the command line with -ei.\n");
2735 GMX_LOG(mdlog.info)
2736 .asParagraph()
2737 .appendText(
2738 "Activating essential dynamics via files passed to "
2739 "gmx mdrun -edi will change in future to be files passed to "
2740 "gmx grompp and the related .mdp options may change also.");
2742 /* Open input and output files, allocate space for ED data structure */
2743 auto edHandle = ed_open(mtop->natoms, oh, ediFileName, edoFileName, startingBehavior, oenv, cr);
2744 auto ed = edHandle->getLegacyED();
2745 GMX_RELEASE_ASSERT(constr != nullptr, "Must have valid constraints object");
2746 constr->saveEdsamPointer(ed);
2748 /* Needed for initializing radacc radius in do_edsam */
2749 ed->bFirst = TRUE;
2751 /* The input file is read by the master and the edi structures are
2752 * initialized here. Input is stored in ed->edpar. Then the edi
2753 * structures are transferred to the other nodes */
2754 if (MASTER(cr))
2756 /* Initialization for every ED/flooding group. Flooding uses one edi group per
2757 * flooding vector, Essential dynamics can be applied to more than one structure
2758 * as well, but will be done in the order given in the edi file, so
2759 * expect different results for different order of edi file concatenation! */
2760 for (auto& edi : ed->edpar)
2762 init_edi(mtop, &edi);
2763 init_flood(&edi, ed, ir->delta_t);
2767 /* The master does the work here. The other nodes get the positions
2768 * not before dd_partition_system which is called after init_edsam */
2769 if (MASTER(cr))
2771 edsamhistory_t* EDstate = oh->edsamHistory.get();
2773 if (!EDstate->bFromCpt)
2775 /* Remove PBC, make molecule(s) subject to ED whole. */
2776 snew(x_pbc, mtop->natoms);
2777 copy_rvecn(globalState->x.rvec_array(), x_pbc, 0, mtop->natoms);
2778 do_pbc_first_mtop(nullptr, ir->pbcType, globalState->box, mtop, x_pbc);
2780 /* Reset pointer to first ED data set which contains the actual ED data */
2781 auto edi = ed->edpar.begin();
2782 GMX_ASSERT(!ed->edpar.empty(), "Essential dynamics parameters is undefined");
2784 /* Loop over all ED/flooding data sets (usually only one, though) */
2785 for (int nr_edi = 1; nr_edi <= EDstate->nED; nr_edi++)
2787 /* For multiple ED groups we use the output frequency that was specified
2788 * in the first set */
2789 if (nr_edi > 1)
2791 edi->outfrq = ed->edpar.begin()->outfrq;
2794 /* Extract the initial reference and average positions. When starting
2795 * from .cpt, these have already been read into sref.x_old
2796 * in init_edsamstate() */
2797 if (!EDstate->bFromCpt)
2799 /* If this is the first run (i.e. no checkpoint present) we assume
2800 * that the starting positions give us the correct PBC representation */
2801 for (i = 0; i < edi->sref.nr; i++)
2803 copy_rvec(x_pbc[edi->sref.anrs[i]], edi->sref.x_old[i]);
2806 for (i = 0; i < edi->sav.nr; i++)
2808 copy_rvec(x_pbc[edi->sav.anrs[i]], edi->sav.x_old[i]);
2812 /* Now we have the PBC-correct start positions of the reference and
2813 average structure. We copy that over to dummy arrays on which we
2814 can apply fitting to print out the RMSD. We srenew the memory since
2815 the size of the buffers is likely different for every ED group */
2816 srenew(xfit, edi->sref.nr);
2817 srenew(xstart, edi->sav.nr);
2818 if (edi->bRefEqAv)
2820 /* Reference indices are the same as average indices */
2821 ref_x_old = edi->sav.x_old;
2823 else
2825 ref_x_old = edi->sref.x_old;
2827 copy_rvecn(ref_x_old, xfit, 0, edi->sref.nr);
2828 copy_rvecn(edi->sav.x_old, xstart, 0, edi->sav.nr);
2830 /* Make the fit to the REFERENCE structure, get translation and rotation */
2831 fit_to_reference(xfit, fit_transvec, fit_rotmat, &(*edi));
2833 /* Output how well we fit to the reference at the start */
2834 translate_and_rotate(xfit, edi->sref.nr, fit_transvec, fit_rotmat);
2835 fprintf(stderr, "ED: Initial RMSD from reference after fit = %f nm",
2836 rmsd_from_structure(xfit, &edi->sref));
2837 if (EDstate->nED > 1)
2839 fprintf(stderr, " (ED group %c)", get_EDgroupChar(nr_edi, EDstate->nED));
2841 fprintf(stderr, "\n");
2843 /* Now apply the translation and rotation to the atoms on which ED sampling will be performed */
2844 translate_and_rotate(xstart, edi->sav.nr, fit_transvec, fit_rotmat);
2846 /* calculate initial projections */
2847 project(xstart, &(*edi));
2849 /* For the target and origin structure both a reference (fit) and an
2850 * average structure can be provided in make_edi. If both structures
2851 * are the same, make_edi only stores one of them in the .edi file.
2852 * If they differ, first the fit and then the average structure is stored
2853 * in star (or sor), thus the number of entries in star/sor is
2854 * (n_fit + n_av) with n_fit the size of the fitting group and n_av
2855 * the size of the average group. */
2857 /* process target structure, if required */
2858 if (edi->star.nr > 0)
2860 fprintf(stderr, "ED: Fitting target structure to reference structure\n");
2862 /* get translation & rotation for fit of target structure to reference structure */
2863 fit_to_reference(edi->star.x, fit_transvec, fit_rotmat, &(*edi));
2864 /* do the fit */
2865 translate_and_rotate(edi->star.x, edi->star.nr, fit_transvec, fit_rotmat);
2866 if (edi->star.nr == edi->sav.nr)
2868 avindex = 0;
2870 else /* edi->star.nr = edi->sref.nr + edi->sav.nr */
2872 /* The last sav.nr indices of the target structure correspond to
2873 * the average structure, which must be projected */
2874 avindex = edi->star.nr - edi->sav.nr;
2876 rad_project(*edi, &edi->star.x[avindex], &edi->vecs.radcon);
2878 else
2880 rad_project(*edi, xstart, &edi->vecs.radcon);
2883 /* process structure that will serve as origin of expansion circle */
2884 if ((EssentialDynamicsType::Flooding == ed->eEDtype) && (!edi->flood.bConstForce))
2886 fprintf(stderr,
2887 "ED: Setting center of flooding potential (0 = average structure)\n");
2890 if (edi->sori.nr > 0)
2892 fprintf(stderr, "ED: Fitting origin structure to reference structure\n");
2894 /* fit this structure to reference structure */
2895 fit_to_reference(edi->sori.x, fit_transvec, fit_rotmat, &(*edi));
2896 /* do the fit */
2897 translate_and_rotate(edi->sori.x, edi->sori.nr, fit_transvec, fit_rotmat);
2898 if (edi->sori.nr == edi->sav.nr)
2900 avindex = 0;
2902 else /* edi->sori.nr = edi->sref.nr + edi->sav.nr */
2904 /* For the projection, we need the last sav.nr indices of sori */
2905 avindex = edi->sori.nr - edi->sav.nr;
2908 rad_project(*edi, &edi->sori.x[avindex], &edi->vecs.radacc);
2909 rad_project(*edi, &edi->sori.x[avindex], &edi->vecs.radfix);
2910 if ((EssentialDynamicsType::Flooding == ed->eEDtype) && (!edi->flood.bConstForce))
2912 fprintf(stderr,
2913 "ED: The ORIGIN structure will define the flooding potential "
2914 "center.\n");
2915 /* Set center of flooding potential to the ORIGIN structure */
2916 rad_project(*edi, &edi->sori.x[avindex], &edi->flood.vecs);
2917 /* We already know that no (moving) reference position was provided,
2918 * therefore we can overwrite refproj[0]*/
2919 copyEvecReference(&edi->flood.vecs, edi->flood.initialReferenceProjection);
2922 else /* No origin structure given */
2924 rad_project(*edi, xstart, &edi->vecs.radacc);
2925 rad_project(*edi, xstart, &edi->vecs.radfix);
2926 if ((EssentialDynamicsType::Flooding == ed->eEDtype) && (!edi->flood.bConstForce))
2928 if (edi->flood.bHarmonic)
2930 fprintf(stderr,
2931 "ED: A (possibly changing) ref. projection will define the "
2932 "flooding potential center.\n");
2933 for (i = 0; i < edi->flood.vecs.neig; i++)
2935 edi->flood.vecs.refproj[i] = edi->flood.initialReferenceProjection[i];
2938 else
2940 fprintf(stderr,
2941 "ED: The AVERAGE structure will define the flooding potential "
2942 "center.\n");
2943 /* Set center of flooding potential to the center of the covariance matrix,
2944 * i.e. the average structure, i.e. zero in the projected system */
2945 for (i = 0; i < edi->flood.vecs.neig; i++)
2947 edi->flood.vecs.refproj[i] = 0.0;
2952 /* For convenience, output the center of the flooding potential for the eigenvectors */
2953 if ((EssentialDynamicsType::Flooding == ed->eEDtype) && (!edi->flood.bConstForce))
2955 for (i = 0; i < edi->flood.vecs.neig; i++)
2957 fprintf(stdout, "ED: EV %d flooding potential center: %11.4e",
2958 edi->flood.vecs.ieig[i], edi->flood.vecs.refproj[i]);
2959 if (edi->flood.bHarmonic)
2961 fprintf(stdout, " (adding %11.4e/timestep)",
2962 edi->flood.referenceProjectionSlope[i]);
2964 fprintf(stdout, "\n");
2968 /* set starting projections for linsam */
2969 rad_project(*edi, xstart, &edi->vecs.linacc);
2970 rad_project(*edi, xstart, &edi->vecs.linfix);
2972 /* Prepare for the next edi data set: */
2973 ++edi;
2975 /* Cleaning up on the master node: */
2976 if (!EDstate->bFromCpt)
2978 sfree(x_pbc);
2980 sfree(xfit);
2981 sfree(xstart);
2983 } /* end of MASTER only section */
2985 if (PAR(cr))
2987 /* Broadcast the essential dynamics / flooding data to all nodes */
2988 broadcast_ed_data(cr, ed);
2990 else
2992 /* In the single-CPU case, point the local atom numbers pointers to the global
2993 * one, so that we can use the same notation in serial and parallel case: */
2994 /* Loop over all ED data sets (usually only one, though) */
2995 for (auto edi = ed->edpar.begin(); edi != ed->edpar.end(); ++edi)
2997 edi->sref.anrs_loc = edi->sref.anrs;
2998 edi->sav.anrs_loc = edi->sav.anrs;
2999 edi->star.anrs_loc = edi->star.anrs;
3000 edi->sori.anrs_loc = edi->sori.anrs;
3001 /* For the same reason as above, make a dummy c_ind array: */
3002 snew(edi->sav.c_ind, edi->sav.nr);
3003 /* Initialize the array */
3004 for (i = 0; i < edi->sav.nr; i++)
3006 edi->sav.c_ind[i] = i;
3008 /* In the general case we will need a different-sized array for the reference indices: */
3009 if (!edi->bRefEqAv)
3011 snew(edi->sref.c_ind, edi->sref.nr);
3012 for (i = 0; i < edi->sref.nr; i++)
3014 edi->sref.c_ind[i] = i;
3017 /* Point to the very same array in case of other structures: */
3018 edi->star.c_ind = edi->sav.c_ind;
3019 edi->sori.c_ind = edi->sav.c_ind;
3020 /* In the serial case, the local number of atoms is the global one: */
3021 edi->sref.nr_loc = edi->sref.nr;
3022 edi->sav.nr_loc = edi->sav.nr;
3023 edi->star.nr_loc = edi->star.nr;
3024 edi->sori.nr_loc = edi->sori.nr;
3028 /* Allocate space for ED buffer variables */
3029 /* Again, loop over ED data sets */
3030 for (auto edi = ed->edpar.begin(); edi != ed->edpar.end(); ++edi)
3032 /* Allocate space for ED buffer variables */
3033 snew_bc(MASTER(cr), edi->buf, 1); /* MASTER has already allocated edi->buf in init_edi() */
3034 snew(edi->buf->do_edsam, 1);
3036 /* Space for collective ED buffer variables */
3038 /* Collective positions of atoms with the average indices */
3039 snew(edi->buf->do_edsam->xcoll, edi->sav.nr);
3040 snew(edi->buf->do_edsam->shifts_xcoll, edi->sav.nr); /* buffer for xcoll shifts */
3041 snew(edi->buf->do_edsam->extra_shifts_xcoll, edi->sav.nr);
3042 /* Collective positions of atoms with the reference indices */
3043 if (!edi->bRefEqAv)
3045 snew(edi->buf->do_edsam->xc_ref, edi->sref.nr);
3046 snew(edi->buf->do_edsam->shifts_xc_ref, edi->sref.nr); /* To store the shifts in */
3047 snew(edi->buf->do_edsam->extra_shifts_xc_ref, edi->sref.nr);
3050 /* Get memory for flooding forces */
3051 snew(edi->flood.forces_cartesian, edi->sav.nr);
3054 /* Flush the edo file so that the user can check some things
3055 * when the simulation has started */
3056 if (ed->edo)
3058 fflush(ed->edo);
3061 return edHandle;
3065 void do_edsam(const t_inputrec* ir, int64_t step, const t_commrec* cr, rvec xs[], rvec v[], const matrix box, gmx_edsam* ed)
3067 int i, edinr, iupdate = 500;
3068 matrix rotmat; /* rotation matrix */
3069 rvec transvec; /* translation vector */
3070 rvec dv, dx, x_unsh; /* tmp vectors for velocity, distance, unshifted x coordinate */
3071 real dt_1; /* 1/dt */
3072 struct t_do_edsam* buf;
3073 real rmsdev = -1; /* RMSD from reference structure prior to applying the constraints */
3074 gmx_bool bSuppress = FALSE; /* Write .xvg output file on master? */
3077 /* Check if ED sampling has to be performed */
3078 if (ed->eEDtype == EssentialDynamicsType::None)
3080 return;
3083 dt_1 = 1.0 / ir->delta_t;
3085 /* Loop over all ED groups (usually one) */
3086 edinr = 0;
3087 for (auto& edi : ed->edpar)
3089 edinr++;
3090 if (bNeedDoEdsam(edi))
3093 buf = edi.buf->do_edsam;
3095 if (ed->bFirst)
3097 /* initialize radacc radius for slope criterion */
3098 buf->oldrad = calc_radius(edi.vecs.radacc);
3101 /* Copy the positions into buf->xc* arrays and after ED
3102 * feed back corrections to the official positions */
3104 /* Broadcast the ED positions such that every node has all of them
3105 * Every node contributes its local positions xs and stores it in
3106 * the collective buf->xcoll array. Note that for edinr > 1
3107 * xs could already have been modified by an earlier ED */
3109 communicate_group_positions(cr, buf->xcoll, buf->shifts_xcoll, buf->extra_shifts_xcoll,
3110 PAR(cr) ? buf->bUpdateShifts : TRUE, xs, edi.sav.nr, edi.sav.nr_loc,
3111 edi.sav.anrs_loc, edi.sav.c_ind, edi.sav.x_old, box);
3113 /* Only assembly reference positions if their indices differ from the average ones */
3114 if (!edi.bRefEqAv)
3116 communicate_group_positions(
3117 cr, buf->xc_ref, buf->shifts_xc_ref, buf->extra_shifts_xc_ref,
3118 PAR(cr) ? buf->bUpdateShifts : TRUE, xs, edi.sref.nr, edi.sref.nr_loc,
3119 edi.sref.anrs_loc, edi.sref.c_ind, edi.sref.x_old, box);
3122 /* If bUpdateShifts was TRUE then the shifts have just been updated in communicate_group_positions.
3123 * We do not need to update the shifts until the next NS step. Note that dd_make_local_ed_indices
3124 * set bUpdateShifts=TRUE in the parallel case. */
3125 buf->bUpdateShifts = FALSE;
3127 /* Now all nodes have all of the ED positions in edi.sav->xcoll,
3128 * as well as the indices in edi.sav.anrs */
3130 /* Fit the reference indices to the reference structure */
3131 if (edi.bRefEqAv)
3133 fit_to_reference(buf->xcoll, transvec, rotmat, &edi);
3135 else
3137 fit_to_reference(buf->xc_ref, transvec, rotmat, &edi);
3140 /* Now apply the translation and rotation to the ED structure */
3141 translate_and_rotate(buf->xcoll, edi.sav.nr, transvec, rotmat);
3143 /* Find out how well we fit to the reference (just for output steps) */
3144 if (do_per_step(step, edi.outfrq) && MASTER(cr))
3146 if (edi.bRefEqAv)
3148 /* Indices of reference and average structures are identical,
3149 * thus we can calculate the rmsd to SREF using xcoll */
3150 rmsdev = rmsd_from_structure(buf->xcoll, &edi.sref);
3152 else
3154 /* We have to translate & rotate the reference atoms first */
3155 translate_and_rotate(buf->xc_ref, edi.sref.nr, transvec, rotmat);
3156 rmsdev = rmsd_from_structure(buf->xc_ref, &edi.sref);
3160 /* update radsam references, when required */
3161 if (do_per_step(step, edi.maxedsteps) && step >= edi.presteps)
3163 project(buf->xcoll, &edi);
3164 rad_project(edi, buf->xcoll, &edi.vecs.radacc);
3165 rad_project(edi, buf->xcoll, &edi.vecs.radfix);
3166 buf->oldrad = -1.e5;
3169 /* update radacc references, when required */
3170 if (do_per_step(step, iupdate) && step >= edi.presteps)
3172 edi.vecs.radacc.radius = calc_radius(edi.vecs.radacc);
3173 if (edi.vecs.radacc.radius - buf->oldrad < edi.slope)
3175 project(buf->xcoll, &edi);
3176 rad_project(edi, buf->xcoll, &edi.vecs.radacc);
3177 buf->oldrad = 0.0;
3179 else
3181 buf->oldrad = edi.vecs.radacc.radius;
3185 /* apply the constraints */
3186 if (step >= edi.presteps && ed_constraints(ed->eEDtype, edi))
3188 /* ED constraints should be applied already in the first MD step
3189 * (which is step 0), therefore we pass step+1 to the routine */
3190 ed_apply_constraints(buf->xcoll, &edi, step + 1 - ir->init_step);
3193 /* write to edo, when required */
3194 if (do_per_step(step, edi.outfrq))
3196 project(buf->xcoll, &edi);
3197 if (MASTER(cr) && !bSuppress)
3199 write_edo(edi, ed->edo, rmsdev);
3203 /* Copy back the positions unless monitoring only */
3204 if (ed_constraints(ed->eEDtype, edi))
3206 /* remove fitting */
3207 rmfit(edi.sav.nr, buf->xcoll, transvec, rotmat);
3209 /* Copy the ED corrected positions into the coordinate array */
3210 /* Each node copies its local part. In the serial case, nat_loc is the
3211 * total number of ED atoms */
3212 for (i = 0; i < edi.sav.nr_loc; i++)
3214 /* Unshift local ED coordinate and store in x_unsh */
3215 ed_unshift_single_coord(box, buf->xcoll[edi.sav.c_ind[i]],
3216 buf->shifts_xcoll[edi.sav.c_ind[i]], x_unsh);
3218 /* dx is the ED correction to the positions: */
3219 rvec_sub(x_unsh, xs[edi.sav.anrs_loc[i]], dx);
3221 if (v != nullptr)
3223 /* dv is the ED correction to the velocity: */
3224 svmul(dt_1, dx, dv);
3225 /* apply the velocity correction: */
3226 rvec_inc(v[edi.sav.anrs_loc[i]], dv);
3228 /* Finally apply the position correction due to ED: */
3229 copy_rvec(x_unsh, xs[edi.sav.anrs_loc[i]]);
3232 } /* END of if ( bNeedDoEdsam(edi) ) */
3234 /* Prepare for the next ED group */
3236 } /* END of loop over ED groups */
3238 ed->bFirst = FALSE;