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,2018,2019, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
39 #include "pull_rotation.h"
50 #include "gromacs/commandline/filenm.h"
51 #include "gromacs/domdec/dlbtiming.h"
52 #include "gromacs/domdec/domdec_struct.h"
53 #include "gromacs/domdec/ga2la.h"
54 #include "gromacs/domdec/localatomset.h"
55 #include "gromacs/domdec/localatomsetmanager.h"
56 #include "gromacs/fileio/gmxfio.h"
57 #include "gromacs/fileio/xvgr.h"
58 #include "gromacs/gmxlib/network.h"
59 #include "gromacs/linearalgebra/nrjac.h"
60 #include "gromacs/math/functions.h"
61 #include "gromacs/math/utilities.h"
62 #include "gromacs/math/vec.h"
63 #include "gromacs/mdlib/groupcoord.h"
64 #include "gromacs/mdlib/stat.h"
65 #include "gromacs/mdrunutility/handlerestart.h"
66 #include "gromacs/mdtypes/commrec.h"
67 #include "gromacs/mdtypes/inputrec.h"
68 #include "gromacs/mdtypes/md_enums.h"
69 #include "gromacs/mdtypes/mdrunoptions.h"
70 #include "gromacs/mdtypes/state.h"
71 #include "gromacs/pbcutil/pbc.h"
72 #include "gromacs/timing/cyclecounter.h"
73 #include "gromacs/timing/wallcycle.h"
74 #include "gromacs/topology/mtop_lookup.h"
75 #include "gromacs/topology/mtop_util.h"
76 #include "gromacs/utility/basedefinitions.h"
77 #include "gromacs/utility/fatalerror.h"
78 #include "gromacs/utility/pleasecite.h"
79 #include "gromacs/utility/smalloc.h"
81 static char const* RotStr = { "Enforced rotation:" };
83 /* Set the minimum weight for the determination of the slab centers */
84 #define WEIGHT_MIN (10 * GMX_FLOAT_MIN)
86 //! Helper structure for sorting positions along rotation vector
87 struct sort_along_vec_t
89 //! Projection of xc on the rotation vector
97 //! Reference position
102 //! Enforced rotation / flexible: determine the angle of each slab
105 //! Number of atoms belonging to this slab
107 /*! \brief The positions belonging to this slab.
109 * In general, this should be all positions of the whole
110 * rotation group, but we leave those away that have a small
113 //! Same for reference
115 //! The weight for each atom
120 //! Helper structure for potential fitting
123 /*! \brief Set of angles for which the potential is calculated.
125 * The optimum fit is determined as the angle for with the
126 * potential is minimal. */
128 //! Potential for the different angles
130 //! Rotation matrix corresponding to the angles
135 //! Enforced rotation data for a single rotation group
138 //! Input parameters for this group
139 const t_rotgrp* rotg = nullptr;
140 //! Index of this group within the set of groups
142 //! Rotation angle in degrees
146 //! The atoms subject to enforced rotation
147 std::unique_ptr<gmx::LocalAtomSet> atomSet;
149 //! The normalized rotation vector
151 //! Rotation potential for this rotation group
153 //! Array to store the forces on the local atoms resulting from enforced rotation potential
156 /* Collective coordinates for the whole rotation group */
157 //! Length of each x_rotref vector after x_rotref has been put into origin
159 //! Center of the rotation group positions, may be mass weighted
161 //! Center of the rotation group reference positions
163 //! Current (collective) positions
165 //! Current (collective) shifts
167 //! Extra shifts since last DD step
169 //! Old (collective) positions
171 //! Normalized form of the current positions
173 //! Reference positions (sorted in the same order as xc when sorted)
175 //! Where is a position found after sorting?
177 //! Collective masses
179 //! Collective masses sorted
181 //! one over the total mass of the rotation group
184 //! Torque in the direction of rotation vector
186 //! Actual angle of the whole rotation group
188 /* Fixed rotation only */
189 //! Weights for angle determination
191 //! Local reference coords, correctly rotated
193 //! Local current coords, correct PBC image
195 //! Masses of the current local atoms
198 /* Flexible rotation only */
199 //! For this many slabs memory is allocated
201 //! Lowermost slab for that the calculation needs to be performed at a given time step
203 //! Uppermost slab ...
205 //! First slab for which ref. center is stored
209 //! Slab buffer region around reference slabs
211 //! First relevant atom for a slab
213 //! Last relevant atom for a slab
215 //! Gaussian-weighted slab center
217 //! Gaussian-weighted slab center for the reference positions
218 rvec* slab_center_ref;
219 //! Sum of gaussian weights in a slab
221 //! Torque T = r x f for each slab. torque_v = m.v = angular momentum in the direction of v
223 //! 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 //! Precalculated gaussians for a single atom
227 //! Tells to which slab each precalculated gaussian belongs
229 //! Inner sum of the flexible2 potential per slab; this is precalculated for optimization reasons
230 rvec* slab_innersumvec;
231 //! Holds atom positions and gaussian weights of atoms belonging to a slab
232 gmx_slabdata* slab_data;
234 /* For potential fits with varying angle: */
235 //! Used for fit type 'potential'
236 gmx_potfit* PotAngleFit;
240 //! Enforced rotation data for all groups
243 //! Input parameters.
244 const t_rot* rot = nullptr;
245 //! Output period for main rotation outfile
247 //! Output period for per-slab data
249 //! Output file for rotation data
250 FILE* out_rot = nullptr;
251 //! Output file for torque data
252 FILE* out_torque = nullptr;
253 //! Output file for slab angles for flexible type
254 FILE* out_angles = nullptr;
255 //! Output file for slab centers
256 FILE* out_slabs = nullptr;
257 //! Allocation size of buf
259 //! Coordinate buffer variable for sorting
260 rvec* xbuf = nullptr;
261 //! Masses buffer variable for sorting
262 real* mbuf = nullptr;
263 //! Buffer variable needed for position sorting
264 sort_along_vec_t* data = nullptr;
266 real* mpi_inbuf = nullptr;
268 real* mpi_outbuf = nullptr;
269 //! Allocation size of in & outbuf
271 //! If true, append output files
272 gmx_bool restartWithAppending = false;
273 //! Used to skip first output when appending to avoid duplicate entries in rotation outfiles
274 gmx_bool bOut = false;
275 //! Stores working data per group
276 std::vector<gmx_enfrotgrp> enfrotgrp;
280 gmx_enfrot::~gmx_enfrot()
284 gmx_fio_fclose(out_rot);
288 gmx_fio_fclose(out_slabs);
292 gmx_fio_fclose(out_angles);
296 gmx_fio_fclose(out_torque);
303 extern template LocalAtomSet LocalAtomSetManager::add<void, void>(ArrayRef<const int> globalAtomIndex);
305 class EnforcedRotation::Impl
308 gmx_enfrot enforcedRotation_;
311 EnforcedRotation::EnforcedRotation() : impl_(new Impl) {}
313 EnforcedRotation::~EnforcedRotation() = default;
315 gmx_enfrot* EnforcedRotation::getLegacyEnfrot()
317 return &impl_->enforcedRotation_;
322 /* Activate output of forces for correctness checks */
323 /* #define PRINT_FORCES */
325 # define PRINT_FORCE_J \
326 fprintf(stderr, "f%d = %15.8f %15.8f %15.8f\n", erg->xc_ref_ind[j], erg->f_rot_loc[j][XX], \
327 erg->f_rot_loc[j][YY], erg->f_rot_loc[j][ZZ]);
328 # define PRINT_POT_TAU \
332 "potential = %15.8f\n" \
333 "torque = %15.8f\n", \
334 erg->V, erg->torque_v); \
337 # define PRINT_FORCE_J
338 # define PRINT_POT_TAU
341 /* Shortcuts for often used queries */
343 (((rg)->eType == erotgFLEX) || ((rg)->eType == erotgFLEXT) || ((rg)->eType == erotgFLEX2) \
344 || ((rg)->eType == erotgFLEX2T))
346 (((rg)->eType == erotgFLEX) || ((rg)->eType == erotgFLEXT) || ((rg)->eType == erotgFLEX2) \
347 || ((rg)->eType == erotgFLEX2T) || ((rg)->eType == erotgRMPF) || ((rg)->eType == erotgRM2PF))
350 /* Does any of the rotation groups use slab decomposition? */
351 static gmx_bool HaveFlexibleGroups(const t_rot* rot)
353 for (int g = 0; g < rot->ngrp; g++)
355 if (ISFLEX(&rot->grp[g]))
365 /* Is for any group the fit angle determined by finding the minimum of the
366 * rotation potential? */
367 static gmx_bool HavePotFitGroups(const t_rot* rot)
369 for (int g = 0; g < rot->ngrp; g++)
371 if (erotgFitPOT == rot->grp[g].eFittype)
381 static double** allocate_square_matrix(int dim)
384 double** mat = nullptr;
388 for (i = 0; i < dim; i++)
397 static void free_square_matrix(double** mat, int dim)
402 for (i = 0; i < dim; i++)
410 /* Return the angle for which the potential is minimal */
411 static real get_fitangle(const gmx_enfrotgrp* erg)
414 real fitangle = -999.9;
415 real pot_min = GMX_FLOAT_MAX;
419 fit = erg->PotAngleFit;
421 for (i = 0; i < erg->rotg->PotAngle_nstep; i++)
423 if (fit->V[i] < pot_min)
426 fitangle = fit->degangle[i];
434 /* Reduce potential angle fit data for this group at this time step? */
435 static inline gmx_bool bPotAngle(const gmx_enfrot* er, const t_rotgrp* rotg, int64_t step)
437 return ((erotgFitPOT == rotg->eFittype)
438 && (do_per_step(step, er->nstsout) || do_per_step(step, er->nstrout)));
441 /* Reduce slab torqe data for this group at this time step? */
442 static inline gmx_bool bSlabTau(const gmx_enfrot* er, const t_rotgrp* rotg, int64_t step)
444 return ((ISFLEX(rotg)) && do_per_step(step, er->nstsout));
447 /* Output rotation energy, torques, etc. for each rotation group */
448 static void reduce_output(const t_commrec* cr, gmx_enfrot* er, real t, int64_t step)
450 int i, islab, nslabs = 0;
451 int count; /* MPI element counter */
455 /* Fill the MPI buffer with stuff to reduce. If items are added for reduction
456 * here, the MPI buffer size has to be enlarged also in calc_mpi_bufsize() */
460 for (auto& ergRef : er->enfrotgrp)
462 gmx_enfrotgrp* erg = &ergRef;
463 const t_rotgrp* rotg = erg->rotg;
464 nslabs = erg->slab_last - erg->slab_first + 1;
465 er->mpi_inbuf[count++] = erg->V;
466 er->mpi_inbuf[count++] = erg->torque_v;
467 er->mpi_inbuf[count++] = erg->angle_v;
468 er->mpi_inbuf[count++] =
469 erg->weight_v; /* weights are not needed for flex types, but this is just a single value */
471 if (bPotAngle(er, rotg, step))
473 for (i = 0; i < rotg->PotAngle_nstep; i++)
475 er->mpi_inbuf[count++] = erg->PotAngleFit->V[i];
478 if (bSlabTau(er, rotg, step))
480 for (i = 0; i < nslabs; i++)
482 er->mpi_inbuf[count++] = erg->slab_torque_v[i];
486 if (count > er->mpi_bufsize)
488 gmx_fatal(FARGS, "%s MPI buffer overflow, please report this error.", RotStr);
492 MPI_Reduce(er->mpi_inbuf, er->mpi_outbuf, count, GMX_MPI_REAL, MPI_SUM, MASTERRANK(cr),
493 cr->mpi_comm_mygroup);
496 /* Copy back the reduced data from the buffer on the master */
500 for (auto& ergRef : er->enfrotgrp)
502 gmx_enfrotgrp* erg = &ergRef;
503 const t_rotgrp* rotg = erg->rotg;
504 nslabs = erg->slab_last - erg->slab_first + 1;
505 erg->V = er->mpi_outbuf[count++];
506 erg->torque_v = er->mpi_outbuf[count++];
507 erg->angle_v = er->mpi_outbuf[count++];
508 erg->weight_v = er->mpi_outbuf[count++];
510 if (bPotAngle(er, rotg, step))
512 for (int i = 0; i < rotg->PotAngle_nstep; i++)
514 erg->PotAngleFit->V[i] = er->mpi_outbuf[count++];
517 if (bSlabTau(er, rotg, step))
519 for (int i = 0; i < nslabs; i++)
521 erg->slab_torque_v[i] = er->mpi_outbuf[count++];
531 /* Angle and torque for each rotation group */
532 for (auto& ergRef : er->enfrotgrp)
534 gmx_enfrotgrp* erg = &ergRef;
535 const t_rotgrp* rotg = erg->rotg;
536 bFlex = ISFLEX(rotg);
538 /* Output to main rotation output file: */
539 if (do_per_step(step, er->nstrout))
541 if (erotgFitPOT == rotg->eFittype)
543 fitangle = get_fitangle(erg);
549 fitangle = erg->angle_v; /* RMSD fit angle */
553 fitangle = (erg->angle_v / erg->weight_v) * 180.0 * M_1_PI;
556 fprintf(er->out_rot, "%12.4f", fitangle);
557 fprintf(er->out_rot, "%12.3e", erg->torque_v);
558 fprintf(er->out_rot, "%12.3e", erg->V);
561 if (do_per_step(step, er->nstsout))
563 /* Output to torque log file: */
566 fprintf(er->out_torque, "%12.3e%6d", t, erg->groupIndex);
567 for (int i = erg->slab_first; i <= erg->slab_last; i++)
569 islab = i - erg->slab_first; /* slab index */
570 /* Only output if enough weight is in slab */
571 if (erg->slab_weights[islab] > rotg->min_gaussian)
573 fprintf(er->out_torque, "%6d%12.3e", i, erg->slab_torque_v[islab]);
576 fprintf(er->out_torque, "\n");
579 /* Output to angles log file: */
580 if (erotgFitPOT == rotg->eFittype)
582 fprintf(er->out_angles, "%12.3e%6d%12.4f", t, erg->groupIndex, erg->degangle);
583 /* Output energies at a set of angles around the reference angle */
584 for (int i = 0; i < rotg->PotAngle_nstep; i++)
586 fprintf(er->out_angles, "%12.3e", erg->PotAngleFit->V[i]);
588 fprintf(er->out_angles, "\n");
592 if (do_per_step(step, er->nstrout))
594 fprintf(er->out_rot, "\n");
600 /* Add the forces from enforced rotation potential to the local forces.
601 * Should be called after the SR forces have been evaluated */
602 real add_rot_forces(gmx_enfrot* er, rvec f[], const t_commrec* cr, int64_t step, real t)
604 real Vrot = 0.0; /* If more than one rotation group is present, Vrot
605 assembles the local parts from all groups */
607 /* Loop over enforced rotation groups (usually 1, though)
608 * Apply the forces from rotation potentials */
609 for (auto& ergRef : er->enfrotgrp)
611 gmx_enfrotgrp* erg = &ergRef;
612 Vrot += erg->V; /* add the local parts from the nodes */
613 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
614 for (gmx::index l = 0; l < localRotationGroupIndex.ssize(); l++)
616 /* Get the right index of the local force */
617 int ii = localRotationGroupIndex[l];
619 rvec_inc(f[ii], erg->f_rot_loc[l]);
623 /* Reduce energy,torque, angles etc. to get the sum values (per rotation group)
624 * on the master and output these values to file. */
625 if ((do_per_step(step, er->nstrout) || do_per_step(step, er->nstsout)) && er->bOut)
627 reduce_output(cr, er, t, step);
630 /* When appending, er->bOut is FALSE the first time to avoid duplicate entries */
639 /* The Gaussian norm is chosen such that the sum of the gaussian functions
640 * over the slabs is approximately 1.0 everywhere */
641 #define GAUSS_NORM 0.569917543430618
644 /* Calculate the maximum beta that leads to a gaussian larger min_gaussian,
645 * also does some checks
647 static double calc_beta_max(real min_gaussian, real slab_dist)
653 /* Actually the next two checks are already made in grompp */
656 gmx_fatal(FARGS, "Slab distance of flexible rotation groups must be >=0 !");
658 if (min_gaussian <= 0)
660 gmx_fatal(FARGS, "Cutoff value for Gaussian must be > 0. (You requested %f)", min_gaussian);
663 /* Define the sigma value */
664 sigma = 0.7 * slab_dist;
666 /* Calculate the argument for the logarithm and check that the log() result is negative or 0 */
667 arg = min_gaussian / GAUSS_NORM;
670 gmx_fatal(FARGS, "min_gaussian of flexible rotation groups must be <%g", GAUSS_NORM);
673 return std::sqrt(-2.0 * sigma * sigma * log(min_gaussian / GAUSS_NORM));
677 static inline real calc_beta(rvec curr_x, const gmx_enfrotgrp* erg, int n)
679 return iprod(curr_x, erg->vec) - erg->rotg->slab_dist * n;
683 static inline real gaussian_weight(rvec curr_x, const gmx_enfrotgrp* erg, int n)
685 const real norm = GAUSS_NORM;
689 /* Define the sigma value */
690 sigma = 0.7 * erg->rotg->slab_dist;
691 /* Calculate the Gaussian value of slab n for position curr_x */
692 return norm * exp(-0.5 * gmx::square(calc_beta(curr_x, erg, n) / sigma));
696 /* Returns the weight in a single slab, also calculates the Gaussian- and mass-
697 * weighted sum of positions for that slab */
698 static real get_slab_weight(int j, const gmx_enfrotgrp* erg, rvec xc[], const real mc[], rvec* x_weighted_sum)
700 rvec curr_x; /* The position of an atom */
701 rvec curr_x_weighted; /* The gaussian-weighted position */
702 real gaussian; /* A single gaussian weight */
703 real wgauss; /* gaussian times current mass */
704 real slabweight = 0.0; /* The sum of weights in the slab */
706 clear_rvec(*x_weighted_sum);
708 /* Loop over all atoms in the rotation group */
709 for (int i = 0; i < erg->rotg->nat; i++)
711 copy_rvec(xc[i], curr_x);
712 gaussian = gaussian_weight(curr_x, erg, j);
713 wgauss = gaussian * mc[i];
714 svmul(wgauss, curr_x, curr_x_weighted);
715 rvec_add(*x_weighted_sum, curr_x_weighted, *x_weighted_sum);
716 slabweight += wgauss;
717 } /* END of loop over rotation group atoms */
723 static void get_slab_centers(gmx_enfrotgrp* erg, /* Enforced rotation group working data */
724 rvec* xc, /* The rotation group positions; will
725 typically be enfrotgrp->xc, but at first call
726 it is enfrotgrp->xc_ref */
727 real* mc, /* The masses of the rotation group atoms */
728 real time, /* Used for output only */
729 FILE* out_slabs, /* For outputting center per slab information */
730 gmx_bool bOutStep, /* Is this an output step? */
731 gmx_bool bReference) /* If this routine is called from
732 init_rot_group we need to store
733 the reference slab centers */
735 /* Loop over slabs */
736 for (int j = erg->slab_first; j <= erg->slab_last; j++)
738 int slabIndex = j - erg->slab_first;
739 erg->slab_weights[slabIndex] = get_slab_weight(j, erg, xc, mc, &erg->slab_center[slabIndex]);
741 /* We can do the calculations ONLY if there is weight in the slab! */
742 if (erg->slab_weights[slabIndex] > WEIGHT_MIN)
744 svmul(1.0 / erg->slab_weights[slabIndex], erg->slab_center[slabIndex],
745 erg->slab_center[slabIndex]);
749 /* We need to check this here, since we divide through slab_weights
750 * in the flexible low-level routines! */
751 gmx_fatal(FARGS, "Not enough weight in slab %d. Slab center cannot be determined!", j);
754 /* At first time step: save the centers of the reference structure */
757 copy_rvec(erg->slab_center[slabIndex], erg->slab_center_ref[slabIndex]);
759 } /* END of loop over slabs */
761 /* Output on the master */
762 if ((nullptr != out_slabs) && bOutStep)
764 fprintf(out_slabs, "%12.3e%6d", time, erg->groupIndex);
765 for (int j = erg->slab_first; j <= erg->slab_last; j++)
767 int slabIndex = j - erg->slab_first;
768 fprintf(out_slabs, "%6d%12.3e%12.3e%12.3e", j, erg->slab_center[slabIndex][XX],
769 erg->slab_center[slabIndex][YY], erg->slab_center[slabIndex][ZZ]);
771 fprintf(out_slabs, "\n");
776 static void calc_rotmat(const rvec vec,
777 real degangle, /* Angle alpha of rotation at time t in degrees */
778 matrix rotmat) /* Rotation matrix */
780 real radangle; /* Rotation angle in radians */
781 real cosa; /* cosine alpha */
782 real sina; /* sine alpha */
783 real OMcosa; /* 1 - cos(alpha) */
784 real dumxy, dumxz, dumyz; /* save computations */
785 rvec rot_vec; /* Rotate around rot_vec ... */
788 radangle = degangle * M_PI / 180.0;
789 copy_rvec(vec, rot_vec);
791 /* Precompute some variables: */
792 cosa = cos(radangle);
793 sina = sin(radangle);
795 dumxy = rot_vec[XX] * rot_vec[YY] * OMcosa;
796 dumxz = rot_vec[XX] * rot_vec[ZZ] * OMcosa;
797 dumyz = rot_vec[YY] * rot_vec[ZZ] * OMcosa;
799 /* Construct the rotation matrix for this rotation group: */
801 rotmat[XX][XX] = cosa + rot_vec[XX] * rot_vec[XX] * OMcosa;
802 rotmat[YY][XX] = dumxy + rot_vec[ZZ] * sina;
803 rotmat[ZZ][XX] = dumxz - rot_vec[YY] * sina;
805 rotmat[XX][YY] = dumxy - rot_vec[ZZ] * sina;
806 rotmat[YY][YY] = cosa + rot_vec[YY] * rot_vec[YY] * OMcosa;
807 rotmat[ZZ][YY] = dumyz + rot_vec[XX] * sina;
809 rotmat[XX][ZZ] = dumxz + rot_vec[YY] * sina;
810 rotmat[YY][ZZ] = dumyz - rot_vec[XX] * sina;
811 rotmat[ZZ][ZZ] = cosa + rot_vec[ZZ] * rot_vec[ZZ] * OMcosa;
816 for (iii = 0; iii < 3; iii++)
818 for (jjj = 0; jjj < 3; jjj++)
820 fprintf(stderr, " %10.8f ", rotmat[iii][jjj]);
822 fprintf(stderr, "\n");
828 /* Calculates torque on the rotation axis tau = position x force */
829 static inline real torque(const rvec rotvec, /* rotation vector; MUST be normalized! */
830 rvec force, /* force */
831 rvec x, /* position of atom on which the force acts */
832 rvec pivot) /* pivot point of rotation axis */
837 /* Subtract offset */
838 rvec_sub(x, pivot, vectmp);
840 /* position x force */
841 cprod(vectmp, force, tau);
843 /* Return the part of the torque which is parallel to the rotation vector */
844 return iprod(tau, rotvec);
848 /* Right-aligned output of value with standard width */
849 static void print_aligned(FILE* fp, char const* str)
851 fprintf(fp, "%12s", str);
855 /* Right-aligned output of value with standard short width */
856 static void print_aligned_short(FILE* fp, char const* str)
858 fprintf(fp, "%6s", str);
862 static FILE* open_output_file(const char* fn, int steps, const char what[])
867 fp = gmx_ffopen(fn, "w");
869 fprintf(fp, "# Output of %s is written in intervals of %d time step%s.\n#\n", what, steps,
870 steps > 1 ? "s" : "");
876 /* Open output file for slab center data. Call on master only */
877 static FILE* open_slab_out(const char* fn, gmx_enfrot* er)
881 if (er->restartWithAppending)
883 fp = gmx_fio_fopen(fn, "a");
887 fp = open_output_file(fn, er->nstsout, "gaussian weighted slab centers");
889 for (auto& ergRef : er->enfrotgrp)
891 gmx_enfrotgrp* erg = &ergRef;
892 if (ISFLEX(erg->rotg))
894 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm, %s.\n", erg->groupIndex,
895 erotg_names[erg->rotg->eType], erg->rotg->slab_dist,
896 erg->rotg->bMassW ? "centers of mass" : "geometrical centers");
900 fprintf(fp, "# Reference centers are listed first (t=-1).\n");
901 fprintf(fp, "# The following columns have the syntax:\n");
903 print_aligned_short(fp, "t");
904 print_aligned_short(fp, "grp");
905 /* Print legend for the first two entries only ... */
906 for (int i = 0; i < 2; i++)
908 print_aligned_short(fp, "slab");
909 print_aligned(fp, "X center");
910 print_aligned(fp, "Y center");
911 print_aligned(fp, "Z center");
913 fprintf(fp, " ...\n");
921 /* Adds 'buf' to 'str' */
922 static void add_to_string(char** str, char* buf)
927 len = strlen(*str) + strlen(buf) + 1;
933 static void add_to_string_aligned(char** str, char* buf)
935 char buf_aligned[STRLEN];
937 sprintf(buf_aligned, "%12s", buf);
938 add_to_string(str, buf_aligned);
942 /* Open output file and print some general information about the rotation groups.
943 * Call on master only */
944 static FILE* open_rot_out(const char* fn, const gmx_output_env_t* oenv, gmx_enfrot* er)
948 const char** setname;
949 char buf[50], buf2[75];
951 char* LegendStr = nullptr;
952 const t_rot* rot = er->rot;
954 if (er->restartWithAppending)
956 fp = gmx_fio_fopen(fn, "a");
960 fp = xvgropen(fn, "Rotation angles and energy", "Time (ps)",
961 "angles (degrees) and energies (kJ/mol)", oenv);
963 "# Output of enforced rotation data is written in intervals of %d time "
965 er->nstrout, er->nstrout > 1 ? "s" : "");
967 "# The scalar tau is the torque (kJ/mol) in the direction of the rotation vector "
969 fprintf(fp, "# To obtain the vectorial torque, multiply tau with the group's rot-vec.\n");
971 "# For flexible groups, tau(t,n) from all slabs n have been summed in a single "
972 "value tau(t) here.\n");
973 fprintf(fp, "# The torques tau(t,n) are found in the rottorque.log (-rt) output file\n");
975 for (int g = 0; g < rot->ngrp; g++)
977 const t_rotgrp* rotg = &rot->grp[g];
978 const gmx_enfrotgrp* erg = &er->enfrotgrp[g];
979 bFlex = ISFLEX(rotg);
982 fprintf(fp, "# ROTATION GROUP %d, potential type '%s':\n", g, erotg_names[rotg->eType]);
983 fprintf(fp, "# rot-massw%d %s\n", g, yesno_names[rotg->bMassW]);
984 fprintf(fp, "# rot-vec%d %12.5e %12.5e %12.5e\n", g, erg->vec[XX],
985 erg->vec[YY], erg->vec[ZZ]);
986 fprintf(fp, "# rot-rate%d %12.5e degrees/ps\n", g, rotg->rate);
987 fprintf(fp, "# rot-k%d %12.5e kJ/(mol*nm^2)\n", g, rotg->k);
988 if (rotg->eType == erotgISO || rotg->eType == erotgPM || rotg->eType == erotgRM
989 || rotg->eType == erotgRM2)
991 fprintf(fp, "# rot-pivot%d %12.5e %12.5e %12.5e nm\n", g, rotg->pivot[XX],
992 rotg->pivot[YY], rotg->pivot[ZZ]);
997 fprintf(fp, "# rot-slab-distance%d %f nm\n", g, rotg->slab_dist);
998 fprintf(fp, "# rot-min-gaussian%d %12.5e\n", g, rotg->min_gaussian);
1001 /* Output the centers of the rotation groups for the pivot-free potentials */
1002 if ((rotg->eType == erotgISOPF) || (rotg->eType == erotgPMPF) || (rotg->eType == erotgRMPF)
1003 || (rotg->eType == erotgRM2PF || (rotg->eType == erotgFLEXT) || (rotg->eType == erotgFLEX2T)))
1005 fprintf(fp, "# ref. grp. %d center %12.5e %12.5e %12.5e\n", g,
1006 erg->xc_ref_center[XX], erg->xc_ref_center[YY], erg->xc_ref_center[ZZ]);
1008 fprintf(fp, "# grp. %d init.center %12.5e %12.5e %12.5e\n", g, erg->xc_center[XX],
1009 erg->xc_center[YY], erg->xc_center[ZZ]);
1012 if ((rotg->eType == erotgRM2) || (rotg->eType == erotgFLEX2) || (rotg->eType == erotgFLEX2T))
1014 fprintf(fp, "# rot-eps%d %12.5e nm^2\n", g, rotg->eps);
1016 if (erotgFitPOT == rotg->eFittype)
1020 "# theta_fit%d is determined by first evaluating the potential for %d "
1021 "angles around theta_ref%d.\n",
1022 g, rotg->PotAngle_nstep, g);
1024 "# The fit angle is the one with the smallest potential. It is given as "
1027 "# from the reference angle, i.e. if theta_ref=X and theta_fit=Y, then the "
1030 "# minimal value of the potential is X+Y. Angular resolution is %g "
1032 rotg->PotAngle_step);
1036 /* Print a nice legend */
1038 LegendStr[0] = '\0';
1039 sprintf(buf, "# %6s", "time");
1040 add_to_string_aligned(&LegendStr, buf);
1043 snew(setname, 4 * rot->ngrp);
1045 for (int g = 0; g < rot->ngrp; g++)
1047 sprintf(buf, "theta_ref%d", g);
1048 add_to_string_aligned(&LegendStr, buf);
1050 sprintf(buf2, "%s (degrees)", buf);
1051 setname[nsets] = gmx_strdup(buf2);
1054 for (int g = 0; g < rot->ngrp; g++)
1056 const t_rotgrp* rotg = &rot->grp[g];
1057 bFlex = ISFLEX(rotg);
1059 /* For flexible axis rotation we use RMSD fitting to determine the
1060 * actual angle of the rotation group */
1061 if (bFlex || erotgFitPOT == rotg->eFittype)
1063 sprintf(buf, "theta_fit%d", g);
1067 sprintf(buf, "theta_av%d", g);
1069 add_to_string_aligned(&LegendStr, buf);
1070 sprintf(buf2, "%s (degrees)", buf);
1071 setname[nsets] = gmx_strdup(buf2);
1074 sprintf(buf, "tau%d", g);
1075 add_to_string_aligned(&LegendStr, buf);
1076 sprintf(buf2, "%s (kJ/mol)", buf);
1077 setname[nsets] = gmx_strdup(buf2);
1080 sprintf(buf, "energy%d", g);
1081 add_to_string_aligned(&LegendStr, buf);
1082 sprintf(buf2, "%s (kJ/mol)", buf);
1083 setname[nsets] = gmx_strdup(buf2);
1090 xvgr_legend(fp, nsets, setname, oenv);
1094 fprintf(fp, "#\n# Legend for the following data columns:\n");
1095 fprintf(fp, "%s\n", LegendStr);
1105 /* Call on master only */
1106 static FILE* open_angles_out(const char* fn, gmx_enfrot* er)
1110 const t_rot* rot = er->rot;
1112 if (er->restartWithAppending)
1114 fp = gmx_fio_fopen(fn, "a");
1118 /* Open output file and write some information about it's structure: */
1119 fp = open_output_file(fn, er->nstsout, "rotation group angles");
1120 fprintf(fp, "# All angles given in degrees, time in ps.\n");
1121 for (int g = 0; g < rot->ngrp; g++)
1123 const t_rotgrp* rotg = &rot->grp[g];
1124 const gmx_enfrotgrp* erg = &er->enfrotgrp[g];
1126 /* Output for this group happens only if potential type is flexible or
1127 * if fit type is potential! */
1128 if (ISFLEX(rotg) || (erotgFitPOT == rotg->eFittype))
1132 sprintf(buf, " slab distance %f nm, ", rotg->slab_dist);
1139 fprintf(fp, "#\n# ROTATION GROUP %d '%s',%s fit type '%s'.\n", g,
1140 erotg_names[rotg->eType], buf, erotg_fitnames[rotg->eFittype]);
1142 /* Special type of fitting using the potential minimum. This is
1143 * done for the whole group only, not for the individual slabs. */
1144 if (erotgFitPOT == rotg->eFittype)
1147 "# To obtain theta_fit%d, the potential is evaluated for %d angles "
1148 "around theta_ref%d\n",
1149 g, rotg->PotAngle_nstep, g);
1151 "# The fit angle in the rotation standard outfile is the one with "
1152 "minimal energy E(theta_fit) [kJ/mol].\n");
1156 fprintf(fp, "# Legend for the group %d data columns:\n", g);
1158 print_aligned_short(fp, "time");
1159 print_aligned_short(fp, "grp");
1160 print_aligned(fp, "theta_ref");
1162 if (erotgFitPOT == rotg->eFittype)
1164 /* Output the set of angles around the reference angle */
1165 for (int i = 0; i < rotg->PotAngle_nstep; i++)
1167 sprintf(buf, "E(%g)", erg->PotAngleFit->degangle[i]);
1168 print_aligned(fp, buf);
1173 /* Output fit angle for each slab */
1174 print_aligned_short(fp, "slab");
1175 print_aligned_short(fp, "atoms");
1176 print_aligned(fp, "theta_fit");
1177 print_aligned_short(fp, "slab");
1178 print_aligned_short(fp, "atoms");
1179 print_aligned(fp, "theta_fit");
1180 fprintf(fp, " ...");
1192 /* Open torque output file and write some information about it's structure.
1193 * Call on master only */
1194 static FILE* open_torque_out(const char* fn, gmx_enfrot* er)
1197 const t_rot* rot = er->rot;
1199 if (er->restartWithAppending)
1201 fp = gmx_fio_fopen(fn, "a");
1205 fp = open_output_file(fn, er->nstsout, "torques");
1207 for (int g = 0; g < rot->ngrp; g++)
1209 const t_rotgrp* rotg = &rot->grp[g];
1210 const gmx_enfrotgrp* erg = &er->enfrotgrp[g];
1213 fprintf(fp, "# Rotation group %d (%s), slab distance %f nm.\n", g,
1214 erotg_names[rotg->eType], rotg->slab_dist);
1216 "# The scalar tau is the torque (kJ/mol) in the direction of the rotation "
1218 fprintf(fp, "# To obtain the vectorial torque, multiply tau with\n");
1219 fprintf(fp, "# rot-vec%d %10.3e %10.3e %10.3e\n", g, erg->vec[XX],
1220 erg->vec[YY], erg->vec[ZZ]);
1224 fprintf(fp, "# Legend for the following data columns: (tau=torque for that slab):\n");
1226 print_aligned_short(fp, "t");
1227 print_aligned_short(fp, "grp");
1228 print_aligned_short(fp, "slab");
1229 print_aligned(fp, "tau");
1230 print_aligned_short(fp, "slab");
1231 print_aligned(fp, "tau");
1232 fprintf(fp, " ...\n");
1240 static void swap_val(double* vec, int i, int j)
1242 double tmp = vec[j];
1250 static void swap_col(double** mat, int i, int j)
1252 double tmp[3] = { mat[0][j], mat[1][j], mat[2][j] };
1255 mat[0][j] = mat[0][i];
1256 mat[1][j] = mat[1][i];
1257 mat[2][j] = mat[2][i];
1265 /* Eigenvectors are stored in columns of eigen_vec */
1266 static void diagonalize_symmetric(double** matrix, double** eigen_vec, double eigenval[3])
1271 jacobi(matrix, 3, eigenval, eigen_vec, &n_rot);
1273 /* sort in ascending order */
1274 if (eigenval[0] > eigenval[1])
1276 swap_val(eigenval, 0, 1);
1277 swap_col(eigen_vec, 0, 1);
1279 if (eigenval[1] > eigenval[2])
1281 swap_val(eigenval, 1, 2);
1282 swap_col(eigen_vec, 1, 2);
1284 if (eigenval[0] > eigenval[1])
1286 swap_val(eigenval, 0, 1);
1287 swap_col(eigen_vec, 0, 1);
1292 static void align_with_z(rvec* s, /* Structure to align */
1297 rvec zet = { 0.0, 0.0, 1.0 };
1298 rvec rot_axis = { 0.0, 0.0, 0.0 };
1299 rvec* rotated_str = nullptr;
1305 snew(rotated_str, natoms);
1307 /* Normalize the axis */
1308 ooanorm = 1.0 / norm(axis);
1309 svmul(ooanorm, axis, axis);
1311 /* Calculate the angle for the fitting procedure */
1312 cprod(axis, zet, rot_axis);
1313 angle = acos(axis[2]);
1319 /* Calculate the rotation matrix */
1320 calc_rotmat(rot_axis, angle * 180.0 / M_PI, rotmat);
1322 /* Apply the rotation matrix to s */
1323 for (i = 0; i < natoms; i++)
1325 for (j = 0; j < 3; j++)
1327 for (k = 0; k < 3; k++)
1329 rotated_str[i][j] += rotmat[j][k] * s[i][k];
1334 /* Rewrite the rotated structure to s */
1335 for (i = 0; i < natoms; i++)
1337 for (j = 0; j < 3; j++)
1339 s[i][j] = rotated_str[i][j];
1347 static void calc_correl_matrix(rvec* Xstr, rvec* Ystr, double** Rmat, int natoms)
1352 for (i = 0; i < 3; i++)
1354 for (j = 0; j < 3; j++)
1360 for (i = 0; i < 3; i++)
1362 for (j = 0; j < 3; j++)
1364 for (k = 0; k < natoms; k++)
1366 Rmat[i][j] += Ystr[k][i] * Xstr[k][j];
1373 static void weigh_coords(rvec* str, real* weight, int natoms)
1378 for (i = 0; i < natoms; i++)
1380 for (j = 0; j < 3; j++)
1382 str[i][j] *= std::sqrt(weight[i]);
1388 static real opt_angle_analytic(rvec* ref_s,
1397 rvec* ref_s_1 = nullptr;
1398 rvec* act_s_1 = nullptr;
1400 double **Rmat, **RtR, **eigvec;
1402 double V[3][3], WS[3][3];
1403 double rot_matrix[3][3];
1407 /* Do not change the original coordinates */
1408 snew(ref_s_1, natoms);
1409 snew(act_s_1, natoms);
1410 for (i = 0; i < natoms; i++)
1412 copy_rvec(ref_s[i], ref_s_1[i]);
1413 copy_rvec(act_s[i], act_s_1[i]);
1416 /* Translate the structures to the origin */
1417 shift[XX] = -ref_com[XX];
1418 shift[YY] = -ref_com[YY];
1419 shift[ZZ] = -ref_com[ZZ];
1420 translate_x(ref_s_1, natoms, shift);
1422 shift[XX] = -act_com[XX];
1423 shift[YY] = -act_com[YY];
1424 shift[ZZ] = -act_com[ZZ];
1425 translate_x(act_s_1, natoms, shift);
1427 /* Align rotation axis with z */
1428 align_with_z(ref_s_1, natoms, axis);
1429 align_with_z(act_s_1, natoms, axis);
1431 /* Correlation matrix */
1432 Rmat = allocate_square_matrix(3);
1434 for (i = 0; i < natoms; i++)
1436 ref_s_1[i][2] = 0.0;
1437 act_s_1[i][2] = 0.0;
1440 /* Weight positions with sqrt(weight) */
1441 if (nullptr != weight)
1443 weigh_coords(ref_s_1, weight, natoms);
1444 weigh_coords(act_s_1, weight, natoms);
1447 /* Calculate correlation matrices R=YXt (X=ref_s; Y=act_s) */
1448 calc_correl_matrix(ref_s_1, act_s_1, Rmat, natoms);
1451 RtR = allocate_square_matrix(3);
1452 for (i = 0; i < 3; i++)
1454 for (j = 0; j < 3; j++)
1456 for (k = 0; k < 3; k++)
1458 RtR[i][j] += Rmat[k][i] * Rmat[k][j];
1462 /* Diagonalize RtR */
1464 for (i = 0; i < 3; i++)
1469 diagonalize_symmetric(RtR, eigvec, eigval);
1470 swap_col(eigvec, 0, 1);
1471 swap_col(eigvec, 1, 2);
1472 swap_val(eigval, 0, 1);
1473 swap_val(eigval, 1, 2);
1476 for (i = 0; i < 3; i++)
1478 for (j = 0; j < 3; j++)
1485 for (i = 0; i < 2; i++)
1487 for (j = 0; j < 2; j++)
1489 WS[i][j] = eigvec[i][j] / std::sqrt(eigval[j]);
1493 for (i = 0; i < 3; i++)
1495 for (j = 0; j < 3; j++)
1497 for (k = 0; k < 3; k++)
1499 V[i][j] += Rmat[i][k] * WS[k][j];
1503 free_square_matrix(Rmat, 3);
1505 /* Calculate optimal rotation matrix */
1506 for (i = 0; i < 3; i++)
1508 for (j = 0; j < 3; j++)
1510 rot_matrix[i][j] = 0.0;
1514 for (i = 0; i < 3; i++)
1516 for (j = 0; j < 3; j++)
1518 for (k = 0; k < 3; k++)
1520 rot_matrix[i][j] += eigvec[i][k] * V[j][k];
1524 rot_matrix[2][2] = 1.0;
1526 /* In some cases abs(rot_matrix[0][0]) can be slighly larger
1527 * than unity due to numerical inacurracies. To be able to calculate
1528 * the acos function, we put these values back in range. */
1529 if (rot_matrix[0][0] > 1.0)
1531 rot_matrix[0][0] = 1.0;
1533 else if (rot_matrix[0][0] < -1.0)
1535 rot_matrix[0][0] = -1.0;
1538 /* Determine the optimal rotation angle: */
1539 opt_angle = (-1.0) * acos(rot_matrix[0][0]) * 180.0 / M_PI;
1540 if (rot_matrix[0][1] < 0.0)
1542 opt_angle = (-1.0) * opt_angle;
1545 /* Give back some memory */
1546 free_square_matrix(RtR, 3);
1549 for (i = 0; i < 3; i++)
1555 return static_cast<real>(opt_angle);
1559 /* Determine angle of the group by RMSD fit to the reference */
1560 /* Not parallelized, call this routine only on the master */
1561 static real flex_fit_angle(gmx_enfrotgrp* erg)
1563 rvec* fitcoords = nullptr;
1564 rvec center; /* Center of positions passed to the fit routine */
1565 real fitangle; /* Angle of the rotation group derived by fitting */
1569 /* Get the center of the rotation group.
1570 * Note, again, erg->xc has been sorted in do_flexible */
1571 get_center(erg->xc, erg->mc_sorted, erg->rotg->nat, center);
1573 /* === Determine the optimal fit angle for the rotation group === */
1574 if (erg->rotg->eFittype == erotgFitNORM)
1576 /* Normalize every position to it's reference length */
1577 for (int i = 0; i < erg->rotg->nat; i++)
1579 /* Put the center of the positions into the origin */
1580 rvec_sub(erg->xc[i], center, coord);
1581 /* Determine the scaling factor for the length: */
1582 scal = erg->xc_ref_length[erg->xc_sortind[i]] / norm(coord);
1583 /* Get position, multiply with the scaling factor and save */
1584 svmul(scal, coord, erg->xc_norm[i]);
1586 fitcoords = erg->xc_norm;
1590 fitcoords = erg->xc;
1592 /* From the point of view of the current positions, the reference has rotated
1593 * backwards. Since we output the angle relative to the fixed reference,
1594 * we need the minus sign. */
1595 fitangle = -opt_angle_analytic(erg->xc_ref_sorted, fitcoords, erg->mc_sorted, erg->rotg->nat,
1596 erg->xc_ref_center, center, erg->vec);
1602 /* Determine actual angle of each slab by RMSD fit to the reference */
1603 /* Not parallelized, call this routine only on the master */
1604 static void flex_fit_angle_perslab(gmx_enfrotgrp* erg, double t, real degangle, FILE* fp)
1607 rvec act_center; /* Center of actual positions that are passed to the fit routine */
1608 rvec ref_center; /* Same for the reference positions */
1609 real fitangle; /* Angle of a slab derived from an RMSD fit to
1610 * the reference structure at t=0 */
1612 real OOm_av; /* 1/average_mass of a rotation group atom */
1613 real m_rel; /* Relative mass of a rotation group atom */
1616 /* Average mass of a rotation group atom: */
1617 OOm_av = erg->invmass * erg->rotg->nat;
1619 /**********************************/
1620 /* First collect the data we need */
1621 /**********************************/
1623 /* Collect the data for the individual slabs */
1624 for (int n = erg->slab_first; n <= erg->slab_last; n++)
1626 int slabIndex = n - erg->slab_first; /* slab index */
1627 sd = &(erg->slab_data[slabIndex]);
1628 sd->nat = erg->lastatom[slabIndex] - erg->firstatom[slabIndex] + 1;
1631 /* Loop over the relevant atoms in the slab */
1632 for (int l = erg->firstatom[slabIndex]; l <= erg->lastatom[slabIndex]; l++)
1634 /* Current position of this atom: x[ii][XX/YY/ZZ] */
1635 copy_rvec(erg->xc[l], curr_x);
1637 /* The (unrotated) reference position of this atom is copied to ref_x.
1638 * Beware, the xc coords have been sorted in do_flexible */
1639 copy_rvec(erg->xc_ref_sorted[l], ref_x);
1641 /* Save data for doing angular RMSD fit later */
1642 /* Save the current atom position */
1643 copy_rvec(curr_x, sd->x[ind]);
1644 /* Save the corresponding reference position */
1645 copy_rvec(ref_x, sd->ref[ind]);
1647 /* Maybe also mass-weighting was requested. If yes, additionally
1648 * multiply the weights with the relative mass of the atom. If not,
1649 * multiply with unity. */
1650 m_rel = erg->mc_sorted[l] * OOm_av;
1652 /* Save the weight for this atom in this slab */
1653 sd->weight[ind] = gaussian_weight(curr_x, erg, n) * m_rel;
1655 /* Next atom in this slab */
1660 /******************************/
1661 /* Now do the fit calculation */
1662 /******************************/
1664 fprintf(fp, "%12.3e%6d%12.3f", t, erg->groupIndex, degangle);
1666 /* === Now do RMSD fitting for each slab === */
1667 /* We require at least SLAB_MIN_ATOMS in a slab, such that the fit makes sense. */
1668 #define SLAB_MIN_ATOMS 4
1670 for (int n = erg->slab_first; n <= erg->slab_last; n++)
1672 int slabIndex = n - erg->slab_first; /* slab index */
1673 sd = &(erg->slab_data[slabIndex]);
1674 if (sd->nat >= SLAB_MIN_ATOMS)
1676 /* Get the center of the slabs reference and current positions */
1677 get_center(sd->ref, sd->weight, sd->nat, ref_center);
1678 get_center(sd->x, sd->weight, sd->nat, act_center);
1679 if (erg->rotg->eFittype == erotgFitNORM)
1681 /* Normalize every position to it's reference length
1682 * prior to performing the fit */
1683 for (int i = 0; i < sd->nat; i++) /* Center */
1685 rvec_dec(sd->ref[i], ref_center);
1686 rvec_dec(sd->x[i], act_center);
1687 /* Normalize x_i such that it gets the same length as ref_i */
1688 svmul(norm(sd->ref[i]) / norm(sd->x[i]), sd->x[i], sd->x[i]);
1690 /* We already subtracted the centers */
1691 clear_rvec(ref_center);
1692 clear_rvec(act_center);
1694 fitangle = -opt_angle_analytic(sd->ref, sd->x, sd->weight, sd->nat, ref_center,
1695 act_center, erg->vec);
1696 fprintf(fp, "%6d%6d%12.3f", n, sd->nat, fitangle);
1701 #undef SLAB_MIN_ATOMS
1705 /* Shift x with is */
1706 static inline void shift_single_coord(const matrix box, rvec x, const ivec is)
1717 x[XX] += tx * box[XX][XX] + ty * box[YY][XX] + tz * box[ZZ][XX];
1718 x[YY] += ty * box[YY][YY] + tz * box[ZZ][YY];
1719 x[ZZ] += tz * box[ZZ][ZZ];
1723 x[XX] += tx * box[XX][XX];
1724 x[YY] += ty * box[YY][YY];
1725 x[ZZ] += tz * box[ZZ][ZZ];
1730 /* Determine the 'home' slab of this atom which is the
1731 * slab with the highest Gaussian weight of all */
1732 static inline int get_homeslab(rvec curr_x, /* The position for which the home slab shall be determined */
1733 const rvec rotvec, /* The rotation vector */
1734 real slabdist) /* The slab distance */
1739 /* The distance of the atom to the coordinate center (where the
1740 * slab with index 0) is */
1741 dist = iprod(rotvec, curr_x);
1743 return gmx::roundToInt(dist / slabdist);
1747 /* For a local atom determine the relevant slabs, i.e. slabs in
1748 * which the gaussian is larger than min_gaussian
1750 static int get_single_atom_gaussians(rvec curr_x, gmx_enfrotgrp* erg)
1753 /* Determine the 'home' slab of this atom: */
1754 int homeslab = get_homeslab(curr_x, erg->vec, erg->rotg->slab_dist);
1756 /* First determine the weight in the atoms home slab: */
1757 real g = gaussian_weight(curr_x, erg, homeslab);
1759 erg->gn_atom[count] = g;
1760 erg->gn_slabind[count] = homeslab;
1764 /* Determine the max slab */
1765 int slab = homeslab;
1766 while (g > erg->rotg->min_gaussian)
1769 g = gaussian_weight(curr_x, erg, slab);
1770 erg->gn_slabind[count] = slab;
1771 erg->gn_atom[count] = g;
1776 /* Determine the min slab */
1781 g = gaussian_weight(curr_x, erg, slab);
1782 erg->gn_slabind[count] = slab;
1783 erg->gn_atom[count] = g;
1785 } while (g > erg->rotg->min_gaussian);
1792 static void flex2_precalc_inner_sum(const gmx_enfrotgrp* erg)
1794 rvec xi; /* positions in the i-sum */
1795 rvec xcn, ycn; /* the current and the reference slab centers */
1798 rvec rin; /* Helper variables */
1801 real OOpsii, OOpsiistar;
1802 real sin_rin; /* s_ii.r_ii */
1803 rvec s_in, tmpvec, tmpvec2;
1804 real mi, wi; /* Mass-weighting of the positions */
1808 N_M = erg->rotg->nat * erg->invmass;
1810 /* Loop over all slabs that contain something */
1811 for (int n = erg->slab_first; n <= erg->slab_last; n++)
1813 int slabIndex = n - erg->slab_first; /* slab index */
1815 /* The current center of this slab is saved in xcn: */
1816 copy_rvec(erg->slab_center[slabIndex], xcn);
1817 /* ... and the reference center in ycn: */
1818 copy_rvec(erg->slab_center_ref[slabIndex + erg->slab_buffer], ycn);
1820 /*** D. Calculate the whole inner sum used for second and third sum */
1821 /* For slab n, we need to loop over all atoms i again. Since we sorted
1822 * the atoms with respect to the rotation vector, we know that it is sufficient
1823 * to calculate from firstatom to lastatom only. All other contributions will
1825 clear_rvec(innersumvec);
1826 for (int i = erg->firstatom[slabIndex]; i <= erg->lastatom[slabIndex]; i++)
1828 /* Coordinate xi of this atom */
1829 copy_rvec(erg->xc[i], xi);
1832 gaussian_xi = gaussian_weight(xi, erg, n);
1833 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1837 copy_rvec(erg->xc_ref_sorted[i], yi0); /* Reference position yi0 */
1838 rvec_sub(yi0, ycn, tmpvec2); /* tmpvec2 = yi0 - ycn */
1839 mvmul(erg->rotmat, tmpvec2, rin); /* rin = Omega.(yi0 - ycn) */
1841 /* Calculate psi_i* and sin */
1842 rvec_sub(xi, xcn, tmpvec2); /* tmpvec2 = xi - xcn */
1844 /* In rare cases, when an atom position coincides with a slab center
1845 * (tmpvec2 == 0) we cannot compute the vector product for s_in.
1846 * However, since the atom is located directly on the pivot, this
1847 * slab's contribution to the force on that atom will be zero
1848 * anyway. Therefore, we continue with the next atom. */
1849 if (gmx_numzero(norm(tmpvec2))) /* 0 == norm(xi - xcn) */
1854 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xi - xcn) */
1855 OOpsiistar = norm2(tmpvec) + erg->rotg->eps; /* OOpsii* = 1/psii* = |v x (xi-xcn)|^2 + eps */
1856 OOpsii = norm(tmpvec); /* OOpsii = 1 / psii = |v x (xi - xcn)| */
1858 /* * v x (xi - xcn) */
1859 unitv(tmpvec, s_in); /* sin = ---------------- */
1860 /* |v x (xi - xcn)| */
1862 sin_rin = iprod(s_in, rin); /* sin_rin = sin . rin */
1864 /* Now the whole sum */
1865 fac = OOpsii / OOpsiistar;
1866 svmul(fac, rin, tmpvec);
1867 fac2 = fac * fac * OOpsii;
1868 svmul(fac2 * sin_rin, s_in, tmpvec2);
1869 rvec_dec(tmpvec, tmpvec2);
1871 svmul(wi * gaussian_xi * sin_rin, tmpvec, tmpvec2);
1873 rvec_inc(innersumvec, tmpvec2);
1874 } /* now we have the inner sum, used both for sum2 and sum3 */
1876 /* Save it to be used in do_flex2_lowlevel */
1877 copy_rvec(innersumvec, erg->slab_innersumvec[slabIndex]);
1878 } /* END of loop over slabs */
1882 static void flex_precalc_inner_sum(const gmx_enfrotgrp* erg)
1884 rvec xi; /* position */
1885 rvec xcn, ycn; /* the current and the reference slab centers */
1886 rvec qin, rin; /* q_i^n and r_i^n */
1889 rvec innersumvec; /* Inner part of sum_n2 */
1890 real gaussian_xi; /* Gaussian weight gn(xi) */
1891 real mi, wi; /* Mass-weighting of the positions */
1894 N_M = erg->rotg->nat * erg->invmass;
1896 /* Loop over all slabs that contain something */
1897 for (int n = erg->slab_first; n <= erg->slab_last; n++)
1899 int slabIndex = n - erg->slab_first; /* slab index */
1901 /* The current center of this slab is saved in xcn: */
1902 copy_rvec(erg->slab_center[slabIndex], xcn);
1903 /* ... and the reference center in ycn: */
1904 copy_rvec(erg->slab_center_ref[slabIndex + erg->slab_buffer], ycn);
1906 /* For slab n, we need to loop over all atoms i again. Since we sorted
1907 * the atoms with respect to the rotation vector, we know that it is sufficient
1908 * to calculate from firstatom to lastatom only. All other contributions will
1910 clear_rvec(innersumvec);
1911 for (int i = erg->firstatom[slabIndex]; i <= erg->lastatom[slabIndex]; i++)
1913 /* Coordinate xi of this atom */
1914 copy_rvec(erg->xc[i], xi);
1917 gaussian_xi = gaussian_weight(xi, erg, n);
1918 mi = erg->mc_sorted[i]; /* need the sorted mass here */
1921 /* Calculate rin and qin */
1922 rvec_sub(erg->xc_ref_sorted[i], ycn, tmpvec); /* tmpvec = yi0-ycn */
1924 /* In rare cases, when an atom position coincides with a slab center
1925 * (tmpvec == 0) we cannot compute the vector product for qin.
1926 * However, since the atom is located directly on the pivot, this
1927 * slab's contribution to the force on that atom will be zero
1928 * anyway. Therefore, we continue with the next atom. */
1929 if (gmx_numzero(norm(tmpvec))) /* 0 == norm(yi0 - ycn) */
1934 mvmul(erg->rotmat, tmpvec, rin); /* rin = Omega.(yi0 - ycn) */
1935 cprod(erg->vec, rin, tmpvec); /* tmpvec = v x Omega*(yi0-ycn) */
1937 /* * v x Omega*(yi0-ycn) */
1938 unitv(tmpvec, qin); /* qin = --------------------- */
1939 /* |v x Omega*(yi0-ycn)| */
1942 rvec_sub(xi, xcn, tmpvec); /* tmpvec = xi-xcn */
1943 bin = iprod(qin, tmpvec); /* bin = qin*(xi-xcn) */
1945 svmul(wi * gaussian_xi * bin, qin, tmpvec);
1947 /* Add this contribution to the inner sum: */
1948 rvec_add(innersumvec, tmpvec, innersumvec);
1949 } /* now we have the inner sum vector S^n for this slab */
1950 /* Save it to be used in do_flex_lowlevel */
1951 copy_rvec(innersumvec, erg->slab_innersumvec[slabIndex]);
1956 static real do_flex2_lowlevel(gmx_enfrotgrp* erg,
1957 real sigma, /* The Gaussian width sigma */
1959 gmx_bool bOutstepRot,
1960 gmx_bool bOutstepSlab,
1963 int count, ii, iigrp;
1964 rvec xj; /* position in the i-sum */
1965 rvec yj0; /* the reference position in the j-sum */
1966 rvec xcn, ycn; /* the current and the reference slab centers */
1967 real V; /* This node's part of the rotation pot. energy */
1968 real gaussian_xj; /* Gaussian weight */
1971 real numerator, fit_numerator;
1972 rvec rjn, fit_rjn; /* Helper variables */
1975 real OOpsij, OOpsijstar;
1976 real OOsigma2; /* 1/(sigma^2) */
1979 rvec sjn, tmpvec, tmpvec2, yj0_ycn;
1980 rvec sum1vec_part, sum1vec, sum2vec_part, sum2vec, sum3vec, sum4vec, innersumvec;
1982 real mj, wj; /* Mass-weighting of the positions */
1984 real Wjn; /* g_n(x_j) m_j / Mjn */
1985 gmx_bool bCalcPotFit;
1987 /* To calculate the torque per slab */
1988 rvec slab_force; /* Single force from slab n on one atom */
1989 rvec slab_sum1vec_part;
1990 real slab_sum3part, slab_sum4part;
1991 rvec slab_sum1vec, slab_sum2vec, slab_sum3vec, slab_sum4vec;
1993 /* Pre-calculate the inner sums, so that we do not have to calculate
1994 * them again for every atom */
1995 flex2_precalc_inner_sum(erg);
1997 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
1999 /********************************************************/
2000 /* Main loop over all local atoms of the rotation group */
2001 /********************************************************/
2002 N_M = erg->rotg->nat * erg->invmass;
2004 OOsigma2 = 1.0 / (sigma * sigma);
2005 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
2006 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2008 for (gmx::index j = 0; j < localRotationGroupIndex.ssize(); j++)
2010 /* Local index of a rotation group atom */
2011 ii = localRotationGroupIndex[j];
2012 /* Position of this atom in the collective array */
2013 iigrp = collectiveRotationGroupIndex[j];
2014 /* Mass-weighting */
2015 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2018 /* Current position of this atom: x[ii][XX/YY/ZZ]
2019 * Note that erg->xc_center contains the center of mass in case the flex2-t
2020 * potential was chosen. For the flex2 potential erg->xc_center must be
2022 rvec_sub(x[ii], erg->xc_center, xj);
2024 /* Shift this atom such that it is near its reference */
2025 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2027 /* Determine the slabs to loop over, i.e. the ones with contributions
2028 * larger than min_gaussian */
2029 count = get_single_atom_gaussians(xj, erg);
2031 clear_rvec(sum1vec_part);
2032 clear_rvec(sum2vec_part);
2035 /* Loop over the relevant slabs for this atom */
2036 for (int ic = 0; ic < count; ic++)
2038 int n = erg->gn_slabind[ic];
2040 /* Get the precomputed Gaussian value of curr_slab for curr_x */
2041 gaussian_xj = erg->gn_atom[ic];
2043 int slabIndex = n - erg->slab_first; /* slab index */
2045 /* The (unrotated) reference position of this atom is copied to yj0: */
2046 copy_rvec(erg->rotg->x_ref[iigrp], yj0);
2048 beta = calc_beta(xj, erg, n);
2050 /* The current center of this slab is saved in xcn: */
2051 copy_rvec(erg->slab_center[slabIndex], xcn);
2052 /* ... and the reference center in ycn: */
2053 copy_rvec(erg->slab_center_ref[slabIndex + erg->slab_buffer], ycn);
2055 rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
2058 mvmul(erg->rotmat, yj0_ycn, rjn); /* rjn = Omega.(yj0 - ycn) */
2060 /* Subtract the slab center from xj */
2061 rvec_sub(xj, xcn, tmpvec2); /* tmpvec2 = xj - xcn */
2063 /* In rare cases, when an atom position coincides with a slab center
2064 * (tmpvec2 == 0) we cannot compute the vector product for sjn.
2065 * However, since the atom is located directly on the pivot, this
2066 * slab's contribution to the force on that atom will be zero
2067 * anyway. Therefore, we directly move on to the next slab. */
2068 if (gmx_numzero(norm(tmpvec2))) /* 0 == norm(xj - xcn) */
2074 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec = v x (xj - xcn) */
2076 OOpsijstar = norm2(tmpvec) + erg->rotg->eps; /* OOpsij* = 1/psij* = |v x (xj-xcn)|^2 + eps */
2078 numerator = gmx::square(iprod(tmpvec, rjn));
2080 /*********************************/
2081 /* Add to the rotation potential */
2082 /*********************************/
2083 V += 0.5 * erg->rotg->k * wj * gaussian_xj * numerator / OOpsijstar;
2085 /* If requested, also calculate the potential for a set of angles
2086 * near the current reference angle */
2089 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2091 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, fit_rjn);
2092 fit_numerator = gmx::square(iprod(tmpvec, fit_rjn));
2093 erg->PotAngleFit->V[ifit] +=
2094 0.5 * erg->rotg->k * wj * gaussian_xj * fit_numerator / OOpsijstar;
2098 /*************************************/
2099 /* Now calculate the force on atom j */
2100 /*************************************/
2102 OOpsij = norm(tmpvec); /* OOpsij = 1 / psij = |v x (xj - xcn)| */
2104 /* * v x (xj - xcn) */
2105 unitv(tmpvec, sjn); /* sjn = ---------------- */
2106 /* |v x (xj - xcn)| */
2108 sjn_rjn = iprod(sjn, rjn); /* sjn_rjn = sjn . rjn */
2111 /*** A. Calculate the first of the four sum terms: ****************/
2112 fac = OOpsij / OOpsijstar;
2113 svmul(fac, rjn, tmpvec);
2114 fac2 = fac * fac * OOpsij;
2115 svmul(fac2 * sjn_rjn, sjn, tmpvec2);
2116 rvec_dec(tmpvec, tmpvec2);
2117 fac2 = wj * gaussian_xj; /* also needed for sum4 */
2118 svmul(fac2 * sjn_rjn, tmpvec, slab_sum1vec_part);
2119 /********************/
2120 /*** Add to sum1: ***/
2121 /********************/
2122 rvec_inc(sum1vec_part, slab_sum1vec_part); /* sum1 still needs to vector multiplied with v */
2124 /*** B. Calculate the forth of the four sum terms: ****************/
2125 betasigpsi = beta * OOsigma2 * OOpsij; /* this is also needed for sum3 */
2126 /********************/
2127 /*** Add to sum4: ***/
2128 /********************/
2129 slab_sum4part = fac2 * betasigpsi * fac * sjn_rjn
2130 * sjn_rjn; /* Note that fac is still valid from above */
2131 sum4 += slab_sum4part;
2133 /*** C. Calculate Wjn for second and third sum */
2134 /* Note that we can safely divide by slab_weights since we check in
2135 * get_slab_centers that it is non-zero. */
2136 Wjn = gaussian_xj * mj / erg->slab_weights[slabIndex];
2138 /* We already have precalculated the inner sum for slab n */
2139 copy_rvec(erg->slab_innersumvec[slabIndex], innersumvec);
2141 /* Weigh the inner sum vector with Wjn */
2142 svmul(Wjn, innersumvec, innersumvec);
2144 /*** E. Calculate the second of the four sum terms: */
2145 /********************/
2146 /*** Add to sum2: ***/
2147 /********************/
2148 rvec_inc(sum2vec_part, innersumvec); /* sum2 still needs to be vector crossproduct'ed with v */
2150 /*** F. Calculate the third of the four sum terms: */
2151 slab_sum3part = betasigpsi * iprod(sjn, innersumvec);
2152 sum3 += slab_sum3part; /* still needs to be multiplied with v */
2154 /*** G. Calculate the torque on the local slab's axis: */
2158 cprod(slab_sum1vec_part, erg->vec, slab_sum1vec);
2160 cprod(innersumvec, erg->vec, slab_sum2vec);
2162 svmul(slab_sum3part, erg->vec, slab_sum3vec);
2164 svmul(slab_sum4part, erg->vec, slab_sum4vec);
2166 /* The force on atom ii from slab n only: */
2167 for (int m = 0; m < DIM; m++)
2169 slab_force[m] = erg->rotg->k
2170 * (-slab_sum1vec[m] + slab_sum2vec[m] - slab_sum3vec[m]
2171 + 0.5 * slab_sum4vec[m]);
2174 erg->slab_torque_v[slabIndex] += torque(erg->vec, slab_force, xj, xcn);
2176 } /* END of loop over slabs */
2178 /* Construct the four individual parts of the vector sum: */
2179 cprod(sum1vec_part, erg->vec, sum1vec); /* sum1vec = { } x v */
2180 cprod(sum2vec_part, erg->vec, sum2vec); /* sum2vec = { } x v */
2181 svmul(sum3, erg->vec, sum3vec); /* sum3vec = { } . v */
2182 svmul(sum4, erg->vec, sum4vec); /* sum4vec = { } . v */
2184 /* Store the additional force so that it can be added to the force
2185 * array after the normal forces have been evaluated */
2186 for (int m = 0; m < DIM; m++)
2188 erg->f_rot_loc[j][m] =
2189 erg->rotg->k * (-sum1vec[m] + sum2vec[m] - sum3vec[m] + 0.5 * sum4vec[m]);
2193 fprintf(stderr, "sum1: %15.8f %15.8f %15.8f\n", -erg->rotg->k * sum1vec[XX],
2194 -erg->rotg->k * sum1vec[YY], -erg->rotg->k * sum1vec[ZZ]);
2195 fprintf(stderr, "sum2: %15.8f %15.8f %15.8f\n", erg->rotg->k * sum2vec[XX],
2196 erg->rotg->k * sum2vec[YY], erg->rotg->k * sum2vec[ZZ]);
2197 fprintf(stderr, "sum3: %15.8f %15.8f %15.8f\n", -erg->rotg->k * sum3vec[XX],
2198 -erg->rotg->k * sum3vec[YY], -erg->rotg->k * sum3vec[ZZ]);
2199 fprintf(stderr, "sum4: %15.8f %15.8f %15.8f\n", 0.5 * erg->rotg->k * sum4vec[XX],
2200 0.5 * erg->rotg->k * sum4vec[YY], 0.5 * erg->rotg->k * sum4vec[ZZ]);
2205 } /* END of loop over local atoms */
2211 static real do_flex_lowlevel(gmx_enfrotgrp* erg,
2212 real sigma, /* The Gaussian width sigma */
2214 gmx_bool bOutstepRot,
2215 gmx_bool bOutstepSlab,
2219 rvec xj, yj0; /* current and reference position */
2220 rvec xcn, ycn; /* the current and the reference slab centers */
2221 rvec yj0_ycn; /* yj0 - ycn */
2222 rvec xj_xcn; /* xj - xcn */
2223 rvec qjn, fit_qjn; /* q_i^n */
2224 rvec sum_n1, sum_n2; /* Two contributions to the rotation force */
2225 rvec innersumvec; /* Inner part of sum_n2 */
2227 rvec force_n; /* Single force from slab n on one atom */
2228 rvec force_n1, force_n2; /* First and second part of force_n */
2229 rvec tmpvec, tmpvec2, tmp_f; /* Helper variables */
2230 real V; /* The rotation potential energy */
2231 real OOsigma2; /* 1/(sigma^2) */
2232 real beta; /* beta_n(xj) */
2233 real bjn, fit_bjn; /* b_j^n */
2234 real gaussian_xj; /* Gaussian weight gn(xj) */
2235 real betan_xj_sigma2;
2236 real mj, wj; /* Mass-weighting of the positions */
2238 gmx_bool bCalcPotFit;
2240 /* Pre-calculate the inner sums, so that we do not have to calculate
2241 * them again for every atom */
2242 flex_precalc_inner_sum(erg);
2244 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
2246 /********************************************************/
2247 /* Main loop over all local atoms of the rotation group */
2248 /********************************************************/
2249 OOsigma2 = 1.0 / (sigma * sigma);
2250 N_M = erg->rotg->nat * erg->invmass;
2252 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
2253 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2255 for (gmx::index j = 0; j < localRotationGroupIndex.ssize(); j++)
2257 /* Local index of a rotation group atom */
2258 int ii = localRotationGroupIndex[j];
2259 /* Position of this atom in the collective array */
2260 iigrp = collectiveRotationGroupIndex[j];
2261 /* Mass-weighting */
2262 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2265 /* Current position of this atom: x[ii][XX/YY/ZZ]
2266 * Note that erg->xc_center contains the center of mass in case the flex-t
2267 * potential was chosen. For the flex potential erg->xc_center must be
2269 rvec_sub(x[ii], erg->xc_center, xj);
2271 /* Shift this atom such that it is near its reference */
2272 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2274 /* Determine the slabs to loop over, i.e. the ones with contributions
2275 * larger than min_gaussian */
2276 count = get_single_atom_gaussians(xj, erg);
2281 /* Loop over the relevant slabs for this atom */
2282 for (int ic = 0; ic < count; ic++)
2284 int n = erg->gn_slabind[ic];
2286 /* Get the precomputed Gaussian for xj in slab n */
2287 gaussian_xj = erg->gn_atom[ic];
2289 int slabIndex = n - erg->slab_first; /* slab index */
2291 /* The (unrotated) reference position of this atom is saved in yj0: */
2292 copy_rvec(erg->rotg->x_ref[iigrp], yj0);
2294 beta = calc_beta(xj, erg, n);
2296 /* The current center of this slab is saved in xcn: */
2297 copy_rvec(erg->slab_center[slabIndex], xcn);
2298 /* ... and the reference center in ycn: */
2299 copy_rvec(erg->slab_center_ref[slabIndex + erg->slab_buffer], ycn);
2301 rvec_sub(yj0, ycn, yj0_ycn); /* yj0_ycn = yj0 - ycn */
2303 /* In rare cases, when an atom position coincides with a reference slab
2304 * center (yj0_ycn == 0) we cannot compute the normal vector qjn.
2305 * However, since the atom is located directly on the pivot, this
2306 * slab's contribution to the force on that atom will be zero
2307 * anyway. Therefore, we directly move on to the next slab. */
2308 if (gmx_numzero(norm(yj0_ycn))) /* 0 == norm(yj0 - ycn) */
2314 mvmul(erg->rotmat, yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2316 /* Subtract the slab center from xj */
2317 rvec_sub(xj, xcn, xj_xcn); /* xj_xcn = xj - xcn */
2320 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2322 /* * v x Omega.(yj0-ycn) */
2323 unitv(tmpvec, qjn); /* qjn = --------------------- */
2324 /* |v x Omega.(yj0-ycn)| */
2326 bjn = iprod(qjn, xj_xcn); /* bjn = qjn * (xj - xcn) */
2328 /*********************************/
2329 /* Add to the rotation potential */
2330 /*********************************/
2331 V += 0.5 * erg->rotg->k * wj * gaussian_xj * gmx::square(bjn);
2333 /* If requested, also calculate the potential for a set of angles
2334 * near the current reference angle */
2337 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2339 /* As above calculate Omega.(yj0-ycn), now for the other angles */
2340 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_ycn, tmpvec2); /* tmpvec2= Omega.(yj0-ycn) */
2341 /* As above calculate qjn */
2342 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec= v x Omega.(yj0-ycn) */
2343 /* * v x Omega.(yj0-ycn) */
2344 unitv(tmpvec, fit_qjn); /* fit_qjn = --------------------- */
2345 /* |v x Omega.(yj0-ycn)| */
2346 fit_bjn = iprod(fit_qjn, xj_xcn); /* fit_bjn = fit_qjn * (xj - xcn) */
2347 /* Add to the rotation potential for this angle */
2348 erg->PotAngleFit->V[ifit] +=
2349 0.5 * erg->rotg->k * wj * gaussian_xj * gmx::square(fit_bjn);
2353 /****************************************************************/
2354 /* sum_n1 will typically be the main contribution to the force: */
2355 /****************************************************************/
2356 betan_xj_sigma2 = beta * OOsigma2; /* beta_n(xj)/sigma^2 */
2358 /* The next lines calculate
2359 * qjn - (bjn*beta(xj)/(2sigma^2))v */
2360 svmul(bjn * 0.5 * betan_xj_sigma2, erg->vec, tmpvec2);
2361 rvec_sub(qjn, tmpvec2, tmpvec);
2363 /* Multiply with gn(xj)*bjn: */
2364 svmul(gaussian_xj * bjn, tmpvec, tmpvec2);
2367 rvec_inc(sum_n1, tmpvec2);
2369 /* We already have precalculated the Sn term for slab n */
2370 copy_rvec(erg->slab_innersumvec[slabIndex], s_n);
2372 svmul(betan_xj_sigma2 * iprod(s_n, xj_xcn), erg->vec, tmpvec); /* tmpvec = ---------- s_n (xj-xcn) */
2375 rvec_sub(s_n, tmpvec, innersumvec);
2377 /* We can safely divide by slab_weights since we check in get_slab_centers
2378 * that it is non-zero. */
2379 svmul(gaussian_xj / erg->slab_weights[slabIndex], innersumvec, innersumvec);
2381 rvec_add(sum_n2, innersumvec, sum_n2);
2383 /* Calculate the torque: */
2386 /* The force on atom ii from slab n only: */
2387 svmul(-erg->rotg->k * wj, tmpvec2, force_n1); /* part 1 */
2388 svmul(erg->rotg->k * mj, innersumvec, force_n2); /* part 2 */
2389 rvec_add(force_n1, force_n2, force_n);
2390 erg->slab_torque_v[slabIndex] += torque(erg->vec, force_n, xj, xcn);
2392 } /* END of loop over slabs */
2394 /* Put both contributions together: */
2395 svmul(wj, sum_n1, sum_n1);
2396 svmul(mj, sum_n2, sum_n2);
2397 rvec_sub(sum_n2, sum_n1, tmp_f); /* F = -grad V */
2399 /* Store the additional force so that it can be added to the force
2400 * array after the normal forces have been evaluated */
2401 for (int m = 0; m < DIM; m++)
2403 erg->f_rot_loc[j][m] = erg->rotg->k * tmp_f[m];
2408 } /* END of loop over local atoms */
2413 static void sort_collective_coordinates(gmx_enfrotgrp* erg,
2414 sort_along_vec_t* data) /* Buffer for sorting the positions */
2416 /* The projection of the position vector on the rotation vector is
2417 * the relevant value for sorting. Fill the 'data' structure */
2418 for (int i = 0; i < erg->rotg->nat; i++)
2420 data[i].xcproj = iprod(erg->xc[i], erg->vec); /* sort criterium */
2421 data[i].m = erg->mc[i];
2423 copy_rvec(erg->xc[i], data[i].x);
2424 copy_rvec(erg->rotg->x_ref[i], data[i].x_ref);
2426 /* Sort the 'data' structure */
2427 std::sort(data, data + erg->rotg->nat, [](const sort_along_vec_t& a, const sort_along_vec_t& b) {
2428 return a.xcproj < b.xcproj;
2431 /* Copy back the sorted values */
2432 for (int i = 0; i < erg->rotg->nat; i++)
2434 copy_rvec(data[i].x, erg->xc[i]);
2435 copy_rvec(data[i].x_ref, erg->xc_ref_sorted[i]);
2436 erg->mc_sorted[i] = data[i].m;
2437 erg->xc_sortind[i] = data[i].ind;
2442 /* For each slab, get the first and the last index of the sorted atom
2444 static void get_firstlast_atom_per_slab(const gmx_enfrotgrp* erg)
2448 /* Find the first atom that needs to enter the calculation for each slab */
2449 int n = erg->slab_first; /* slab */
2450 int i = 0; /* start with the first atom */
2453 /* Find the first atom that significantly contributes to this slab */
2454 do /* move forward in position until a large enough beta is found */
2456 beta = calc_beta(erg->xc[i], erg, n);
2458 } while ((beta < -erg->max_beta) && (i < erg->rotg->nat));
2460 int slabIndex = n - erg->slab_first; /* slab index */
2461 erg->firstatom[slabIndex] = i;
2462 /* Proceed to the next slab */
2464 } while (n <= erg->slab_last);
2466 /* Find the last atom for each slab */
2467 n = erg->slab_last; /* start with last slab */
2468 i = erg->rotg->nat - 1; /* start with the last atom */
2471 do /* move backward in position until a large enough beta is found */
2473 beta = calc_beta(erg->xc[i], erg, n);
2475 } while ((beta > erg->max_beta) && (i > -1));
2477 int slabIndex = n - erg->slab_first; /* slab index */
2478 erg->lastatom[slabIndex] = i;
2479 /* Proceed to the next slab */
2481 } while (n >= erg->slab_first);
2485 /* Determine the very first and very last slab that needs to be considered
2486 * For the first slab that needs to be considered, we have to find the smallest
2489 * x_first * v - n*Delta_x <= beta_max
2491 * slab index n, slab distance Delta_x, rotation vector v. For the last slab we
2492 * have to find the largest n that obeys
2494 * x_last * v - n*Delta_x >= -beta_max
2497 static inline int get_first_slab(const gmx_enfrotgrp* erg,
2498 rvec firstatom) /* First atom after sorting along the rotation vector v */
2500 /* Find the first slab for the first atom */
2501 return static_cast<int>(ceil(
2502 static_cast<double>((iprod(firstatom, erg->vec) - erg->max_beta) / erg->rotg->slab_dist)));
2506 static inline int get_last_slab(const gmx_enfrotgrp* erg, rvec lastatom) /* Last atom along v */
2508 /* Find the last slab for the last atom */
2509 return static_cast<int>(floor(
2510 static_cast<double>((iprod(lastatom, erg->vec) + erg->max_beta) / erg->rotg->slab_dist)));
2514 static void get_firstlast_slab_check(gmx_enfrotgrp* erg, /* The rotation group (data only accessible in this file) */
2515 rvec firstatom, /* First atom after sorting along the rotation vector v */
2516 rvec lastatom) /* Last atom along v */
2518 erg->slab_first = get_first_slab(erg, firstatom);
2519 erg->slab_last = get_last_slab(erg, lastatom);
2521 /* Calculate the slab buffer size, which changes when slab_first changes */
2522 erg->slab_buffer = erg->slab_first - erg->slab_first_ref;
2524 /* Check whether we have reference data to compare against */
2525 if (erg->slab_first < erg->slab_first_ref)
2527 gmx_fatal(FARGS, "%s No reference data for first slab (n=%d), unable to proceed.", RotStr,
2531 /* Check whether we have reference data to compare against */
2532 if (erg->slab_last > erg->slab_last_ref)
2534 gmx_fatal(FARGS, "%s No reference data for last slab (n=%d), unable to proceed.", RotStr,
2540 /* Enforced rotation with a flexible axis */
2541 static void do_flexible(gmx_bool bMaster,
2542 gmx_enfrot* enfrot, /* Other rotation data */
2544 rvec x[], /* The local positions */
2546 double t, /* Time in picoseconds */
2547 gmx_bool bOutstepRot, /* Output to main rotation output file */
2548 gmx_bool bOutstepSlab) /* Output per-slab data */
2551 real sigma; /* The Gaussian width sigma */
2553 /* Define the sigma value */
2554 sigma = 0.7 * erg->rotg->slab_dist;
2556 /* Sort the collective coordinates erg->xc along the rotation vector. This is
2557 * an optimization for the inner loop. */
2558 sort_collective_coordinates(erg, enfrot->data);
2560 /* Determine the first relevant slab for the first atom and the last
2561 * relevant slab for the last atom */
2562 get_firstlast_slab_check(erg, erg->xc[0], erg->xc[erg->rotg->nat - 1]);
2564 /* Determine for each slab depending on the min_gaussian cutoff criterium,
2565 * a first and a last atom index inbetween stuff needs to be calculated */
2566 get_firstlast_atom_per_slab(erg);
2568 /* Determine the gaussian-weighted center of positions for all slabs */
2569 get_slab_centers(erg, erg->xc, erg->mc_sorted, t, enfrot->out_slabs, bOutstepSlab, FALSE);
2571 /* Clear the torque per slab from last time step: */
2572 nslabs = erg->slab_last - erg->slab_first + 1;
2573 for (int l = 0; l < nslabs; l++)
2575 erg->slab_torque_v[l] = 0.0;
2578 /* Call the rotational forces kernel */
2579 if (erg->rotg->eType == erotgFLEX || erg->rotg->eType == erotgFLEXT)
2581 erg->V = do_flex_lowlevel(erg, sigma, x, bOutstepRot, bOutstepSlab, box);
2583 else if (erg->rotg->eType == erotgFLEX2 || erg->rotg->eType == erotgFLEX2T)
2585 erg->V = do_flex2_lowlevel(erg, sigma, x, bOutstepRot, bOutstepSlab, box);
2589 gmx_fatal(FARGS, "Unknown flexible rotation type");
2592 /* Determine angle by RMSD fit to the reference - Let's hope this */
2593 /* only happens once in a while, since this is not parallelized! */
2594 if (bMaster && (erotgFitPOT != erg->rotg->eFittype))
2598 /* Fit angle of the whole rotation group */
2599 erg->angle_v = flex_fit_angle(erg);
2603 /* Fit angle of each slab */
2604 flex_fit_angle_perslab(erg, t, erg->degangle, enfrot->out_angles);
2608 /* Lump together the torques from all slabs: */
2609 erg->torque_v = 0.0;
2610 for (int l = 0; l < nslabs; l++)
2612 erg->torque_v += erg->slab_torque_v[l];
2617 /* Calculate the angle between reference and actual rotation group atom,
2618 * both projected into a plane perpendicular to the rotation vector: */
2619 static void angle(const gmx_enfrotgrp* erg,
2623 real* weight) /* atoms near the rotation axis should count less than atoms far away */
2625 rvec xp, xrp; /* current and reference positions projected on a plane perpendicular to pg->vec */
2629 /* Project x_ref and x into a plane through the origin perpendicular to rot_vec: */
2630 /* Project x_ref: xrp = x_ref - (vec * x_ref) * vec */
2631 svmul(iprod(erg->vec, x_ref), erg->vec, dum);
2632 rvec_sub(x_ref, dum, xrp);
2633 /* Project x_act: */
2634 svmul(iprod(erg->vec, x_act), erg->vec, dum);
2635 rvec_sub(x_act, dum, xp);
2637 /* Retrieve information about which vector precedes. gmx_angle always
2638 * returns a positive angle. */
2639 cprod(xp, xrp, dum); /* if reference precedes, this is pointing into the same direction as vec */
2641 if (iprod(erg->vec, dum) >= 0)
2643 *alpha = -gmx_angle(xrp, xp);
2647 *alpha = +gmx_angle(xrp, xp);
2650 /* Also return the weight */
2655 /* Project first vector onto a plane perpendicular to the second vector
2657 * Note that v must be of unit length.
2659 static inline void project_onto_plane(rvec dr, const rvec v)
2664 svmul(iprod(dr, v), v, tmp); /* tmp = (dr.v)v */
2665 rvec_dec(dr, tmp); /* dr = dr - (dr.v)v */
2669 /* Fixed rotation: The rotation reference group rotates around the v axis. */
2670 /* The atoms of the actual rotation group are attached with imaginary */
2671 /* springs to the reference atoms. */
2672 static void do_fixed(gmx_enfrotgrp* erg,
2673 gmx_bool bOutstepRot, /* Output to main rotation output file */
2674 gmx_bool bOutstepSlab) /* Output per-slab data */
2677 rvec tmp_f; /* Force */
2678 real alpha; /* a single angle between an actual and a reference position */
2679 real weight; /* single weight for a single angle */
2680 rvec xi_xc; /* xi - xc */
2681 gmx_bool bCalcPotFit;
2684 /* for mass weighting: */
2685 real wi; /* Mass-weighting of the positions */
2687 real k_wi; /* k times wi */
2691 bProject = (erg->rotg->eType == erotgPM) || (erg->rotg->eType == erotgPMPF);
2692 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
2694 N_M = erg->rotg->nat * erg->invmass;
2695 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2696 /* Each process calculates the forces on its local atoms */
2697 for (size_t j = 0; j < erg->atomSet->numAtomsLocal(); j++)
2699 /* Calculate (x_i-x_c) resp. (x_i-u) */
2700 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xi_xc);
2702 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2703 rvec_sub(erg->xr_loc[j], xi_xc, dr);
2707 project_onto_plane(dr, erg->vec);
2710 /* Mass-weighting */
2711 wi = N_M * erg->m_loc[j];
2713 /* Store the additional force so that it can be added to the force
2714 * array after the normal forces have been evaluated */
2715 k_wi = erg->rotg->k * wi;
2716 for (int m = 0; m < DIM; m++)
2718 tmp_f[m] = k_wi * dr[m];
2719 erg->f_rot_loc[j][m] = tmp_f[m];
2720 erg->V += 0.5 * k_wi * gmx::square(dr[m]);
2723 /* If requested, also calculate the potential for a set of angles
2724 * near the current reference angle */
2727 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2729 /* Index of this rotation group atom with respect to the whole rotation group */
2730 int jj = collectiveRotationGroupIndex[j];
2732 /* Rotate with the alternative angle. Like rotate_local_reference(),
2733 * just for a single local atom */
2734 mvmul(erg->PotAngleFit->rotmat[ifit], erg->rotg->x_ref[jj],
2735 fit_xr_loc); /* fit_xr_loc = Omega*(y_i-y_c) */
2737 /* Calculate Omega*(y_i-y_c)-(x_i-x_c) */
2738 rvec_sub(fit_xr_loc, xi_xc, dr);
2742 project_onto_plane(dr, erg->vec);
2745 /* Add to the rotation potential for this angle: */
2746 erg->PotAngleFit->V[ifit] += 0.5 * k_wi * norm2(dr);
2752 /* Add to the torque of this rotation group */
2753 erg->torque_v += torque(erg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2755 /* Calculate the angle between reference and actual rotation group atom. */
2756 angle(erg, xi_xc, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2757 erg->angle_v += alpha * weight;
2758 erg->weight_v += weight;
2760 /* If you want enforced rotation to contribute to the virial,
2761 * activate the following lines:
2764 Add the rotation contribution to the virial
2765 for(j=0; j<DIM; j++)
2767 vir[j][m] += 0.5*f[ii][j]*dr[m];
2773 } /* end of loop over local rotation group atoms */
2777 /* Calculate the radial motion potential and forces */
2778 static void do_radial_motion(gmx_enfrotgrp* erg,
2779 gmx_bool bOutstepRot, /* Output to main rotation output file */
2780 gmx_bool bOutstepSlab) /* Output per-slab data */
2782 rvec tmp_f; /* Force */
2783 real alpha; /* a single angle between an actual and a reference position */
2784 real weight; /* single weight for a single angle */
2785 rvec xj_u; /* xj - u */
2786 rvec tmpvec, fit_tmpvec;
2787 real fac, fac2, sum = 0.0;
2789 gmx_bool bCalcPotFit;
2791 /* For mass weighting: */
2792 real wj; /* Mass-weighting of the positions */
2795 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
2797 N_M = erg->rotg->nat * erg->invmass;
2798 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2799 /* Each process calculates the forces on its local atoms */
2800 for (size_t j = 0; j < erg->atomSet->numAtomsLocal(); j++)
2802 /* Calculate (xj-u) */
2803 rvec_sub(erg->x_loc_pbc[j], erg->xc_center, xj_u); /* xj_u = xj-u */
2805 /* Calculate Omega.(yj0-u) */
2806 cprod(erg->vec, erg->xr_loc[j], tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2808 /* * v x Omega.(yj0-u) */
2809 unitv(tmpvec, pj); /* pj = --------------------- */
2810 /* | v x Omega.(yj0-u) | */
2812 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2815 /* Mass-weighting */
2816 wj = N_M * erg->m_loc[j];
2818 /* Store the additional force so that it can be added to the force
2819 * array after the normal forces have been evaluated */
2820 svmul(-erg->rotg->k * wj * fac, pj, tmp_f);
2821 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2824 /* If requested, also calculate the potential for a set of angles
2825 * near the current reference angle */
2828 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2830 /* Index of this rotation group atom with respect to the whole rotation group */
2831 int jj = collectiveRotationGroupIndex[j];
2833 /* Rotate with the alternative angle. Like rotate_local_reference(),
2834 * just for a single local atom */
2835 mvmul(erg->PotAngleFit->rotmat[ifit], erg->rotg->x_ref[jj],
2836 fit_tmpvec); /* fit_tmpvec = Omega*(yj0-u) */
2838 /* Calculate Omega.(yj0-u) */
2839 cprod(erg->vec, fit_tmpvec, tmpvec); /* tmpvec = v x Omega.(yj0-u) */
2840 /* * v x Omega.(yj0-u) */
2841 unitv(tmpvec, pj); /* pj = --------------------- */
2842 /* | v x Omega.(yj0-u) | */
2844 fac = iprod(pj, xj_u); /* fac = pj.(xj-u) */
2847 /* Add to the rotation potential for this angle: */
2848 erg->PotAngleFit->V[ifit] += 0.5 * erg->rotg->k * wj * fac2;
2854 /* Add to the torque of this rotation group */
2855 erg->torque_v += torque(erg->vec, tmp_f, erg->x_loc_pbc[j], erg->xc_center);
2857 /* Calculate the angle between reference and actual rotation group atom. */
2858 angle(erg, xj_u, erg->xr_loc[j], &alpha, &weight); /* angle in rad, weighted */
2859 erg->angle_v += alpha * weight;
2860 erg->weight_v += weight;
2865 } /* end of loop over local rotation group atoms */
2866 erg->V = 0.5 * erg->rotg->k * sum;
2870 /* Calculate the radial motion pivot-free potential and forces */
2871 static void do_radial_motion_pf(gmx_enfrotgrp* erg,
2872 rvec x[], /* The positions */
2873 const matrix box, /* The simulation box */
2874 gmx_bool bOutstepRot, /* Output to main rotation output file */
2875 gmx_bool bOutstepSlab) /* Output per-slab data */
2877 rvec xj; /* Current position */
2878 rvec xj_xc; /* xj - xc */
2879 rvec yj0_yc0; /* yj0 - yc0 */
2880 rvec tmp_f; /* Force */
2881 real alpha; /* a single angle between an actual and a reference position */
2882 real weight; /* single weight for a single angle */
2883 rvec tmpvec, tmpvec2;
2884 rvec innersumvec; /* Precalculation of the inner sum */
2886 real fac, fac2, V = 0.0;
2888 gmx_bool bCalcPotFit;
2890 /* For mass weighting: */
2891 real mj, wi, wj; /* Mass-weighting of the positions */
2894 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
2896 N_M = erg->rotg->nat * erg->invmass;
2898 /* Get the current center of the rotation group: */
2899 get_center(erg->xc, erg->mc, erg->rotg->nat, erg->xc_center);
2901 /* Precalculate Sum_i [ wi qi.(xi-xc) qi ] which is needed for every single j */
2902 clear_rvec(innersumvec);
2903 for (int i = 0; i < erg->rotg->nat; i++)
2905 /* Mass-weighting */
2906 wi = N_M * erg->mc[i];
2908 /* Calculate qi. Note that xc_ref_center has already been subtracted from
2909 * x_ref in init_rot_group.*/
2910 mvmul(erg->rotmat, erg->rotg->x_ref[i], tmpvec); /* tmpvec = Omega.(yi0-yc0) */
2912 cprod(erg->vec, tmpvec, tmpvec2); /* tmpvec2 = v x Omega.(yi0-yc0) */
2914 /* * v x Omega.(yi0-yc0) */
2915 unitv(tmpvec2, qi); /* qi = ----------------------- */
2916 /* | v x Omega.(yi0-yc0) | */
2918 rvec_sub(erg->xc[i], erg->xc_center, tmpvec); /* tmpvec = xi-xc */
2920 svmul(wi * iprod(qi, tmpvec), qi, tmpvec2);
2922 rvec_inc(innersumvec, tmpvec2);
2924 svmul(erg->rotg->k * erg->invmass, innersumvec, innersumveckM);
2926 /* Each process calculates the forces on its local atoms */
2927 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
2928 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
2929 for (gmx::index j = 0; j < localRotationGroupIndex.ssize(); j++)
2931 /* Local index of a rotation group atom */
2932 int ii = localRotationGroupIndex[j];
2933 /* Position of this atom in the collective array */
2934 int iigrp = collectiveRotationGroupIndex[j];
2935 /* Mass-weighting */
2936 mj = erg->mc[iigrp]; /* need the unsorted mass here */
2939 /* Current position of this atom: x[ii][XX/YY/ZZ] */
2940 copy_rvec(x[ii], xj);
2942 /* Shift this atom such that it is near its reference */
2943 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
2945 /* The (unrotated) reference position is yj0. yc0 has already
2946 * been subtracted in init_rot_group */
2947 copy_rvec(erg->rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
2949 /* Calculate Omega.(yj0-yc0) */
2950 mvmul(erg->rotmat, yj0_yc0, tmpvec2); /* tmpvec2 = Omega.(yj0 - yc0) */
2952 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
2954 /* * v x Omega.(yj0-yc0) */
2955 unitv(tmpvec, qj); /* qj = ----------------------- */
2956 /* | v x Omega.(yj0-yc0) | */
2958 /* Calculate (xj-xc) */
2959 rvec_sub(xj, erg->xc_center, xj_xc); /* xj_xc = xj-xc */
2961 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
2964 /* Store the additional force so that it can be added to the force
2965 * array after the normal forces have been evaluated */
2966 svmul(-erg->rotg->k * wj * fac, qj, tmp_f); /* part 1 of force */
2967 svmul(mj, innersumveckM, tmpvec); /* part 2 of force */
2968 rvec_inc(tmp_f, tmpvec);
2969 copy_rvec(tmp_f, erg->f_rot_loc[j]);
2972 /* If requested, also calculate the potential for a set of angles
2973 * near the current reference angle */
2976 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
2978 /* Rotate with the alternative angle. Like rotate_local_reference(),
2979 * just for a single local atom */
2980 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, tmpvec2); /* tmpvec2 = Omega*(yj0-yc0) */
2982 /* Calculate Omega.(yj0-u) */
2983 cprod(erg->vec, tmpvec2, tmpvec); /* tmpvec = v x Omega.(yj0-yc0) */
2984 /* * v x Omega.(yj0-yc0) */
2985 unitv(tmpvec, qj); /* qj = ----------------------- */
2986 /* | v x Omega.(yj0-yc0) | */
2988 fac = iprod(qj, xj_xc); /* fac = qj.(xj-xc) */
2991 /* Add to the rotation potential for this angle: */
2992 erg->PotAngleFit->V[ifit] += 0.5 * erg->rotg->k * wj * fac2;
2998 /* Add to the torque of this rotation group */
2999 erg->torque_v += torque(erg->vec, tmp_f, xj, erg->xc_center);
3001 /* Calculate the angle between reference and actual rotation group atom. */
3002 angle(erg, xj_xc, yj0_yc0, &alpha, &weight); /* angle in rad, weighted */
3003 erg->angle_v += alpha * weight;
3004 erg->weight_v += weight;
3009 } /* end of loop over local rotation group atoms */
3010 erg->V = 0.5 * erg->rotg->k * V;
3014 /* Precalculate the inner sum for the radial motion 2 forces */
3015 static void radial_motion2_precalc_inner_sum(const gmx_enfrotgrp* erg, rvec innersumvec)
3018 rvec xi_xc; /* xj - xc */
3019 rvec tmpvec, tmpvec2;
3023 rvec v_xi_xc; /* v x (xj - u) */
3024 real psii, psiistar;
3025 real wi; /* Mass-weighting of the positions */
3029 N_M = erg->rotg->nat * erg->invmass;
3031 /* Loop over the collective set of positions */
3033 for (i = 0; i < erg->rotg->nat; i++)
3035 /* Mass-weighting */
3036 wi = N_M * erg->mc[i];
3038 rvec_sub(erg->xc[i], erg->xc_center, xi_xc); /* xi_xc = xi-xc */
3040 /* Calculate ri. Note that xc_ref_center has already been subtracted from
3041 * x_ref in init_rot_group.*/
3042 mvmul(erg->rotmat, erg->rotg->x_ref[i], ri); /* ri = Omega.(yi0-yc0) */
3044 cprod(erg->vec, xi_xc, v_xi_xc); /* v_xi_xc = v x (xi-u) */
3046 fac = norm2(v_xi_xc);
3048 psiistar = 1.0 / (fac + erg->rotg->eps); /* psiistar = --------------------- */
3049 /* |v x (xi-xc)|^2 + eps */
3051 psii = gmx::invsqrt(fac); /* 1 */
3052 /* psii = ------------- */
3055 svmul(psii, v_xi_xc, si); /* si = psii * (v x (xi-xc) ) */
3057 siri = iprod(si, ri); /* siri = si.ri */
3059 svmul(psiistar / psii, ri, tmpvec);
3060 svmul(psiistar * psiistar / (psii * psii * psii) * siri, si, tmpvec2);
3061 rvec_dec(tmpvec, tmpvec2);
3062 cprod(tmpvec, erg->vec, tmpvec2);
3064 svmul(wi * siri, tmpvec2, tmpvec);
3066 rvec_inc(sumvec, tmpvec);
3068 svmul(erg->rotg->k * erg->invmass, sumvec, innersumvec);
3072 /* Calculate the radial motion 2 potential and forces */
3073 static void do_radial_motion2(gmx_enfrotgrp* erg,
3074 rvec x[], /* The positions */
3075 const matrix box, /* The simulation box */
3076 gmx_bool bOutstepRot, /* Output to main rotation output file */
3077 gmx_bool bOutstepSlab) /* Output per-slab data */
3079 rvec xj; /* Position */
3080 real alpha; /* a single angle between an actual and a reference position */
3081 real weight; /* single weight for a single angle */
3082 rvec xj_u; /* xj - u */
3083 rvec yj0_yc0; /* yj0 -yc0 */
3084 rvec tmpvec, tmpvec2;
3085 real fac, fit_fac, fac2, Vpart = 0.0;
3086 rvec rj, fit_rj, sj;
3088 rvec v_xj_u; /* v x (xj - u) */
3089 real psij, psijstar;
3090 real mj, wj; /* For mass-weighting of the positions */
3094 gmx_bool bCalcPotFit;
3096 bPF = erg->rotg->eType == erotgRM2PF;
3097 bCalcPotFit = (bOutstepRot || bOutstepSlab) && (erotgFitPOT == erg->rotg->eFittype);
3099 clear_rvec(yj0_yc0); /* Make the compiler happy */
3101 clear_rvec(innersumvec);
3104 /* For the pivot-free variant we have to use the current center of
3105 * mass of the rotation group instead of the pivot u */
3106 get_center(erg->xc, erg->mc, erg->rotg->nat, erg->xc_center);
3108 /* Also, we precalculate the second term of the forces that is identical
3109 * (up to the weight factor mj) for all forces */
3110 radial_motion2_precalc_inner_sum(erg, innersumvec);
3113 N_M = erg->rotg->nat * erg->invmass;
3115 /* Each process calculates the forces on its local atoms */
3116 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
3117 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
3118 for (gmx::index j = 0; j < localRotationGroupIndex.ssize(); j++)
3122 /* Local index of a rotation group atom */
3123 int ii = localRotationGroupIndex[j];
3124 /* Position of this atom in the collective array */
3125 int iigrp = collectiveRotationGroupIndex[j];
3126 /* Mass-weighting */
3127 mj = erg->mc[iigrp];
3129 /* Current position of this atom: x[ii] */
3130 copy_rvec(x[ii], xj);
3132 /* Shift this atom such that it is near its reference */
3133 shift_single_coord(box, xj, erg->xc_shifts[iigrp]);
3135 /* The (unrotated) reference position is yj0. yc0 has already
3136 * been subtracted in init_rot_group */
3137 copy_rvec(erg->rotg->x_ref[iigrp], yj0_yc0); /* yj0_yc0 = yj0 - yc0 */
3139 /* Calculate Omega.(yj0-yc0) */
3140 mvmul(erg->rotmat, yj0_yc0, rj); /* rj = Omega.(yj0-yc0) */
3145 copy_rvec(erg->x_loc_pbc[j], xj);
3146 copy_rvec(erg->xr_loc[j], rj); /* rj = Omega.(yj0-u) */
3148 /* Mass-weighting */
3151 /* Calculate (xj-u) resp. (xj-xc) */
3152 rvec_sub(xj, erg->xc_center, xj_u); /* xj_u = xj-u */
3154 cprod(erg->vec, xj_u, v_xj_u); /* v_xj_u = v x (xj-u) */
3156 fac = norm2(v_xj_u);
3158 psijstar = 1.0 / (fac + erg->rotg->eps); /* psistar = -------------------- */
3159 /* * |v x (xj-u)|^2 + eps */
3161 psij = gmx::invsqrt(fac); /* 1 */
3162 /* psij = ------------ */
3165 svmul(psij, v_xj_u, sj); /* sj = psij * (v x (xj-u) ) */
3167 fac = iprod(v_xj_u, rj); /* fac = (v x (xj-u)).rj */
3170 sjrj = iprod(sj, rj); /* sjrj = sj.rj */
3172 svmul(psijstar / psij, rj, tmpvec);
3173 svmul(psijstar * psijstar / (psij * psij * psij) * sjrj, sj, tmpvec2);
3174 rvec_dec(tmpvec, tmpvec2);
3175 cprod(tmpvec, erg->vec, tmpvec2);
3177 /* Store the additional force so that it can be added to the force
3178 * array after the normal forces have been evaluated */
3179 svmul(-erg->rotg->k * wj * sjrj, tmpvec2, tmpvec);
3180 svmul(mj, innersumvec, tmpvec2); /* This is != 0 only for the pivot-free variant */
3182 rvec_add(tmpvec2, tmpvec, erg->f_rot_loc[j]);
3183 Vpart += wj * psijstar * fac2;
3185 /* If requested, also calculate the potential for a set of angles
3186 * near the current reference angle */
3189 for (int ifit = 0; ifit < erg->rotg->PotAngle_nstep; ifit++)
3193 mvmul(erg->PotAngleFit->rotmat[ifit], yj0_yc0, fit_rj); /* fit_rj = Omega.(yj0-yc0) */
3197 /* Position of this atom in the collective array */
3198 int iigrp = collectiveRotationGroupIndex[j];
3199 /* Rotate with the alternative angle. Like rotate_local_reference(),
3200 * just for a single local atom */
3201 mvmul(erg->PotAngleFit->rotmat[ifit], erg->rotg->x_ref[iigrp],
3202 fit_rj); /* fit_rj = Omega*(yj0-u) */
3204 fit_fac = iprod(v_xj_u, fit_rj); /* fac = (v x (xj-u)).fit_rj */
3205 /* Add to the rotation potential for this angle: */
3206 erg->PotAngleFit->V[ifit] += 0.5 * erg->rotg->k * wj * psijstar * fit_fac * fit_fac;
3212 /* Add to the torque of this rotation group */
3213 erg->torque_v += torque(erg->vec, erg->f_rot_loc[j], xj, erg->xc_center);
3215 /* Calculate the angle between reference and actual rotation group atom. */
3216 angle(erg, xj_u, rj, &alpha, &weight); /* angle in rad, weighted */
3217 erg->angle_v += alpha * weight;
3218 erg->weight_v += weight;
3223 } /* end of loop over local rotation group atoms */
3224 erg->V = 0.5 * erg->rotg->k * Vpart;
3228 /* Determine the smallest and largest position vector (with respect to the
3229 * rotation vector) for the reference group */
3230 static void get_firstlast_atom_ref(const gmx_enfrotgrp* erg, int* firstindex, int* lastindex)
3233 real xcproj; /* The projection of a reference position on the
3235 real minproj, maxproj; /* Smallest and largest projection on v */
3237 /* Start with some value */
3238 minproj = iprod(erg->rotg->x_ref[0], erg->vec);
3241 /* This is just to ensure that it still works if all the atoms of the
3242 * reference structure are situated in a plane perpendicular to the rotation
3245 *lastindex = erg->rotg->nat - 1;
3247 /* Loop over all atoms of the reference group,
3248 * project them on the rotation vector to find the extremes */
3249 for (i = 0; i < erg->rotg->nat; i++)
3251 xcproj = iprod(erg->rotg->x_ref[i], erg->vec);
3252 if (xcproj < minproj)
3257 if (xcproj > maxproj)
3266 /* Allocate memory for the slabs */
3267 static void allocate_slabs(gmx_enfrotgrp* erg, FILE* fplog, gmx_bool bVerbose)
3269 /* More slabs than are defined for the reference are never needed */
3270 int nslabs = erg->slab_last_ref - erg->slab_first_ref + 1;
3272 /* Remember how many we allocated */
3273 erg->nslabs_alloc = nslabs;
3275 if ((nullptr != fplog) && bVerbose)
3277 fprintf(fplog, "%s allocating memory to store data for %d slabs (rotation group %d).\n",
3278 RotStr, nslabs, erg->groupIndex);
3280 snew(erg->slab_center, nslabs);
3281 snew(erg->slab_center_ref, nslabs);
3282 snew(erg->slab_weights, nslabs);
3283 snew(erg->slab_torque_v, nslabs);
3284 snew(erg->slab_data, nslabs);
3285 snew(erg->gn_atom, nslabs);
3286 snew(erg->gn_slabind, nslabs);
3287 snew(erg->slab_innersumvec, nslabs);
3288 for (int i = 0; i < nslabs; i++)
3290 snew(erg->slab_data[i].x, erg->rotg->nat);
3291 snew(erg->slab_data[i].ref, erg->rotg->nat);
3292 snew(erg->slab_data[i].weight, erg->rotg->nat);
3294 snew(erg->xc_ref_sorted, erg->rotg->nat);
3295 snew(erg->xc_sortind, erg->rotg->nat);
3296 snew(erg->firstatom, nslabs);
3297 snew(erg->lastatom, nslabs);
3301 /* From the extreme positions of the reference group, determine the first
3302 * and last slab of the reference. We can never have more slabs in the real
3303 * simulation than calculated here for the reference.
3305 static void get_firstlast_slab_ref(gmx_enfrotgrp* erg, real mc[], int ref_firstindex, int ref_lastindex)
3309 int first = get_first_slab(erg, erg->rotg->x_ref[ref_firstindex]);
3310 int last = get_last_slab(erg, erg->rotg->x_ref[ref_lastindex]);
3312 while (get_slab_weight(first, erg, erg->rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3316 erg->slab_first_ref = first + 1;
3317 while (get_slab_weight(last, erg, erg->rotg->x_ref, mc, &dummy) > WEIGHT_MIN)
3321 erg->slab_last_ref = last - 1;
3325 /* Special version of copy_rvec:
3326 * During the copy procedure of xcurr to b, the correct PBC image is chosen
3327 * such that the copied vector ends up near its reference position xref */
3328 static inline void copy_correct_pbc_image(const rvec xcurr, /* copy vector xcurr ... */
3329 rvec b, /* ... to b ... */
3330 const rvec xref, /* choosing the PBC image such that b ends up near xref */
3339 /* Shortest PBC distance between the atom and its reference */
3340 rvec_sub(xcurr, xref, dx);
3342 /* Determine the shift for this atom */
3344 for (m = npbcdim - 1; m >= 0; m--)
3346 while (dx[m] < -0.5 * box[m][m])
3348 for (d = 0; d < DIM; d++)
3354 while (dx[m] >= 0.5 * box[m][m])
3356 for (d = 0; d < DIM; d++)
3364 /* Apply the shift to the position */
3365 copy_rvec(xcurr, b);
3366 shift_single_coord(box, b, shift);
3370 static void init_rot_group(FILE* fplog,
3371 const t_commrec* cr,
3379 gmx_bool bOutputCenters)
3381 rvec coord, xref, *xdum;
3382 gmx_bool bFlex, bColl;
3383 int ref_firstindex, ref_lastindex;
3384 real mass, totalmass;
3387 const t_rotgrp* rotg = erg->rotg;
3390 /* Do we have a flexible axis? */
3391 bFlex = ISFLEX(rotg);
3392 /* Do we use a global set of coordinates? */
3393 bColl = ISCOLL(rotg);
3395 /* Allocate space for collective coordinates if needed */
3398 snew(erg->xc, erg->rotg->nat);
3399 snew(erg->xc_shifts, erg->rotg->nat);
3400 snew(erg->xc_eshifts, erg->rotg->nat);
3401 snew(erg->xc_old, erg->rotg->nat);
3403 if (erg->rotg->eFittype == erotgFitNORM)
3405 snew(erg->xc_ref_length, erg->rotg->nat); /* in case fit type NORM is chosen */
3406 snew(erg->xc_norm, erg->rotg->nat);
3411 snew(erg->xr_loc, erg->rotg->nat);
3412 snew(erg->x_loc_pbc, erg->rotg->nat);
3415 copy_rvec(erg->rotg->inputVec, erg->vec);
3416 snew(erg->f_rot_loc, erg->rotg->nat);
3418 /* Make space for the calculation of the potential at other angles (used
3419 * for fitting only) */
3420 if (erotgFitPOT == erg->rotg->eFittype)
3422 snew(erg->PotAngleFit, 1);
3423 snew(erg->PotAngleFit->degangle, erg->rotg->PotAngle_nstep);
3424 snew(erg->PotAngleFit->V, erg->rotg->PotAngle_nstep);
3425 snew(erg->PotAngleFit->rotmat, erg->rotg->PotAngle_nstep);
3427 /* Get the set of angles around the reference angle */
3428 start = -0.5 * (erg->rotg->PotAngle_nstep - 1) * erg->rotg->PotAngle_step;
3429 for (int i = 0; i < erg->rotg->PotAngle_nstep; i++)
3431 erg->PotAngleFit->degangle[i] = start + i * erg->rotg->PotAngle_step;
3436 erg->PotAngleFit = nullptr;
3439 /* Copy the masses so that the center can be determined. For all types of
3440 * enforced rotation, we store the masses in the erg->mc array. */
3441 snew(erg->mc, erg->rotg->nat);
3444 snew(erg->mc_sorted, erg->rotg->nat);
3448 snew(erg->m_loc, erg->rotg->nat);
3452 for (int i = 0; i < erg->rotg->nat; i++)
3454 if (erg->rotg->bMassW)
3456 mass = mtopGetAtomMass(mtop, erg->rotg->ind[i], &molb);
3465 erg->invmass = 1.0 / totalmass;
3467 /* Set xc_ref_center for any rotation potential */
3468 if ((erg->rotg->eType == erotgISO) || (erg->rotg->eType == erotgPM)
3469 || (erg->rotg->eType == erotgRM) || (erg->rotg->eType == erotgRM2))
3471 /* Set the pivot point for the fixed, stationary-axis potentials. This
3472 * won't change during the simulation */
3473 copy_rvec(erg->rotg->pivot, erg->xc_ref_center);
3474 copy_rvec(erg->rotg->pivot, erg->xc_center);
3478 /* Center of the reference positions */
3479 get_center(erg->rotg->x_ref, erg->mc, erg->rotg->nat, erg->xc_ref_center);
3481 /* Center of the actual positions */
3484 snew(xdum, erg->rotg->nat);
3485 for (int i = 0; i < erg->rotg->nat; i++)
3487 int ii = erg->rotg->ind[i];
3488 copy_rvec(x[ii], xdum[i]);
3490 get_center(xdum, erg->mc, erg->rotg->nat, erg->xc_center);
3496 gmx_bcast(sizeof(erg->xc_center), erg->xc_center, cr);
3503 /* Save the original (whole) set of positions in xc_old such that at later
3504 * steps the rotation group can always be made whole again. If the simulation is
3505 * restarted, we compute the starting reference positions (given the time)
3506 * and assume that the correct PBC image of each position is the one nearest
3507 * to the current reference */
3510 /* Calculate the rotation matrix for this angle: */
3511 t_start = ir->init_t + ir->init_step * ir->delta_t;
3512 erg->degangle = erg->rotg->rate * t_start;
3513 calc_rotmat(erg->vec, erg->degangle, erg->rotmat);
3515 for (int i = 0; i < erg->rotg->nat; i++)
3517 int ii = erg->rotg->ind[i];
3519 /* Subtract pivot, rotate, and add pivot again. This will yield the
3520 * reference position for time t */
3521 rvec_sub(erg->rotg->x_ref[i], erg->xc_ref_center, coord);
3522 mvmul(erg->rotmat, coord, xref);
3523 rvec_inc(xref, erg->xc_ref_center);
3525 copy_correct_pbc_image(x[ii], erg->xc_old[i], xref, box, 3);
3531 gmx_bcast(erg->rotg->nat * sizeof(erg->xc_old[0]), erg->xc_old, cr);
3536 if ((erg->rotg->eType != erotgFLEX) && (erg->rotg->eType != erotgFLEX2))
3538 /* Put the reference positions into origin: */
3539 for (int i = 0; i < erg->rotg->nat; i++)
3541 rvec_dec(erg->rotg->x_ref[i], erg->xc_ref_center);
3545 /* Enforced rotation with flexible axis */
3548 /* Calculate maximum beta value from minimum gaussian (performance opt.) */
3549 erg->max_beta = calc_beta_max(erg->rotg->min_gaussian, erg->rotg->slab_dist);
3551 /* Determine the smallest and largest coordinate with respect to the rotation vector */
3552 get_firstlast_atom_ref(erg, &ref_firstindex, &ref_lastindex);
3554 /* From the extreme positions of the reference group, determine the first
3555 * and last slab of the reference. */
3556 get_firstlast_slab_ref(erg, erg->mc, ref_firstindex, ref_lastindex);
3558 /* Allocate memory for the slabs */
3559 allocate_slabs(erg, fplog, bVerbose);
3561 /* Flexible rotation: determine the reference centers for the rest of the simulation */
3562 erg->slab_first = erg->slab_first_ref;
3563 erg->slab_last = erg->slab_last_ref;
3564 get_slab_centers(erg, erg->rotg->x_ref, erg->mc, -1, out_slabs, bOutputCenters, TRUE);
3566 /* Length of each x_rotref vector from center (needed if fit routine NORM is chosen): */
3567 if (erg->rotg->eFittype == erotgFitNORM)
3569 for (int i = 0; i < erg->rotg->nat; i++)
3571 rvec_sub(erg->rotg->x_ref[i], erg->xc_ref_center, coord);
3572 erg->xc_ref_length[i] = norm(coord);
3578 /* Calculate the size of the MPI buffer needed in reduce_output() */
3579 static int calc_mpi_bufsize(const gmx_enfrot* er)
3582 int count_total = 0;
3583 for (int g = 0; g < er->rot->ngrp; g++)
3585 const t_rotgrp* rotg = &er->rot->grp[g];
3586 const gmx_enfrotgrp* erg = &er->enfrotgrp[g];
3588 /* Count the items that are transferred for this group: */
3589 int count_group = 4; /* V, torque, angle, weight */
3591 /* Add the maximum number of slabs for flexible groups */
3594 count_group += erg->slab_last_ref - erg->slab_first_ref + 1;
3597 /* Add space for the potentials at different angles: */
3598 if (erotgFitPOT == erg->rotg->eFittype)
3600 count_group += erg->rotg->PotAngle_nstep;
3603 /* Add to the total number: */
3604 count_total += count_group;
3611 std::unique_ptr<gmx::EnforcedRotation> init_rot(FILE* fplog,
3614 const t_filenm fnm[],
3615 const t_commrec* cr,
3616 gmx::LocalAtomSetManager* atomSets,
3617 const t_state* globalState,
3619 const gmx_output_env_t* oenv,
3620 const gmx::MdrunOptions& mdrunOptions,
3621 const gmx::StartingBehavior startingBehavior)
3623 int nat_max = 0; /* Size of biggest rotation group */
3624 rvec* x_pbc = nullptr; /* Space for the pbc-correct atom positions */
3626 if (MASTER(cr) && mdrunOptions.verbose)
3628 fprintf(stdout, "%s Initializing ...\n", RotStr);
3631 auto enforcedRotation = std::make_unique<gmx::EnforcedRotation>();
3632 gmx_enfrot* er = enforcedRotation->getLegacyEnfrot();
3633 // TODO When this module implements IMdpOptions, the ownership will become more clear.
3635 er->restartWithAppending = (startingBehavior == gmx::StartingBehavior::RestartWithAppending);
3637 /* When appending, skip first output to avoid duplicate entries in the data files */
3638 if (er->restartWithAppending)
3647 if (MASTER(cr) && er->bOut)
3649 please_cite(fplog, "Kutzner2011");
3652 /* Output every step for reruns */
3653 if (mdrunOptions.rerun)
3655 if (nullptr != fplog)
3657 fprintf(fplog, "%s rerun - will write rotation output every available step.\n", RotStr);
3664 er->nstrout = er->rot->nstrout;
3665 er->nstsout = er->rot->nstsout;
3668 er->out_slabs = nullptr;
3669 if (MASTER(cr) && HaveFlexibleGroups(er->rot))
3671 er->out_slabs = open_slab_out(opt2fn("-rs", nfile, fnm), er);
3676 /* Remove pbc, make molecule whole.
3677 * When ir->bContinuation=TRUE this has already been done, but ok. */
3678 snew(x_pbc, mtop->natoms);
3679 copy_rvecn(globalState->x.rvec_array(), x_pbc, 0, mtop->natoms);
3680 do_pbc_first_mtop(nullptr, ir->ePBC, globalState->box, mtop, x_pbc);
3681 /* All molecules will be whole now, but not necessarily in the home box.
3682 * Additionally, if a rotation group consists of more than one molecule
3683 * (e.g. two strands of DNA), each one of them can end up in a different
3684 * periodic box. This is taken care of in init_rot_group. */
3687 /* Allocate space for the per-rotation-group data: */
3688 er->enfrotgrp.resize(er->rot->ngrp);
3690 for (auto& ergRef : er->enfrotgrp)
3692 gmx_enfrotgrp* erg = &ergRef;
3693 erg->rotg = &er->rot->grp[groupIndex];
3694 erg->atomSet = std::make_unique<gmx::LocalAtomSet>(
3695 atomSets->add({ erg->rotg->ind, erg->rotg->ind + erg->rotg->nat }));
3696 erg->groupIndex = groupIndex;
3698 if (nullptr != fplog)
3700 fprintf(fplog, "%s group %d type '%s'\n", RotStr, groupIndex, erotg_names[erg->rotg->eType]);
3703 if (erg->rotg->nat > 0)
3705 nat_max = std::max(nat_max, erg->rotg->nat);
3707 init_rot_group(fplog, cr, erg, x_pbc, mtop, mdrunOptions.verbose, er->out_slabs,
3708 MASTER(cr) ? globalState->box : nullptr, ir,
3709 !er->restartWithAppending); /* Do not output the reference centers
3710 * again if we are appending */
3715 /* Allocate space for enforced rotation buffer variables */
3716 er->bufsize = nat_max;
3717 snew(er->data, nat_max);
3718 snew(er->xbuf, nat_max);
3719 snew(er->mbuf, nat_max);
3721 /* Buffers for MPI reducing torques, angles, weights (for each group), and V */
3724 er->mpi_bufsize = calc_mpi_bufsize(er) + 100; /* larger to catch errors */
3725 snew(er->mpi_inbuf, er->mpi_bufsize);
3726 snew(er->mpi_outbuf, er->mpi_bufsize);
3730 er->mpi_bufsize = 0;
3731 er->mpi_inbuf = nullptr;
3732 er->mpi_outbuf = nullptr;
3735 /* Only do I/O on the MASTER */
3736 er->out_angles = nullptr;
3737 er->out_rot = nullptr;
3738 er->out_torque = nullptr;
3741 er->out_rot = open_rot_out(opt2fn("-ro", nfile, fnm), oenv, er);
3743 if (er->nstsout > 0)
3745 if (HaveFlexibleGroups(er->rot) || HavePotFitGroups(er->rot))
3747 er->out_angles = open_angles_out(opt2fn("-ra", nfile, fnm), er);
3749 if (HaveFlexibleGroups(er->rot))
3751 er->out_torque = open_torque_out(opt2fn("-rt", nfile, fnm), er);
3757 return enforcedRotation;
3760 /* Rotate the local reference positions and store them in
3761 * erg->xr_loc[0...(nat_loc-1)]
3763 * Note that we already subtracted u or y_c from the reference positions
3764 * in init_rot_group().
3766 static void rotate_local_reference(gmx_enfrotgrp* erg)
3768 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
3769 for (size_t i = 0; i < erg->atomSet->numAtomsLocal(); i++)
3771 /* Index of this rotation group atom with respect to the whole rotation group */
3772 int ii = collectiveRotationGroupIndex[i];
3774 mvmul(erg->rotmat, erg->rotg->x_ref[ii], erg->xr_loc[i]);
3779 /* Select the PBC representation for each local x position and store that
3780 * for later usage. We assume the right PBC image of an x is the one nearest to
3781 * its rotated reference */
3782 static void choose_pbc_image(rvec x[], gmx_enfrotgrp* erg, const matrix box, int npbcdim)
3784 const auto& localRotationGroupIndex = erg->atomSet->localIndex();
3785 for (gmx::index i = 0; i < localRotationGroupIndex.ssize(); i++)
3787 /* Index of a rotation group atom */
3788 int ii = localRotationGroupIndex[i];
3790 /* Get the correctly rotated reference position. The pivot was already
3791 * subtracted in init_rot_group() from the reference positions. Also,
3792 * the reference positions have already been rotated in
3793 * rotate_local_reference(). For the current reference position we thus
3794 * only need to add the pivot again. */
3796 copy_rvec(erg->xr_loc[i], xref);
3797 rvec_inc(xref, erg->xc_ref_center);
3799 copy_correct_pbc_image(x[ii], erg->x_loc_pbc[i], xref, box, npbcdim);
3804 void do_rotation(const t_commrec* cr, gmx_enfrot* er, const matrix box, rvec x[], real t, int64_t step, gmx_bool bNS)
3806 gmx_bool outstep_slab, outstep_rot;
3809 gmx_potfit* fit = nullptr; /* For fit type 'potential' determine the fit
3810 angle via the potential minimum */
3816 /* When to output in main rotation output file */
3817 outstep_rot = do_per_step(step, er->nstrout) && er->bOut;
3818 /* When to output per-slab data */
3819 outstep_slab = do_per_step(step, er->nstsout) && er->bOut;
3821 /* Output time into rotation output file */
3822 if (outstep_rot && MASTER(cr))
3824 fprintf(er->out_rot, "%12.3e", t);
3827 /**************************************************************************/
3828 /* First do ALL the communication! */
3829 for (auto& ergRef : er->enfrotgrp)
3831 gmx_enfrotgrp* erg = &ergRef;
3832 const t_rotgrp* rotg = erg->rotg;
3834 /* Do we use a collective (global) set of coordinates? */
3835 bColl = ISCOLL(rotg);
3837 /* Calculate the rotation matrix for this angle: */
3838 erg->degangle = rotg->rate * t;
3839 calc_rotmat(erg->vec, erg->degangle, erg->rotmat);
3843 /* Transfer the rotation group's positions such that every node has
3844 * all of them. Every node contributes its local positions x and stores
3845 * it in the collective erg->xc array. */
3846 communicate_group_positions(cr, erg->xc, erg->xc_shifts, erg->xc_eshifts, bNS, x,
3847 rotg->nat, erg->atomSet->numAtomsLocal(),
3848 erg->atomSet->localIndex().data(),
3849 erg->atomSet->collectiveIndex().data(), erg->xc_old, box);
3853 /* Fill the local masses array;
3854 * this array changes in DD/neighborsearching steps */
3857 const auto& collectiveRotationGroupIndex = erg->atomSet->collectiveIndex();
3858 for (gmx::index i = 0; i < collectiveRotationGroupIndex.ssize(); i++)
3860 /* Index of local atom w.r.t. the collective rotation group */
3861 int ii = collectiveRotationGroupIndex[i];
3862 erg->m_loc[i] = erg->mc[ii];
3866 /* Calculate Omega*(y_i-y_c) for the local positions */
3867 rotate_local_reference(erg);
3869 /* Choose the nearest PBC images of the group atoms with respect
3870 * to the rotated reference positions */
3871 choose_pbc_image(x, erg, box, 3);
3873 /* Get the center of the rotation group */
3874 if ((rotg->eType == erotgISOPF) || (rotg->eType == erotgPMPF))
3876 get_center_comm(cr, erg->x_loc_pbc, erg->m_loc, erg->atomSet->numAtomsLocal(),
3877 rotg->nat, erg->xc_center);
3881 } /* End of loop over rotation groups */
3883 /**************************************************************************/
3884 /* Done communicating, we can start to count cycles for the load balancing now ... */
3885 if (DOMAINDECOMP(cr))
3887 ddReopenBalanceRegionCpu(cr->dd);
3894 for (auto& ergRef : er->enfrotgrp)
3896 gmx_enfrotgrp* erg = &ergRef;
3897 const t_rotgrp* rotg = erg->rotg;
3899 if (outstep_rot && MASTER(cr))
3901 fprintf(er->out_rot, "%12.4f", erg->degangle);
3904 /* Calculate angles and rotation matrices for potential fitting: */
3905 if ((outstep_rot || outstep_slab) && (erotgFitPOT == rotg->eFittype))
3907 fit = erg->PotAngleFit;
3908 for (int i = 0; i < rotg->PotAngle_nstep; i++)
3910 calc_rotmat(erg->vec, erg->degangle + fit->degangle[i], fit->rotmat[i]);
3912 /* Clear value from last step */
3913 erg->PotAngleFit->V[i] = 0.0;
3917 /* Clear values from last time step */
3919 erg->torque_v = 0.0;
3921 erg->weight_v = 0.0;
3923 switch (rotg->eType)
3928 case erotgPMPF: do_fixed(erg, outstep_rot, outstep_slab); break;
3929 case erotgRM: do_radial_motion(erg, outstep_rot, outstep_slab); break;
3930 case erotgRMPF: do_radial_motion_pf(erg, x, box, outstep_rot, outstep_slab); break;
3932 case erotgRM2PF: do_radial_motion2(erg, x, box, outstep_rot, outstep_slab); break;
3935 /* Subtract the center of the rotation group from the collective positions array
3936 * Also store the center in erg->xc_center since it needs to be subtracted
3937 * in the low level routines from the local coordinates as well */
3938 get_center(erg->xc, erg->mc, rotg->nat, erg->xc_center);
3939 svmul(-1.0, erg->xc_center, transvec);
3940 translate_x(erg->xc, rotg->nat, transvec);
3941 do_flexible(MASTER(cr), er, erg, x, box, t, outstep_rot, outstep_slab);
3945 /* Do NOT subtract the center of mass in the low level routines! */
3946 clear_rvec(erg->xc_center);
3947 do_flexible(MASTER(cr), er, erg, x, box, t, outstep_rot, outstep_slab);
3949 default: gmx_fatal(FARGS, "No such rotation potential.");
3956 fprintf(stderr, "%s calculation (step %d) took %g seconds.\n", RotStr, step, MPI_Wtime() - t0);