Split lines with many copyright years
[gromacs.git] / src / gromacs / pulling / pull_rotation.cpp
blob1ae8f45edda3309d9d9efcfaa2a90ca19e9b2a83
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-2008, 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 "pull_rotation.h"
42 #include "config.h"
44 #include <cstdio>
45 #include <cstdlib>
46 #include <cstring>
48 #include <algorithm>
49 #include <memory>
51 #include "gromacs/commandline/filenm.h"
52 #include "gromacs/domdec/dlbtiming.h"
53 #include "gromacs/domdec/domdec_struct.h"
54 #include "gromacs/domdec/ga2la.h"
55 #include "gromacs/domdec/localatomset.h"
56 #include "gromacs/domdec/localatomsetmanager.h"
57 #include "gromacs/fileio/gmxfio.h"
58 #include "gromacs/fileio/xvgr.h"
59 #include "gromacs/gmxlib/network.h"
60 #include "gromacs/linearalgebra/nrjac.h"
61 #include "gromacs/math/functions.h"
62 #include "gromacs/math/utilities.h"
63 #include "gromacs/math/vec.h"
64 #include "gromacs/mdlib/groupcoord.h"
65 #include "gromacs/mdlib/stat.h"
66 #include "gromacs/mdrunutility/handlerestart.h"
67 #include "gromacs/mdtypes/commrec.h"
68 #include "gromacs/mdtypes/inputrec.h"
69 #include "gromacs/mdtypes/md_enums.h"
70 #include "gromacs/mdtypes/mdrunoptions.h"
71 #include "gromacs/mdtypes/state.h"
72 #include "gromacs/pbcutil/pbc.h"
73 #include "gromacs/timing/cyclecounter.h"
74 #include "gromacs/timing/wallcycle.h"
75 #include "gromacs/topology/mtop_lookup.h"
76 #include "gromacs/topology/mtop_util.h"
77 #include "gromacs/utility/basedefinitions.h"
78 #include "gromacs/utility/fatalerror.h"
79 #include "gromacs/utility/pleasecite.h"
80 #include "gromacs/utility/smalloc.h"
82 static char const* RotStr = { "Enforced rotation:" };
84 /* Set the minimum weight for the determination of the slab centers */
85 #define WEIGHT_MIN (10 * GMX_FLOAT_MIN)
87 //! Helper structure for sorting positions along rotation vector
88 struct sort_along_vec_t
90 //! Projection of xc on the rotation vector
91 real xcproj;
92 //! Index of xc
93 int ind;
94 //! Mass
95 real m;
96 //! Position
97 rvec x;
98 //! Reference position
99 rvec x_ref;
103 //! Enforced rotation / flexible: determine the angle of each slab
104 struct gmx_slabdata
106 //! Number of atoms belonging to this slab
107 int nat;
108 /*! \brief The positions belonging to this slab.
110 * In general, this should be all positions of the whole
111 * rotation group, but we leave those away that have a small
112 * enough weight. */
113 rvec* x;
114 //! Same for reference
115 rvec* ref;
116 //! The weight for each atom
117 real* weight;
121 //! Helper structure for potential fitting
122 struct gmx_potfit
124 /*! \brief Set of angles for which the potential is calculated.
126 * The optimum fit is determined as the angle for with the
127 * potential is minimal. */
128 real* degangle;
129 //! Potential for the different angles
130 real* V;
131 //! Rotation matrix corresponding to the angles
132 matrix* rotmat;
136 //! Enforced rotation data for a single rotation group
137 struct gmx_enfrotgrp
139 //! Input parameters for this group
140 const t_rotgrp* rotg = nullptr;
141 //! Index of this group within the set of groups
142 int groupIndex;
143 //! Rotation angle in degrees
144 real degangle;
145 //! Rotation matrix
146 matrix rotmat;
147 //! The atoms subject to enforced rotation
148 std::unique_ptr<gmx::LocalAtomSet> atomSet;
150 //! The normalized rotation vector
151 rvec vec;
152 //! Rotation potential for this rotation group
153 real V;
154 //! Array to store the forces on the local atoms resulting from enforced rotation potential
155 rvec* f_rot_loc;
157 /* Collective coordinates for the whole rotation group */
158 //! Length of each x_rotref vector after x_rotref has been put into origin
159 real* xc_ref_length;
160 //! Center of the rotation group positions, may be mass weighted
161 rvec xc_center;
162 //! Center of the rotation group reference positions
163 rvec xc_ref_center;
164 //! Current (collective) positions
165 rvec* xc;
166 //! Current (collective) shifts
167 ivec* xc_shifts;
168 //! Extra shifts since last DD step
169 ivec* xc_eshifts;
170 //! Old (collective) positions
171 rvec* xc_old;
172 //! Normalized form of the current positions
173 rvec* xc_norm;
174 //! Reference positions (sorted in the same order as xc when sorted)
175 rvec* xc_ref_sorted;
176 //! Where is a position found after sorting?
177 int* xc_sortind;
178 //! Collective masses
179 real* mc;
180 //! Collective masses sorted
181 real* mc_sorted;
182 //! one over the total mass of the rotation group
183 real invmass;
185 //! Torque in the direction of rotation vector
186 real torque_v;
187 //! Actual angle of the whole rotation group
188 real angle_v;
189 /* Fixed rotation only */
190 //! Weights for angle determination
191 real weight_v;
192 //! Local reference coords, correctly rotated
193 rvec* xr_loc;
194 //! Local current coords, correct PBC image
195 rvec* x_loc_pbc;
196 //! Masses of the current local atoms
197 real* m_loc;
199 /* Flexible rotation only */
200 //! For this many slabs memory is allocated
201 int nslabs_alloc;
202 //! Lowermost slab for that the calculation needs to be performed at a given time step
203 int slab_first;
204 //! Uppermost slab ...
205 int slab_last;
206 //! First slab for which ref. center is stored
207 int slab_first_ref;
208 //! Last ...
209 int slab_last_ref;
210 //! Slab buffer region around reference slabs
211 int slab_buffer;
212 //! First relevant atom for a slab
213 int* firstatom;
214 //! Last relevant atom for a slab
215 int* lastatom;
216 //! Gaussian-weighted slab center
217 rvec* slab_center;
218 //! Gaussian-weighted slab center for the reference positions
219 rvec* slab_center_ref;
220 //! Sum of gaussian weights in a slab
221 real* slab_weights;
222 //! Torque T = r x f for each slab. torque_v = m.v = angular momentum in the direction of v
223 real* slab_torque_v;
224 //! min_gaussian from t_rotgrp is the minimum value the gaussian must have so that the force is actually evaluated. max_beta is just another way to put it
225 real max_beta;
226 //! Precalculated gaussians for a single atom
227 real* gn_atom;
228 //! Tells to which slab each precalculated gaussian belongs
229 int* gn_slabind;
230 //! Inner sum of the flexible2 potential per slab; this is precalculated for optimization reasons
231 rvec* slab_innersumvec;
232 //! Holds atom positions and gaussian weights of atoms belonging to a slab
233 gmx_slabdata* slab_data;
235 /* For potential fits with varying angle: */
236 //! Used for fit type 'potential'
237 gmx_potfit* PotAngleFit;
241 //! Enforced rotation data for all groups
242 struct gmx_enfrot
244 //! Input parameters.
245 const t_rot* rot = nullptr;
246 //! Output period for main rotation outfile
247 int nstrout;
248 //! Output period for per-slab data
249 int nstsout;
250 //! Output file for rotation data
251 FILE* out_rot = nullptr;
252 //! Output file for torque data
253 FILE* out_torque = nullptr;
254 //! Output file for slab angles for flexible type
255 FILE* out_angles = nullptr;
256 //! Output file for slab centers
257 FILE* out_slabs = nullptr;
258 //! Allocation size of buf
259 int bufsize = 0;
260 //! Coordinate buffer variable for sorting
261 rvec* xbuf = nullptr;
262 //! Masses buffer variable for sorting
263 real* mbuf = nullptr;
264 //! Buffer variable needed for position sorting
265 sort_along_vec_t* data = nullptr;
266 //! MPI buffer
267 real* mpi_inbuf = nullptr;
268 //! MPI buffer
269 real* mpi_outbuf = nullptr;
270 //! Allocation size of in & outbuf
271 int mpi_bufsize = 0;
272 //! If true, append output files
273 gmx_bool restartWithAppending = false;
274 //! Used to skip first output when appending to avoid duplicate entries in rotation outfiles
275 gmx_bool bOut = false;
276 //! Stores working data per group
277 std::vector<gmx_enfrotgrp> enfrotgrp;
278 ~gmx_enfrot();
281 gmx_enfrot::~gmx_enfrot()
283 if (out_rot)
285 gmx_fio_fclose(out_rot);
287 if (out_slabs)
289 gmx_fio_fclose(out_slabs);
291 if (out_angles)
293 gmx_fio_fclose(out_angles);
295 if (out_torque)
297 gmx_fio_fclose(out_torque);
301 namespace gmx
304 extern template LocalAtomSet LocalAtomSetManager::add<void, void>(ArrayRef<const int> globalAtomIndex);
306 class EnforcedRotation::Impl
308 public:
309 gmx_enfrot enforcedRotation_;
312 EnforcedRotation::EnforcedRotation() : impl_(new Impl) {}
314 EnforcedRotation::~EnforcedRotation() = default;
316 gmx_enfrot* EnforcedRotation::getLegacyEnfrot()
318 return &impl_->enforcedRotation_;
321 } // namespace gmx
323 /* Activate output of forces for correctness checks */
324 /* #define PRINT_FORCES */
325 #ifdef PRINT_FORCES
326 # define PRINT_FORCE_J \
327 fprintf(stderr, "f%d = %15.8f %15.8f %15.8f\n", erg->xc_ref_ind[j], erg->f_rot_loc[j][XX], \
328 erg->f_rot_loc[j][YY], erg->f_rot_loc[j][ZZ]);
329 # define PRINT_POT_TAU \
330 if (MASTER(cr)) \
332 fprintf(stderr, \
333 "potential = %15.8f\n" \
334 "torque = %15.8f\n", \
335 erg->V, erg->torque_v); \
337 #else
338 # define PRINT_FORCE_J
339 # define PRINT_POT_TAU
340 #endif
342 /* Shortcuts for often used queries */
343 #define ISFLEX(rg) \
344 (((rg)->eType == erotgFLEX) || ((rg)->eType == erotgFLEXT) || ((rg)->eType == erotgFLEX2) \
345 || ((rg)->eType == erotgFLEX2T))
346 #define ISCOLL(rg) \
347 (((rg)->eType == erotgFLEX) || ((rg)->eType == erotgFLEXT) || ((rg)->eType == erotgFLEX2) \
348 || ((rg)->eType == erotgFLEX2T) || ((rg)->eType == erotgRMPF) || ((rg)->eType == erotgRM2PF))
351 /* Does any of the rotation groups use slab decomposition? */
352 static gmx_bool HaveFlexibleGroups(const t_rot* rot)
354 for (int g = 0; g < rot->ngrp; g++)
356 if (ISFLEX(&rot->grp[g]))
358 return TRUE;
362 return FALSE;
366 /* Is for any group the fit angle determined by finding the minimum of the
367 * rotation potential? */
368 static gmx_bool HavePotFitGroups(const t_rot* rot)
370 for (int g = 0; g < rot->ngrp; g++)
372 if (erotgFitPOT == rot->grp[g].eFittype)
374 return TRUE;
378 return FALSE;
382 static double** allocate_square_matrix(int dim)
384 int i;
385 double** mat = nullptr;
388 snew(mat, dim);
389 for (i = 0; i < dim; i++)
391 snew(mat[i], dim);
394 return mat;
398 static void free_square_matrix(double** mat, int dim)
400 int i;
403 for (i = 0; i < dim; i++)
405 sfree(mat[i]);
407 sfree(mat);
411 /* Return the angle for which the potential is minimal */
412 static real get_fitangle(const gmx_enfrotgrp* erg)
414 int i;
415 real fitangle = -999.9;
416 real pot_min = GMX_FLOAT_MAX;
417 gmx_potfit* fit;
420 fit = erg->PotAngleFit;
422 for (i = 0; i < erg->rotg->PotAngle_nstep; i++)
424 if (fit->V[i] < pot_min)
426 pot_min = fit->V[i];
427 fitangle = fit->degangle[i];
431 return fitangle;
435 /* Reduce potential angle fit data for this group at this time step? */
436 static inline gmx_bool bPotAngle(const gmx_enfrot* er, const t_rotgrp* rotg, int64_t step)
438 return ((erotgFitPOT == rotg->eFittype)
439 && (do_per_step(step, er->nstsout) || do_per_step(step, er->nstrout)));
442 /* Reduce slab torqe data for this group at this time step? */
443 static inline gmx_bool bSlabTau(const gmx_enfrot* er, const t_rotgrp* rotg, int64_t step)
445 return ((ISFLEX(rotg)) && do_per_step(step, er->nstsout));
448 /* Output rotation energy, torques, etc. for each rotation group */
449 static void reduce_output(const t_commrec* cr, gmx_enfrot* er, real t, int64_t step)
451 int i, islab, nslabs = 0;
452 int count; /* MPI element counter */
453 real fitangle;
454 gmx_bool bFlex;
456 /* Fill the MPI buffer with stuff to reduce. If items are added for reduction
457 * here, the MPI buffer size has to be enlarged also in calc_mpi_bufsize() */
458 if (PAR(cr))
460 count = 0;
461 for (auto& ergRef : er->enfrotgrp)
463 gmx_enfrotgrp* erg = &ergRef;
464 const t_rotgrp* rotg = erg->rotg;
465 nslabs = erg->slab_last - erg->slab_first + 1;
466 er->mpi_inbuf[count++] = erg->V;
467 er->mpi_inbuf[count++] = erg->torque_v;
468 er->mpi_inbuf[count++] = erg->angle_v;
469 er->mpi_inbuf[count++] =
470 erg->weight_v; /* weights are not needed for flex types, but this is just a single value */
472 if (bPotAngle(er, rotg, step))
474 for (i = 0; i < rotg->PotAngle_nstep; i++)
476 er->mpi_inbuf[count++] = erg->PotAngleFit->V[i];
479 if (bSlabTau(er, rotg, step))
481 for (i = 0; i < nslabs; i++)
483 er->mpi_inbuf[count++] = erg->slab_torque_v[i];
487 if (count > er->mpi_bufsize)
489 gmx_fatal(FARGS, "%s MPI buffer overflow, please report this error.", RotStr);
492 #if GMX_MPI
493 MPI_Reduce(er->mpi_inbuf, er->mpi_outbuf, count, GMX_MPI_REAL, MPI_SUM, MASTERRANK(cr),
494 cr->mpi_comm_mygroup);
495 #endif
497 /* Copy back the reduced data from the buffer on the master */
498 if (MASTER(cr))
500 count = 0;
501 for (auto& ergRef : er->enfrotgrp)
503 gmx_enfrotgrp* erg = &ergRef;
504 const t_rotgrp* rotg = erg->rotg;
505 nslabs = erg->slab_last - erg->slab_first + 1;
506 erg->V = er->mpi_outbuf[count++];
507 erg->torque_v = er->mpi_outbuf[count++];
508 erg->angle_v = er->mpi_outbuf[count++];
509 erg->weight_v = er->mpi_outbuf[count++];
511 if (bPotAngle(er, rotg, step))
513 for (int i = 0; i < rotg->PotAngle_nstep; i++)
515 erg->PotAngleFit->V[i] = er->mpi_outbuf[count++];
518 if (bSlabTau(er, rotg, step))
520 for (int i = 0; i < nslabs; i++)
522 erg->slab_torque_v[i] = er->mpi_outbuf[count++];
529 /* Output */
530 if (MASTER(cr))
532 /* Angle and torque for each rotation group */
533 for (auto& ergRef : er->enfrotgrp)
535 gmx_enfrotgrp* erg = &ergRef;
536 const t_rotgrp* rotg = erg->rotg;
537 bFlex = ISFLEX(rotg);
539 /* Output to main rotation output file: */
540 if (do_per_step(step, er->nstrout))
542 if (erotgFitPOT == rotg->eFittype)
544 fitangle = get_fitangle(erg);
546 else
548 if (bFlex)
550 fitangle = erg->angle_v; /* RMSD fit angle */
552 else
554 fitangle = (erg->angle_v / erg->weight_v) * 180.0 * M_1_PI;
557 fprintf(er->out_rot, "%12.4f", fitangle);
558 fprintf(er->out_rot, "%12.3e", erg->torque_v);
559 fprintf(er->out_rot, "%12.3e", erg->V);
562 if (do_per_step(step, er->nstsout))
564 /* Output to torque log file: */
565 if (bFlex)
567 fprintf(er->out_torque, "%12.3e%6d", t, erg->groupIndex);
568 for (int i = erg->slab_first; i <= erg->slab_last; i++)
570 islab = i - erg->slab_first; /* slab index */
571 /* Only output if enough weight is in slab */
572 if (erg->slab_weights[islab] > rotg->min_gaussian)
574 fprintf(er->out_torque, "%6d%12.3e", i, erg->slab_torque_v[islab]);
577 fprintf(er->out_torque, "\n");
580 /* Output to angles log file: */
581 if (erotgFitPOT == rotg->eFittype)
583 fprintf(er->out_angles, "%12.3e%6d%12.4f", t, erg->groupIndex, erg->degangle);
584 /* Output energies at a set of angles around the reference angle */
585 for (int i = 0; i < rotg->PotAngle_nstep; i++)
587 fprintf(er->out_angles, "%12.3e", erg->PotAngleFit->V[i]);
589 fprintf(er->out_angles, "\n");
593 if (do_per_step(step, er->nstrout))
595 fprintf(er->out_rot, "\n");
601 /* Add the forces from enforced rotation potential to the local forces.
602 * Should be called after the SR forces have been evaluated */
603 real add_rot_forces(gmx_enfrot* er, rvec f[], const t_commrec* cr, int64_t step, real t)
605 real Vrot = 0.0; /* If more than one rotation group is present, Vrot
606 assembles the local parts from all groups */
608 /* Loop over enforced rotation groups (usually 1, though)
609 * Apply the forces from rotation potentials */
610 for (auto& ergRef : er->enfrotgrp)
612 gmx_enfrotgrp* erg = &ergRef;
613 Vrot += erg->V; /* add the local parts from the nodes */
614 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
615 for (gmx::index l = 0; l < localRotationGroupIndex.ssize(); l++)
617 /* Get the right index of the local force */
618 int ii = localRotationGroupIndex[l];
619 /* Add */
620 rvec_inc(f[ii], erg->f_rot_loc[l]);
624 /* Reduce energy,torque, angles etc. to get the sum values (per rotation group)
625 * on the master and output these values to file. */
626 if ((do_per_step(step, er->nstrout) || do_per_step(step, er->nstsout)) && er->bOut)
628 reduce_output(cr, er, t, step);
631 /* When appending, er->bOut is FALSE the first time to avoid duplicate entries */
632 er->bOut = TRUE;
634 PRINT_POT_TAU
636 return Vrot;
640 /* The Gaussian norm is chosen such that the sum of the gaussian functions
641 * over the slabs is approximately 1.0 everywhere */
642 #define GAUSS_NORM 0.569917543430618
645 /* Calculate the maximum beta that leads to a gaussian larger min_gaussian,
646 * also does some checks
648 static double calc_beta_max(real min_gaussian, real slab_dist)
650 double sigma;
651 double arg;
654 /* Actually the next two checks are already made in grompp */
655 if (slab_dist <= 0)
657 gmx_fatal(FARGS, "Slab distance of flexible rotation groups must be >=0 !");
659 if (min_gaussian <= 0)
661 gmx_fatal(FARGS, "Cutoff value for Gaussian must be > 0. (You requested %f)", min_gaussian);
664 /* Define the sigma value */
665 sigma = 0.7 * slab_dist;
667 /* Calculate the argument for the logarithm and check that the log() result is negative or 0 */
668 arg = min_gaussian / GAUSS_NORM;
669 if (arg > 1.0)
671 gmx_fatal(FARGS, "min_gaussian of flexible rotation groups must be <%g", GAUSS_NORM);
674 return std::sqrt(-2.0 * sigma * sigma * log(min_gaussian / GAUSS_NORM));
678 static inline real calc_beta(rvec curr_x, const gmx_enfrotgrp* erg, int n)
680 return iprod(curr_x, erg->vec) - erg->rotg->slab_dist * n;
684 static inline real gaussian_weight(rvec curr_x, const gmx_enfrotgrp* erg, int n)
686 const real norm = GAUSS_NORM;
687 real sigma;
690 /* Define the sigma value */
691 sigma = 0.7 * erg->rotg->slab_dist;
692 /* Calculate the Gaussian value of slab n for position curr_x */
693 return norm * exp(-0.5 * gmx::square(calc_beta(curr_x, erg, n) / sigma));
697 /* Returns the weight in a single slab, also calculates the Gaussian- and mass-
698 * weighted sum of positions for that slab */
699 static real get_slab_weight(int j, const gmx_enfrotgrp* erg, rvec xc[], const real mc[], rvec* x_weighted_sum)
701 rvec curr_x; /* The position of an atom */
702 rvec curr_x_weighted; /* The gaussian-weighted position */
703 real gaussian; /* A single gaussian weight */
704 real wgauss; /* gaussian times current mass */
705 real slabweight = 0.0; /* The sum of weights in the slab */
707 clear_rvec(*x_weighted_sum);
709 /* Loop over all atoms in the rotation group */
710 for (int i = 0; i < erg->rotg->nat; i++)
712 copy_rvec(xc[i], curr_x);
713 gaussian = gaussian_weight(curr_x, erg, j);
714 wgauss = gaussian * mc[i];
715 svmul(wgauss, curr_x, curr_x_weighted);
716 rvec_add(*x_weighted_sum, curr_x_weighted, *x_weighted_sum);
717 slabweight += wgauss;
718 } /* END of loop over rotation group atoms */
720 return slabweight;
724 static void get_slab_centers(gmx_enfrotgrp* erg, /* Enforced rotation group working data */
725 rvec* xc, /* The rotation group positions; will
726 typically be enfrotgrp->xc, but at first call
727 it is enfrotgrp->xc_ref */
728 real* mc, /* The masses of the rotation group atoms */
729 real time, /* Used for output only */
730 FILE* out_slabs, /* For outputting center per slab information */
731 gmx_bool bOutStep, /* Is this an output step? */
732 gmx_bool bReference) /* If this routine is called from
733 init_rot_group we need to store
734 the reference slab centers */
736 /* Loop over slabs */
737 for (int j = erg->slab_first; j <= erg->slab_last; j++)
739 int slabIndex = j - erg->slab_first;
740 erg->slab_weights[slabIndex] = get_slab_weight(j, erg, xc, mc, &erg->slab_center[slabIndex]);
742 /* We can do the calculations ONLY if there is weight in the slab! */
743 if (erg->slab_weights[slabIndex] > WEIGHT_MIN)
745 svmul(1.0 / erg->slab_weights[slabIndex], erg->slab_center[slabIndex],
746 erg->slab_center[slabIndex]);
748 else
750 /* We need to check this here, since we divide through slab_weights
751 * in the flexible low-level routines! */
752 gmx_fatal(FARGS, "Not enough weight in slab %d. Slab center cannot be determined!", j);
755 /* At first time step: save the centers of the reference structure */
756 if (bReference)
758 copy_rvec(erg->slab_center[slabIndex], erg->slab_center_ref[slabIndex]);
760 } /* END of loop over slabs */
762 /* Output on the master */
763 if ((nullptr != out_slabs) && bOutStep)
765 fprintf(out_slabs, "%12.3e%6d", time, erg->groupIndex);
766 for (int j = erg->slab_first; j <= erg->slab_last; j++)
768 int slabIndex = j - erg->slab_first;
769 fprintf(out_slabs, "%6d%12.3e%12.3e%12.3e", j, erg->slab_center[slabIndex][XX],
770 erg->slab_center[slabIndex][YY], erg->slab_center[slabIndex][ZZ]);
772 fprintf(out_slabs, "\n");
777 static void calc_rotmat(const rvec vec,
778 real degangle, /* Angle alpha of rotation at time t in degrees */
779 matrix rotmat) /* Rotation matrix */
781 real radangle; /* Rotation angle in radians */
782 real cosa; /* cosine alpha */
783 real sina; /* sine alpha */
784 real OMcosa; /* 1 - cos(alpha) */
785 real dumxy, dumxz, dumyz; /* save computations */
786 rvec rot_vec; /* Rotate around rot_vec ... */
789 radangle = degangle * M_PI / 180.0;
790 copy_rvec(vec, rot_vec);
792 /* Precompute some variables: */
793 cosa = cos(radangle);
794 sina = sin(radangle);
795 OMcosa = 1.0 - cosa;
796 dumxy = rot_vec[XX] * rot_vec[YY] * OMcosa;
797 dumxz = rot_vec[XX] * rot_vec[ZZ] * OMcosa;
798 dumyz = rot_vec[YY] * rot_vec[ZZ] * OMcosa;
800 /* Construct the rotation matrix for this rotation group: */
801 /* 1st column: */
802 rotmat[XX][XX] = cosa + rot_vec[XX] * rot_vec[XX] * OMcosa;
803 rotmat[YY][XX] = dumxy + rot_vec[ZZ] * sina;
804 rotmat[ZZ][XX] = dumxz - rot_vec[YY] * sina;
805 /* 2nd column: */
806 rotmat[XX][YY] = dumxy - rot_vec[ZZ] * sina;
807 rotmat[YY][YY] = cosa + rot_vec[YY] * rot_vec[YY] * OMcosa;
808 rotmat[ZZ][YY] = dumyz + rot_vec[XX] * sina;
809 /* 3rd column: */
810 rotmat[XX][ZZ] = dumxz + rot_vec[YY] * sina;
811 rotmat[YY][ZZ] = dumyz - rot_vec[XX] * sina;
812 rotmat[ZZ][ZZ] = cosa + rot_vec[ZZ] * rot_vec[ZZ] * OMcosa;
814 #ifdef PRINTMATRIX
815 int iii, jjj;
817 for (iii = 0; iii < 3; iii++)
819 for (jjj = 0; jjj < 3; jjj++)
821 fprintf(stderr, " %10.8f ", rotmat[iii][jjj]);
823 fprintf(stderr, "\n");
825 #endif
829 /* Calculates torque on the rotation axis tau = position x force */
830 static inline real torque(const rvec rotvec, /* rotation vector; MUST be normalized! */
831 rvec force, /* force */
832 rvec x, /* position of atom on which the force acts */
833 rvec pivot) /* pivot point of rotation axis */
835 rvec vectmp, tau;
838 /* Subtract offset */
839 rvec_sub(x, pivot, vectmp);
841 /* position x force */
842 cprod(vectmp, force, tau);
844 /* Return the part of the torque which is parallel to the rotation vector */
845 return iprod(tau, rotvec);
849 /* Right-aligned output of value with standard width */
850 static void print_aligned(FILE* fp, char const* str)
852 fprintf(fp, "%12s", str);
856 /* Right-aligned output of value with standard short width */
857 static void print_aligned_short(FILE* fp, char const* str)
859 fprintf(fp, "%6s", str);
863 static FILE* open_output_file(const char* fn, int steps, const char what[])
865 FILE* fp;
868 fp = gmx_ffopen(fn, "w");
870 fprintf(fp, "# Output of %s is written in intervals of %d time step%s.\n#\n", what, steps,
871 steps > 1 ? "s" : "");
873 return fp;
877 /* Open output file for slab center data. Call on master only */
878 static FILE* open_slab_out(const char* fn, gmx_enfrot* er)
880 FILE* fp;
882 if (er->restartWithAppending)
884 fp = gmx_fio_fopen(fn, "a");
886 else
888 fp = open_output_file(fn, er->nstsout, "gaussian weighted slab centers");
890 for (auto& ergRef : er->enfrotgrp)
892 gmx_enfrotgrp* erg = &ergRef;
893 if (ISFLEX(erg->rotg))
895 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm, %s.\n", erg->groupIndex,
896 erotg_names[erg->rotg->eType], erg->rotg->slab_dist,
897 erg->rotg->bMassW ? "centers of mass" : "geometrical centers");
901 fprintf(fp, "# Reference centers are listed first (t=-1).\n");
902 fprintf(fp, "# The following columns have the syntax:\n");
903 fprintf(fp, "# ");
904 print_aligned_short(fp, "t");
905 print_aligned_short(fp, "grp");
906 /* Print legend for the first two entries only ... */
907 for (int i = 0; i < 2; i++)
909 print_aligned_short(fp, "slab");
910 print_aligned(fp, "X center");
911 print_aligned(fp, "Y center");
912 print_aligned(fp, "Z center");
914 fprintf(fp, " ...\n");
915 fflush(fp);
918 return fp;
922 /* Adds 'buf' to 'str' */
923 static void add_to_string(char** str, char* buf)
925 int len;
928 len = strlen(*str) + strlen(buf) + 1;
929 srenew(*str, len);
930 strcat(*str, buf);
934 static void add_to_string_aligned(char** str, char* buf)
936 char buf_aligned[STRLEN];
938 sprintf(buf_aligned, "%12s", buf);
939 add_to_string(str, buf_aligned);
943 /* Open output file and print some general information about the rotation groups.
944 * Call on master only */
945 static FILE* open_rot_out(const char* fn, const gmx_output_env_t* oenv, gmx_enfrot* er)
947 FILE* fp;
948 int nsets;
949 const char** setname;
950 char buf[50], buf2[75];
951 gmx_bool bFlex;
952 char* LegendStr = nullptr;
953 const t_rot* rot = er->rot;
955 if (er->restartWithAppending)
957 fp = gmx_fio_fopen(fn, "a");
959 else
961 fp = xvgropen(fn, "Rotation angles and energy", "Time (ps)",
962 "angles (degrees) and energies (kJ/mol)", oenv);
963 fprintf(fp,
964 "# Output of enforced rotation data is written in intervals of %d time "
965 "step%s.\n#\n",
966 er->nstrout, er->nstrout > 1 ? "s" : "");
967 fprintf(fp,
968 "# The scalar tau is the torque (kJ/mol) in the direction of the rotation vector "
969 "v.\n");
970 fprintf(fp, "# To obtain the vectorial torque, multiply tau with the group's rot-vec.\n");
971 fprintf(fp,
972 "# For flexible groups, tau(t,n) from all slabs n have been summed in a single "
973 "value tau(t) here.\n");
974 fprintf(fp, "# The torques tau(t,n) are found in the rottorque.log (-rt) output file\n");
976 for (int g = 0; g < rot->ngrp; g++)
978 const t_rotgrp* rotg = &rot->grp[g];
979 const gmx_enfrotgrp* erg = &er->enfrotgrp[g];
980 bFlex = ISFLEX(rotg);
982 fprintf(fp, "#\n");
983 fprintf(fp, "# ROTATION GROUP %d, potential type '%s':\n", g, erotg_names[rotg->eType]);
984 fprintf(fp, "# rot-massw%d %s\n", g, yesno_names[rotg->bMassW]);
985 fprintf(fp, "# rot-vec%d %12.5e %12.5e %12.5e\n", g, erg->vec[XX],
986 erg->vec[YY], erg->vec[ZZ]);
987 fprintf(fp, "# rot-rate%d %12.5e degrees/ps\n", g, rotg->rate);
988 fprintf(fp, "# rot-k%d %12.5e kJ/(mol*nm^2)\n", g, rotg->k);
989 if (rotg->eType == erotgISO || rotg->eType == erotgPM || rotg->eType == erotgRM
990 || rotg->eType == erotgRM2)
992 fprintf(fp, "# rot-pivot%d %12.5e %12.5e %12.5e nm\n", g, rotg->pivot[XX],
993 rotg->pivot[YY], rotg->pivot[ZZ]);
996 if (bFlex)
998 fprintf(fp, "# rot-slab-distance%d %f nm\n", g, rotg->slab_dist);
999 fprintf(fp, "# rot-min-gaussian%d %12.5e\n", g, rotg->min_gaussian);
1002 /* Output the centers of the rotation groups for the pivot-free potentials */
1003 if ((rotg->eType == erotgISOPF) || (rotg->eType == erotgPMPF) || (rotg->eType == erotgRMPF)
1004 || (rotg->eType == erotgRM2PF || (rotg->eType == erotgFLEXT) || (rotg->eType == erotgFLEX2T)))
1006 fprintf(fp, "# ref. grp. %d center %12.5e %12.5e %12.5e\n", g,
1007 erg->xc_ref_center[XX], erg->xc_ref_center[YY], erg->xc_ref_center[ZZ]);
1009 fprintf(fp, "# grp. %d init.center %12.5e %12.5e %12.5e\n", g, erg->xc_center[XX],
1010 erg->xc_center[YY], erg->xc_center[ZZ]);
1013 if ((rotg->eType == erotgRM2) || (rotg->eType == erotgFLEX2) || (rotg->eType == erotgFLEX2T))
1015 fprintf(fp, "# rot-eps%d %12.5e nm^2\n", g, rotg->eps);
1017 if (erotgFitPOT == rotg->eFittype)
1019 fprintf(fp, "#\n");
1020 fprintf(fp,
1021 "# theta_fit%d is determined by first evaluating the potential for %d "
1022 "angles around theta_ref%d.\n",
1023 g, rotg->PotAngle_nstep, g);
1024 fprintf(fp,
1025 "# The fit angle is the one with the smallest potential. It is given as "
1026 "the deviation\n");
1027 fprintf(fp,
1028 "# from the reference angle, i.e. if theta_ref=X and theta_fit=Y, then the "
1029 "angle with\n");
1030 fprintf(fp,
1031 "# minimal value of the potential is X+Y. Angular resolution is %g "
1032 "degrees.\n",
1033 rotg->PotAngle_step);
1037 /* Print a nice legend */
1038 snew(LegendStr, 1);
1039 LegendStr[0] = '\0';
1040 sprintf(buf, "# %6s", "time");
1041 add_to_string_aligned(&LegendStr, buf);
1043 nsets = 0;
1044 snew(setname, 4 * rot->ngrp);
1046 for (int g = 0; g < rot->ngrp; g++)
1048 sprintf(buf, "theta_ref%d", g);
1049 add_to_string_aligned(&LegendStr, buf);
1051 sprintf(buf2, "%s (degrees)", buf);
1052 setname[nsets] = gmx_strdup(buf2);
1053 nsets++;
1055 for (int g = 0; g < rot->ngrp; g++)
1057 const t_rotgrp* rotg = &rot->grp[g];
1058 bFlex = ISFLEX(rotg);
1060 /* For flexible axis rotation we use RMSD fitting to determine the
1061 * actual angle of the rotation group */
1062 if (bFlex || erotgFitPOT == rotg->eFittype)
1064 sprintf(buf, "theta_fit%d", g);
1066 else
1068 sprintf(buf, "theta_av%d", g);
1070 add_to_string_aligned(&LegendStr, buf);
1071 sprintf(buf2, "%s (degrees)", buf);
1072 setname[nsets] = gmx_strdup(buf2);
1073 nsets++;
1075 sprintf(buf, "tau%d", g);
1076 add_to_string_aligned(&LegendStr, buf);
1077 sprintf(buf2, "%s (kJ/mol)", buf);
1078 setname[nsets] = gmx_strdup(buf2);
1079 nsets++;
1081 sprintf(buf, "energy%d", g);
1082 add_to_string_aligned(&LegendStr, buf);
1083 sprintf(buf2, "%s (kJ/mol)", buf);
1084 setname[nsets] = gmx_strdup(buf2);
1085 nsets++;
1087 fprintf(fp, "#\n");
1089 if (nsets > 1)
1091 xvgr_legend(fp, nsets, setname, oenv);
1093 sfree(setname);
1095 fprintf(fp, "#\n# Legend for the following data columns:\n");
1096 fprintf(fp, "%s\n", LegendStr);
1097 sfree(LegendStr);
1099 fflush(fp);
1102 return fp;
1106 /* Call on master only */
1107 static FILE* open_angles_out(const char* fn, gmx_enfrot* er)
1109 FILE* fp;
1110 char buf[100];
1111 const t_rot* rot = er->rot;
1113 if (er->restartWithAppending)
1115 fp = gmx_fio_fopen(fn, "a");
1117 else
1119 /* Open output file and write some information about it's structure: */
1120 fp = open_output_file(fn, er->nstsout, "rotation group angles");
1121 fprintf(fp, "# All angles given in degrees, time in ps.\n");
1122 for (int g = 0; g < rot->ngrp; g++)
1124 const t_rotgrp* rotg = &rot->grp[g];
1125 const gmx_enfrotgrp* erg = &er->enfrotgrp[g];
1127 /* Output for this group happens only if potential type is flexible or
1128 * if fit type is potential! */
1129 if (ISFLEX(rotg) || (erotgFitPOT == rotg->eFittype))
1131 if (ISFLEX(rotg))
1133 sprintf(buf, " slab distance %f nm, ", rotg->slab_dist);
1135 else
1137 buf[0] = '\0';
1140 fprintf(fp, "#\n# ROTATION GROUP %d '%s',%s fit type '%s'.\n", g,
1141 erotg_names[rotg->eType], buf, erotg_fitnames[rotg->eFittype]);
1143 /* Special type of fitting using the potential minimum. This is
1144 * done for the whole group only, not for the individual slabs. */
1145 if (erotgFitPOT == rotg->eFittype)
1147 fprintf(fp,
1148 "# To obtain theta_fit%d, the potential is evaluated for %d angles "
1149 "around theta_ref%d\n",
1150 g, rotg->PotAngle_nstep, g);
1151 fprintf(fp,
1152 "# The fit angle in the rotation standard outfile is the one with "
1153 "minimal energy E(theta_fit) [kJ/mol].\n");
1154 fprintf(fp, "#\n");
1157 fprintf(fp, "# Legend for the group %d data columns:\n", g);
1158 fprintf(fp, "# ");
1159 print_aligned_short(fp, "time");
1160 print_aligned_short(fp, "grp");
1161 print_aligned(fp, "theta_ref");
1163 if (erotgFitPOT == rotg->eFittype)
1165 /* Output the set of angles around the reference angle */
1166 for (int i = 0; i < rotg->PotAngle_nstep; i++)
1168 sprintf(buf, "E(%g)", erg->PotAngleFit->degangle[i]);
1169 print_aligned(fp, buf);
1172 else
1174 /* Output fit angle for each slab */
1175 print_aligned_short(fp, "slab");
1176 print_aligned_short(fp, "atoms");
1177 print_aligned(fp, "theta_fit");
1178 print_aligned_short(fp, "slab");
1179 print_aligned_short(fp, "atoms");
1180 print_aligned(fp, "theta_fit");
1181 fprintf(fp, " ...");
1183 fprintf(fp, "\n");
1186 fflush(fp);
1189 return fp;
1193 /* Open torque output file and write some information about it's structure.
1194 * Call on master only */
1195 static FILE* open_torque_out(const char* fn, gmx_enfrot* er)
1197 FILE* fp;
1198 const t_rot* rot = er->rot;
1200 if (er->restartWithAppending)
1202 fp = gmx_fio_fopen(fn, "a");
1204 else
1206 fp = open_output_file(fn, er->nstsout, "torques");
1208 for (int g = 0; g < rot->ngrp; g++)
1210 const t_rotgrp* rotg = &rot->grp[g];
1211 const gmx_enfrotgrp* erg = &er->enfrotgrp[g];
1212 if (ISFLEX(rotg))
1214 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm.\n", g,
1215 erotg_names[rotg->eType], rotg->slab_dist);
1216 fprintf(fp,
1217 "# The scalar tau is the torque (kJ/mol) in the direction of the rotation "
1218 "vector.\n");
1219 fprintf(fp, "# To obtain the vectorial torque, multiply tau with\n");
1220 fprintf(fp, "# rot-vec%d %10.3e %10.3e %10.3e\n", g, erg->vec[XX],
1221 erg->vec[YY], erg->vec[ZZ]);
1222 fprintf(fp, "#\n");
1225 fprintf(fp, "# Legend for the following data columns: (tau=torque for that slab):\n");
1226 fprintf(fp, "# ");
1227 print_aligned_short(fp, "t");
1228 print_aligned_short(fp, "grp");
1229 print_aligned_short(fp, "slab");
1230 print_aligned(fp, "tau");
1231 print_aligned_short(fp, "slab");
1232 print_aligned(fp, "tau");
1233 fprintf(fp, " ...\n");
1234 fflush(fp);
1237 return fp;
1241 static void swap_val(double* vec, int i, int j)
1243 double tmp = vec[j];
1246 vec[j] = vec[i];
1247 vec[i] = tmp;
1251 static void swap_col(double** mat, int i, int j)
1253 double tmp[3] = { mat[0][j], mat[1][j], mat[2][j] };
1256 mat[0][j] = mat[0][i];
1257 mat[1][j] = mat[1][i];
1258 mat[2][j] = mat[2][i];
1260 mat[0][i] = tmp[0];
1261 mat[1][i] = tmp[1];
1262 mat[2][i] = tmp[2];
1266 /* Eigenvectors are stored in columns of eigen_vec */
1267 static void diagonalize_symmetric(double** matrix, double** eigen_vec, double eigenval[3])
1269 int n_rot;
1272 jacobi(matrix, 3, eigenval, eigen_vec, &n_rot);
1274 /* sort in ascending order */
1275 if (eigenval[0] > eigenval[1])
1277 swap_val(eigenval, 0, 1);
1278 swap_col(eigen_vec, 0, 1);
1280 if (eigenval[1] > eigenval[2])
1282 swap_val(eigenval, 1, 2);
1283 swap_col(eigen_vec, 1, 2);
1285 if (eigenval[0] > eigenval[1])
1287 swap_val(eigenval, 0, 1);
1288 swap_col(eigen_vec, 0, 1);
1293 static void align_with_z(rvec* s, /* Structure to align */
1294 int natoms,
1295 rvec axis)
1297 int i, j, k;
1298 rvec zet = { 0.0, 0.0, 1.0 };
1299 rvec rot_axis = { 0.0, 0.0, 0.0 };
1300 rvec* rotated_str = nullptr;
1301 real ooanorm;
1302 real angle;
1303 matrix rotmat;
1306 snew(rotated_str, natoms);
1308 /* Normalize the axis */
1309 ooanorm = 1.0 / norm(axis);
1310 svmul(ooanorm, axis, axis);
1312 /* Calculate the angle for the fitting procedure */
1313 cprod(axis, zet, rot_axis);
1314 angle = acos(axis[2]);
1315 if (angle < 0.0)
1317 angle += M_PI;
1320 /* Calculate the rotation matrix */
1321 calc_rotmat(rot_axis, angle * 180.0 / M_PI, rotmat);
1323 /* Apply the rotation matrix to s */
1324 for (i = 0; i < natoms; i++)
1326 for (j = 0; j < 3; j++)
1328 for (k = 0; k < 3; k++)
1330 rotated_str[i][j] += rotmat[j][k] * s[i][k];
1335 /* Rewrite the rotated structure to s */
1336 for (i = 0; i < natoms; i++)
1338 for (j = 0; j < 3; j++)
1340 s[i][j] = rotated_str[i][j];
1344 sfree(rotated_str);
1348 static void calc_correl_matrix(rvec* Xstr, rvec* Ystr, double** Rmat, int natoms)
1350 int i, j, k;
1353 for (i = 0; i < 3; i++)
1355 for (j = 0; j < 3; j++)
1357 Rmat[i][j] = 0.0;
1361 for (i = 0; i < 3; i++)
1363 for (j = 0; j < 3; j++)
1365 for (k = 0; k < natoms; k++)
1367 Rmat[i][j] += Ystr[k][i] * Xstr[k][j];
1374 static void weigh_coords(rvec* str, real* weight, int natoms)
1376 int i, j;
1379 for (i = 0; i < natoms; i++)
1381 for (j = 0; j < 3; j++)
1383 str[i][j] *= std::sqrt(weight[i]);
1389 static real opt_angle_analytic(rvec* ref_s,
1390 rvec* act_s,
1391 real* weight,
1392 int natoms,
1393 const rvec ref_com,
1394 const rvec act_com,
1395 rvec axis)
1397 int i, j, k;
1398 rvec* ref_s_1 = nullptr;
1399 rvec* act_s_1 = nullptr;
1400 rvec shift;
1401 double **Rmat, **RtR, **eigvec;
1402 double eigval[3];
1403 double V[3][3], WS[3][3];
1404 double rot_matrix[3][3];
1405 double opt_angle;
1408 /* Do not change the original coordinates */
1409 snew(ref_s_1, natoms);
1410 snew(act_s_1, natoms);
1411 for (i = 0; i < natoms; i++)
1413 copy_rvec(ref_s[i], ref_s_1[i]);
1414 copy_rvec(act_s[i], act_s_1[i]);
1417 /* Translate the structures to the origin */
1418 shift[XX] = -ref_com[XX];
1419 shift[YY] = -ref_com[YY];
1420 shift[ZZ] = -ref_com[ZZ];
1421 translate_x(ref_s_1, natoms, shift);
1423 shift[XX] = -act_com[XX];
1424 shift[YY] = -act_com[YY];
1425 shift[ZZ] = -act_com[ZZ];
1426 translate_x(act_s_1, natoms, shift);
1428 /* Align rotation axis with z */
1429 align_with_z(ref_s_1, natoms, axis);
1430 align_with_z(act_s_1, natoms, axis);
1432 /* Correlation matrix */
1433 Rmat = allocate_square_matrix(3);
1435 for (i = 0; i < natoms; i++)
1437 ref_s_1[i][2] = 0.0;
1438 act_s_1[i][2] = 0.0;
1441 /* Weight positions with sqrt(weight) */
1442 if (nullptr != weight)
1444 weigh_coords(ref_s_1, weight, natoms);
1445 weigh_coords(act_s_1, weight, natoms);
1448 /* Calculate correlation matrices R=YXt (X=ref_s; Y=act_s) */
1449 calc_correl_matrix(ref_s_1, act_s_1, Rmat, natoms);
1451 /* Calculate RtR */
1452 RtR = allocate_square_matrix(3);
1453 for (i = 0; i < 3; i++)
1455 for (j = 0; j < 3; j++)
1457 for (k = 0; k < 3; k++)
1459 RtR[i][j] += Rmat[k][i] * Rmat[k][j];
1463 /* Diagonalize RtR */
1464 snew(eigvec, 3);
1465 for (i = 0; i < 3; i++)
1467 snew(eigvec[i], 3);
1470 diagonalize_symmetric(RtR, eigvec, eigval);
1471 swap_col(eigvec, 0, 1);
1472 swap_col(eigvec, 1, 2);
1473 swap_val(eigval, 0, 1);
1474 swap_val(eigval, 1, 2);
1476 /* Calculate V */
1477 for (i = 0; i < 3; i++)
1479 for (j = 0; j < 3; j++)
1481 V[i][j] = 0.0;
1482 WS[i][j] = 0.0;
1486 for (i = 0; i < 2; i++)
1488 for (j = 0; j < 2; j++)
1490 WS[i][j] = eigvec[i][j] / std::sqrt(eigval[j]);
1494 for (i = 0; i < 3; i++)
1496 for (j = 0; j < 3; j++)
1498 for (k = 0; k < 3; k++)
1500 V[i][j] += Rmat[i][k] * WS[k][j];
1504 free_square_matrix(Rmat, 3);
1506 /* Calculate optimal rotation matrix */
1507 for (i = 0; i < 3; i++)
1509 for (j = 0; j < 3; j++)
1511 rot_matrix[i][j] = 0.0;
1515 for (i = 0; i < 3; i++)
1517 for (j = 0; j < 3; j++)
1519 for (k = 0; k < 3; k++)
1521 rot_matrix[i][j] += eigvec[i][k] * V[j][k];
1525 rot_matrix[2][2] = 1.0;
1527 /* In some cases abs(rot_matrix[0][0]) can be slighly larger
1528 * than unity due to numerical inacurracies. To be able to calculate
1529 * the acos function, we put these values back in range. */
1530 if (rot_matrix[0][0] > 1.0)
1532 rot_matrix[0][0] = 1.0;
1534 else if (rot_matrix[0][0] < -1.0)
1536 rot_matrix[0][0] = -1.0;
1539 /* Determine the optimal rotation angle: */
1540 opt_angle = (-1.0) * acos(rot_matrix[0][0]) * 180.0 / M_PI;
1541 if (rot_matrix[0][1] < 0.0)
1543 opt_angle = (-1.0) * opt_angle;
1546 /* Give back some memory */
1547 free_square_matrix(RtR, 3);
1548 sfree(ref_s_1);
1549 sfree(act_s_1);
1550 for (i = 0; i < 3; i++)
1552 sfree(eigvec[i]);
1554 sfree(eigvec);
1556 return static_cast<real>(opt_angle);
1560 /* Determine angle of the group by RMSD fit to the reference */
1561 /* Not parallelized, call this routine only on the master */
1562 static real flex_fit_angle(gmx_enfrotgrp* erg)
1564 rvec* fitcoords = nullptr;
1565 rvec center; /* Center of positions passed to the fit routine */
1566 real fitangle; /* Angle of the rotation group derived by fitting */
1567 rvec coord;
1568 real scal;
1570 /* Get the center of the rotation group.
1571 * Note, again, erg->xc has been sorted in do_flexible */
1572 get_center(erg->xc, erg->mc_sorted, erg->rotg->nat, center);
1574 /* === Determine the optimal fit angle for the rotation group === */
1575 if (erg->rotg->eFittype == erotgFitNORM)
1577 /* Normalize every position to it's reference length */
1578 for (int i = 0; i < erg->rotg->nat; i++)
1580 /* Put the center of the positions into the origin */
1581 rvec_sub(erg->xc[i], center, coord);
1582 /* Determine the scaling factor for the length: */
1583 scal = erg->xc_ref_length[erg->xc_sortind[i]] / norm(coord);
1584 /* Get position, multiply with the scaling factor and save */
1585 svmul(scal, coord, erg->xc_norm[i]);
1587 fitcoords = erg->xc_norm;
1589 else
1591 fitcoords = erg->xc;
1593 /* From the point of view of the current positions, the reference has rotated
1594 * backwards. Since we output the angle relative to the fixed reference,
1595 * we need the minus sign. */
1596 fitangle = -opt_angle_analytic(erg->xc_ref_sorted, fitcoords, erg->mc_sorted, erg->rotg->nat,
1597 erg->xc_ref_center, center, erg->vec);
1599 return fitangle;
1603 /* Determine actual angle of each slab by RMSD fit to the reference */
1604 /* Not parallelized, call this routine only on the master */
1605 static void flex_fit_angle_perslab(gmx_enfrotgrp* erg, double t, real degangle, FILE* fp)
1607 rvec curr_x, ref_x;
1608 rvec act_center; /* Center of actual positions that are passed to the fit routine */
1609 rvec ref_center; /* Same for the reference positions */
1610 real fitangle; /* Angle of a slab derived from an RMSD fit to
1611 * the reference structure at t=0 */
1612 gmx_slabdata* sd;
1613 real OOm_av; /* 1/average_mass of a rotation group atom */
1614 real m_rel; /* Relative mass of a rotation group atom */
1617 /* Average mass of a rotation group atom: */
1618 OOm_av = erg->invmass * erg->rotg->nat;
1620 /**********************************/
1621 /* First collect the data we need */
1622 /**********************************/
1624 /* Collect the data for the individual slabs */
1625 for (int n = erg->slab_first; n <= erg->slab_last; n++)
1627 int slabIndex = n - erg->slab_first; /* slab index */
1628 sd = &(erg->slab_data[slabIndex]);
1629 sd->nat = erg->lastatom[slabIndex] - erg->firstatom[slabIndex] + 1;
1630 int ind = 0;
1632 /* Loop over the relevant atoms in the slab */
1633 for (int l = erg->firstatom[slabIndex]; l <= erg->lastatom[slabIndex]; l++)
1635 /* Current position of this atom: x[ii][XX/YY/ZZ] */
1636 copy_rvec(erg->xc[l], curr_x);
1638 /* The (unrotated) reference position of this atom is copied to ref_x.
1639 * Beware, the xc coords have been sorted in do_flexible */
1640 copy_rvec(erg->xc_ref_sorted[l], ref_x);
1642 /* Save data for doing angular RMSD fit later */
1643 /* Save the current atom position */
1644 copy_rvec(curr_x, sd->x[ind]);
1645 /* Save the corresponding reference position */
1646 copy_rvec(ref_x, sd->ref[ind]);
1648 /* Maybe also mass-weighting was requested. If yes, additionally
1649 * multiply the weights with the relative mass of the atom. If not,
1650 * multiply with unity. */
1651 m_rel = erg->mc_sorted[l] * OOm_av;
1653 /* Save the weight for this atom in this slab */
1654 sd->weight[ind] = gaussian_weight(curr_x, erg, n) * m_rel;
1656 /* Next atom in this slab */
1657 ind++;
1661 /******************************/
1662 /* Now do the fit calculation */
1663 /******************************/
1665 fprintf(fp, "%12.3e%6d%12.3f", t, erg->groupIndex, degangle);
1667 /* === Now do RMSD fitting for each slab === */
1668 /* We require at least SLAB_MIN_ATOMS in a slab, such that the fit makes sense. */
1669 #define SLAB_MIN_ATOMS 4
1671 for (int n = erg->slab_first; n <= erg->slab_last; n++)
1673 int slabIndex = n - erg->slab_first; /* slab index */
1674 sd = &(erg->slab_data[slabIndex]);
1675 if (sd->nat >= SLAB_MIN_ATOMS)
1677 /* Get the center of the slabs reference and current positions */
1678 get_center(sd->ref, sd->weight, sd->nat, ref_center);
1679 get_center(sd->x, sd->weight, sd->nat, act_center);
1680 if (erg->rotg->eFittype == erotgFitNORM)
1682 /* Normalize every position to it's reference length
1683 * prior to performing the fit */
1684 for (int i = 0; i < sd->nat; i++) /* Center */
1686 rvec_dec(sd->ref[i], ref_center);
1687 rvec_dec(sd->x[i], act_center);
1688 /* Normalize x_i such that it gets the same length as ref_i */
1689 svmul(norm(sd->ref[i]) / norm(sd->x[i]), sd->x[i], sd->x[i]);
1691 /* We already subtracted the centers */
1692 clear_rvec(ref_center);
1693 clear_rvec(act_center);
1695 fitangle = -opt_angle_analytic(sd->ref, sd->x, sd->weight, sd->nat, ref_center,
1696 act_center, erg->vec);
1697 fprintf(fp, "%6d%6d%12.3f", n, sd->nat, fitangle);
1700 fprintf(fp, "\n");
1702 #undef SLAB_MIN_ATOMS
1706 /* Shift x with is */
1707 static inline void shift_single_coord(const matrix box, rvec x, const ivec is)
1709 int tx, ty, tz;
1712 tx = is[XX];
1713 ty = is[YY];
1714 tz = is[ZZ];
1716 if (TRICLINIC(box))
1718 x[XX] += tx * box[XX][XX] + ty * box[YY][XX] + tz * box[ZZ][XX];
1719 x[YY] += ty * box[YY][YY] + tz * box[ZZ][YY];
1720 x[ZZ] += tz * box[ZZ][ZZ];
1722 else
1724 x[XX] += tx * box[XX][XX];
1725 x[YY] += ty * box[YY][YY];
1726 x[ZZ] += tz * box[ZZ][ZZ];
1731 /* Determine the 'home' slab of this atom which is the
1732 * slab with the highest Gaussian weight of all */
1733 static inline int get_homeslab(rvec curr_x, /* The position for which the home slab shall be determined */
1734 const rvec rotvec, /* The rotation vector */
1735 real slabdist) /* The slab distance */
1737 real dist;
1740 /* The distance of the atom to the coordinate center (where the
1741 * slab with index 0) is */
1742 dist = iprod(rotvec, curr_x);
1744 return gmx::roundToInt(dist / slabdist);
1748 /* For a local atom determine the relevant slabs, i.e. slabs in
1749 * which the gaussian is larger than min_gaussian
1751 static int get_single_atom_gaussians(rvec curr_x, gmx_enfrotgrp* erg)
1754 /* Determine the 'home' slab of this atom: */
1755 int homeslab = get_homeslab(curr_x, erg->vec, erg->rotg->slab_dist);
1757 /* First determine the weight in the atoms home slab: */
1758 real g = gaussian_weight(curr_x, erg, homeslab);
1759 int count = 0;
1760 erg->gn_atom[count] = g;
1761 erg->gn_slabind[count] = homeslab;
1762 count++;
1765 /* Determine the max slab */
1766 int slab = homeslab;
1767 while (g > erg->rotg->min_gaussian)
1769 slab++;
1770 g = gaussian_weight(curr_x, erg, slab);
1771 erg->gn_slabind[count] = slab;
1772 erg->gn_atom[count] = g;
1773 count++;
1775 count--;
1777 /* Determine the min slab */
1778 slab = homeslab;
1781 slab--;
1782 g = gaussian_weight(curr_x, erg, slab);
1783 erg->gn_slabind[count] = slab;
1784 erg->gn_atom[count] = g;
1785 count++;
1786 } while (g > erg->rotg->min_gaussian);
1787 count--;
1789 return count;
1793 static void flex2_precalc_inner_sum(const gmx_enfrotgrp* erg)
1795 rvec xi; /* positions in the i-sum */
1796 rvec xcn, ycn; /* the current and the reference slab centers */
1797 real gaussian_xi;
1798 rvec yi0;
1799 rvec rin; /* Helper variables */
1800 real fac, fac2;
1801 rvec innersumvec;
1802 real OOpsii, OOpsiistar;
1803 real sin_rin; /* s_ii.r_ii */
1804 rvec s_in, tmpvec, tmpvec2;
1805 real mi, wi; /* Mass-weighting of the positions */
1806 real N_M; /* N/M */
1809 N_M = erg->rotg->nat * erg->invmass;
1811 /* Loop over all slabs that contain something */
1812 for (int n = erg->slab_first; n <= erg->slab_last; n++)
1814 int slabIndex = n - erg->slab_first; /* slab index */
1816 /* The current center of this slab is saved in xcn: */
1817 copy_rvec(erg->slab_center[slabIndex], xcn);
1818 /* ... and the reference center in ycn: */
1819 copy_rvec(erg->slab_center_ref[slabIndex + erg->slab_buffer], ycn);
1821 /*** D. Calculate the whole inner sum used for second and third sum */
1822 /* For slab n, we need to loop over all atoms i again. Since we sorted
1823 * the atoms with respect to the rotation vector, we know that it is sufficient
1824 * to calculate from firstatom to lastatom only. All other contributions will
1825 * be very small. */
1826 clear_rvec(innersumvec);
1827 for (int i = erg->firstatom[slabIndex]; i <= erg->lastatom[slabIndex]; i++)
1829 /* Coordinate xi of this atom */
1830 copy_rvec(erg->xc[i], xi);
1832 /* The i-weights */
1833 gaussian_xi = gaussian_weight(xi, erg, n);
1834 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1835 wi = N_M * mi;
1837 /* Calculate rin */
1838 copy_rvec(erg->xc_ref_sorted[i], yi0); /* Reference position yi0 */
1839 rvec_sub(yi0, ycn, tmpvec2); /* tmpvec2 = yi0 - ycn */
1840 mvmul(erg->rotmat, tmpvec2, rin); /* rin = Omega.(yi0 - ycn) */
1842 /* Calculate psi_i* and sin */
1843 rvec_sub(xi, xcn, tmpvec2); /* tmpvec2 = xi - xcn */
1845 /* In rare cases, when an atom position coincides with a slab center
1846 * (tmpvec2 == 0) we cannot compute the vector product for s_in.
1847 * However, since the atom is located directly on the pivot, this
1848 * slab's contribution to the force on that atom will be zero
1849 * anyway. Therefore, we continue with the next atom. */
1850 if (gmx_numzero(norm(tmpvec2))) /* 0 == norm(xi - xcn) */
1852 continue;
1855 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xi - xcn) */
1856 OOpsiistar = norm2(tmpvec) + erg->rotg->eps; /* OOpsii* = 1/psii* = |v x (xi-xcn)|^2 + eps */
1857 OOpsii = norm(tmpvec); /* OOpsii = 1 / psii = |v x (xi - xcn)| */
1859 /* * v x (xi - xcn) */
1860 unitv(tmpvec, s_in); /* sin = ---------------- */
1861 /* |v x (xi - xcn)| */
1863 sin_rin = iprod(s_in, rin); /* sin_rin = sin . rin */
1865 /* Now the whole sum */
1866 fac = OOpsii / OOpsiistar;
1867 svmul(fac, rin, tmpvec);
1868 fac2 = fac * fac * OOpsii;
1869 svmul(fac2 * sin_rin, s_in, tmpvec2);
1870 rvec_dec(tmpvec, tmpvec2);
1872 svmul(wi * gaussian_xi * sin_rin, tmpvec, tmpvec2);
1874 rvec_inc(innersumvec, tmpvec2);
1875 } /* now we have the inner sum, used both for sum2 and sum3 */
1877 /* Save it to be used in do_flex2_lowlevel */
1878 copy_rvec(innersumvec, erg->slab_innersumvec[slabIndex]);
1879 } /* END of loop over slabs */
1883 static void flex_precalc_inner_sum(const gmx_enfrotgrp* erg)
1885 rvec xi; /* position */
1886 rvec xcn, ycn; /* the current and the reference slab centers */
1887 rvec qin, rin; /* q_i^n and r_i^n */
1888 real bin;
1889 rvec tmpvec;
1890 rvec innersumvec; /* Inner part of sum_n2 */
1891 real gaussian_xi; /* Gaussian weight gn(xi) */
1892 real mi, wi; /* Mass-weighting of the positions */
1893 real N_M; /* N/M */
1895 N_M = erg->rotg->nat * erg->invmass;
1897 /* Loop over all slabs that contain something */
1898 for (int n = erg->slab_first; n <= erg->slab_last; n++)
1900 int slabIndex = n - erg->slab_first; /* slab index */
1902 /* The current center of this slab is saved in xcn: */
1903 copy_rvec(erg->slab_center[slabIndex], xcn);
1904 /* ... and the reference center in ycn: */
1905 copy_rvec(erg->slab_center_ref[slabIndex + erg->slab_buffer], ycn);
1907 /* For slab n, we need to loop over all atoms i again. Since we sorted
1908 * the atoms with respect to the rotation vector, we know that it is sufficient
1909 * to calculate from firstatom to lastatom only. All other contributions will
1910 * be very small. */
1911 clear_rvec(innersumvec);
1912 for (int i = erg->firstatom[slabIndex]; i <= erg->lastatom[slabIndex]; i++)
1914 /* Coordinate xi of this atom */
1915 copy_rvec(erg->xc[i], xi);
1917 /* The i-weights */
1918 gaussian_xi = gaussian_weight(xi, erg, n);
1919 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1920 wi = N_M * mi;
1922 /* Calculate rin and qin */
1923 rvec_sub(erg->xc_ref_sorted[i], ycn, tmpvec); /* tmpvec = yi0-ycn */
1925 /* In rare cases, when an atom position coincides with a slab center
1926 * (tmpvec == 0) we cannot compute the vector product for qin.
1927 * However, since the atom is located directly on the pivot, this
1928 * slab's contribution to the force on that atom will be zero
1929 * anyway. Therefore, we continue with the next atom. */
1930 if (gmx_numzero(norm(tmpvec))) /* 0 == norm(yi0 - ycn) */
1932 continue;
1935 mvmul(erg->rotmat, tmpvec, rin); /* rin = Omega.(yi0 - ycn) */
1936 cprod(erg->vec, rin, tmpvec); /* tmpvec = v x Omega*(yi0-ycn) */
1938 /* * v x Omega*(yi0-ycn) */
1939 unitv(tmpvec, qin); /* qin = --------------------- */
1940 /* |v x Omega*(yi0-ycn)| */
1942 /* Calculate bin */
1943 rvec_sub(xi, xcn, tmpvec); /* tmpvec = xi-xcn */
1944 bin = iprod(qin, tmpvec); /* bin = qin*(xi-xcn) */
1946 svmul(wi * gaussian_xi * bin, qin, tmpvec);
1948 /* Add this contribution to the inner sum: */
1949 rvec_add(innersumvec, tmpvec, innersumvec);
1950 } /* now we have the inner sum vector S^n for this slab */
1951 /* Save it to be used in do_flex_lowlevel */
1952 copy_rvec(innersumvec, erg->slab_innersumvec[slabIndex]);
1957 static real do_flex2_lowlevel(gmx_enfrotgrp* erg,
1958 real sigma, /* The Gaussian width sigma */
1959 rvec x[],
1960 gmx_bool bOutstepRot,
1961 gmx_bool bOutstepSlab,
1962 const matrix box)
1964 int count, ii, iigrp;
1965 rvec xj; /* position in the i-sum */
1966 rvec yj0; /* the reference position in the j-sum */
1967 rvec xcn, ycn; /* the current and the reference slab centers */
1968 real V; /* This node's part of the rotation pot. energy */
1969 real gaussian_xj; /* Gaussian weight */
1970 real beta;
1972 real numerator, fit_numerator;
1973 rvec rjn, fit_rjn; /* Helper variables */
1974 real fac, fac2;
1976 real OOpsij, OOpsijstar;
1977 real OOsigma2; /* 1/(sigma^2) */
1978 real sjn_rjn;
1979 real betasigpsi;
1980 rvec sjn, tmpvec, tmpvec2, yj0_ycn;
1981 rvec sum1vec_part, sum1vec, sum2vec_part, sum2vec, sum3vec, sum4vec, innersumvec;
1982 real sum3, sum4;
1983 real mj, wj; /* Mass-weighting of the positions */
1984 real N_M; /* N/M */
1985 real Wjn; /* g_n(x_j) m_j / Mjn */
1986 gmx_bool bCalcPotFit;
1988 /* To calculate the torque per slab */
1989 rvec slab_force; /* Single force from slab n on one atom */
1990 rvec slab_sum1vec_part;
1991 real slab_sum3part, slab_sum4part;
1992 rvec slab_sum1vec, slab_sum2vec, slab_sum3vec, slab_sum4vec;
1994 /* Pre-calculate the inner sums, so that we do not have to calculate
1995 * them again for every atom */
1996 flex2_precalc_inner_sum(erg);
1998 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
2000 /********************************************************/
2001 /* Main loop over all local atoms of the rotation group */
2002 /********************************************************/
2003 N_M = erg->rotg->nat * erg->invmass;
2004 V = 0.0;
2005 OOsigma2 = 1.0 / (sigma * sigma);
2006 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
2007 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2009 for (gmx::index j = 0; j < localRotationGroupIndex.ssize(); j++)
2011 /* Local index of a rotation group atom */
2012 ii = localRotationGroupIndex[j];
2013 /* Position of this atom in the collective array */
2014 iigrp = collectiveRotationGroupIndex[j];
2015 /* Mass-weighting */
2016 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2017 wj = N_M * mj;
2019 /* Current position of this atom: x[ii][XX/YY/ZZ]
2020 * Note that erg->xc_center contains the center of mass in case the flex2-t
2021 * potential was chosen. For the flex2 potential erg->xc_center must be
2022 * zero. */
2023 rvec_sub(x[ii], erg->xc_center, xj);
2025 /* Shift this atom such that it is near its reference */
2026 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2028 /* Determine the slabs to loop over, i.e. the ones with contributions
2029 * larger than min_gaussian */
2030 count = get_single_atom_gaussians(xj, erg);
2032 clear_rvec(sum1vec_part);
2033 clear_rvec(sum2vec_part);
2034 sum3 = 0.0;
2035 sum4 = 0.0;
2036 /* Loop over the relevant slabs for this atom */
2037 for (int ic = 0; ic < count; ic++)
2039 int n = erg->gn_slabind[ic];
2041 /* Get the precomputed Gaussian value of curr_slab for curr_x */
2042 gaussian_xj = erg->gn_atom[ic];
2044 int slabIndex = n - erg->slab_first; /* slab index */
2046 /* The (unrotated) reference position of this atom is copied to yj0: */
2047 copy_rvec(erg->rotg->x_ref[iigrp], yj0);
2049 beta = calc_beta(xj, erg, n);
2051 /* The current center of this slab is saved in xcn: */
2052 copy_rvec(erg->slab_center[slabIndex], xcn);
2053 /* ... and the reference center in ycn: */
2054 copy_rvec(erg->slab_center_ref[slabIndex + erg->slab_buffer], ycn);
2056 rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
2058 /* Rotate: */
2059 mvmul(erg->rotmat, yj0_ycn, rjn); /* rjn = Omega.(yj0 - ycn) */
2061 /* Subtract the slab center from xj */
2062 rvec_sub(xj, xcn, tmpvec2); /* tmpvec2 = xj - xcn */
2064 /* In rare cases, when an atom position coincides with a slab center
2065 * (tmpvec2 == 0) we cannot compute the vector product for sjn.
2066 * However, since the atom is located directly on the pivot, this
2067 * slab's contribution to the force on that atom will be zero
2068 * anyway. Therefore, we directly move on to the next slab. */
2069 if (gmx_numzero(norm(tmpvec2))) /* 0 == norm(xj - xcn) */
2071 continue;
2074 /* Calculate sjn */
2075 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xj - xcn) */
2077 OOpsijstar = norm2(tmpvec) + erg->rotg->eps; /* OOpsij* = 1/psij* = |v x (xj-xcn)|^2 + eps */
2079 numerator = gmx::square(iprod(tmpvec, rjn));
2081 /*********************************/
2082 /* Add to the rotation potential */
2083 /*********************************/
2084 V += 0.5 * erg->rotg->k * wj * gaussian_xj * numerator / OOpsijstar;
2086 /* If requested, also calculate the potential for a set of angles
2087 * near the current reference angle */
2088 if (bCalcPotFit)
2090 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2092 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, fit_rjn);
2093 fit_numerator = gmx::square(iprod(tmpvec, fit_rjn));
2094 erg->PotAngleFit->V[ifit] +=
2095 0.5 * erg->rotg->k * wj * gaussian_xj * fit_numerator / OOpsijstar;
2099 /*************************************/
2100 /* Now calculate the force on atom j */
2101 /*************************************/
2103 OOpsij = norm(tmpvec); /* OOpsij = 1 / psij = |v x (xj - xcn)| */
2105 /* * v x (xj - xcn) */
2106 unitv(tmpvec, sjn); /* sjn = ---------------- */
2107 /* |v x (xj - xcn)| */
2109 sjn_rjn = iprod(sjn, rjn); /* sjn_rjn = sjn . rjn */
2112 /*** A. Calculate the first of the four sum terms: ****************/
2113 fac = OOpsij / OOpsijstar;
2114 svmul(fac, rjn, tmpvec);
2115 fac2 = fac * fac * OOpsij;
2116 svmul(fac2 * sjn_rjn, sjn, tmpvec2);
2117 rvec_dec(tmpvec, tmpvec2);
2118 fac2 = wj * gaussian_xj; /* also needed for sum4 */
2119 svmul(fac2 * sjn_rjn, tmpvec, slab_sum1vec_part);
2120 /********************/
2121 /*** Add to sum1: ***/
2122 /********************/
2123 rvec_inc(sum1vec_part, slab_sum1vec_part); /* sum1 still needs to vector multiplied with v */
2125 /*** B. Calculate the forth of the four sum terms: ****************/
2126 betasigpsi = beta * OOsigma2 * OOpsij; /* this is also needed for sum3 */
2127 /********************/
2128 /*** Add to sum4: ***/
2129 /********************/
2130 slab_sum4part = fac2 * betasigpsi * fac * sjn_rjn
2131 * sjn_rjn; /* Note that fac is still valid from above */
2132 sum4 += slab_sum4part;
2134 /*** C. Calculate Wjn for second and third sum */
2135 /* Note that we can safely divide by slab_weights since we check in
2136 * get_slab_centers that it is non-zero. */
2137 Wjn = gaussian_xj * mj / erg->slab_weights[slabIndex];
2139 /* We already have precalculated the inner sum for slab n */
2140 copy_rvec(erg->slab_innersumvec[slabIndex], innersumvec);
2142 /* Weigh the inner sum vector with Wjn */
2143 svmul(Wjn, innersumvec, innersumvec);
2145 /*** E. Calculate the second of the four sum terms: */
2146 /********************/
2147 /*** Add to sum2: ***/
2148 /********************/
2149 rvec_inc(sum2vec_part, innersumvec); /* sum2 still needs to be vector crossproduct'ed with v */
2151 /*** F. Calculate the third of the four sum terms: */
2152 slab_sum3part = betasigpsi * iprod(sjn, innersumvec);
2153 sum3 += slab_sum3part; /* still needs to be multiplied with v */
2155 /*** G. Calculate the torque on the local slab's axis: */
2156 if (bOutstepRot)
2158 /* Sum1 */
2159 cprod(slab_sum1vec_part, erg->vec, slab_sum1vec);
2160 /* Sum2 */
2161 cprod(innersumvec, erg->vec, slab_sum2vec);
2162 /* Sum3 */
2163 svmul(slab_sum3part, erg->vec, slab_sum3vec);
2164 /* Sum4 */
2165 svmul(slab_sum4part, erg->vec, slab_sum4vec);
2167 /* The force on atom ii from slab n only: */
2168 for (int m = 0; m < DIM; m++)
2170 slab_force[m] = erg->rotg->k
2171 * (-slab_sum1vec[m] + slab_sum2vec[m] - slab_sum3vec[m]
2172 + 0.5 * slab_sum4vec[m]);
2175 erg->slab_torque_v[slabIndex] += torque(erg->vec, slab_force, xj, xcn);
2177 } /* END of loop over slabs */
2179 /* Construct the four individual parts of the vector sum: */
2180 cprod(sum1vec_part, erg->vec, sum1vec); /* sum1vec = { } x v */
2181 cprod(sum2vec_part, erg->vec, sum2vec); /* sum2vec = { } x v */
2182 svmul(sum3, erg->vec, sum3vec); /* sum3vec = { } . v */
2183 svmul(sum4, erg->vec, sum4vec); /* sum4vec = { } . v */
2185 /* Store the additional force so that it can be added to the force
2186 * array after the normal forces have been evaluated */
2187 for (int m = 0; m < DIM; m++)
2189 erg->f_rot_loc[j][m] =
2190 erg->rotg->k * (-sum1vec[m] + sum2vec[m] - sum3vec[m] + 0.5 * sum4vec[m]);
2193 #ifdef SUM_PARTS
2194 fprintf(stderr, "sum1: %15.8f %15.8f %15.8f\n", -erg->rotg->k * sum1vec[XX],
2195 -erg->rotg->k * sum1vec[YY], -erg->rotg->k * sum1vec[ZZ]);
2196 fprintf(stderr, "sum2: %15.8f %15.8f %15.8f\n", erg->rotg->k * sum2vec[XX],
2197 erg->rotg->k * sum2vec[YY], erg->rotg->k * sum2vec[ZZ]);
2198 fprintf(stderr, "sum3: %15.8f %15.8f %15.8f\n", -erg->rotg->k * sum3vec[XX],
2199 -erg->rotg->k * sum3vec[YY], -erg->rotg->k * sum3vec[ZZ]);
2200 fprintf(stderr, "sum4: %15.8f %15.8f %15.8f\n", 0.5 * erg->rotg->k * sum4vec[XX],
2201 0.5 * erg->rotg->k * sum4vec[YY], 0.5 * erg->rotg->k * sum4vec[ZZ]);
2202 #endif
2204 PRINT_FORCE_J
2206 } /* END of loop over local atoms */
2208 return V;
2212 static real do_flex_lowlevel(gmx_enfrotgrp* erg,
2213 real sigma, /* The Gaussian width sigma */
2214 rvec x[],
2215 gmx_bool bOutstepRot,
2216 gmx_bool bOutstepSlab,
2217 const matrix box)
2219 int count, iigrp;
2220 rvec xj, yj0; /* current and reference position */
2221 rvec xcn, ycn; /* the current and the reference slab centers */
2222 rvec yj0_ycn; /* yj0 - ycn */
2223 rvec xj_xcn; /* xj - xcn */
2224 rvec qjn, fit_qjn; /* q_i^n */
2225 rvec sum_n1, sum_n2; /* Two contributions to the rotation force */
2226 rvec innersumvec; /* Inner part of sum_n2 */
2227 rvec s_n;
2228 rvec force_n; /* Single force from slab n on one atom */
2229 rvec force_n1, force_n2; /* First and second part of force_n */
2230 rvec tmpvec, tmpvec2, tmp_f; /* Helper variables */
2231 real V; /* The rotation potential energy */
2232 real OOsigma2; /* 1/(sigma^2) */
2233 real beta; /* beta_n(xj) */
2234 real bjn, fit_bjn; /* b_j^n */
2235 real gaussian_xj; /* Gaussian weight gn(xj) */
2236 real betan_xj_sigma2;
2237 real mj, wj; /* Mass-weighting of the positions */
2238 real N_M; /* N/M */
2239 gmx_bool bCalcPotFit;
2241 /* Pre-calculate the inner sums, so that we do not have to calculate
2242 * them again for every atom */
2243 flex_precalc_inner_sum(erg);
2245 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
2247 /********************************************************/
2248 /* Main loop over all local atoms of the rotation group */
2249 /********************************************************/
2250 OOsigma2 = 1.0 / (sigma * sigma);
2251 N_M = erg->rotg->nat * erg->invmass;
2252 V = 0.0;
2253 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
2254 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2256 for (gmx::index j = 0; j < localRotationGroupIndex.ssize(); j++)
2258 /* Local index of a rotation group atom */
2259 int ii = localRotationGroupIndex[j];
2260 /* Position of this atom in the collective array */
2261 iigrp = collectiveRotationGroupIndex[j];
2262 /* Mass-weighting */
2263 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2264 wj = N_M * mj;
2266 /* Current position of this atom: x[ii][XX/YY/ZZ]
2267 * Note that erg->xc_center contains the center of mass in case the flex-t
2268 * potential was chosen. For the flex potential erg->xc_center must be
2269 * zero. */
2270 rvec_sub(x[ii], erg->xc_center, xj);
2272 /* Shift this atom such that it is near its reference */
2273 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2275 /* Determine the slabs to loop over, i.e. the ones with contributions
2276 * larger than min_gaussian */
2277 count = get_single_atom_gaussians(xj, erg);
2279 clear_rvec(sum_n1);
2280 clear_rvec(sum_n2);
2282 /* Loop over the relevant slabs for this atom */
2283 for (int ic = 0; ic < count; ic++)
2285 int n = erg->gn_slabind[ic];
2287 /* Get the precomputed Gaussian for xj in slab n */
2288 gaussian_xj = erg->gn_atom[ic];
2290 int slabIndex = n - erg->slab_first; /* slab index */
2292 /* The (unrotated) reference position of this atom is saved in yj0: */
2293 copy_rvec(erg->rotg->x_ref[iigrp], yj0);
2295 beta = calc_beta(xj, erg, n);
2297 /* The current center of this slab is saved in xcn: */
2298 copy_rvec(erg->slab_center[slabIndex], xcn);
2299 /* ... and the reference center in ycn: */
2300 copy_rvec(erg->slab_center_ref[slabIndex + erg->slab_buffer], ycn);
2302 rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
2304 /* In rare cases, when an atom position coincides with a reference slab
2305 * center (yj0_ycn == 0) we cannot compute the normal vector qjn.
2306 * However, since the atom is located directly on the pivot, this
2307 * slab's contribution to the force on that atom will be zero
2308 * anyway. Therefore, we directly move on to the next slab. */
2309 if (gmx_numzero(norm(yj0_ycn))) /* 0 == norm(yj0 - ycn) */
2311 continue;
2314 /* Rotate: */
2315 mvmul(erg->rotmat, yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2317 /* Subtract the slab center from xj */
2318 rvec_sub(xj, xcn, xj_xcn); /* xj_xcn = xj - xcn */
2320 /* Calculate qjn */
2321 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2323 /* * v x Omega.(yj0-ycn) */
2324 unitv(tmpvec, qjn); /* qjn = --------------------- */
2325 /* |v x Omega.(yj0-ycn)| */
2327 bjn = iprod(qjn, xj_xcn); /* bjn = qjn * (xj - xcn) */
2329 /*********************************/
2330 /* Add to the rotation potential */
2331 /*********************************/
2332 V += 0.5 * erg->rotg->k * wj * gaussian_xj * gmx::square(bjn);
2334 /* If requested, also calculate the potential for a set of angles
2335 * near the current reference angle */
2336 if (bCalcPotFit)
2338 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2340 /* As above calculate Omega.(yj0-ycn), now for the other angles */
2341 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2342 /* As above calculate qjn */
2343 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2344 /* * v x Omega.(yj0-ycn) */
2345 unitv(tmpvec, fit_qjn); /* fit_qjn = --------------------- */
2346 /* |v x Omega.(yj0-ycn)| */
2347 fit_bjn = iprod(fit_qjn, xj_xcn); /* fit_bjn = fit_qjn * (xj - xcn) */
2348 /* Add to the rotation potential for this angle */
2349 erg->PotAngleFit->V[ifit] +=
2350 0.5 * erg->rotg->k * wj * gaussian_xj * gmx::square(fit_bjn);
2354 /****************************************************************/
2355 /* sum_n1 will typically be the main contribution to the force: */
2356 /****************************************************************/
2357 betan_xj_sigma2 = beta * OOsigma2; /* beta_n(xj)/sigma^2 */
2359 /* The next lines calculate
2360 * qjn - (bjn*beta(xj)/(2sigma^2))v */
2361 svmul(bjn * 0.5 * betan_xj_sigma2, erg->vec, tmpvec2);
2362 rvec_sub(qjn, tmpvec2, tmpvec);
2364 /* Multiply with gn(xj)*bjn: */
2365 svmul(gaussian_xj * bjn, tmpvec, tmpvec2);
2367 /* Sum over n: */
2368 rvec_inc(sum_n1, tmpvec2);
2370 /* We already have precalculated the Sn term for slab n */
2371 copy_rvec(erg->slab_innersumvec[slabIndex], s_n);
2372 /* * beta_n(xj) */
2373 svmul(betan_xj_sigma2 * iprod(s_n, xj_xcn), erg->vec, tmpvec); /* tmpvec = ---------- s_n (xj-xcn) */
2374 /* sigma^2 */
2376 rvec_sub(s_n, tmpvec, innersumvec);
2378 /* We can safely divide by slab_weights since we check in get_slab_centers
2379 * that it is non-zero. */
2380 svmul(gaussian_xj / erg->slab_weights[slabIndex], innersumvec, innersumvec);
2382 rvec_add(sum_n2, innersumvec, sum_n2);
2384 /* Calculate the torque: */
2385 if (bOutstepRot)
2387 /* The force on atom ii from slab n only: */
2388 svmul(-erg->rotg->k * wj, tmpvec2, force_n1); /* part 1 */
2389 svmul(erg->rotg->k * mj, innersumvec, force_n2); /* part 2 */
2390 rvec_add(force_n1, force_n2, force_n);
2391 erg->slab_torque_v[slabIndex] += torque(erg->vec, force_n, xj, xcn);
2393 } /* END of loop over slabs */
2395 /* Put both contributions together: */
2396 svmul(wj, sum_n1, sum_n1);
2397 svmul(mj, sum_n2, sum_n2);
2398 rvec_sub(sum_n2, sum_n1, tmp_f); /* F = -grad V */
2400 /* Store the additional force so that it can be added to the force
2401 * array after the normal forces have been evaluated */
2402 for (int m = 0; m < DIM; m++)
2404 erg->f_rot_loc[j][m] = erg->rotg->k * tmp_f[m];
2407 PRINT_FORCE_J
2409 } /* END of loop over local atoms */
2411 return V;
2414 static void sort_collective_coordinates(gmx_enfrotgrp* erg,
2415 sort_along_vec_t* data) /* Buffer for sorting the positions */
2417 /* The projection of the position vector on the rotation vector is
2418 * the relevant value for sorting. Fill the 'data' structure */
2419 for (int i = 0; i < erg->rotg->nat; i++)
2421 data[i].xcproj = iprod(erg->xc[i], erg->vec); /* sort criterium */
2422 data[i].m = erg->mc[i];
2423 data[i].ind = i;
2424 copy_rvec(erg->xc[i], data[i].x);
2425 copy_rvec(erg->rotg->x_ref[i], data[i].x_ref);
2427 /* Sort the 'data' structure */
2428 std::sort(data, data + erg->rotg->nat, [](const sort_along_vec_t& a, const sort_along_vec_t& b) {
2429 return a.xcproj < b.xcproj;
2432 /* Copy back the sorted values */
2433 for (int i = 0; i < erg->rotg->nat; i++)
2435 copy_rvec(data[i].x, erg->xc[i]);
2436 copy_rvec(data[i].x_ref, erg->xc_ref_sorted[i]);
2437 erg->mc_sorted[i] = data[i].m;
2438 erg->xc_sortind[i] = data[i].ind;
2443 /* For each slab, get the first and the last index of the sorted atom
2444 * indices */
2445 static void get_firstlast_atom_per_slab(const gmx_enfrotgrp* erg)
2447 real beta;
2449 /* Find the first atom that needs to enter the calculation for each slab */
2450 int n = erg->slab_first; /* slab */
2451 int i = 0; /* start with the first atom */
2454 /* Find the first atom that significantly contributes to this slab */
2455 do /* move forward in position until a large enough beta is found */
2457 beta = calc_beta(erg->xc[i], erg, n);
2458 i++;
2459 } while ((beta < -erg->max_beta) && (i < erg->rotg->nat));
2460 i--;
2461 int slabIndex = n - erg->slab_first; /* slab index */
2462 erg->firstatom[slabIndex] = i;
2463 /* Proceed to the next slab */
2464 n++;
2465 } while (n <= erg->slab_last);
2467 /* Find the last atom for each slab */
2468 n = erg->slab_last; /* start with last slab */
2469 i = erg->rotg->nat - 1; /* start with the last atom */
2472 do /* move backward in position until a large enough beta is found */
2474 beta = calc_beta(erg->xc[i], erg, n);
2475 i--;
2476 } while ((beta > erg->max_beta) && (i > -1));
2477 i++;
2478 int slabIndex = n - erg->slab_first; /* slab index */
2479 erg->lastatom[slabIndex] = i;
2480 /* Proceed to the next slab */
2481 n--;
2482 } while (n >= erg->slab_first);
2486 /* Determine the very first and very last slab that needs to be considered
2487 * For the first slab that needs to be considered, we have to find the smallest
2488 * n that obeys:
2490 * x_first * v - n*Delta_x <= beta_max
2492 * slab index n, slab distance Delta_x, rotation vector v. For the last slab we
2493 * have to find the largest n that obeys
2495 * x_last * v - n*Delta_x >= -beta_max
2498 static inline int get_first_slab(const gmx_enfrotgrp* erg,
2499 rvec firstatom) /* First atom after sorting along the rotation vector v */
2501 /* Find the first slab for the first atom */
2502 return static_cast<int>(ceil(
2503 static_cast<double>((iprod(firstatom, erg->vec) - erg->max_beta) / erg->rotg->slab_dist)));
2507 static inline int get_last_slab(const gmx_enfrotgrp* erg, rvec lastatom) /* Last atom along v */
2509 /* Find the last slab for the last atom */
2510 return static_cast<int>(floor(
2511 static_cast<double>((iprod(lastatom, erg->vec) + erg->max_beta) / erg->rotg->slab_dist)));
2515 static void get_firstlast_slab_check(gmx_enfrotgrp* erg, /* The rotation group (data only accessible in this file) */
2516 rvec firstatom, /* First atom after sorting along the rotation vector v */
2517 rvec lastatom) /* Last atom along v */
2519 erg->slab_first = get_first_slab(erg, firstatom);
2520 erg->slab_last = get_last_slab(erg, lastatom);
2522 /* Calculate the slab buffer size, which changes when slab_first changes */
2523 erg->slab_buffer = erg->slab_first - erg->slab_first_ref;
2525 /* Check whether we have reference data to compare against */
2526 if (erg->slab_first < erg->slab_first_ref)
2528 gmx_fatal(FARGS, "%s No reference data for first slab (n=%d), unable to proceed.", RotStr,
2529 erg->slab_first);
2532 /* Check whether we have reference data to compare against */
2533 if (erg->slab_last > erg->slab_last_ref)
2535 gmx_fatal(FARGS, "%s No reference data for last slab (n=%d), unable to proceed.", RotStr,
2536 erg->slab_last);
2541 /* Enforced rotation with a flexible axis */
2542 static void do_flexible(gmx_bool bMaster,
2543 gmx_enfrot* enfrot, /* Other rotation data */
2544 gmx_enfrotgrp* erg,
2545 rvec x[], /* The local positions */
2546 const matrix box,
2547 double t, /* Time in picoseconds */
2548 gmx_bool bOutstepRot, /* Output to main rotation output file */
2549 gmx_bool bOutstepSlab) /* Output per-slab data */
2551 int nslabs;
2552 real sigma; /* The Gaussian width sigma */
2554 /* Define the sigma value */
2555 sigma = 0.7 * erg->rotg->slab_dist;
2557 /* Sort the collective coordinates erg->xc along the rotation vector. This is
2558 * an optimization for the inner loop. */
2559 sort_collective_coordinates(erg, enfrot->data);
2561 /* Determine the first relevant slab for the first atom and the last
2562 * relevant slab for the last atom */
2563 get_firstlast_slab_check(erg, erg->xc[0], erg->xc[erg->rotg->nat - 1]);
2565 /* Determine for each slab depending on the min_gaussian cutoff criterium,
2566 * a first and a last atom index inbetween stuff needs to be calculated */
2567 get_firstlast_atom_per_slab(erg);
2569 /* Determine the gaussian-weighted center of positions for all slabs */
2570 get_slab_centers(erg, erg->xc, erg->mc_sorted, t, enfrot->out_slabs, bOutstepSlab, FALSE);
2572 /* Clear the torque per slab from last time step: */
2573 nslabs = erg->slab_last - erg->slab_first + 1;
2574 for (int l = 0; l < nslabs; l++)
2576 erg->slab_torque_v[l] = 0.0;
2579 /* Call the rotational forces kernel */
2580 if (erg->rotg->eType == erotgFLEX || erg->rotg->eType == erotgFLEXT)
2582 erg->V = do_flex_lowlevel(erg, sigma, x, bOutstepRot, bOutstepSlab, box);
2584 else if (erg->rotg->eType == erotgFLEX2 || erg->rotg->eType == erotgFLEX2T)
2586 erg->V = do_flex2_lowlevel(erg, sigma, x, bOutstepRot, bOutstepSlab, box);
2588 else
2590 gmx_fatal(FARGS, "Unknown flexible rotation type");
2593 /* Determine angle by RMSD fit to the reference - Let's hope this */
2594 /* only happens once in a while, since this is not parallelized! */
2595 if (bMaster && (erotgFitPOT != erg->rotg->eFittype))
2597 if (bOutstepRot)
2599 /* Fit angle of the whole rotation group */
2600 erg->angle_v = flex_fit_angle(erg);
2602 if (bOutstepSlab)
2604 /* Fit angle of each slab */
2605 flex_fit_angle_perslab(erg, t, erg->degangle, enfrot->out_angles);
2609 /* Lump together the torques from all slabs: */
2610 erg->torque_v = 0.0;
2611 for (int l = 0; l < nslabs; l++)
2613 erg->torque_v += erg->slab_torque_v[l];
2618 /* Calculate the angle between reference and actual rotation group atom,
2619 * both projected into a plane perpendicular to the rotation vector: */
2620 static void angle(const gmx_enfrotgrp* erg,
2621 rvec x_act,
2622 rvec x_ref,
2623 real* alpha,
2624 real* weight) /* atoms near the rotation axis should count less than atoms far away */
2626 rvec xp, xrp; /* current and reference positions projected on a plane perpendicular to pg->vec */
2627 rvec dum;
2630 /* Project x_ref and x into a plane through the origin perpendicular to rot_vec: */
2631 /* Project x_ref: xrp = x_ref - (vec * x_ref) * vec */
2632 svmul(iprod(erg->vec, x_ref), erg->vec, dum);
2633 rvec_sub(x_ref, dum, xrp);
2634 /* Project x_act: */
2635 svmul(iprod(erg->vec, x_act), erg->vec, dum);
2636 rvec_sub(x_act, dum, xp);
2638 /* Retrieve information about which vector precedes. gmx_angle always
2639 * returns a positive angle. */
2640 cprod(xp, xrp, dum); /* if reference precedes, this is pointing into the same direction as vec */
2642 if (iprod(erg->vec, dum) >= 0)
2644 *alpha = -gmx_angle(xrp, xp);
2646 else
2648 *alpha = +gmx_angle(xrp, xp);
2651 /* Also return the weight */
2652 *weight = norm(xp);
2656 /* Project first vector onto a plane perpendicular to the second vector
2657 * dr = dr - (dr.v)v
2658 * Note that v must be of unit length.
2660 static inline void project_onto_plane(rvec dr, const rvec v)
2662 rvec tmp;
2665 svmul(iprod(dr, v), v, tmp); /* tmp = (dr.v)v */
2666 rvec_dec(dr, tmp); /* dr = dr - (dr.v)v */
2670 /* Fixed rotation: The rotation reference group rotates around the v axis. */
2671 /* The atoms of the actual rotation group are attached with imaginary */
2672 /* springs to the reference atoms. */
2673 static void do_fixed(gmx_enfrotgrp* erg,
2674 gmx_bool bOutstepRot, /* Output to main rotation output file */
2675 gmx_bool bOutstepSlab) /* Output per-slab data */
2677 rvec dr;
2678 rvec tmp_f; /* Force */
2679 real alpha; /* a single angle between an actual and a reference position */
2680 real weight; /* single weight for a single angle */
2681 rvec xi_xc; /* xi - xc */
2682 gmx_bool bCalcPotFit;
2683 rvec fit_xr_loc;
2685 /* for mass weighting: */
2686 real wi; /* Mass-weighting of the positions */
2687 real N_M; /* N/M */
2688 real k_wi; /* k times wi */
2690 gmx_bool bProject;
2692 bProject = (erg->rotg->eType == erotgPM) || (erg->rotg->eType == erotgPMPF);
2693 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
2695 N_M = erg->rotg->nat * erg->invmass;
2696 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2697 /* Each process calculates the forces on its local atoms */
2698 for (size_t j = 0; j < erg->atomSet->numAtomsLocal(); j++)
2700 /* Calculate (x_i-x_c) resp. (x_i-u) */
2701 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xi_xc);
2703 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2704 rvec_sub(erg->xr_loc[j], xi_xc, dr);
2706 if (bProject)
2708 project_onto_plane(dr, erg->vec);
2711 /* Mass-weighting */
2712 wi = N_M * erg->m_loc[j];
2714 /* Store the additional force so that it can be added to the force
2715 * array after the normal forces have been evaluated */
2716 k_wi = erg->rotg->k * wi;
2717 for (int m = 0; m < DIM; m++)
2719 tmp_f[m] = k_wi * dr[m];
2720 erg->f_rot_loc[j][m] = tmp_f[m];
2721 erg->V += 0.5 * k_wi * gmx::square(dr[m]);
2724 /* If requested, also calculate the potential for a set of angles
2725 * near the current reference angle */
2726 if (bCalcPotFit)
2728 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2730 /* Index of this rotation group atom with respect to the whole rotation group */
2731 int jj = collectiveRotationGroupIndex[j];
2733 /* Rotate with the alternative angle. Like rotate_local_reference(),
2734 * just for a single local atom */
2735 mvmul(erg->PotAngleFit->rotmat[ifit], erg->rotg->x_ref[jj],
2736 fit_xr_loc); /* fit_xr_loc = Omega*(y_i-y_c) */
2738 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2739 rvec_sub(fit_xr_loc, xi_xc, dr);
2741 if (bProject)
2743 project_onto_plane(dr, erg->vec);
2746 /* Add to the rotation potential for this angle: */
2747 erg->PotAngleFit->V[ifit] += 0.5 * k_wi * norm2(dr);
2751 if (bOutstepRot)
2753 /* Add to the torque of this rotation group */
2754 erg->torque_v += torque(erg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2756 /* Calculate the angle between reference and actual rotation group atom. */
2757 angle(erg, xi_xc, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2758 erg->angle_v += alpha * weight;
2759 erg->weight_v += weight;
2761 /* If you want enforced rotation to contribute to the virial,
2762 * activate the following lines:
2763 if (MASTER(cr))
2765 Add the rotation contribution to the virial
2766 for(j=0; j<DIM; j++)
2767 for(m=0;m<DIM;m++)
2768 vir[j][m] += 0.5*f[ii][j]*dr[m];
2772 PRINT_FORCE_J
2774 } /* end of loop over local rotation group atoms */
2778 /* Calculate the radial motion potential and forces */
2779 static void do_radial_motion(gmx_enfrotgrp* erg,
2780 gmx_bool bOutstepRot, /* Output to main rotation output file */
2781 gmx_bool bOutstepSlab) /* Output per-slab data */
2783 rvec tmp_f; /* Force */
2784 real alpha; /* a single angle between an actual and a reference position */
2785 real weight; /* single weight for a single angle */
2786 rvec xj_u; /* xj - u */
2787 rvec tmpvec, fit_tmpvec;
2788 real fac, fac2, sum = 0.0;
2789 rvec pj;
2790 gmx_bool bCalcPotFit;
2792 /* For mass weighting: */
2793 real wj; /* Mass-weighting of the positions */
2794 real N_M; /* N/M */
2796 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
2798 N_M = erg->rotg->nat * erg->invmass;
2799 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2800 /* Each process calculates the forces on its local atoms */
2801 for (size_t j = 0; j < erg->atomSet->numAtomsLocal(); j++)
2803 /* Calculate (xj-u) */
2804 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xj_u); /* xj_u = xj-u */
2806 /* Calculate Omega.(yj0-u) */
2807 cprod(erg->vec, erg->xr_loc[j], tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2809 /* * v x Omega.(yj0-u) */
2810 unitv(tmpvec, pj); /* pj = --------------------- */
2811 /* | v x Omega.(yj0-u) | */
2813 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2814 fac2 = fac * fac;
2816 /* Mass-weighting */
2817 wj = N_M * erg->m_loc[j];
2819 /* Store the additional force so that it can be added to the force
2820 * array after the normal forces have been evaluated */
2821 svmul(-erg->rotg->k * wj * fac, pj, tmp_f);
2822 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2823 sum += wj * fac2;
2825 /* If requested, also calculate the potential for a set of angles
2826 * near the current reference angle */
2827 if (bCalcPotFit)
2829 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2831 /* Index of this rotation group atom with respect to the whole rotation group */
2832 int jj = collectiveRotationGroupIndex[j];
2834 /* Rotate with the alternative angle. Like rotate_local_reference(),
2835 * just for a single local atom */
2836 mvmul(erg->PotAngleFit->rotmat[ifit], erg->rotg->x_ref[jj],
2837 fit_tmpvec); /* fit_tmpvec = Omega*(yj0-u) */
2839 /* Calculate Omega.(yj0-u) */
2840 cprod(erg->vec, fit_tmpvec, tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2841 /* * v x Omega.(yj0-u) */
2842 unitv(tmpvec, pj); /* pj = --------------------- */
2843 /* | v x Omega.(yj0-u) | */
2845 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2846 fac2 = fac * fac;
2848 /* Add to the rotation potential for this angle: */
2849 erg->PotAngleFit->V[ifit] += 0.5 * erg->rotg->k * wj * fac2;
2853 if (bOutstepRot)
2855 /* Add to the torque of this rotation group */
2856 erg->torque_v += torque(erg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2858 /* Calculate the angle between reference and actual rotation group atom. */
2859 angle(erg, xj_u, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2860 erg->angle_v += alpha * weight;
2861 erg->weight_v += weight;
2864 PRINT_FORCE_J
2866 } /* end of loop over local rotation group atoms */
2867 erg->V = 0.5 * erg->rotg->k * sum;
2871 /* Calculate the radial motion pivot-free potential and forces */
2872 static void do_radial_motion_pf(gmx_enfrotgrp* erg,
2873 rvec x[], /* The positions */
2874 const matrix box, /* The simulation box */
2875 gmx_bool bOutstepRot, /* Output to main rotation output file */
2876 gmx_bool bOutstepSlab) /* Output per-slab data */
2878 rvec xj; /* Current position */
2879 rvec xj_xc; /* xj - xc */
2880 rvec yj0_yc0; /* yj0 - yc0 */
2881 rvec tmp_f; /* Force */
2882 real alpha; /* a single angle between an actual and a reference position */
2883 real weight; /* single weight for a single angle */
2884 rvec tmpvec, tmpvec2;
2885 rvec innersumvec; /* Precalculation of the inner sum */
2886 rvec innersumveckM;
2887 real fac, fac2, V = 0.0;
2888 rvec qi, qj;
2889 gmx_bool bCalcPotFit;
2891 /* For mass weighting: */
2892 real mj, wi, wj; /* Mass-weighting of the positions */
2893 real N_M; /* N/M */
2895 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
2897 N_M = erg->rotg->nat * erg->invmass;
2899 /* Get the current center of the rotation group: */
2900 get_center(erg->xc, erg->mc, erg->rotg->nat, erg->xc_center);
2902 /* Precalculate Sum_i [ wi qi.(xi-xc) qi ] which is needed for every single j */
2903 clear_rvec(innersumvec);
2904 for (int i = 0; i < erg->rotg->nat; i++)
2906 /* Mass-weighting */
2907 wi = N_M * erg->mc[i];
2909 /* Calculate qi. Note that xc_ref_center has already been subtracted from
2910 * x_ref in init_rot_group.*/
2911 mvmul(erg->rotmat, erg->rotg->x_ref[i], tmpvec); /* tmpvec = Omega.(yi0-yc0) */
2913 cprod(erg->vec, tmpvec, tmpvec2); /* tmpvec2 = v x Omega.(yi0-yc0) */
2915 /* * v x Omega.(yi0-yc0) */
2916 unitv(tmpvec2, qi); /* qi = ----------------------- */
2917 /* | v x Omega.(yi0-yc0) | */
2919 rvec_sub(erg->xc[i], erg->xc_center, tmpvec); /* tmpvec = xi-xc */
2921 svmul(wi * iprod(qi, tmpvec), qi, tmpvec2);
2923 rvec_inc(innersumvec, tmpvec2);
2925 svmul(erg->rotg->k * erg->invmass, innersumvec, innersumveckM);
2927 /* Each process calculates the forces on its local atoms */
2928 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
2929 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2930 for (gmx::index j = 0; j < localRotationGroupIndex.ssize(); j++)
2932 /* Local index of a rotation group atom */
2933 int ii = localRotationGroupIndex[j];
2934 /* Position of this atom in the collective array */
2935 int iigrp = collectiveRotationGroupIndex[j];
2936 /* Mass-weighting */
2937 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2938 wj = N_M * mj;
2940 /* Current position of this atom: x[ii][XX/YY/ZZ] */
2941 copy_rvec(x[ii], xj);
2943 /* Shift this atom such that it is near its reference */
2944 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2946 /* The (unrotated) reference position is yj0. yc0 has already
2947 * been subtracted in init_rot_group */
2948 copy_rvec(erg->rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
2950 /* Calculate Omega.(yj0-yc0) */
2951 mvmul(erg->rotmat, yj0_yc0, tmpvec2); /* tmpvec2 = Omega.(yj0 - yc0) */
2953 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
2955 /* * v x Omega.(yj0-yc0) */
2956 unitv(tmpvec, qj); /* qj = ----------------------- */
2957 /* | v x Omega.(yj0-yc0) | */
2959 /* Calculate (xj-xc) */
2960 rvec_sub(xj, erg->xc_center, xj_xc); /* xj_xc = xj-xc */
2962 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
2963 fac2 = fac * fac;
2965 /* Store the additional force so that it can be added to the force
2966 * array after the normal forces have been evaluated */
2967 svmul(-erg->rotg->k * wj * fac, qj, tmp_f); /* part 1 of force */
2968 svmul(mj, innersumveckM, tmpvec); /* part 2 of force */
2969 rvec_inc(tmp_f, tmpvec);
2970 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2971 V += wj * fac2;
2973 /* If requested, also calculate the potential for a set of angles
2974 * near the current reference angle */
2975 if (bCalcPotFit)
2977 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2979 /* Rotate with the alternative angle. Like rotate_local_reference(),
2980 * just for a single local atom */
2981 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, tmpvec2); /* tmpvec2 = Omega*(yj0-yc0) */
2983 /* Calculate Omega.(yj0-u) */
2984 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
2985 /* * v x Omega.(yj0-yc0) */
2986 unitv(tmpvec, qj); /* qj = ----------------------- */
2987 /* | v x Omega.(yj0-yc0) | */
2989 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
2990 fac2 = fac * fac;
2992 /* Add to the rotation potential for this angle: */
2993 erg->PotAngleFit->V[ifit] += 0.5 * erg->rotg->k * wj * fac2;
2997 if (bOutstepRot)
2999 /* Add to the torque of this rotation group */
3000 erg->torque_v += torque(erg->vec, tmp_f, xj, erg->xc_center);
3002 /* Calculate the angle between reference and actual rotation group atom. */
3003 angle(erg, xj_xc, yj0_yc0, &alpha, &weight); /* angle in rad, weighted */
3004 erg->angle_v += alpha * weight;
3005 erg->weight_v += weight;
3008 PRINT_FORCE_J
3010 } /* end of loop over local rotation group atoms */
3011 erg->V = 0.5 * erg->rotg->k * V;
3015 /* Precalculate the inner sum for the radial motion 2 forces */
3016 static void radial_motion2_precalc_inner_sum(const gmx_enfrotgrp* erg, rvec innersumvec)
3018 int i;
3019 rvec xi_xc; /* xj - xc */
3020 rvec tmpvec, tmpvec2;
3021 real fac;
3022 rvec ri, si;
3023 real siri;
3024 rvec v_xi_xc; /* v x (xj - u) */
3025 real psii, psiistar;
3026 real wi; /* Mass-weighting of the positions */
3027 real N_M; /* N/M */
3028 rvec sumvec;
3030 N_M = erg->rotg->nat * erg->invmass;
3032 /* Loop over the collective set of positions */
3033 clear_rvec(sumvec);
3034 for (i = 0; i < erg->rotg->nat; i++)
3036 /* Mass-weighting */
3037 wi = N_M * erg->mc[i];
3039 rvec_sub(erg->xc[i], erg->xc_center, xi_xc); /* xi_xc = xi-xc */
3041 /* Calculate ri. Note that xc_ref_center has already been subtracted from
3042 * x_ref in init_rot_group.*/
3043 mvmul(erg->rotmat, erg->rotg->x_ref[i], ri); /* ri = Omega.(yi0-yc0) */
3045 cprod(erg->vec, xi_xc, v_xi_xc); /* v_xi_xc = v x (xi-u) */
3047 fac = norm2(v_xi_xc);
3048 /* * 1 */
3049 psiistar = 1.0 / (fac + erg->rotg->eps); /* psiistar = --------------------- */
3050 /* |v x (xi-xc)|^2 + eps */
3052 psii = gmx::invsqrt(fac); /* 1 */
3053 /* psii = ------------- */
3054 /* |v x (xi-xc)| */
3056 svmul(psii, v_xi_xc, si); /* si = psii * (v x (xi-xc) ) */
3058 siri = iprod(si, ri); /* siri = si.ri */
3060 svmul(psiistar / psii, ri, tmpvec);
3061 svmul(psiistar * psiistar / (psii * psii * psii) * siri, si, tmpvec2);
3062 rvec_dec(tmpvec, tmpvec2);
3063 cprod(tmpvec, erg->vec, tmpvec2);
3065 svmul(wi * siri, tmpvec2, tmpvec);
3067 rvec_inc(sumvec, tmpvec);
3069 svmul(erg->rotg->k * erg->invmass, sumvec, innersumvec);
3073 /* Calculate the radial motion 2 potential and forces */
3074 static void do_radial_motion2(gmx_enfrotgrp* erg,
3075 rvec x[], /* The positions */
3076 const matrix box, /* The simulation box */
3077 gmx_bool bOutstepRot, /* Output to main rotation output file */
3078 gmx_bool bOutstepSlab) /* Output per-slab data */
3080 rvec xj; /* Position */
3081 real alpha; /* a single angle between an actual and a reference position */
3082 real weight; /* single weight for a single angle */
3083 rvec xj_u; /* xj - u */
3084 rvec yj0_yc0; /* yj0 -yc0 */
3085 rvec tmpvec, tmpvec2;
3086 real fac, fit_fac, fac2, Vpart = 0.0;
3087 rvec rj, fit_rj, sj;
3088 real sjrj;
3089 rvec v_xj_u; /* v x (xj - u) */
3090 real psij, psijstar;
3091 real mj, wj; /* For mass-weighting of the positions */
3092 real N_M; /* N/M */
3093 gmx_bool bPF;
3094 rvec innersumvec;
3095 gmx_bool bCalcPotFit;
3097 bPF = erg->rotg->eType == erotgRM2PF;
3098 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
3100 clear_rvec(yj0_yc0); /* Make the compiler happy */
3102 clear_rvec(innersumvec);
3103 if (bPF)
3105 /* For the pivot-free variant we have to use the current center of
3106 * mass of the rotation group instead of the pivot u */
3107 get_center(erg->xc, erg->mc, erg->rotg->nat, erg->xc_center);
3109 /* Also, we precalculate the second term of the forces that is identical
3110 * (up to the weight factor mj) for all forces */
3111 radial_motion2_precalc_inner_sum(erg, innersumvec);
3114 N_M = erg->rotg->nat * erg->invmass;
3116 /* Each process calculates the forces on its local atoms */
3117 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
3118 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
3119 for (gmx::index j = 0; j < localRotationGroupIndex.ssize(); j++)
3121 if (bPF)
3123 /* Local index of a rotation group atom */
3124 int ii = localRotationGroupIndex[j];
3125 /* Position of this atom in the collective array */
3126 int iigrp = collectiveRotationGroupIndex[j];
3127 /* Mass-weighting */
3128 mj = erg->mc[iigrp];
3130 /* Current position of this atom: x[ii] */
3131 copy_rvec(x[ii], xj);
3133 /* Shift this atom such that it is near its reference */
3134 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
3136 /* The (unrotated) reference position is yj0. yc0 has already
3137 * been subtracted in init_rot_group */
3138 copy_rvec(erg->rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
3140 /* Calculate Omega.(yj0-yc0) */
3141 mvmul(erg->rotmat, yj0_yc0, rj); /* rj = Omega.(yj0-yc0) */
3143 else
3145 mj = erg->m_loc[j];
3146 copy_rvec(erg->x_loc_pbc[j], xj);
3147 copy_rvec(erg->xr_loc[j], rj); /* rj = Omega.(yj0-u) */
3149 /* Mass-weighting */
3150 wj = N_M * mj;
3152 /* Calculate (xj-u) resp. (xj-xc) */
3153 rvec_sub(xj, erg->xc_center, xj_u); /* xj_u = xj-u */
3155 cprod(erg->vec, xj_u, v_xj_u); /* v_xj_u = v x (xj-u) */
3157 fac = norm2(v_xj_u);
3158 /* * 1 */
3159 psijstar = 1.0 / (fac + erg->rotg->eps); /* psistar = -------------------- */
3160 /* * |v x (xj-u)|^2 + eps */
3162 psij = gmx::invsqrt(fac); /* 1 */
3163 /* psij = ------------ */
3164 /* |v x (xj-u)| */
3166 svmul(psij, v_xj_u, sj); /* sj = psij * (v x (xj-u) ) */
3168 fac = iprod(v_xj_u, rj); /* fac = (v x (xj-u)).rj */
3169 fac2 = fac * fac;
3171 sjrj = iprod(sj, rj); /* sjrj = sj.rj */
3173 svmul(psijstar / psij, rj, tmpvec);
3174 svmul(psijstar * psijstar / (psij * psij * psij) * sjrj, sj, tmpvec2);
3175 rvec_dec(tmpvec, tmpvec2);
3176 cprod(tmpvec, erg->vec, tmpvec2);
3178 /* Store the additional force so that it can be added to the force
3179 * array after the normal forces have been evaluated */
3180 svmul(-erg->rotg->k * wj * sjrj, tmpvec2, tmpvec);
3181 svmul(mj, innersumvec, tmpvec2); /* This is != 0 only for the pivot-free variant */
3183 rvec_add(tmpvec2, tmpvec, erg->f_rot_loc[j]);
3184 Vpart += wj * psijstar * fac2;
3186 /* If requested, also calculate the potential for a set of angles
3187 * near the current reference angle */
3188 if (bCalcPotFit)
3190 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
3192 if (bPF)
3194 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, fit_rj); /* fit_rj = Omega.(yj0-yc0) */
3196 else
3198 /* Position of this atom in the collective array */
3199 int iigrp = collectiveRotationGroupIndex[j];
3200 /* Rotate with the alternative angle. Like rotate_local_reference(),
3201 * just for a single local atom */
3202 mvmul(erg->PotAngleFit->rotmat[ifit], erg->rotg->x_ref[iigrp],
3203 fit_rj); /* fit_rj = Omega*(yj0-u) */
3205 fit_fac = iprod(v_xj_u, fit_rj); /* fac = (v x (xj-u)).fit_rj */
3206 /* Add to the rotation potential for this angle: */
3207 erg->PotAngleFit->V[ifit] += 0.5 * erg->rotg->k * wj * psijstar * fit_fac * fit_fac;
3211 if (bOutstepRot)
3213 /* Add to the torque of this rotation group */
3214 erg->torque_v += torque(erg->vec, erg->f_rot_loc[j], xj, erg->xc_center);
3216 /* Calculate the angle between reference and actual rotation group atom. */
3217 angle(erg, xj_u, rj, &alpha, &weight); /* angle in rad, weighted */
3218 erg->angle_v += alpha * weight;
3219 erg->weight_v += weight;
3222 PRINT_FORCE_J
3224 } /* end of loop over local rotation group atoms */
3225 erg->V = 0.5 * erg->rotg->k * Vpart;
3229 /* Determine the smallest and largest position vector (with respect to the
3230 * rotation vector) for the reference group */
3231 static void get_firstlast_atom_ref(const gmx_enfrotgrp* erg, int* firstindex, int* lastindex)
3233 int i;
3234 real xcproj; /* The projection of a reference position on the
3235 rotation vector */
3236 real minproj, maxproj; /* Smallest and largest projection on v */
3238 /* Start with some value */
3239 minproj = iprod(erg->rotg->x_ref[0], erg->vec);
3240 maxproj = minproj;
3242 /* This is just to ensure that it still works if all the atoms of the
3243 * reference structure are situated in a plane perpendicular to the rotation
3244 * vector */
3245 *firstindex = 0;
3246 *lastindex = erg->rotg->nat - 1;
3248 /* Loop over all atoms of the reference group,
3249 * project them on the rotation vector to find the extremes */
3250 for (i = 0; i < erg->rotg->nat; i++)
3252 xcproj = iprod(erg->rotg->x_ref[i], erg->vec);
3253 if (xcproj < minproj)
3255 minproj = xcproj;
3256 *firstindex = i;
3258 if (xcproj > maxproj)
3260 maxproj = xcproj;
3261 *lastindex = i;
3267 /* Allocate memory for the slabs */
3268 static void allocate_slabs(gmx_enfrotgrp* erg, FILE* fplog, gmx_bool bVerbose)
3270 /* More slabs than are defined for the reference are never needed */
3271 int nslabs = erg->slab_last_ref - erg->slab_first_ref + 1;
3273 /* Remember how many we allocated */
3274 erg->nslabs_alloc = nslabs;
3276 if ((nullptr != fplog) && bVerbose)
3278 fprintf(fplog, "%s allocating memory to store data for %d slabs (rotation group %d).\n",
3279 RotStr, nslabs, erg->groupIndex);
3281 snew(erg->slab_center, nslabs);
3282 snew(erg->slab_center_ref, nslabs);
3283 snew(erg->slab_weights, nslabs);
3284 snew(erg->slab_torque_v, nslabs);
3285 snew(erg->slab_data, nslabs);
3286 snew(erg->gn_atom, nslabs);
3287 snew(erg->gn_slabind, nslabs);
3288 snew(erg->slab_innersumvec, nslabs);
3289 for (int i = 0; i < nslabs; i++)
3291 snew(erg->slab_data[i].x, erg->rotg->nat);
3292 snew(erg->slab_data[i].ref, erg->rotg->nat);
3293 snew(erg->slab_data[i].weight, erg->rotg->nat);
3295 snew(erg->xc_ref_sorted, erg->rotg->nat);
3296 snew(erg->xc_sortind, erg->rotg->nat);
3297 snew(erg->firstatom, nslabs);
3298 snew(erg->lastatom, nslabs);
3302 /* From the extreme positions of the reference group, determine the first
3303 * and last slab of the reference. We can never have more slabs in the real
3304 * simulation than calculated here for the reference.
3306 static void get_firstlast_slab_ref(gmx_enfrotgrp* erg, real mc[], int ref_firstindex, int ref_lastindex)
3308 rvec dummy;
3310 int first = get_first_slab(erg, erg->rotg->x_ref[ref_firstindex]);
3311 int last = get_last_slab(erg, erg->rotg->x_ref[ref_lastindex]);
3313 while (get_slab_weight(first, erg, erg->rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3315 first--;
3317 erg->slab_first_ref = first + 1;
3318 while (get_slab_weight(last, erg, erg->rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3320 last++;
3322 erg->slab_last_ref = last - 1;
3326 /* Special version of copy_rvec:
3327 * During the copy procedure of xcurr to b, the correct PBC image is chosen
3328 * such that the copied vector ends up near its reference position xref */
3329 static inline void copy_correct_pbc_image(const rvec xcurr, /* copy vector xcurr ... */
3330 rvec b, /* ... to b ... */
3331 const rvec xref, /* choosing the PBC image such that b ends up near xref */
3332 const matrix box,
3333 int npbcdim)
3335 rvec dx;
3336 int d, m;
3337 ivec shift;
3340 /* Shortest PBC distance between the atom and its reference */
3341 rvec_sub(xcurr, xref, dx);
3343 /* Determine the shift for this atom */
3344 clear_ivec(shift);
3345 for (m = npbcdim - 1; m >= 0; m--)
3347 while (dx[m] < -0.5 * box[m][m])
3349 for (d = 0; d < DIM; d++)
3351 dx[d] += box[m][d];
3353 shift[m]++;
3355 while (dx[m] >= 0.5 * box[m][m])
3357 for (d = 0; d < DIM; d++)
3359 dx[d] -= box[m][d];
3361 shift[m]--;
3365 /* Apply the shift to the position */
3366 copy_rvec(xcurr, b);
3367 shift_single_coord(box, b, shift);
3371 static void init_rot_group(FILE* fplog,
3372 const t_commrec* cr,
3373 gmx_enfrotgrp* erg,
3374 rvec* x,
3375 gmx_mtop_t* mtop,
3376 gmx_bool bVerbose,
3377 FILE* out_slabs,
3378 const matrix box,
3379 t_inputrec* ir,
3380 gmx_bool bOutputCenters)
3382 rvec coord, xref, *xdum;
3383 gmx_bool bFlex, bColl;
3384 int ref_firstindex, ref_lastindex;
3385 real mass, totalmass;
3386 real start = 0.0;
3387 double t_start;
3388 const t_rotgrp* rotg = erg->rotg;
3391 /* Do we have a flexible axis? */
3392 bFlex = ISFLEX(rotg);
3393 /* Do we use a global set of coordinates? */
3394 bColl = ISCOLL(rotg);
3396 /* Allocate space for collective coordinates if needed */
3397 if (bColl)
3399 snew(erg->xc, erg->rotg->nat);
3400 snew(erg->xc_shifts, erg->rotg->nat);
3401 snew(erg->xc_eshifts, erg->rotg->nat);
3402 snew(erg->xc_old, erg->rotg->nat);
3404 if (erg->rotg->eFittype == erotgFitNORM)
3406 snew(erg->xc_ref_length, erg->rotg->nat); /* in case fit type NORM is chosen */
3407 snew(erg->xc_norm, erg->rotg->nat);
3410 else
3412 snew(erg->xr_loc, erg->rotg->nat);
3413 snew(erg->x_loc_pbc, erg->rotg->nat);
3416 copy_rvec(erg->rotg->inputVec, erg->vec);
3417 snew(erg->f_rot_loc, erg->rotg->nat);
3419 /* Make space for the calculation of the potential at other angles (used
3420 * for fitting only) */
3421 if (erotgFitPOT == erg->rotg->eFittype)
3423 snew(erg->PotAngleFit, 1);
3424 snew(erg->PotAngleFit->degangle, erg->rotg->PotAngle_nstep);
3425 snew(erg->PotAngleFit->V, erg->rotg->PotAngle_nstep);
3426 snew(erg->PotAngleFit->rotmat, erg->rotg->PotAngle_nstep);
3428 /* Get the set of angles around the reference angle */
3429 start = -0.5 * (erg->rotg->PotAngle_nstep - 1) * erg->rotg->PotAngle_step;
3430 for (int i = 0; i < erg->rotg->PotAngle_nstep; i++)
3432 erg->PotAngleFit->degangle[i] = start + i * erg->rotg->PotAngle_step;
3435 else
3437 erg->PotAngleFit = nullptr;
3440 /* Copy the masses so that the center can be determined. For all types of
3441 * enforced rotation, we store the masses in the erg->mc array. */
3442 snew(erg->mc, erg->rotg->nat);
3443 if (bFlex)
3445 snew(erg->mc_sorted, erg->rotg->nat);
3447 if (!bColl)
3449 snew(erg->m_loc, erg->rotg->nat);
3451 totalmass = 0.0;
3452 int molb = 0;
3453 for (int i = 0; i < erg->rotg->nat; i++)
3455 if (erg->rotg->bMassW)
3457 mass = mtopGetAtomMass(mtop, erg->rotg->ind[i], &molb);
3459 else
3461 mass = 1.0;
3463 erg->mc[i] = mass;
3464 totalmass += mass;
3466 erg->invmass = 1.0 / totalmass;
3468 /* Set xc_ref_center for any rotation potential */
3469 if ((erg->rotg->eType == erotgISO) || (erg->rotg->eType == erotgPM)
3470 || (erg->rotg->eType == erotgRM) || (erg->rotg->eType == erotgRM2))
3472 /* Set the pivot point for the fixed, stationary-axis potentials. This
3473 * won't change during the simulation */
3474 copy_rvec(erg->rotg->pivot, erg->xc_ref_center);
3475 copy_rvec(erg->rotg->pivot, erg->xc_center);
3477 else
3479 /* Center of the reference positions */
3480 get_center(erg->rotg->x_ref, erg->mc, erg->rotg->nat, erg->xc_ref_center);
3482 /* Center of the actual positions */
3483 if (MASTER(cr))
3485 snew(xdum, erg->rotg->nat);
3486 for (int i = 0; i < erg->rotg->nat; i++)
3488 int ii = erg->rotg->ind[i];
3489 copy_rvec(x[ii], xdum[i]);
3491 get_center(xdum, erg->mc, erg->rotg->nat, erg->xc_center);
3492 sfree(xdum);
3494 #if GMX_MPI
3495 if (PAR(cr))
3497 gmx_bcast(sizeof(erg->xc_center), erg->xc_center, cr);
3499 #endif
3502 if (bColl)
3504 /* Save the original (whole) set of positions in xc_old such that at later
3505 * steps the rotation group can always be made whole again. If the simulation is
3506 * restarted, we compute the starting reference positions (given the time)
3507 * and assume that the correct PBC image of each position is the one nearest
3508 * to the current reference */
3509 if (MASTER(cr))
3511 /* Calculate the rotation matrix for this angle: */
3512 t_start = ir->init_t + ir->init_step * ir->delta_t;
3513 erg->degangle = erg->rotg->rate * t_start;
3514 calc_rotmat(erg->vec, erg->degangle, erg->rotmat);
3516 for (int i = 0; i < erg->rotg->nat; i++)
3518 int ii = erg->rotg->ind[i];
3520 /* Subtract pivot, rotate, and add pivot again. This will yield the
3521 * reference position for time t */
3522 rvec_sub(erg->rotg->x_ref[i], erg->xc_ref_center, coord);
3523 mvmul(erg->rotmat, coord, xref);
3524 rvec_inc(xref, erg->xc_ref_center);
3526 copy_correct_pbc_image(x[ii], erg->xc_old[i], xref, box, 3);
3529 #if GMX_MPI
3530 if (PAR(cr))
3532 gmx_bcast(erg->rotg->nat * sizeof(erg->xc_old[0]), erg->xc_old, cr);
3534 #endif
3537 if ((erg->rotg->eType != erotgFLEX) && (erg->rotg->eType != erotgFLEX2))
3539 /* Put the reference positions into origin: */
3540 for (int i = 0; i < erg->rotg->nat; i++)
3542 rvec_dec(erg->rotg->x_ref[i], erg->xc_ref_center);
3546 /* Enforced rotation with flexible axis */
3547 if (bFlex)
3549 /* Calculate maximum beta value from minimum gaussian (performance opt.) */
3550 erg->max_beta = calc_beta_max(erg->rotg->min_gaussian, erg->rotg->slab_dist);
3552 /* Determine the smallest and largest coordinate with respect to the rotation vector */
3553 get_firstlast_atom_ref(erg, &ref_firstindex, &ref_lastindex);
3555 /* From the extreme positions of the reference group, determine the first
3556 * and last slab of the reference. */
3557 get_firstlast_slab_ref(erg, erg->mc, ref_firstindex, ref_lastindex);
3559 /* Allocate memory for the slabs */
3560 allocate_slabs(erg, fplog, bVerbose);
3562 /* Flexible rotation: determine the reference centers for the rest of the simulation */
3563 erg->slab_first = erg->slab_first_ref;
3564 erg->slab_last = erg->slab_last_ref;
3565 get_slab_centers(erg, erg->rotg->x_ref, erg->mc, -1, out_slabs, bOutputCenters, TRUE);
3567 /* Length of each x_rotref vector from center (needed if fit routine NORM is chosen): */
3568 if (erg->rotg->eFittype == erotgFitNORM)
3570 for (int i = 0; i < erg->rotg->nat; i++)
3572 rvec_sub(erg->rotg->x_ref[i], erg->xc_ref_center, coord);
3573 erg->xc_ref_length[i] = norm(coord);
3579 /* Calculate the size of the MPI buffer needed in reduce_output() */
3580 static int calc_mpi_bufsize(const gmx_enfrot* er)
3583 int count_total = 0;
3584 for (int g = 0; g < er->rot->ngrp; g++)
3586 const t_rotgrp* rotg = &er->rot->grp[g];
3587 const gmx_enfrotgrp* erg = &er->enfrotgrp[g];
3589 /* Count the items that are transferred for this group: */
3590 int count_group = 4; /* V, torque, angle, weight */
3592 /* Add the maximum number of slabs for flexible groups */
3593 if (ISFLEX(rotg))
3595 count_group += erg->slab_last_ref - erg->slab_first_ref + 1;
3598 /* Add space for the potentials at different angles: */
3599 if (erotgFitPOT == erg->rotg->eFittype)
3601 count_group += erg->rotg->PotAngle_nstep;
3604 /* Add to the total number: */
3605 count_total += count_group;
3608 return count_total;
3612 std::unique_ptr<gmx::EnforcedRotation> init_rot(FILE* fplog,
3613 t_inputrec* ir,
3614 int nfile,
3615 const t_filenm fnm[],
3616 const t_commrec* cr,
3617 gmx::LocalAtomSetManager* atomSets,
3618 const t_state* globalState,
3619 gmx_mtop_t* mtop,
3620 const gmx_output_env_t* oenv,
3621 const gmx::MdrunOptions& mdrunOptions,
3622 const gmx::StartingBehavior startingBehavior)
3624 int nat_max = 0; /* Size of biggest rotation group */
3625 rvec* x_pbc = nullptr; /* Space for the pbc-correct atom positions */
3627 if (MASTER(cr) && mdrunOptions.verbose)
3629 fprintf(stdout, "%s Initializing ...\n", RotStr);
3632 auto enforcedRotation = std::make_unique<gmx::EnforcedRotation>();
3633 gmx_enfrot* er = enforcedRotation->getLegacyEnfrot();
3634 // TODO When this module implements IMdpOptions, the ownership will become more clear.
3635 er->rot = ir->rot;
3636 er->restartWithAppending = (startingBehavior == gmx::StartingBehavior::RestartWithAppending);
3638 /* When appending, skip first output to avoid duplicate entries in the data files */
3639 if (er->restartWithAppending)
3641 er->bOut = FALSE;
3643 else
3645 er->bOut = TRUE;
3648 if (MASTER(cr) && er->bOut)
3650 please_cite(fplog, "Kutzner2011");
3653 /* Output every step for reruns */
3654 if (mdrunOptions.rerun)
3656 if (nullptr != fplog)
3658 fprintf(fplog, "%s rerun - will write rotation output every available step.\n", RotStr);
3660 er->nstrout = 1;
3661 er->nstsout = 1;
3663 else
3665 er->nstrout = er->rot->nstrout;
3666 er->nstsout = er->rot->nstsout;
3669 er->out_slabs = nullptr;
3670 if (MASTER(cr) && HaveFlexibleGroups(er->rot))
3672 er->out_slabs = open_slab_out(opt2fn("-rs", nfile, fnm), er);
3675 if (MASTER(cr))
3677 /* Remove pbc, make molecule whole.
3678 * When ir->bContinuation=TRUE this has already been done, but ok. */
3679 snew(x_pbc, mtop->natoms);
3680 copy_rvecn(globalState->x.rvec_array(), x_pbc, 0, mtop->natoms);
3681 do_pbc_first_mtop(nullptr, ir->ePBC, globalState->box, mtop, x_pbc);
3682 /* All molecules will be whole now, but not necessarily in the home box.
3683 * Additionally, if a rotation group consists of more than one molecule
3684 * (e.g. two strands of DNA), each one of them can end up in a different
3685 * periodic box. This is taken care of in init_rot_group. */
3688 /* Allocate space for the per-rotation-group data: */
3689 er->enfrotgrp.resize(er->rot->ngrp);
3690 int groupIndex = 0;
3691 for (auto& ergRef : er->enfrotgrp)
3693 gmx_enfrotgrp* erg = &ergRef;
3694 erg->rotg = &er->rot->grp[groupIndex];
3695 erg->atomSet = std::make_unique<gmx::LocalAtomSet>(
3696 atomSets->add({ erg->rotg->ind, erg->rotg->ind + erg->rotg->nat }));
3697 erg->groupIndex = groupIndex;
3699 if (nullptr != fplog)
3701 fprintf(fplog, "%s group %d type '%s'\n", RotStr, groupIndex, erotg_names[erg->rotg->eType]);
3704 if (erg->rotg->nat > 0)
3706 nat_max = std::max(nat_max, erg->rotg->nat);
3708 init_rot_group(fplog, cr, erg, x_pbc, mtop, mdrunOptions.verbose, er->out_slabs,
3709 MASTER(cr) ? globalState->box : nullptr, ir,
3710 !er->restartWithAppending); /* Do not output the reference centers
3711 * again if we are appending */
3713 ++groupIndex;
3716 /* Allocate space for enforced rotation buffer variables */
3717 er->bufsize = nat_max;
3718 snew(er->data, nat_max);
3719 snew(er->xbuf, nat_max);
3720 snew(er->mbuf, nat_max);
3722 /* Buffers for MPI reducing torques, angles, weights (for each group), and V */
3723 if (PAR(cr))
3725 er->mpi_bufsize = calc_mpi_bufsize(er) + 100; /* larger to catch errors */
3726 snew(er->mpi_inbuf, er->mpi_bufsize);
3727 snew(er->mpi_outbuf, er->mpi_bufsize);
3729 else
3731 er->mpi_bufsize = 0;
3732 er->mpi_inbuf = nullptr;
3733 er->mpi_outbuf = nullptr;
3736 /* Only do I/O on the MASTER */
3737 er->out_angles = nullptr;
3738 er->out_rot = nullptr;
3739 er->out_torque = nullptr;
3740 if (MASTER(cr))
3742 er->out_rot = open_rot_out(opt2fn("-ro", nfile, fnm), oenv, er);
3744 if (er->nstsout > 0)
3746 if (HaveFlexibleGroups(er->rot) || HavePotFitGroups(er->rot))
3748 er->out_angles = open_angles_out(opt2fn("-ra", nfile, fnm), er);
3750 if (HaveFlexibleGroups(er->rot))
3752 er->out_torque = open_torque_out(opt2fn("-rt", nfile, fnm), er);
3756 sfree(x_pbc);
3758 return enforcedRotation;
3761 /* Rotate the local reference positions and store them in
3762 * erg->xr_loc[0...(nat_loc-1)]
3764 * Note that we already subtracted u or y_c from the reference positions
3765 * in init_rot_group().
3767 static void rotate_local_reference(gmx_enfrotgrp* erg)
3769 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
3770 for (size_t i = 0; i < erg->atomSet->numAtomsLocal(); i++)
3772 /* Index of this rotation group atom with respect to the whole rotation group */
3773 int ii = collectiveRotationGroupIndex[i];
3774 /* Rotate */
3775 mvmul(erg->rotmat, erg->rotg->x_ref[ii], erg->xr_loc[i]);
3780 /* Select the PBC representation for each local x position and store that
3781 * for later usage. We assume the right PBC image of an x is the one nearest to
3782 * its rotated reference */
3783 static void choose_pbc_image(rvec x[], gmx_enfrotgrp* erg, const matrix box, int npbcdim)
3785 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
3786 for (gmx::index i = 0; i < localRotationGroupIndex.ssize(); i++)
3788 /* Index of a rotation group atom */
3789 int ii = localRotationGroupIndex[i];
3791 /* Get the correctly rotated reference position. The pivot was already
3792 * subtracted in init_rot_group() from the reference positions. Also,
3793 * the reference positions have already been rotated in
3794 * rotate_local_reference(). For the current reference position we thus
3795 * only need to add the pivot again. */
3796 rvec xref;
3797 copy_rvec(erg->xr_loc[i], xref);
3798 rvec_inc(xref, erg->xc_ref_center);
3800 copy_correct_pbc_image(x[ii], erg->x_loc_pbc[i], xref, box, npbcdim);
3805 void do_rotation(const t_commrec* cr, gmx_enfrot* er, const matrix box, rvec x[], real t, int64_t step, gmx_bool bNS)
3807 gmx_bool outstep_slab, outstep_rot;
3808 gmx_bool bColl;
3809 rvec transvec;
3810 gmx_potfit* fit = nullptr; /* For fit type 'potential' determine the fit
3811 angle via the potential minimum */
3813 #ifdef TAKETIME
3814 double t0;
3815 #endif
3817 /* When to output in main rotation output file */
3818 outstep_rot = do_per_step(step, er->nstrout) && er->bOut;
3819 /* When to output per-slab data */
3820 outstep_slab = do_per_step(step, er->nstsout) && er->bOut;
3822 /* Output time into rotation output file */
3823 if (outstep_rot && MASTER(cr))
3825 fprintf(er->out_rot, "%12.3e", t);
3828 /**************************************************************************/
3829 /* First do ALL the communication! */
3830 for (auto& ergRef : er->enfrotgrp)
3832 gmx_enfrotgrp* erg = &ergRef;
3833 const t_rotgrp* rotg = erg->rotg;
3835 /* Do we use a collective (global) set of coordinates? */
3836 bColl = ISCOLL(rotg);
3838 /* Calculate the rotation matrix for this angle: */
3839 erg->degangle = rotg->rate * t;
3840 calc_rotmat(erg->vec, erg->degangle, erg->rotmat);
3842 if (bColl)
3844 /* Transfer the rotation group's positions such that every node has
3845 * all of them. Every node contributes its local positions x and stores
3846 * it in the collective erg->xc array. */
3847 communicate_group_positions(cr, erg->xc, erg->xc_shifts, erg->xc_eshifts, bNS, x,
3848 rotg->nat, erg->atomSet->numAtomsLocal(),
3849 erg->atomSet->localIndex().data(),
3850 erg->atomSet->collectiveIndex().data(), erg->xc_old, box);
3852 else
3854 /* Fill the local masses array;
3855 * this array changes in DD/neighborsearching steps */
3856 if (bNS)
3858 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
3859 for (gmx::index i = 0; i < collectiveRotationGroupIndex.ssize(); i++)
3861 /* Index of local atom w.r.t. the collective rotation group */
3862 int ii = collectiveRotationGroupIndex[i];
3863 erg->m_loc[i] = erg->mc[ii];
3867 /* Calculate Omega*(y_i-y_c) for the local positions */
3868 rotate_local_reference(erg);
3870 /* Choose the nearest PBC images of the group atoms with respect
3871 * to the rotated reference positions */
3872 choose_pbc_image(x, erg, box, 3);
3874 /* Get the center of the rotation group */
3875 if ((rotg->eType == erotgISOPF) || (rotg->eType == erotgPMPF))
3877 get_center_comm(cr, erg->x_loc_pbc, erg->m_loc, erg->atomSet->numAtomsLocal(),
3878 rotg->nat, erg->xc_center);
3882 } /* End of loop over rotation groups */
3884 /**************************************************************************/
3885 /* Done communicating, we can start to count cycles for the load balancing now ... */
3886 if (DOMAINDECOMP(cr))
3888 ddReopenBalanceRegionCpu(cr->dd);
3891 #ifdef TAKETIME
3892 t0 = MPI_Wtime();
3893 #endif
3895 for (auto& ergRef : er->enfrotgrp)
3897 gmx_enfrotgrp* erg = &ergRef;
3898 const t_rotgrp* rotg = erg->rotg;
3900 if (outstep_rot && MASTER(cr))
3902 fprintf(er->out_rot, "%12.4f", erg->degangle);
3905 /* Calculate angles and rotation matrices for potential fitting: */
3906 if ((outstep_rot || outstep_slab) && (erotgFitPOT == rotg->eFittype))
3908 fit = erg->PotAngleFit;
3909 for (int i = 0; i < rotg->PotAngle_nstep; i++)
3911 calc_rotmat(erg->vec, erg->degangle + fit->degangle[i], fit->rotmat[i]);
3913 /* Clear value from last step */
3914 erg->PotAngleFit->V[i] = 0.0;
3918 /* Clear values from last time step */
3919 erg->V = 0.0;
3920 erg->torque_v = 0.0;
3921 erg->angle_v = 0.0;
3922 erg->weight_v = 0.0;
3924 switch (rotg->eType)
3926 case erotgISO:
3927 case erotgISOPF:
3928 case erotgPM:
3929 case erotgPMPF: do_fixed(erg, outstep_rot, outstep_slab); break;
3930 case erotgRM: do_radial_motion(erg, outstep_rot, outstep_slab); break;
3931 case erotgRMPF: do_radial_motion_pf(erg, x, box, outstep_rot, outstep_slab); break;
3932 case erotgRM2:
3933 case erotgRM2PF: do_radial_motion2(erg, x, box, outstep_rot, outstep_slab); break;
3934 case erotgFLEXT:
3935 case erotgFLEX2T:
3936 /* Subtract the center of the rotation group from the collective positions array
3937 * Also store the center in erg->xc_center since it needs to be subtracted
3938 * in the low level routines from the local coordinates as well */
3939 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
3940 svmul(-1.0, erg->xc_center, transvec);
3941 translate_x(erg->xc, rotg->nat, transvec);
3942 do_flexible(MASTER(cr), er, erg, x, box, t, outstep_rot, outstep_slab);
3943 break;
3944 case erotgFLEX:
3945 case erotgFLEX2:
3946 /* Do NOT subtract the center of mass in the low level routines! */
3947 clear_rvec(erg->xc_center);
3948 do_flexible(MASTER(cr), er, erg, x, box, t, outstep_rot, outstep_slab);
3949 break;
3950 default: gmx_fatal(FARGS, "No such rotation potential.");
3954 #ifdef TAKETIME
3955 if (MASTER(cr))
3957 fprintf(stderr, "%s calculation (step %d) took %g seconds.\n", RotStr, step, MPI_Wtime() - t0);
3959 #endif